Actual source code: sbaijfact4.c

  1: #define PETSCMAT_DLL

 3:  #include src/mat/impls/sbaij/seq/sbaij.h
 4:  #include src/inline/ilu.h

  6: /*
  7:       Version for when blocks are 3 by 3 Using natural ordering
  8: */
 11: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B)
 12: {
 13:   Mat            C = *B;
 14:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
 16:   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
 17:   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
 18:   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
 19:   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;

 22: 
 23:   /* initialization */
 24:   PetscMalloc(9*mbs*sizeof(MatScalar),&rtmp);
 25:   PetscMemzero(rtmp,9*mbs*sizeof(MatScalar));
 26:   PetscMalloc(2*mbs*sizeof(PetscInt),&il);
 27:   jl   = il + mbs;
 28:   for (i=0; i<mbs; i++) {
 29:     jl[i] = mbs; il[0] = 0;
 30:   }
 31:   PetscMalloc(18*sizeof(MatScalar),&dk);
 32:   uik  = dk + 9;

 34:   ai   = a->i; aj = a->j; aa = a->a;

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

 39:     /*initialize k-th row with elements nonzero in row k of A */
 40:     jmin = ai[k]; jmax = ai[k+1];
 41:     if (jmin < jmax) {
 42:       ap = aa + jmin*9;
 43:       for (j = jmin; j < jmax; j++){
 44:         vj = aj[j];         /* block col. index */
 45:         rtmp_ptr = rtmp + vj*9;
 46:         for (i=0; i<9; i++) *rtmp_ptr++ = *ap++;
 47:       }
 48:     }

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

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

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

 60:       /* uik = -inv(Di)*U_bar(i,k) */
 61:       diag = ba + i*9;
 62:       u    = ba + ili*9;

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

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

 72:       uik[6] = -(diag[0]*u[6] + diag[3]*u[7] + diag[6]*u[8]);
 73:       uik[7] = -(diag[1]*u[6] + diag[4]*u[7] + diag[7]*u[8]);
 74:       uik[8] = -(diag[2]*u[6] + diag[5]*u[7] + diag[8]*u[8]);

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

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

 85:       dk[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8];
 86:       dk[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8];
 87:       dk[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8];

 89:       PetscLogFlops(27*4);

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

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

105:           rtmp_ptr[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5];
106:           rtmp_ptr[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5];
107:           rtmp_ptr[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5];

109:           rtmp_ptr[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8];
110:           rtmp_ptr[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8];
111:           rtmp_ptr[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8];
112:         }
113:         PetscLogFlops(2*27*(jmax-jmin));
114: 
115:         /* ... add i to row list for next nonzero entry */
116:         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
117:         j     = bj[jmin];
118:         jl[i] = jl[j]; jl[j] = i; /* update jl */
119:       }
120:       i = nexti;
121:     }

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

125:     /* invert diagonal block */
126:     diag = ba+k*9;
127:     PetscMemcpy(diag,dk,9*sizeof(MatScalar));
128:     Kernel_A_gets_inverse_A_3(diag);
129: 
130:     jmin = bi[k]; jmax = bi[k+1];
131:     if (jmin < jmax) {
132:       for (j=jmin; j<jmax; j++){
133:          vj = bj[j];           /* block col. index of U */
134:          u   = ba + j*9;
135:          rtmp_ptr = rtmp + vj*9;
136:          for (k1=0; k1<9; k1++){
137:            *u++        = *rtmp_ptr;
138:            *rtmp_ptr++ = 0.0;
139:          }
140:       }
141: 
142:       /* ... add k to row list for first nonzero entry in k-th row */
143:       il[k] = jmin;
144:       i     = bj[jmin];
145:       jl[k] = jl[i]; jl[i] = k;
146:     }
147:   }

149:   PetscFree(rtmp);
150:   PetscFree(il);
151:   PetscFree(dk);

153:   C->factor    = FACTOR_CHOLESKY;
154:   C->assembled = PETSC_TRUE;
155:   C->preallocated = PETSC_TRUE;
156:   PetscLogFlops(1.3333*27*b->mbs); /* from inverting diagonal blocks */
157:   return(0);
158: }