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