1 #define PETSCMAT_DLL 2 3 #include "src/mat/impls/sbaij/seq/sbaij.h" 4 #include "src/inline/ilu.h" 5 6 /* Version for when blocks are 3 by 3 */ 7 #undef __FUNCT__ 8 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_3" 9 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat A,MatFactorInfo *info,Mat *B) 10 { 11 Mat C = *B; 12 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 13 IS perm = b->row; 14 PetscErrorCode ierr; 15 PetscInt *perm_ptr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 16 PetscInt *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 17 MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 18 MatScalar *u,*diag,*rtmp,*rtmp_ptr; 19 PetscReal shift = info->shiftinblocks; 20 21 PetscFunctionBegin; 22 23 /* initialization */ 24 ierr = PetscMalloc(9*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 25 ierr = PetscMemzero(rtmp,9*mbs*sizeof(MatScalar));CHKERRQ(ierr); 26 ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); 27 jl = il + mbs; 28 for (i=0; i<mbs; i++) { 29 jl[i] = mbs; il[0] = 0; 30 } 31 ierr = PetscMalloc(18*sizeof(MatScalar),&dk);CHKERRQ(ierr); 32 uik = dk + 9; 33 ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr); 34 35 /* check permutation */ 36 if (!a->permute){ 37 ai = a->i; aj = a->j; aa = a->a; 38 } else { 39 ai = a->inew; aj = a->jnew; 40 ierr = PetscMalloc(9*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); 41 ierr = PetscMemcpy(aa,a->a,9*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); 42 ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr); 43 ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr); 44 45 for (i=0; i<mbs; i++){ 46 jmin = ai[i]; jmax = ai[i+1]; 47 for (j=jmin; j<jmax; j++){ 48 while (a2anew[j] != j){ 49 k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k; 50 for (k1=0; k1<9; k1++){ 51 dk[k1] = aa[k*9+k1]; 52 aa[k*9+k1] = aa[j*9+k1]; 53 aa[j*9+k1] = dk[k1]; 54 } 55 } 56 /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */ 57 if (i > aj[j]){ 58 /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */ 59 ap = aa + j*9; /* ptr to the beginning of j-th block of aa */ 60 for (k=0; k<9; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */ 61 for (k=0; k<3; k++){ /* j-th block of aa <- dk^T */ 62 for (k1=0; k1<3; k1++) *ap++ = dk[k + 3*k1]; 63 } 64 } 65 } 66 } 67 ierr = PetscFree(a2anew);CHKERRQ(ierr); 68 } 69 70 /* for each row k */ 71 for (k = 0; k<mbs; k++){ 72 73 /*initialize k-th row with elements nonzero in row perm(k) of A */ 74 jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1]; 75 if (jmin < jmax) { 76 ap = aa + jmin*9; 77 for (j = jmin; j < jmax; j++){ 78 vj = perm_ptr[aj[j]]; /* block col. index */ 79 rtmp_ptr = rtmp + vj*9; 80 for (i=0; i<9; i++) *rtmp_ptr++ = *ap++; 81 } 82 } 83 84 /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 85 ierr = PetscMemcpy(dk,rtmp+k*9,9*sizeof(MatScalar));CHKERRQ(ierr); 86 i = jl[k]; /* first row to be added to k_th row */ 87 88 while (i < mbs){ 89 nexti = jl[i]; /* next row to be added to k_th row */ 90 91 /* compute multiplier */ 92 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 93 94 /* uik = -inv(Di)*U_bar(i,k) */ 95 diag = ba + i*9; 96 u = ba + ili*9; 97 98 uik[0] = -(diag[0]*u[0] + diag[3]*u[1] + diag[6]*u[2]); 99 uik[1] = -(diag[1]*u[0] + diag[4]*u[1] + diag[7]*u[2]); 100 uik[2] = -(diag[2]*u[0] + diag[5]*u[1] + diag[8]*u[2]); 101 102 uik[3] = -(diag[0]*u[3] + diag[3]*u[4] + diag[6]*u[5]); 103 uik[4] = -(diag[1]*u[3] + diag[4]*u[4] + diag[7]*u[5]); 104 uik[5] = -(diag[2]*u[3] + diag[5]*u[4] + diag[8]*u[5]); 105 106 uik[6] = -(diag[0]*u[6] + diag[3]*u[7] + diag[6]*u[8]); 107 uik[7] = -(diag[1]*u[6] + diag[4]*u[7] + diag[7]*u[8]); 108 uik[8] = -(diag[2]*u[6] + diag[5]*u[7] + diag[8]*u[8]); 109 110 /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 111 dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2]; 112 dk[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2]; 113 dk[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2]; 114 115 dk[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5]; 116 dk[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5]; 117 dk[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5]; 118 119 dk[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8]; 120 dk[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8]; 121 dk[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8]; 122 123 ierr = PetscLogFlops(27*4);CHKERRQ(ierr); 124 125 /* update -U(i,k) */ 126 ierr = PetscMemcpy(ba+ili*9,uik,9*sizeof(MatScalar));CHKERRQ(ierr); 127 128 /* add multiple of row i to k-th row ... */ 129 jmin = ili + 1; jmax = bi[i+1]; 130 if (jmin < jmax){ 131 for (j=jmin; j<jmax; j++) { 132 /* rtmp += -U(i,k)^T * U_bar(i,j) */ 133 rtmp_ptr = rtmp + bj[j]*9; 134 u = ba + j*9; 135 rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2]; 136 rtmp_ptr[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2]; 137 rtmp_ptr[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2]; 138 139 rtmp_ptr[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5]; 140 rtmp_ptr[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5]; 141 rtmp_ptr[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5]; 142 143 rtmp_ptr[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8]; 144 rtmp_ptr[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8]; 145 rtmp_ptr[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8]; 146 } 147 ierr = PetscLogFlops(2*27*(jmax-jmin));CHKERRQ(ierr); 148 149 /* ... add i to row list for next nonzero entry */ 150 il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 151 j = bj[jmin]; 152 jl[i] = jl[j]; jl[j] = i; /* update jl */ 153 } 154 i = nexti; 155 } 156 157 /* save nonzero entries in k-th row of U ... */ 158 159 /* invert diagonal block */ 160 diag = ba+k*9; 161 ierr = PetscMemcpy(diag,dk,9*sizeof(MatScalar));CHKERRQ(ierr); 162 ierr = Kernel_A_gets_inverse_A_3(diag,shift);CHKERRQ(ierr); 163 164 jmin = bi[k]; jmax = bi[k+1]; 165 if (jmin < jmax) { 166 for (j=jmin; j<jmax; j++){ 167 vj = bj[j]; /* block col. index of U */ 168 u = ba + j*9; 169 rtmp_ptr = rtmp + vj*9; 170 for (k1=0; k1<9; k1++){ 171 *u++ = *rtmp_ptr; 172 *rtmp_ptr++ = 0.0; 173 } 174 } 175 176 /* ... add k to row list for first nonzero entry in k-th row */ 177 il[k] = jmin; 178 i = bj[jmin]; 179 jl[k] = jl[i]; jl[i] = k; 180 } 181 } 182 183 ierr = PetscFree(rtmp);CHKERRQ(ierr); 184 ierr = PetscFree(il);CHKERRQ(ierr); 185 ierr = PetscFree(dk);CHKERRQ(ierr); 186 if (a->permute){ 187 ierr = PetscFree(aa);CHKERRQ(ierr); 188 } 189 190 ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr); 191 C->ops->solve = MatSolve_SeqSBAIJ_3; 192 C->ops->solvetranspose = MatSolve_SeqSBAIJ_3; 193 C->assembled = PETSC_TRUE; 194 C->preallocated = PETSC_TRUE; 195 ierr = PetscLogFlops(1.3333*27*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 196 PetscFunctionReturn(0); 197 } 198