Actual source code: sbaijfact4.c

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

  4: /*
  5:       Version for when blocks are 3 by 3 Using natural ordering
  6: */
  7: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering(Mat C, Mat A, const MatFactorInfo *info)
  8: {
  9:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data, *b = (Mat_SeqSBAIJ *)C->data;
 10:   PetscInt      i, j, mbs = a->mbs, *bi = b->i, *bj = b->j;
 11:   PetscInt     *ai, *aj, 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:   PetscReal     shift = info->shiftamount;
 15:   PetscBool     allowzeropivot, zeropivotdetected;

 17:   PetscFunctionBegin;
 18:   /* initialization */
 19:   allowzeropivot = PetscNot(A->erroriffailure);
 20:   PetscCall(PetscCalloc1(9 * 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(9, &dk, 9, &uik));
 26:   ai = a->i;
 27:   aj = a->j;
 28:   aa = a->a;

 30:   /* for each row k */
 31:   for (k = 0; k < mbs; k++) {
 32:     /*initialize k-th row with elements nonzero in row k of A */
 33:     jmin = ai[k];
 34:     jmax = ai[k + 1];
 35:     if (jmin < jmax) {
 36:       ap = aa + jmin * 9;
 37:       for (j = jmin; j < jmax; j++) {
 38:         vj       = aj[j]; /* block col. index */
 39:         rtmp_ptr = rtmp + vj * 9;
 40:         for (i = 0; i < 9; i++) *rtmp_ptr++ = *ap++;
 41:       }
 42:     }

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

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

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

 54:       /* uik = -inv(Di)*U_bar(i,k) */
 55:       diag = ba + i * 9;
 56:       u    = ba + ili * 9;

 58:       uik[0] = -(diag[0] * u[0] + diag[3] * u[1] + diag[6] * u[2]);
 59:       uik[1] = -(diag[1] * u[0] + diag[4] * u[1] + diag[7] * u[2]);
 60:       uik[2] = -(diag[2] * u[0] + diag[5] * u[1] + diag[8] * u[2]);

 62:       uik[3] = -(diag[0] * u[3] + diag[3] * u[4] + diag[6] * u[5]);
 63:       uik[4] = -(diag[1] * u[3] + diag[4] * u[4] + diag[7] * u[5]);
 64:       uik[5] = -(diag[2] * u[3] + diag[5] * u[4] + diag[8] * u[5]);

 66:       uik[6] = -(diag[0] * u[6] + diag[3] * u[7] + diag[6] * u[8]);
 67:       uik[7] = -(diag[1] * u[6] + diag[4] * u[7] + diag[7] * u[8]);
 68:       uik[8] = -(diag[2] * u[6] + diag[5] * u[7] + diag[8] * u[8]);

 70:       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
 71:       dk[0] += uik[0] * u[0] + uik[1] * u[1] + uik[2] * u[2];
 72:       dk[1] += uik[3] * u[0] + uik[4] * u[1] + uik[5] * u[2];
 73:       dk[2] += uik[6] * u[0] + uik[7] * u[1] + uik[8] * u[2];

 75:       dk[3] += uik[0] * u[3] + uik[1] * u[4] + uik[2] * u[5];
 76:       dk[4] += uik[3] * u[3] + uik[4] * u[4] + uik[5] * u[5];
 77:       dk[5] += uik[6] * u[3] + uik[7] * u[4] + uik[8] * u[5];

 79:       dk[6] += uik[0] * u[6] + uik[1] * u[7] + uik[2] * u[8];
 80:       dk[7] += uik[3] * u[6] + uik[4] * u[7] + uik[5] * u[8];
 81:       dk[8] += uik[6] * u[6] + uik[7] * u[7] + uik[8] * u[8];

 83:       PetscCall(PetscLogFlops(27.0 * 4.0));

 85:       /* update -U(i,k) */
 86:       PetscCall(PetscArraycpy(ba + ili * 9, uik, 9));

 88:       /* add multiple of row i to k-th row ... */
 89:       jmin = ili + 1;
 90:       jmax = bi[i + 1];
 91:       if (jmin < jmax) {
 92:         for (j = jmin; j < jmax; j++) {
 93:           /* rtmp += -U(i,k)^T * U_bar(i,j) */
 94:           rtmp_ptr = rtmp + bj[j] * 9;
 95:           u        = ba + j * 9;
 96:           rtmp_ptr[0] += uik[0] * u[0] + uik[1] * u[1] + uik[2] * u[2];
 97:           rtmp_ptr[1] += uik[3] * u[0] + uik[4] * u[1] + uik[5] * u[2];
 98:           rtmp_ptr[2] += uik[6] * u[0] + uik[7] * u[1] + uik[8] * u[2];

100:           rtmp_ptr[3] += uik[0] * u[3] + uik[1] * u[4] + uik[2] * u[5];
101:           rtmp_ptr[4] += uik[3] * u[3] + uik[4] * u[4] + uik[5] * u[5];
102:           rtmp_ptr[5] += uik[6] * u[3] + uik[7] * u[4] + uik[8] * u[5];

104:           rtmp_ptr[6] += uik[0] * u[6] + uik[1] * u[7] + uik[2] * u[8];
105:           rtmp_ptr[7] += uik[3] * u[6] + uik[4] * u[7] + uik[5] * u[8];
106:           rtmp_ptr[8] += uik[6] * u[6] + uik[7] * u[7] + uik[8] * u[8];
107:         }
108:         PetscCall(PetscLogFlops(2.0 * 27.0 * (jmax - jmin)));

110:         /* ... add i to row list for next nonzero entry */
111:         il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
112:         j     = bj[jmin];
113:         jl[i] = jl[j];
114:         jl[j] = i; /* update jl */
115:       }
116:       i = nexti;
117:     }

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

121:     /* invert diagonal block */
122:     diag = ba + k * 9;
123:     PetscCall(PetscArraycpy(diag, dk, 9));
124:     PetscCall(PetscKernel_A_gets_inverse_A_3(diag, shift, allowzeropivot, &zeropivotdetected));
125:     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;

127:     jmin = bi[k];
128:     jmax = bi[k + 1];
129:     if (jmin < jmax) {
130:       for (j = jmin; j < jmax; j++) {
131:         vj       = bj[j]; /* block col. index of U */
132:         u        = ba + j * 9;
133:         rtmp_ptr = rtmp + vj * 9;
134:         for (k1 = 0; k1 < 9; k1++) {
135:           *u++        = *rtmp_ptr;
136:           *rtmp_ptr++ = 0.0;
137:         }
138:       }

140:       /* ... add k to row list for first nonzero entry in k-th row */
141:       il[k] = jmin;
142:       i     = bj[jmin];
143:       jl[k] = jl[i];
144:       jl[i] = k;
145:     }
146:   }

148:   PetscCall(PetscFree(rtmp));
149:   PetscCall(PetscFree2(il, jl));
150:   PetscCall(PetscFree2(dk, uik));

152:   C->ops->solve          = MatSolve_SeqSBAIJ_3_NaturalOrdering_inplace;
153:   C->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering_inplace;
154:   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace;
155:   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace;

157:   C->assembled    = PETSC_TRUE;
158:   C->preallocated = PETSC_TRUE;

160:   PetscCall(PetscLogFlops(1.3333 * 27 * b->mbs)); /* from inverting diagonal blocks */
161:   PetscFunctionReturn(PETSC_SUCCESS);
162: }