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