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