Actual source code: baijsolvtrannat6.c

  1: #include <../src/mat/impls/baij/seq/baij.h>

  3: PetscErrorCode MatSolveTranspose_SeqBAIJ_6_NaturalOrdering_inplace(Mat A, Vec bb, Vec xx)
  4: {
  5:   Mat_SeqBAIJ     *a    = (Mat_SeqBAIJ *)A->data;
  6:   const PetscInt  *diag = a->diag, n = a->mbs, *vi, *ai = a->i, *aj = a->j;
  7:   PetscInt         i, nz, idx, idt, oidx;
  8:   const MatScalar *aa = a->a, *v;
  9:   PetscScalar      s1, s2, s3, s4, s5, s6, x1, x2, x3, x4, x5, x6, *x;

 11:   PetscFunctionBegin;
 12:   PetscCall(VecCopy(bb, xx));
 13:   PetscCall(VecGetArray(xx, &x));

 15:   /* forward solve the U^T */
 16:   idx = 0;
 17:   for (i = 0; i < n; i++) {
 18:     v = aa + 36 * diag[i];
 19:     /* multiply by the inverse of the block diagonal */
 20:     x1 = x[idx];
 21:     x2 = x[1 + idx];
 22:     x3 = x[2 + idx];
 23:     x4 = x[3 + idx];
 24:     x5 = x[4 + idx];
 25:     x6 = x[5 + idx];
 26:     s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4 + v[4] * x5 + v[5] * x6;
 27:     s2 = v[6] * x1 + v[7] * x2 + v[8] * x3 + v[9] * x4 + v[10] * x5 + v[11] * x6;
 28:     s3 = v[12] * x1 + v[13] * x2 + v[14] * x3 + v[15] * x4 + v[16] * x5 + v[17] * x6;
 29:     s4 = v[18] * x1 + v[19] * x2 + v[20] * x3 + v[21] * x4 + v[22] * x5 + v[23] * x6;
 30:     s5 = v[24] * x1 + v[25] * x2 + v[26] * x3 + v[27] * x4 + v[28] * x5 + v[29] * x6;
 31:     s6 = v[30] * x1 + v[31] * x2 + v[32] * x3 + v[33] * x4 + v[34] * x5 + v[35] * x6;
 32:     v += 36;

 34:     vi = aj + diag[i] + 1;
 35:     nz = ai[i + 1] - diag[i] - 1;
 36:     while (nz--) {
 37:       oidx = 6 * (*vi++);
 38:       x[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6;
 39:       x[oidx + 1] -= v[6] * s1 + v[7] * s2 + v[8] * s3 + v[9] * s4 + v[10] * s5 + v[11] * s6;
 40:       x[oidx + 2] -= v[12] * s1 + v[13] * s2 + v[14] * s3 + v[15] * s4 + v[16] * s5 + v[17] * s6;
 41:       x[oidx + 3] -= v[18] * s1 + v[19] * s2 + v[20] * s3 + v[21] * s4 + v[22] * s5 + v[23] * s6;
 42:       x[oidx + 4] -= v[24] * s1 + v[25] * s2 + v[26] * s3 + v[27] * s4 + v[28] * s5 + v[29] * s6;
 43:       x[oidx + 5] -= v[30] * s1 + v[31] * s2 + v[32] * s3 + v[33] * s4 + v[34] * s5 + v[35] * s6;
 44:       v += 36;
 45:     }
 46:     x[idx]     = s1;
 47:     x[1 + idx] = s2;
 48:     x[2 + idx] = s3;
 49:     x[3 + idx] = s4;
 50:     x[4 + idx] = s5;
 51:     x[5 + idx] = s6;
 52:     idx += 6;
 53:   }
 54:   /* backward solve the L^T */
 55:   for (i = n - 1; i >= 0; i--) {
 56:     v   = aa + 36 * diag[i] - 36;
 57:     vi  = aj + diag[i] - 1;
 58:     nz  = diag[i] - ai[i];
 59:     idt = 6 * i;
 60:     s1  = x[idt];
 61:     s2  = x[1 + idt];
 62:     s3  = x[2 + idt];
 63:     s4  = x[3 + idt];
 64:     s5  = x[4 + idt];
 65:     s6  = x[5 + idt];
 66:     while (nz--) {
 67:       idx = 6 * (*vi--);
 68:       x[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6;
 69:       x[idx + 1] -= v[6] * s1 + v[7] * s2 + v[8] * s3 + v[9] * s4 + v[10] * s5 + v[11] * s6;
 70:       x[idx + 2] -= v[12] * s1 + v[13] * s2 + v[14] * s3 + v[15] * s4 + v[16] * s5 + v[17] * s6;
 71:       x[idx + 3] -= v[18] * s1 + v[19] * s2 + v[20] * s3 + v[21] * s4 + v[22] * s5 + v[23] * s6;
 72:       x[idx + 4] -= v[24] * s1 + v[25] * s2 + v[26] * s3 + v[27] * s4 + v[28] * s5 + v[29] * s6;
 73:       x[idx + 5] -= v[30] * s1 + v[31] * s2 + v[32] * s3 + v[33] * s4 + v[34] * s5 + v[35] * s6;
 74:       v -= 36;
 75:     }
 76:   }
 77:   PetscCall(VecRestoreArray(xx, &x));
 78:   PetscCall(PetscLogFlops(2.0 * 36 * (a->nz) - 6.0 * A->cmap->n));
 79:   PetscFunctionReturn(PETSC_SUCCESS);
 80: }

 82: PetscErrorCode MatSolveTranspose_SeqBAIJ_6_NaturalOrdering(Mat A, Vec bb, Vec xx)
 83: {
 84:   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *)A->data;
 85:   const PetscInt   n = a->mbs, *vi, *ai = a->i, *aj = a->j, *diag = a->diag;
 86:   PetscInt         nz, idx, idt, j, i, oidx;
 87:   const PetscInt   bs = A->rmap->bs, bs2 = a->bs2;
 88:   const MatScalar *aa = a->a, *v;
 89:   PetscScalar      s1, s2, s3, s4, s5, s6, x1, x2, x3, x4, x5, x6, *x;

 91:   PetscFunctionBegin;
 92:   PetscCall(VecCopy(bb, xx));
 93:   PetscCall(VecGetArray(xx, &x));

 95:   /* forward solve the U^T */
 96:   idx = 0;
 97:   for (i = 0; i < n; i++) {
 98:     v = aa + bs2 * diag[i];
 99:     /* multiply by the inverse of the block diagonal */
100:     x1 = x[idx];
101:     x2 = x[1 + idx];
102:     x3 = x[2 + idx];
103:     x4 = x[3 + idx];
104:     x5 = x[4 + idx];
105:     x6 = x[5 + idx];
106:     s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4 + v[4] * x5 + v[5] * x6;
107:     s2 = v[6] * x1 + v[7] * x2 + v[8] * x3 + v[9] * x4 + v[10] * x5 + v[11] * x6;
108:     s3 = v[12] * x1 + v[13] * x2 + v[14] * x3 + v[15] * x4 + v[16] * x5 + v[17] * x6;
109:     s4 = v[18] * x1 + v[19] * x2 + v[20] * x3 + v[21] * x4 + v[22] * x5 + v[23] * x6;
110:     s5 = v[24] * x1 + v[25] * x2 + v[26] * x3 + v[27] * x4 + v[28] * x5 + v[29] * x6;
111:     s6 = v[30] * x1 + v[31] * x2 + v[32] * x3 + v[33] * x4 + v[34] * x5 + v[35] * x6;
112:     v -= bs2;

114:     vi = aj + diag[i] - 1;
115:     nz = diag[i] - diag[i + 1] - 1;
116:     for (j = 0; j > -nz; j--) {
117:       oidx = bs * vi[j];
118:       x[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6;
119:       x[oidx + 1] -= v[6] * s1 + v[7] * s2 + v[8] * s3 + v[9] * s4 + v[10] * s5 + v[11] * s6;
120:       x[oidx + 2] -= v[12] * s1 + v[13] * s2 + v[14] * s3 + v[15] * s4 + v[16] * s5 + v[17] * s6;
121:       x[oidx + 3] -= v[18] * s1 + v[19] * s2 + v[20] * s3 + v[21] * s4 + v[22] * s5 + v[23] * s6;
122:       x[oidx + 4] -= v[24] * s1 + v[25] * s2 + v[26] * s3 + v[27] * s4 + v[28] * s5 + v[29] * s6;
123:       x[oidx + 5] -= v[30] * s1 + v[31] * s2 + v[32] * s3 + v[33] * s4 + v[34] * s5 + v[35] * s6;
124:       v -= bs2;
125:     }
126:     x[idx]     = s1;
127:     x[1 + idx] = s2;
128:     x[2 + idx] = s3;
129:     x[3 + idx] = s4;
130:     x[4 + idx] = s5;
131:     x[5 + idx] = s6;
132:     idx += bs;
133:   }
134:   /* backward solve the L^T */
135:   for (i = n - 1; i >= 0; i--) {
136:     v   = aa + bs2 * ai[i];
137:     vi  = aj + ai[i];
138:     nz  = ai[i + 1] - ai[i];
139:     idt = bs * i;
140:     s1  = x[idt];
141:     s2  = x[1 + idt];
142:     s3  = x[2 + idt];
143:     s4  = x[3 + idt];
144:     s5  = x[4 + idt];
145:     s6  = x[5 + idt];
146:     for (j = 0; j < nz; j++) {
147:       idx = bs * vi[j];
148:       x[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6;
149:       x[idx + 1] -= v[6] * s1 + v[7] * s2 + v[8] * s3 + v[9] * s4 + v[10] * s5 + v[11] * s6;
150:       x[idx + 2] -= v[12] * s1 + v[13] * s2 + v[14] * s3 + v[15] * s4 + v[16] * s5 + v[17] * s6;
151:       x[idx + 3] -= v[18] * s1 + v[19] * s2 + v[20] * s3 + v[21] * s4 + v[22] * s5 + v[23] * s6;
152:       x[idx + 4] -= v[24] * s1 + v[25] * s2 + v[26] * s3 + v[27] * s4 + v[28] * s5 + v[29] * s6;
153:       x[idx + 5] -= v[30] * s1 + v[31] * s2 + v[32] * s3 + v[33] * s4 + v[34] * s5 + v[35] * s6;
154:       v += bs2;
155:     }
156:   }
157:   PetscCall(VecRestoreArray(xx, &x));
158:   PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n));
159:   PetscFunctionReturn(PETSC_SUCCESS);
160: }