Actual source code: sbaijfact6.c


  2: #include <../src/mat/impls/sbaij/seq/sbaij.h>
  3: #include <petsc/private/kernels/blockinvert.h>

  5: /* Version for when blocks are 4 by 4  */
  6: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat C, Mat A, const MatFactorInfo *info)
  7: {
  8:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ *)A->data, *b = (Mat_SeqSBAIJ *)C->data;
  9:   IS              perm = b->row;
 10:   const PetscInt *ai, *aj, *perm_ptr, mbs = a->mbs, *bi = b->i, *bj = b->j;
 11:   PetscInt        i, j, *a2anew, k, k1, jmin, jmax, *jl, *il, vj, nexti, ili;
 12:   MatScalar      *ba = b->a, *aa, *ap, *dk, *uik;
 13:   MatScalar      *u, *diag, *rtmp, *rtmp_ptr;
 14:   PetscBool       pivotinblocks = b->pivotinblocks;
 15:   PetscReal       shift         = info->shiftamount;
 16:   PetscBool       allowzeropivot, zeropivotdetected = PETSC_FALSE;

 18:   PetscFunctionBegin;
 19:   /* initialization */
 20:   allowzeropivot = PetscNot(A->erroriffailure);
 21:   PetscCall(PetscCalloc1(16 * mbs, &rtmp));
 22:   PetscCall(PetscMalloc2(mbs, &il, mbs, &jl));
 23:   il[0] = 0;
 24:   for (i = 0; i < mbs; i++) jl[i] = mbs;

 26:   PetscCall(PetscMalloc2(16, &dk, 16, &uik));
 27:   PetscCall(ISGetIndices(perm, &perm_ptr));

 29:   /* check permutation */
 30:   if (!a->permute) {
 31:     ai = a->i;
 32:     aj = a->j;
 33:     aa = a->a;
 34:   } else {
 35:     ai = a->inew;
 36:     aj = a->jnew;
 37:     PetscCall(PetscMalloc1(16 * ai[mbs], &aa));
 38:     PetscCall(PetscArraycpy(aa, a->a, 16 * ai[mbs]));
 39:     PetscCall(PetscMalloc1(ai[mbs], &a2anew));
 40:     PetscCall(PetscArraycpy(a2anew, a->a2anew, ai[mbs]));

 42:     for (i = 0; i < mbs; i++) {
 43:       jmin = ai[i];
 44:       jmax = ai[i + 1];
 45:       for (j = jmin; j < jmax; j++) {
 46:         while (a2anew[j] != j) {
 47:           k         = a2anew[j];
 48:           a2anew[j] = a2anew[k];
 49:           a2anew[k] = k;
 50:           for (k1 = 0; k1 < 16; k1++) {
 51:             dk[k1]          = aa[k * 16 + k1];
 52:             aa[k * 16 + k1] = aa[j * 16 + k1];
 53:             aa[j * 16 + k1] = dk[k1];
 54:           }
 55:         }
 56:         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
 57:         if (i > aj[j]) {
 58:           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
 59:           ap = aa + j * 16;                       /* ptr to the beginning of j-th block of aa */
 60:           for (k = 0; k < 16; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
 61:           for (k = 0; k < 4; k++) {               /* j-th block of aa <- dk^T */
 62:             for (k1 = 0; k1 < 4; k1++) *ap++ = dk[k + 4 * k1];
 63:           }
 64:         }
 65:       }
 66:     }
 67:     PetscCall(PetscFree(a2anew));
 68:   }

 70:   /* for each row k */
 71:   for (k = 0; k < mbs; k++) {
 72:     /*initialize k-th row with elements nonzero in row perm(k) of A */
 73:     jmin = ai[perm_ptr[k]];
 74:     jmax = ai[perm_ptr[k] + 1];
 75:     if (jmin < jmax) {
 76:       ap = aa + jmin * 16;
 77:       for (j = jmin; j < jmax; j++) {
 78:         vj       = perm_ptr[aj[j]]; /* block col. index */
 79:         rtmp_ptr = rtmp + vj * 16;
 80:         for (i = 0; i < 16; i++) *rtmp_ptr++ = *ap++;
 81:       }
 82:     }

 84:     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
 85:     PetscCall(PetscArraycpy(dk, rtmp + k * 16, 16));
 86:     i = jl[k]; /* first row to be added to k_th row  */

 88:     while (i < mbs) {
 89:       nexti = jl[i]; /* next row to be added to k_th row */

 91:       /* compute multiplier */
 92:       ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */

 94:       /* uik = -inv(Di)*U_bar(i,k) */
 95:       diag = ba + i * 16;
 96:       u    = ba + ili * 16;

 98:       uik[0] = -(diag[0] * u[0] + diag[4] * u[1] + diag[8] * u[2] + diag[12] * u[3]);
 99:       uik[1] = -(diag[1] * u[0] + diag[5] * u[1] + diag[9] * u[2] + diag[13] * u[3]);
100:       uik[2] = -(diag[2] * u[0] + diag[6] * u[1] + diag[10] * u[2] + diag[14] * u[3]);
101:       uik[3] = -(diag[3] * u[0] + diag[7] * u[1] + diag[11] * u[2] + diag[15] * u[3]);

103:       uik[4] = -(diag[0] * u[4] + diag[4] * u[5] + diag[8] * u[6] + diag[12] * u[7]);
104:       uik[5] = -(diag[1] * u[4] + diag[5] * u[5] + diag[9] * u[6] + diag[13] * u[7]);
105:       uik[6] = -(diag[2] * u[4] + diag[6] * u[5] + diag[10] * u[6] + diag[14] * u[7]);
106:       uik[7] = -(diag[3] * u[4] + diag[7] * u[5] + diag[11] * u[6] + diag[15] * u[7]);

108:       uik[8]  = -(diag[0] * u[8] + diag[4] * u[9] + diag[8] * u[10] + diag[12] * u[11]);
109:       uik[9]  = -(diag[1] * u[8] + diag[5] * u[9] + diag[9] * u[10] + diag[13] * u[11]);
110:       uik[10] = -(diag[2] * u[8] + diag[6] * u[9] + diag[10] * u[10] + diag[14] * u[11]);
111:       uik[11] = -(diag[3] * u[8] + diag[7] * u[9] + diag[11] * u[10] + diag[15] * u[11]);

113:       uik[12] = -(diag[0] * u[12] + diag[4] * u[13] + diag[8] * u[14] + diag[12] * u[15]);
114:       uik[13] = -(diag[1] * u[12] + diag[5] * u[13] + diag[9] * u[14] + diag[13] * u[15]);
115:       uik[14] = -(diag[2] * u[12] + diag[6] * u[13] + diag[10] * u[14] + diag[14] * u[15]);
116:       uik[15] = -(diag[3] * u[12] + diag[7] * u[13] + diag[11] * u[14] + diag[15] * u[15]);

118:       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
119:       dk[0] += uik[0] * u[0] + uik[1] * u[1] + uik[2] * u[2] + uik[3] * u[3];
120:       dk[1] += uik[4] * u[0] + uik[5] * u[1] + uik[6] * u[2] + uik[7] * u[3];
121:       dk[2] += uik[8] * u[0] + uik[9] * u[1] + uik[10] * u[2] + uik[11] * u[3];
122:       dk[3] += uik[12] * u[0] + uik[13] * u[1] + uik[14] * u[2] + uik[15] * u[3];

124:       dk[4] += uik[0] * u[4] + uik[1] * u[5] + uik[2] * u[6] + uik[3] * u[7];
125:       dk[5] += uik[4] * u[4] + uik[5] * u[5] + uik[6] * u[6] + uik[7] * u[7];
126:       dk[6] += uik[8] * u[4] + uik[9] * u[5] + uik[10] * u[6] + uik[11] * u[7];
127:       dk[7] += uik[12] * u[4] + uik[13] * u[5] + uik[14] * u[6] + uik[15] * u[7];

129:       dk[8] += uik[0] * u[8] + uik[1] * u[9] + uik[2] * u[10] + uik[3] * u[11];
130:       dk[9] += uik[4] * u[8] + uik[5] * u[9] + uik[6] * u[10] + uik[7] * u[11];
131:       dk[10] += uik[8] * u[8] + uik[9] * u[9] + uik[10] * u[10] + uik[11] * u[11];
132:       dk[11] += uik[12] * u[8] + uik[13] * u[9] + uik[14] * u[10] + uik[15] * u[11];

134:       dk[12] += uik[0] * u[12] + uik[1] * u[13] + uik[2] * u[14] + uik[3] * u[15];
135:       dk[13] += uik[4] * u[12] + uik[5] * u[13] + uik[6] * u[14] + uik[7] * u[15];
136:       dk[14] += uik[8] * u[12] + uik[9] * u[13] + uik[10] * u[14] + uik[11] * u[15];
137:       dk[15] += uik[12] * u[12] + uik[13] * u[13] + uik[14] * u[14] + uik[15] * u[15];

139:       PetscCall(PetscLogFlops(64.0 * 4.0));

141:       /* update -U(i,k) */
142:       PetscCall(PetscArraycpy(ba + ili * 16, uik, 16));

144:       /* add multiple of row i to k-th row ... */
145:       jmin = ili + 1;
146:       jmax = bi[i + 1];
147:       if (jmin < jmax) {
148:         for (j = jmin; j < jmax; j++) {
149:           /* rtmp += -U(i,k)^T * U_bar(i,j) */
150:           rtmp_ptr = rtmp + bj[j] * 16;
151:           u        = ba + j * 16;
152:           rtmp_ptr[0] += uik[0] * u[0] + uik[1] * u[1] + uik[2] * u[2] + uik[3] * u[3];
153:           rtmp_ptr[1] += uik[4] * u[0] + uik[5] * u[1] + uik[6] * u[2] + uik[7] * u[3];
154:           rtmp_ptr[2] += uik[8] * u[0] + uik[9] * u[1] + uik[10] * u[2] + uik[11] * u[3];
155:           rtmp_ptr[3] += uik[12] * u[0] + uik[13] * u[1] + uik[14] * u[2] + uik[15] * u[3];

157:           rtmp_ptr[4] += uik[0] * u[4] + uik[1] * u[5] + uik[2] * u[6] + uik[3] * u[7];
158:           rtmp_ptr[5] += uik[4] * u[4] + uik[5] * u[5] + uik[6] * u[6] + uik[7] * u[7];
159:           rtmp_ptr[6] += uik[8] * u[4] + uik[9] * u[5] + uik[10] * u[6] + uik[11] * u[7];
160:           rtmp_ptr[7] += uik[12] * u[4] + uik[13] * u[5] + uik[14] * u[6] + uik[15] * u[7];

162:           rtmp_ptr[8] += uik[0] * u[8] + uik[1] * u[9] + uik[2] * u[10] + uik[3] * u[11];
163:           rtmp_ptr[9] += uik[4] * u[8] + uik[5] * u[9] + uik[6] * u[10] + uik[7] * u[11];
164:           rtmp_ptr[10] += uik[8] * u[8] + uik[9] * u[9] + uik[10] * u[10] + uik[11] * u[11];
165:           rtmp_ptr[11] += uik[12] * u[8] + uik[13] * u[9] + uik[14] * u[10] + uik[15] * u[11];

167:           rtmp_ptr[12] += uik[0] * u[12] + uik[1] * u[13] + uik[2] * u[14] + uik[3] * u[15];
168:           rtmp_ptr[13] += uik[4] * u[12] + uik[5] * u[13] + uik[6] * u[14] + uik[7] * u[15];
169:           rtmp_ptr[14] += uik[8] * u[12] + uik[9] * u[13] + uik[10] * u[14] + uik[11] * u[15];
170:           rtmp_ptr[15] += uik[12] * u[12] + uik[13] * u[13] + uik[14] * u[14] + uik[15] * u[15];
171:         }
172:         PetscCall(PetscLogFlops(2.0 * 64.0 * (jmax - jmin)));

174:         /* ... add i to row list for next nonzero entry */
175:         il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
176:         j     = bj[jmin];
177:         jl[i] = jl[j];
178:         jl[j] = i; /* update jl */
179:       }
180:       i = nexti;
181:     }

183:     /* save nonzero entries in k-th row of U ... */

185:     /* invert diagonal block */
186:     diag = ba + k * 16;
187:     PetscCall(PetscArraycpy(diag, dk, 16));

189:     if (pivotinblocks) {
190:       PetscCall(PetscKernel_A_gets_inverse_A_4(diag, shift, allowzeropivot, &zeropivotdetected));
191:       if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
192:     } else {
193:       PetscCall(PetscKernel_A_gets_inverse_A_4_nopivot(diag));
194:     }

196:     jmin = bi[k];
197:     jmax = bi[k + 1];
198:     if (jmin < jmax) {
199:       for (j = jmin; j < jmax; j++) {
200:         vj       = bj[j]; /* block col. index of U */
201:         u        = ba + j * 16;
202:         rtmp_ptr = rtmp + vj * 16;
203:         for (k1 = 0; k1 < 16; k1++) {
204:           *u++        = *rtmp_ptr;
205:           *rtmp_ptr++ = 0.0;
206:         }
207:       }

209:       /* ... add k to row list for first nonzero entry in k-th row */
210:       il[k] = jmin;
211:       i     = bj[jmin];
212:       jl[k] = jl[i];
213:       jl[i] = k;
214:     }
215:   }

217:   PetscCall(PetscFree(rtmp));
218:   PetscCall(PetscFree2(il, jl));
219:   PetscCall(PetscFree2(dk, uik));
220:   if (a->permute) PetscCall(PetscFree(aa));

222:   PetscCall(ISRestoreIndices(perm, &perm_ptr));

224:   C->ops->solve          = MatSolve_SeqSBAIJ_4_inplace;
225:   C->ops->solvetranspose = MatSolve_SeqSBAIJ_4_inplace;
226:   C->assembled           = PETSC_TRUE;
227:   C->preallocated        = PETSC_TRUE;

229:   PetscCall(PetscLogFlops(1.3333 * 64 * b->mbs)); /* from inverting diagonal blocks */
230:   PetscFunctionReturn(PETSC_SUCCESS);
231: }