Actual source code: sbaijfact7.c


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

  5: /* Version for when blocks are 5 by 5  */
  6: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5(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, ipvt[5];
 12:   MatScalar      *ba = b->a, *aa, *ap, *dk, *uik;
 13:   MatScalar      *u, *d, *rtmp, *rtmp_ptr, work[25];
 14:   PetscReal       shift = info->shiftamount;
 15:   PetscBool       allowzeropivot, zeropivotdetected;

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

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

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

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

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

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

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

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

 93:       /* uik = -inv(Di)*U_bar(i,k) */
 94:       d = ba + i * 25;
 95:       u = ba + ili * 25;

 97:       uik[0] = -(d[0] * u[0] + d[5] * u[1] + d[10] * u[2] + d[15] * u[3] + d[20] * u[4]);
 98:       uik[1] = -(d[1] * u[0] + d[6] * u[1] + d[11] * u[2] + d[16] * u[3] + d[21] * u[4]);
 99:       uik[2] = -(d[2] * u[0] + d[7] * u[1] + d[12] * u[2] + d[17] * u[3] + d[22] * u[4]);
100:       uik[3] = -(d[3] * u[0] + d[8] * u[1] + d[13] * u[2] + d[18] * u[3] + d[23] * u[4]);
101:       uik[4] = -(d[4] * u[0] + d[9] * u[1] + d[14] * u[2] + d[19] * u[3] + d[24] * u[4]);

103:       uik[5] = -(d[0] * u[5] + d[5] * u[6] + d[10] * u[7] + d[15] * u[8] + d[20] * u[9]);
104:       uik[6] = -(d[1] * u[5] + d[6] * u[6] + d[11] * u[7] + d[16] * u[8] + d[21] * u[9]);
105:       uik[7] = -(d[2] * u[5] + d[7] * u[6] + d[12] * u[7] + d[17] * u[8] + d[22] * u[9]);
106:       uik[8] = -(d[3] * u[5] + d[8] * u[6] + d[13] * u[7] + d[18] * u[8] + d[23] * u[9]);
107:       uik[9] = -(d[4] * u[5] + d[9] * u[6] + d[14] * u[7] + d[19] * u[8] + d[24] * u[9]);

109:       uik[10] = -(d[0] * u[10] + d[5] * u[11] + d[10] * u[12] + d[15] * u[13] + d[20] * u[14]);
110:       uik[11] = -(d[1] * u[10] + d[6] * u[11] + d[11] * u[12] + d[16] * u[13] + d[21] * u[14]);
111:       uik[12] = -(d[2] * u[10] + d[7] * u[11] + d[12] * u[12] + d[17] * u[13] + d[22] * u[14]);
112:       uik[13] = -(d[3] * u[10] + d[8] * u[11] + d[13] * u[12] + d[18] * u[13] + d[23] * u[14]);
113:       uik[14] = -(d[4] * u[10] + d[9] * u[11] + d[14] * u[12] + d[19] * u[13] + d[24] * u[14]);

115:       uik[15] = -(d[0] * u[15] + d[5] * u[16] + d[10] * u[17] + d[15] * u[18] + d[20] * u[19]);
116:       uik[16] = -(d[1] * u[15] + d[6] * u[16] + d[11] * u[17] + d[16] * u[18] + d[21] * u[19]);
117:       uik[17] = -(d[2] * u[15] + d[7] * u[16] + d[12] * u[17] + d[17] * u[18] + d[22] * u[19]);
118:       uik[18] = -(d[3] * u[15] + d[8] * u[16] + d[13] * u[17] + d[18] * u[18] + d[23] * u[19]);
119:       uik[19] = -(d[4] * u[15] + d[9] * u[16] + d[14] * u[17] + d[19] * u[18] + d[24] * u[19]);

121:       uik[20] = -(d[0] * u[20] + d[5] * u[21] + d[10] * u[22] + d[15] * u[23] + d[20] * u[24]);
122:       uik[21] = -(d[1] * u[20] + d[6] * u[21] + d[11] * u[22] + d[16] * u[23] + d[21] * u[24]);
123:       uik[22] = -(d[2] * u[20] + d[7] * u[21] + d[12] * u[22] + d[17] * u[23] + d[22] * u[24]);
124:       uik[23] = -(d[3] * u[20] + d[8] * u[21] + d[13] * u[22] + d[18] * u[23] + d[23] * u[24]);
125:       uik[24] = -(d[4] * u[20] + d[9] * u[21] + d[14] * u[22] + d[19] * u[23] + d[24] * u[24]);

127:       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
128:       dk[0] += uik[0] * u[0] + uik[1] * u[1] + uik[2] * u[2] + uik[3] * u[3] + uik[4] * u[4];
129:       dk[1] += uik[5] * u[0] + uik[6] * u[1] + uik[7] * u[2] + uik[8] * u[3] + uik[9] * u[4];
130:       dk[2] += uik[10] * u[0] + uik[11] * u[1] + uik[12] * u[2] + uik[13] * u[3] + uik[14] * u[4];
131:       dk[3] += uik[15] * u[0] + uik[16] * u[1] + uik[17] * u[2] + uik[18] * u[3] + uik[19] * u[4];
132:       dk[4] += uik[20] * u[0] + uik[21] * u[1] + uik[22] * u[2] + uik[23] * u[3] + uik[24] * u[4];

134:       dk[5] += uik[0] * u[5] + uik[1] * u[6] + uik[2] * u[7] + uik[3] * u[8] + uik[4] * u[9];
135:       dk[6] += uik[5] * u[5] + uik[6] * u[6] + uik[7] * u[7] + uik[8] * u[8] + uik[9] * u[9];
136:       dk[7] += uik[10] * u[5] + uik[11] * u[6] + uik[12] * u[7] + uik[13] * u[8] + uik[14] * u[9];
137:       dk[8] += uik[15] * u[5] + uik[16] * u[6] + uik[17] * u[7] + uik[18] * u[8] + uik[19] * u[9];
138:       dk[9] += uik[20] * u[5] + uik[21] * u[6] + uik[22] * u[7] + uik[23] * u[8] + uik[24] * u[9];

140:       dk[10] += uik[0] * u[10] + uik[1] * u[11] + uik[2] * u[12] + uik[3] * u[13] + uik[4] * u[14];
141:       dk[11] += uik[5] * u[10] + uik[6] * u[11] + uik[7] * u[12] + uik[8] * u[13] + uik[9] * u[14];
142:       dk[12] += uik[10] * u[10] + uik[11] * u[11] + uik[12] * u[12] + uik[13] * u[13] + uik[14] * u[14];
143:       dk[13] += uik[15] * u[10] + uik[16] * u[11] + uik[17] * u[12] + uik[18] * u[13] + uik[19] * u[14];
144:       dk[14] += uik[20] * u[10] + uik[21] * u[11] + uik[22] * u[12] + uik[23] * u[13] + uik[24] * u[14];

146:       dk[15] += uik[0] * u[15] + uik[1] * u[16] + uik[2] * u[17] + uik[3] * u[18] + uik[4] * u[19];
147:       dk[16] += uik[5] * u[15] + uik[6] * u[16] + uik[7] * u[17] + uik[8] * u[18] + uik[9] * u[19];
148:       dk[17] += uik[10] * u[15] + uik[11] * u[16] + uik[12] * u[17] + uik[13] * u[18] + uik[14] * u[19];
149:       dk[18] += uik[15] * u[15] + uik[16] * u[16] + uik[17] * u[17] + uik[18] * u[18] + uik[19] * u[19];
150:       dk[19] += uik[20] * u[15] + uik[21] * u[16] + uik[22] * u[17] + uik[23] * u[18] + uik[24] * u[19];

152:       dk[20] += uik[0] * u[20] + uik[1] * u[21] + uik[2] * u[22] + uik[3] * u[23] + uik[4] * u[24];
153:       dk[21] += uik[5] * u[20] + uik[6] * u[21] + uik[7] * u[22] + uik[8] * u[23] + uik[9] * u[24];
154:       dk[22] += uik[10] * u[20] + uik[11] * u[21] + uik[12] * u[22] + uik[13] * u[23] + uik[14] * u[24];
155:       dk[23] += uik[15] * u[20] + uik[16] * u[21] + uik[17] * u[22] + uik[18] * u[23] + uik[19] * u[24];
156:       dk[24] += uik[20] * u[20] + uik[21] * u[21] + uik[22] * u[22] + uik[23] * u[23] + uik[24] * u[24];

158:       PetscCall(PetscLogFlops(125.0 * 4.0));

160:       /* update -U(i,k) */
161:       PetscCall(PetscArraycpy(ba + ili * 25, uik, 25));

163:       /* add multiple of row i to k-th row ... */
164:       jmin = ili + 1;
165:       jmax = bi[i + 1];
166:       if (jmin < jmax) {
167:         for (j = jmin; j < jmax; j++) {
168:           /* rtmp += -U(i,k)^T * U_bar(i,j) */
169:           rtmp_ptr = rtmp + bj[j] * 25;
170:           u        = ba + j * 25;
171:           rtmp_ptr[0] += uik[0] * u[0] + uik[1] * u[1] + uik[2] * u[2] + uik[3] * u[3] + uik[4] * u[4];
172:           rtmp_ptr[1] += uik[5] * u[0] + uik[6] * u[1] + uik[7] * u[2] + uik[8] * u[3] + uik[9] * u[4];
173:           rtmp_ptr[2] += uik[10] * u[0] + uik[11] * u[1] + uik[12] * u[2] + uik[13] * u[3] + uik[14] * u[4];
174:           rtmp_ptr[3] += uik[15] * u[0] + uik[16] * u[1] + uik[17] * u[2] + uik[18] * u[3] + uik[19] * u[4];
175:           rtmp_ptr[4] += uik[20] * u[0] + uik[21] * u[1] + uik[22] * u[2] + uik[23] * u[3] + uik[24] * u[4];

177:           rtmp_ptr[5] += uik[0] * u[5] + uik[1] * u[6] + uik[2] * u[7] + uik[3] * u[8] + uik[4] * u[9];
178:           rtmp_ptr[6] += uik[5] * u[5] + uik[6] * u[6] + uik[7] * u[7] + uik[8] * u[8] + uik[9] * u[9];
179:           rtmp_ptr[7] += uik[10] * u[5] + uik[11] * u[6] + uik[12] * u[7] + uik[13] * u[8] + uik[14] * u[9];
180:           rtmp_ptr[8] += uik[15] * u[5] + uik[16] * u[6] + uik[17] * u[7] + uik[18] * u[8] + uik[19] * u[9];
181:           rtmp_ptr[9] += uik[20] * u[5] + uik[21] * u[6] + uik[22] * u[7] + uik[23] * u[8] + uik[24] * u[9];

183:           rtmp_ptr[10] += uik[0] * u[10] + uik[1] * u[11] + uik[2] * u[12] + uik[3] * u[13] + uik[4] * u[14];
184:           rtmp_ptr[11] += uik[5] * u[10] + uik[6] * u[11] + uik[7] * u[12] + uik[8] * u[13] + uik[9] * u[14];
185:           rtmp_ptr[12] += uik[10] * u[10] + uik[11] * u[11] + uik[12] * u[12] + uik[13] * u[13] + uik[14] * u[14];
186:           rtmp_ptr[13] += uik[15] * u[10] + uik[16] * u[11] + uik[17] * u[12] + uik[18] * u[13] + uik[19] * u[14];
187:           rtmp_ptr[14] += uik[20] * u[10] + uik[21] * u[11] + uik[22] * u[12] + uik[23] * u[13] + uik[24] * u[14];

189:           rtmp_ptr[15] += uik[0] * u[15] + uik[1] * u[16] + uik[2] * u[17] + uik[3] * u[18] + uik[4] * u[19];
190:           rtmp_ptr[16] += uik[5] * u[15] + uik[6] * u[16] + uik[7] * u[17] + uik[8] * u[18] + uik[9] * u[19];
191:           rtmp_ptr[17] += uik[10] * u[15] + uik[11] * u[16] + uik[12] * u[17] + uik[13] * u[18] + uik[14] * u[19];
192:           rtmp_ptr[18] += uik[15] * u[15] + uik[16] * u[16] + uik[17] * u[17] + uik[18] * u[18] + uik[19] * u[19];
193:           rtmp_ptr[19] += uik[20] * u[15] + uik[21] * u[16] + uik[22] * u[17] + uik[23] * u[18] + uik[24] * u[19];

195:           rtmp_ptr[20] += uik[0] * u[20] + uik[1] * u[21] + uik[2] * u[22] + uik[3] * u[23] + uik[4] * u[24];
196:           rtmp_ptr[21] += uik[5] * u[20] + uik[6] * u[21] + uik[7] * u[22] + uik[8] * u[23] + uik[9] * u[24];
197:           rtmp_ptr[22] += uik[10] * u[20] + uik[11] * u[21] + uik[12] * u[22] + uik[13] * u[23] + uik[14] * u[24];
198:           rtmp_ptr[23] += uik[15] * u[20] + uik[16] * u[21] + uik[17] * u[22] + uik[18] * u[23] + uik[19] * u[24];
199:           rtmp_ptr[24] += uik[20] * u[20] + uik[21] * u[21] + uik[22] * u[22] + uik[23] * u[23] + uik[24] * u[24];
200:         }
201:         PetscCall(PetscLogFlops(2.0 * 125.0 * (jmax - jmin)));

203:         /* ... add i to row list for next nonzero entry */
204:         il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
205:         j     = bj[jmin];
206:         jl[i] = jl[j];
207:         jl[j] = i; /* update jl */
208:       }
209:       i = nexti;
210:     }

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

