Actual source code: sbaijfact4.c


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

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

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

 25:   PetscMalloc2(9,&dk,9,&uik);
 26:   ai   = a->i; aj = a->j; aa = a->a;

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

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

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

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

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

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

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

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

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

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

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

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

 81:       PetscLogFlops(27.0*4.0);

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

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

 97:           rtmp_ptr[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5];
 98:           rtmp_ptr[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5];
 99:           rtmp_ptr[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5];

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

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

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

117:     /* invert diagonal block */
118:     diag = ba+k*9;
119:     PetscArraycpy(diag,dk,9);
120:     PetscKernel_A_gets_inverse_A_3(diag,shift,allowzeropivot,&zeropivotdetected);
121:     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;

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

135:       /* ... add k to row list for first nonzero entry in k-th row */
136:       il[k] = jmin;
137:       i     = bj[jmin];
138:       jl[k] = jl[i]; jl[i] = k;
139:     }
140:   }

142:   PetscFree(rtmp);
143:   PetscFree2(il,jl);
144:   PetscFree2(dk,uik);

146:   C->ops->solve          = MatSolve_SeqSBAIJ_3_NaturalOrdering_inplace;
147:   C->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering_inplace;
148:   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace;
149:   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace;

151:   C->assembled    = PETSC_TRUE;
152:   C->preallocated = PETSC_TRUE;

154:   PetscLogFlops(1.3333*27*b->mbs); /* from inverting diagonal blocks */
155:   return 0;
156: }