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

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

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

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

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

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

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

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

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

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

 91:       uik[0] = -(diag[0]*u[0] + diag[4]*u[1] + diag[8]*u[2] + diag[12]*u[3]);
 92:       uik[1] = -(diag[1]*u[0] + diag[5]*u[1] + diag[9]*u[2] + diag[13]*u[3]);
 93:       uik[2] = -(diag[2]*u[0] + diag[6]*u[1] + diag[10]*u[2]+ diag[14]*u[3]);
 94:       uik[3] = -(diag[3]*u[0] + diag[7]*u[1] + diag[11]*u[2]+ diag[15]*u[3]);

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

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

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

111:       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
112:       dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2] + uik[3]*u[3];
113:       dk[1] += uik[4]*u[0] + uik[5]*u[1] + uik[6]*u[2] + uik[7]*u[3];
114:       dk[2] += uik[8]*u[0] + uik[9]*u[1] + uik[10]*u[2]+ uik[11]*u[3];
115:       dk[3] += uik[12]*u[0]+ uik[13]*u[1]+ uik[14]*u[2]+ uik[15]*u[3];

117:       dk[4] += uik[0]*u[4] + uik[1]*u[5] + uik[2]*u[6] + uik[3]*u[7];
118:       dk[5] += uik[4]*u[4] + uik[5]*u[5] + uik[6]*u[6] + uik[7]*u[7];
119:       dk[6] += uik[8]*u[4] + uik[9]*u[5] + uik[10]*u[6]+ uik[11]*u[7];
120:       dk[7] += uik[12]*u[4]+ uik[13]*u[5]+ uik[14]*u[6]+ uik[15]*u[7];

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

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

132:       PetscLogFlops(64.0*4.0);

134:       /* update -U(i,k) */
135:       PetscArraycpy(ba+ili*16,uik,16);

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

149:           rtmp_ptr[4] += uik[0]*u[4] + uik[1]*u[5] + uik[2]*u[6] + uik[3]*u[7];
150:           rtmp_ptr[5] += uik[4]*u[4] + uik[5]*u[5] + uik[6]*u[6] + uik[7]*u[7];
151:           rtmp_ptr[6] += uik[8]*u[4] + uik[9]*u[5] + uik[10]*u[6]+ uik[11]*u[7];
152:           rtmp_ptr[7] += uik[12]*u[4]+ uik[13]*u[5]+ uik[14]*u[6]+ uik[15]*u[7];

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

159:           rtmp_ptr[12]+= uik[0]*u[12] + uik[1]*u[13] + uik[2]*u[14] + uik[3]*u[15];
160:           rtmp_ptr[13]+= uik[4]*u[12] + uik[5]*u[13] + uik[6]*u[14] + uik[7]*u[15];
161:           rtmp_ptr[14]+= uik[8]*u[12] + uik[9]*u[13] + uik[10]*u[14]+ uik[11]*u[15];
162:           rtmp_ptr[15]+= uik[12]*u[12]+ uik[13]*u[13]+ uik[14]*u[14]+ uik[15]*u[15];
163:         }
164:         PetscLogFlops(2.0*64.0*(jmax-jmin));

166:         /* ... add i to row list for next nonzero entry */
167:         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
168:         j     = bj[jmin];
169:         jl[i] = jl[j]; jl[j] = i; /* update jl */
170:       }
171:       i = nexti;
172:     }

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

176:     /* invert diagonal block */
177:     diag = ba+k*16;
178:     PetscArraycpy(diag,dk,16);

180:     if (pivotinblocks) {
181:       PetscKernel_A_gets_inverse_A_4(diag,shift, allowzeropivot,&zeropivotdetected);
182:       if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
183:     } else {
184:       PetscKernel_A_gets_inverse_A_4_nopivot(diag);
185:     }

187:     jmin = bi[k]; jmax = bi[k+1];
188:     if (jmin < jmax) {
189:       for (j=jmin; j<jmax; j++) {
190:         vj       = bj[j];      /* block col. index of U */
191:         u        = ba + j*16;
192:         rtmp_ptr = rtmp + vj*16;
193:         for (k1=0; k1<16; k1++) {
194:           *u++        = *rtmp_ptr;
195:           *rtmp_ptr++ = 0.0;
196:         }
197:       }

199:       /* ... add k to row list for first nonzero entry in k-th row */
200:       il[k] = jmin;
201:       i     = bj[jmin];
202:       jl[k] = jl[i]; jl[i] = k;
203:     }
204:   }

206:   PetscFree(rtmp);
207:   PetscFree2(il,jl);
208:   PetscFree2(dk,uik);
209:   if (a->permute) {
210:     PetscFree(aa);
211:   }

213:   ISRestoreIndices(perm,&perm_ptr);

215:   C->ops->solve          = MatSolve_SeqSBAIJ_4_inplace;
216:   C->ops->solvetranspose = MatSolve_SeqSBAIJ_4_inplace;
217:   C->assembled           = PETSC_TRUE;
218:   C->preallocated        = PETSC_TRUE;

220:   PetscLogFlops(1.3333*64*b->mbs); /* from inverting diagonal blocks */
221:   return 0;
222: }