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 4 by 4 */ 7 #undef __FUNCT__ 8 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_4" 9 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(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 PetscTruth pivotinblocks = b->pivotinblocks; 20 PetscReal shift = info->shiftinblocks; 21 22 PetscFunctionBegin; 23 /* initialization */ 24 ierr = PetscMalloc(16*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 25 ierr = PetscMemzero(rtmp,16*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(32*sizeof(MatScalar),&dk);CHKERRQ(ierr); 32 uik = dk + 16; 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(16*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); 41 ierr = PetscMemcpy(aa,a->a,16*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<16; k1++){ 51 dk[k1] = aa[k*16+k1]; 52 aa[k*16+k1] = aa[j*16+k1]; 53 aa[j*16+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*16; /* ptr to the beginning of j-th block of aa */ 60 for (k=0; k<16; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */ 61 for (k=0; k<4; k++){ /* j-th block of aa <- dk^T */ 62 for (k1=0; k1<4; k1++) *ap++ = dk[k + 4*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*16; 77 for (j = jmin; j < jmax; j++){ 78 vj = perm_ptr[aj[j]]; /* block col. index */ 79 rtmp_ptr = rtmp + vj*16; 80 for (i=0; i<16; 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*16,16*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*16; 96 u = ba + ili*16; 97 98 uik[0] = -(diag[0]*u[0] + diag[4]*u[1] + diag[8]*u[2] + diag[12]*u[3]); 99 uik[1] = -(diag[1]*u[0] + diag[5]*u[1] + diag[9]*u[2] + diag[13]*u[3]); 100 uik[2] = -(diag[2]*u[0] + diag[6]*u[1] + diag[10]*u[2]+ diag[14]*u[3]); 101 uik[3] = -(diag[3]*u[0] + diag[7]*u[1] + diag[11]*u[2]+ diag[15]*u[3]); 102 103 uik[4] = -(diag[0]*u[4] + diag[4]*u[5] + diag[8]*u[6] + diag[12]*u[7]); 104 uik[5] = -(diag[1]*u[4] + diag[5]*u[5] + diag[9]*u[6] + diag[13]*u[7]); 105 uik[6] = -(diag[2]*u[4] + diag[6]*u[5] + diag[10]*u[6]+ diag[14]*u[7]); 106 uik[7] = -(diag[3]*u[4] + diag[7]*u[5] + diag[11]*u[6]+ diag[15]*u[7]); 107 108 uik[8] = -(diag[0]*u[8] + diag[4]*u[9] + diag[8]*u[10] + diag[12]*u[11]); 109 uik[9] = -(diag[1]*u[8] + diag[5]*u[9] + diag[9]*u[10] + diag[13]*u[11]); 110 uik[10]= -(diag[2]*u[8] + diag[6]*u[9] + diag[10]*u[10]+ diag[14]*u[11]); 111 uik[11]= -(diag[3]*u[8] + diag[7]*u[9] + diag[11]*u[10]+ diag[15]*u[11]); 112 113 uik[12]= -(diag[0]*u[12] + diag[4]*u[13] + diag[8]*u[14] + diag[12]*u[15]); 114 uik[13]= -(diag[1]*u[12] + diag[5]*u[13] + diag[9]*u[14] + diag[13]*u[15]); 115 uik[14]= -(diag[2]*u[12] + diag[6]*u[13] + diag[10]*u[14]+ diag[14]*u[15]); 116 uik[15]= -(diag[3]*u[12] + diag[7]*u[13] + diag[11]*u[14]+ diag[15]*u[15]); 117 118 /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 119 dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2] + uik[3]*u[3]; 120 dk[1] += uik[4]*u[0] + uik[5]*u[1] + uik[6]*u[2] + uik[7]*u[3]; 121 dk[2] += uik[8]*u[0] + uik[9]*u[1] + uik[10]*u[2]+ uik[11]*u[3]; 122 dk[3] += uik[12]*u[0]+ uik[13]*u[1]+ uik[14]*u[2]+ uik[15]*u[3]; 123 124 dk[4] += uik[0]*u[4] + uik[1]*u[5] + uik[2]*u[6] + uik[3]*u[7]; 125 dk[5] += uik[4]*u[4] + uik[5]*u[5] + uik[6]*u[6] + uik[7]*u[7]; 126 dk[6] += uik[8]*u[4] + uik[9]*u[5] + uik[10]*u[6]+ uik[11]*u[7]; 127 dk[7] += uik[12]*u[4]+ uik[13]*u[5]+ uik[14]*u[6]+ uik[15]*u[7]; 128 129 dk[8] += uik[0]*u[8] + uik[1]*u[9] + uik[2]*u[10] + uik[3]*u[11]; 130 dk[9] += uik[4]*u[8] + uik[5]*u[9] + uik[6]*u[10] + uik[7]*u[11]; 131 dk[10]+= uik[8]*u[8] + uik[9]*u[9] + uik[10]*u[10]+ uik[11]*u[11]; 132 dk[11]+= uik[12]*u[8]+ uik[13]*u[9]+ uik[14]*u[10]+ uik[15]*u[11]; 133 134 dk[12]+= uik[0]*u[12] + uik[1]*u[13] + uik[2]*u[14] + uik[3]*u[15]; 135 dk[13]+= uik[4]*u[12] + uik[5]*u[13] + uik[6]*u[14] + uik[7]*u[15]; 136 dk[14]+= uik[8]*u[12] + uik[9]*u[13] + uik[10]*u[14]+ uik[11]*u[15]; 137 dk[15]+= uik[12]*u[12]+ uik[13]*u[13]+ uik[14]*u[14]+ uik[15]*u[15]; 138 139 ierr = PetscLogFlops(64*4);CHKERRQ(ierr); 140 141 /* update -U(i,k) */ 142 ierr = PetscMemcpy(ba+ili*16,uik,16*sizeof(MatScalar));CHKERRQ(ierr); 143 144 /* add multiple of row i to k-th row ... */ 145 jmin = ili + 1; jmax = bi[i+1]; 146 if (jmin < jmax){ 147 for (j=jmin; j<jmax; j++) { 148 /* rtmp += -U(i,k)^T * U_bar(i,j) */ 149 rtmp_ptr = rtmp + bj[j]*16; 150 u = ba + j*16; 151 rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2] + uik[3]*u[3]; 152 rtmp_ptr[1] += uik[4]*u[0] + uik[5]*u[1] + uik[6]*u[2] + uik[7]*u[3]; 153 rtmp_ptr[2] += uik[8]*u[0] + uik[9]*u[1] + uik[10]*u[2]+ uik[11]*u[3]; 154 rtmp_ptr[3] += uik[12]*u[0]+ uik[13]*u[1]+ uik[14]*u[2]+ uik[15]*u[3]; 155 156 rtmp_ptr[4] += uik[0]*u[4] + uik[1]*u[5] + uik[2]*u[6] + uik[3]*u[7]; 157 rtmp_ptr[5] += uik[4]*u[4] + uik[5]*u[5] + uik[6]*u[6] + uik[7]*u[7]; 158 rtmp_ptr[6] += uik[8]*u[4] + uik[9]*u[5] + uik[10]*u[6]+ uik[11]*u[7]; 159 rtmp_ptr[7] += uik[12]*u[4]+ uik[13]*u[5]+ uik[14]*u[6]+ uik[15]*u[7]; 160 161 rtmp_ptr[8] += uik[0]*u[8] + uik[1]*u[9] + uik[2]*u[10] + uik[3]*u[11]; 162 rtmp_ptr[9] += uik[4]*u[8] + uik[5]*u[9] + uik[6]*u[10] + uik[7]*u[11]; 163 rtmp_ptr[10]+= uik[8]*u[8] + uik[9]*u[9] + uik[10]*u[10]+ uik[11]*u[11]; 164 rtmp_ptr[11]+= uik[12]*u[8]+ uik[13]*u[9]+ uik[14]*u[10]+ uik[15]*u[11]; 165 166 rtmp_ptr[12]+= uik[0]*u[12] + uik[1]*u[13] + uik[2]*u[14] + uik[3]*u[15]; 167 rtmp_ptr[13]+= uik[4]*u[12] + uik[5]*u[13] + uik[6]*u[14] + uik[7]*u[15]; 168 rtmp_ptr[14]+= uik[8]*u[12] + uik[9]*u[13] + uik[10]*u[14]+ uik[11]*u[15]; 169 rtmp_ptr[15]+= uik[12]*u[12]+ uik[13]*u[13]+ uik[14]*u[14]+ uik[15]*u[15]; 170 } 171 ierr = PetscLogFlops(2*64*(jmax-jmin));CHKERRQ(ierr); 172 173 /* ... add i to row list for next nonzero entry */ 174 il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 175 j = bj[jmin]; 176 jl[i] = jl[j]; jl[j] = i; /* update jl */ 177 } 178 i = nexti; 179 } 180 181 /* save nonzero entries in k-th row of U ... */ 182 183 /* invert diagonal block */ 184 diag = ba+k*16; 185 ierr = PetscMemcpy(diag,dk,16*sizeof(MatScalar));CHKERRQ(ierr); 186 187 if (pivotinblocks) { 188 ierr = Kernel_A_gets_inverse_A_4(diag,shift);CHKERRQ(ierr); 189 } else { 190 ierr = Kernel_A_gets_inverse_A_4_nopivot(diag);CHKERRQ(ierr); 191 } 192 193 jmin = bi[k]; jmax = bi[k+1]; 194 if (jmin < jmax) { 195 for (j=jmin; j<jmax; j++){ 196 vj = bj[j]; /* block col. index of U */ 197 u = ba + j*16; 198 rtmp_ptr = rtmp + vj*16; 199 for (k1=0; k1<16; k1++){ 200 *u++ = *rtmp_ptr; 201 *rtmp_ptr++ = 0.0; 202 } 203 } 204 205 /* ... add k to row list for first nonzero entry in k-th row */ 206 il[k] = jmin; 207 i = bj[jmin]; 208 jl[k] = jl[i]; jl[i] = k; 209 } 210 } 211 212 ierr = PetscFree(rtmp);CHKERRQ(ierr); 213 ierr = PetscFree(il);CHKERRQ(ierr); 214 ierr = PetscFree(dk);CHKERRQ(ierr); 215 if (a->permute){ 216 ierr = PetscFree(aa);CHKERRQ(ierr); 217 } 218 219 ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr); 220 C->ops->solve = MatSolve_SeqSBAIJ_4; 221 C->ops->solvetranspose = MatSolve_SeqSBAIJ_4; 222 C->assembled = PETSC_TRUE; 223 C->preallocated = PETSC_TRUE; 224 ierr = PetscLogFlops(1.3333*64*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 225 PetscFunctionReturn(0); 226 } 227