Actual source code: sbaijfact3.c

petsc-3.12.5 2020-03-29
Report Typos and Errors

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

  5: /* Version for when blocks are 3 by 3  */
  6: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(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;
 11:   const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j;
 12:   PetscInt       *a2anew,i,j,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
 13:   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
 14:   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
 15:   PetscReal      shift = info->shiftamount;
 16:   PetscBool      allowzeropivot,zeropivotdetected;

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

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

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

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

 64:   /* for each row k */
 65:   for (k = 0; k<mbs; k++) {

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

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

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

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

 88:       /* uik = -inv(Di)*U_bar(i,k) */
 89:       diag = ba + i*9;
 90:       u    = ba + ili*9;

 92:       uik[0] = -(diag[0]*u[0] + diag[3]*u[1] + diag[6]*u[2]);
 93:       uik[1] = -(diag[1]*u[0] + diag[4]*u[1] + diag[7]*u[2]);
 94:       uik[2] = -(diag[2]*u[0] + diag[5]*u[1] + diag[8]*u[2]);

 96:       uik[3] = -(diag[0]*u[3] + diag[3]*u[4] + diag[6]*u[5]);
 97:       uik[4] = -(diag[1]*u[3] + diag[4]*u[4] + diag[7]*u[5]);
 98:       uik[5] = -(diag[2]*u[3] + diag[5]*u[4] + diag[8]*u[5]);

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

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

109:       dk[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5];
110:       dk[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5];
111:       dk[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5];

113:       dk[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8];
114:       dk[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8];
115:       dk[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8];

117:       PetscLogFlops(27.0*4.0);

119:       /* update -U(i,k) */
120:       PetscArraycpy(ba+ili*9,uik,9);

122:       /* add multiple of row i to k-th row ... */
123:       jmin = ili + 1; jmax = bi[i+1];
124:       if (jmin < jmax) {
125:         for (j=jmin; j<jmax; j++) {
126:           /* rtmp += -U(i,k)^T * U_bar(i,j) */
127:           rtmp_ptr     = rtmp + bj[j]*9;
128:           u            = ba + j*9;
129:           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2];
130:           rtmp_ptr[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2];
131:           rtmp_ptr[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2];

133:           rtmp_ptr[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5];
134:           rtmp_ptr[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5];
135:           rtmp_ptr[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5];

137:           rtmp_ptr[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8];
138:           rtmp_ptr[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8];
139:           rtmp_ptr[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8];
140:         }
141:         PetscLogFlops(2.0*27.0*(jmax-jmin));

143:         /* ... add i to row list for next nonzero entry */
144:         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
145:         j     = bj[jmin];
146:         jl[i] = jl[j]; jl[j] = i; /* update jl */
147:       }
148:       i = nexti;
149:     }

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

153:     /* invert diagonal block */
154:     diag = ba+k*9;
155:     PetscArraycpy(diag,dk,9);
156:     PetscKernel_A_gets_inverse_A_3(diag,shift,allowzeropivot,&zeropivotdetected);
157:     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;

159:     jmin = bi[k]; jmax = bi[k+1];
160:     if (jmin < jmax) {
161:       for (j=jmin; j<jmax; j++) {
162:         vj       = bj[j];      /* block col. index of U */
163:         u        = ba + j*9;
164:         rtmp_ptr = rtmp + vj*9;
165:         for (k1=0; k1<9; k1++) {
166:           *u++        = *rtmp_ptr;
167:           *rtmp_ptr++ = 0.0;
168:         }
169:       }

171:       /* ... add k to row list for first nonzero entry in k-th row */
172:       il[k] = jmin;
173:       i     = bj[jmin];
174:       jl[k] = jl[i]; jl[i] = k;
175:     }
176:   }

178:   PetscFree(rtmp);
179:   PetscFree2(il,jl);
180:   PetscFree2(dk,uik);
181:   if (a->permute) {
182:     PetscFree(aa);
183:   }

185:   ISRestoreIndices(perm,&perm_ptr);

187:   C->ops->solve          = MatSolve_SeqSBAIJ_3_inplace;
188:   C->ops->solvetranspose = MatSolve_SeqSBAIJ_3_inplace;
189:   C->assembled           = PETSC_TRUE;
190:   C->preallocated        = PETSC_TRUE;

192:   PetscLogFlops(1.3333*27*b->mbs); /* from inverting diagonal blocks */
193:   return(0);
194: }