1 2 #include <../src/mat/impls/sbaij/seq/sbaij.h> 3 #include <petsc/private/kernels/blockinvert.h> 4 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; 17 18 PetscFunctionBegin; 19 /* initialization */ 20 allowzeropivot = PetscNot(A->erroriffailure); 21 PetscCall(PetscCalloc1(9*mbs,&rtmp)); 22 PetscCall(PetscMalloc2(mbs,&il,mbs,&jl)); 23 il[0] = 0; 24 for (i=0; i<mbs; i++) jl[i] = mbs; 25 26 PetscCall(PetscMalloc2(9,&dk,9,&uik)); 27 ai = a->i; aj = a->j; aa = a->a; 28 29 /* for each row k */ 30 for (k = 0; k<mbs; k++) { 31 32 /*initialize k-th row with elements nonzero in row k of A */ 33 jmin = ai[k]; jmax = ai[k+1]; 34 if (jmin < jmax) { 35 ap = aa + jmin*9; 36 for (j = jmin; j < jmax; j++) { 37 vj = aj[j]; /* block col. index */ 38 rtmp_ptr = rtmp + vj*9; 39 for (i=0; i<9; i++) *rtmp_ptr++ = *ap++; 40 } 41 } 42 43 /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 44 PetscCall(PetscArraycpy(dk,rtmp+k*9,9)); 45 i = jl[k]; /* first row to be added to k_th row */ 46 47 while (i < mbs) { 48 nexti = jl[i]; /* next row to be added to k_th row */ 49 50 /* compute multiplier */ 51 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 52 53 /* uik = -inv(Di)*U_bar(i,k) */ 54 diag = ba + i*9; 55 u = ba + ili*9; 56 57 uik[0] = -(diag[0]*u[0] + diag[3]*u[1] + diag[6]*u[2]); 58 uik[1] = -(diag[1]*u[0] + diag[4]*u[1] + diag[7]*u[2]); 59 uik[2] = -(diag[2]*u[0] + diag[5]*u[1] + diag[8]*u[2]); 60 61 uik[3] = -(diag[0]*u[3] + diag[3]*u[4] + diag[6]*u[5]); 62 uik[4] = -(diag[1]*u[3] + diag[4]*u[4] + diag[7]*u[5]); 63 uik[5] = -(diag[2]*u[3] + diag[5]*u[4] + diag[8]*u[5]); 64 65 uik[6] = -(diag[0]*u[6] + diag[3]*u[7] + diag[6]*u[8]); 66 uik[7] = -(diag[1]*u[6] + diag[4]*u[7] + diag[7]*u[8]); 67 uik[8] = -(diag[2]*u[6] + diag[5]*u[7] + diag[8]*u[8]); 68 69 /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 70 dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2]; 71 dk[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2]; 72 dk[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2]; 73 74 dk[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5]; 75 dk[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5]; 76 dk[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5]; 77 78 dk[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8]; 79 dk[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8]; 80 dk[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8]; 81 82 PetscCall(PetscLogFlops(27.0*4.0)); 83 84 /* update -U(i,k) */ 85 PetscCall(PetscArraycpy(ba+ili*9,uik,9)); 86 87 /* add multiple of row i to k-th row ... */ 88 jmin = ili + 1; jmax = bi[i+1]; 89 if (jmin < jmax) { 90 for (j=jmin; j<jmax; j++) { 91 /* rtmp += -U(i,k)^T * U_bar(i,j) */ 92 rtmp_ptr = rtmp + bj[j]*9; 93 u = ba + j*9; 94 rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2]; 95 rtmp_ptr[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2]; 96 rtmp_ptr[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2]; 97 98 rtmp_ptr[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5]; 99 rtmp_ptr[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5]; 100 rtmp_ptr[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5]; 101 102 rtmp_ptr[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8]; 103 rtmp_ptr[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8]; 104 rtmp_ptr[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8]; 105 } 106 PetscCall(PetscLogFlops(2.0*27.0*(jmax-jmin))); 107 108 /* ... add i to row list for next nonzero entry */ 109 il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 110 j = bj[jmin]; 111 jl[i] = jl[j]; jl[j] = i; /* update jl */ 112 } 113 i = nexti; 114 } 115 116 /* save nonzero entries in k-th row of U ... */ 117 118 /* invert diagonal block */ 119 diag = ba+k*9; 120 PetscCall(PetscArraycpy(diag,dk,9)); 121 PetscCall(PetscKernel_A_gets_inverse_A_3(diag,shift,allowzeropivot,&zeropivotdetected)); 122 if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 123 124 jmin = bi[k]; jmax = bi[k+1]; 125 if (jmin < jmax) { 126 for (j=jmin; j<jmax; j++) { 127 vj = bj[j]; /* block col. index of U */ 128 u = ba + j*9; 129 rtmp_ptr = rtmp + vj*9; 130 for (k1=0; k1<9; k1++) { 131 *u++ = *rtmp_ptr; 132 *rtmp_ptr++ = 0.0; 133 } 134 } 135 136 /* ... add k to row list for first nonzero entry in k-th row */ 137 il[k] = jmin; 138 i = bj[jmin]; 139 jl[k] = jl[i]; jl[i] = k; 140 } 141 } 142 143 PetscCall(PetscFree(rtmp)); 144 PetscCall(PetscFree2(il,jl)); 145 PetscCall(PetscFree2(dk,uik)); 146 147 C->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering_inplace; 148 C->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering_inplace; 149 C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace; 150 C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace; 151 152 C->assembled = PETSC_TRUE; 153 C->preallocated = PETSC_TRUE; 154 155 PetscCall(PetscLogFlops(1.3333*27*b->mbs)); /* from inverting diagonal blocks */ 156 PetscFunctionReturn(0); 157 } 158