Actual source code: sbaijfact3.c

petsc-3.7.3 2016-08-01
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  */
  8: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat C,Mat A,const MatFactorInfo *info)
  9: {
 10:   Mat_SeqSBAIJ   *a   = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
 11:   IS             perm = b->row;
 13:   const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j;
 14:   PetscInt       *a2anew,i,j,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
 15:   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
 16:   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
 17:   PetscReal      shift = info->shiftamount;
 18:   PetscBool      allowzeropivot,zeropivotdetected;

 21:   /* initialization */
 22:   allowzeropivot = PetscNot(A->erroriffailure);
 23:   PetscCalloc1(9*mbs,&rtmp);
 24:   PetscMalloc2(mbs,&il,mbs,&jl);
 25:   il[0] = 0;
 26:   for (i=0; i<mbs; i++) jl[i] = mbs;
 27: 
 28:   PetscMalloc2(9,&dk,9,&uik);
 29:   ISGetIndices(perm,&perm_ptr);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

119:       PetscLogFlops(27.0*4.0);

121:       /* update -U(i,k) */
122:       PetscMemcpy(ba+ili*9,uik,9*sizeof(MatScalar));

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

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

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

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

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

155:     /* invert diagonal block */
156:     diag = ba+k*9;
157:     PetscMemcpy(diag,dk,9*sizeof(MatScalar));
158:     PetscKernel_A_gets_inverse_A_3(diag,shift,allowzeropivot,&zeropivotdetected);
159:     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;

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

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

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

187:   ISRestoreIndices(perm,&perm_ptr);

189:   C->ops->solve          = MatSolve_SeqSBAIJ_3_inplace;
190:   C->ops->solvetranspose = MatSolve_SeqSBAIJ_3_inplace;
191:   C->assembled           = PETSC_TRUE;
192:   C->preallocated        = PETSC_TRUE;

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