214:     /* invert diagonal block */
215:     d = ba + k * 25;
216:     PetscCall(PetscArraycpy(d, dk, 25));
217:     PetscCall(PetscKernel_A_gets_inverse_A_5(d, ipvt, work, shift, allowzeropivot, &zeropivotdetected));
218:     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;

220:     jmin = bi[k];
221:     jmax = bi[k + 1];
222:     if (jmin < jmax) {
223:       for (j = jmin; j < jmax; j++) {
224:         vj       = bj[j]; /* block col. index of U */
225:         u        = ba + j * 25;
226:         rtmp_ptr = rtmp + vj * 25;
227:         for (k1 = 0; k1 < 25; k1++) {
228:           *u++        = *rtmp_ptr;
229:           *rtmp_ptr++ = 0.0;
230:         }
231:       }

233:       /* ... add k to row list for first nonzero entry in k-th row */
234:       il[k] = jmin;
235:       i     = bj[jmin];
236:       jl[k] = jl[i];
237:       jl[i] = k;
238:     }
239:   }

241:   PetscCall(PetscFree(rtmp));
242:   PetscCall(PetscFree2(il, jl));
243:   PetscCall(PetscFree2(dk, uik));
244:   if (a->permute) PetscCall(PetscFree(aa));

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

248:   C->ops->solve          = MatSolve_SeqSBAIJ_5_inplace;
249:   C->ops->solvetranspose = MatSolve_SeqSBAIJ_5_inplace;
250:   C->assembled           = PETSC_TRUE;
251:   C->preallocated        = PETSC_TRUE;

253:   PetscCall(PetscLogFlops(1.3333 * 125 * b->mbs)); /* from inverting diagonal blocks */
254:   PetscFunctionReturn(PETSC_SUCCESS);
255: }