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