1 2 #include <../src/mat/impls/sbaij/seq/sbaij.h> 3 #include <petsc/private/kernels/blockinvert.h> 4 5 /* Version for when blocks are 3 by 3 */ 6 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(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 *a2anew,i,j,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 PetscReal shift = info->shiftamount; 15 PetscBool allowzeropivot,zeropivotdetected; 16 17 PetscFunctionBegin; 18 /* initialization */ 19 allowzeropivot = PetscNot(A->erroriffailure); 20 PetscCall(PetscCalloc1(9*mbs,&rtmp)); 21 PetscCall(PetscMalloc2(mbs,&il,mbs,&jl)); 22 il[0] = 0; 23 for (i=0; i<mbs; i++) jl[i] = mbs; 24 25 PetscCall(PetscMalloc2(9,&dk,9,&uik)); 26 PetscCall(ISGetIndices(perm,&perm_ptr)); 27 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 PetscCall(PetscMalloc1(9*ai[mbs],&aa)); 34 PetscCall(PetscArraycpy(aa,a->a,9*ai[mbs])); 35 PetscCall(PetscMalloc1(ai[mbs],&a2anew)); 36 PetscCall(PetscArraycpy(a2anew,a->a2anew,ai[mbs])); 37 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<9; k1++) { 44 dk[k1] = aa[k*9+k1]; 45 aa[k*9+k1] = aa[j*9+k1]; 46 aa[j*9+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*9; /* ptr to the beginning of j-th block of aa */ 53 for (k=0; k<9; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */ 54 for (k=0; k<3; k++) { /* j-th block of aa <- dk^T */ 55 for (k1=0; k1<3; k1++) *ap++ = dk[k + 3*k1]; 56 } 57 } 58 } 59 } 60 PetscCall(PetscFree(a2anew)); 61 } 62 63 /* for each row k */ 64 for (k = 0; k<mbs; k++) { 65 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*9; 70 for (j = jmin; j < jmax; j++) { 71 vj = perm_ptr[aj[j]]; /* block col. index */ 72 rtmp_ptr = rtmp + vj*9; 73 for (i=0; i<9; i++) *rtmp_ptr++ = *ap++; 74 } 75 } 76 77 /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 78 PetscCall(PetscArraycpy(dk,rtmp+k*9,9)); 79 i = jl[k]; /* first row to be added to k_th row */ 80 81 while (i < mbs) { 82 nexti = jl[i]; /* next row to be added to k_th row */ 83 84 /* compute multiplier */ 85 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 86 87 /* uik = -inv(Di)*U_bar(i,k) */ 88 diag = ba + i*9; 89 u = ba + ili*9; 90 91 uik[0] = -(diag[0]*u[0] + diag[3]*u[1] + diag[6]*u[2]); 92 uik[1] = -(diag[1]*u[0] + diag[4]*u[1] + diag[7]*u[2]); 93 uik[2] = -(diag[2]*u[0] + diag[5]*u[1] + diag[8]*u[2]); 94 95 uik[3] = -(diag[0]*u[3] + diag[3]*u[4] + diag[6]*u[5]); 96 uik[4] = -(diag[1]*u[3] + diag[4]*u[4] + diag[7]*u[5]); 97 uik[5] = -(diag[2]*u[3] + diag[5]*u[4] + diag[8]*u[5]); 98 99 uik[6] = -(diag[0]*u[6] + diag[3]*u[7] + diag[6]*u[8]); 100 uik[7] = -(diag[1]*u[6] + diag[4]*u[7] + diag[7]*u[8]); 101 uik[8] = -(diag[2]*u[6] + diag[5]*u[7] + diag[8]*u[8]); 102 103 /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 104 dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2]; 105 dk[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2]; 106 dk[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2]; 107 108 dk[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5]; 109 dk[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5]; 110 dk[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5]; 111 112 dk[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8]; 113 dk[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8]; 114 dk[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8]; 115 116 PetscCall(PetscLogFlops(27.0*4.0)); 117 118 /* update -U(i,k) */ 119 PetscCall(PetscArraycpy(ba+ili*9,uik,9)); 120 121 /* add multiple of row i to k-th row ... */ 122 jmin = ili + 1; jmax = bi[i+1]; 123 if (jmin < jmax) { 124 for (j=jmin; j<jmax; j++) { 125 /* rtmp += -U(i,k)^T * U_bar(i,j) */ 126 rtmp_ptr = rtmp + bj[j]*9; 127 u = ba + j*9; 128 rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2]; 129 rtmp_ptr[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2]; 130 rtmp_ptr[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2]; 131 132 rtmp_ptr[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5]; 133 rtmp_ptr[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5]; 134 rtmp_ptr[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5]; 135 136 rtmp_ptr[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8]; 137 rtmp_ptr[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8]; 138 rtmp_ptr[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8]; 139 } 140 PetscCall(PetscLogFlops(2.0*27.0*(jmax-jmin))); 141 142 /* ... add i to row list for next nonzero entry */ 143 il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 144 j = bj[jmin]; 145 jl[i] = jl[j]; jl[j] = i; /* update jl */ 146 } 147 i = nexti; 148 } 149 150 /* save nonzero entries in k-th row of U ... */ 151 152 /* invert diagonal block */ 153 diag = ba+k*9; 154 PetscCall(PetscArraycpy(diag,dk,9)); 155 PetscCall(PetscKernel_A_gets_inverse_A_3(diag,shift,allowzeropivot,&zeropivotdetected)); 156 if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 157 158 jmin = bi[k]; jmax = bi[k+1]; 159 if (jmin < jmax) { 160 for (j=jmin; j<jmax; j++) { 161 vj = bj[j]; /* block col. index of U */ 162 u = ba + j*9; 163 rtmp_ptr = rtmp + vj*9; 164 for (k1=0; k1<9; k1++) { 165 *u++ = *rtmp_ptr; 166 *rtmp_ptr++ = 0.0; 167 } 168 } 169 170 /* ... add k to row list for first nonzero entry in k-th row */ 171 il[k] = jmin; 172 i = bj[jmin]; 173 jl[k] = jl[i]; jl[i] = k; 174 } 175 } 176 177 PetscCall(PetscFree(rtmp)); 178 PetscCall(PetscFree2(il,jl)); 179 PetscCall(PetscFree2(dk,uik)); 180 if (a->permute) { 181 PetscCall(PetscFree(aa)); 182 } 183 184 PetscCall(ISRestoreIndices(perm,&perm_ptr)); 185 186 C->ops->solve = MatSolve_SeqSBAIJ_3_inplace; 187 C->ops->solvetranspose = MatSolve_SeqSBAIJ_3_inplace; 188 C->assembled = PETSC_TRUE; 189 C->preallocated = PETSC_TRUE; 190 191 PetscCall(PetscLogFlops(1.3333*27*b->mbs)); /* from inverting diagonal blocks */ 192 PetscFunctionReturn(0); 193 } 194