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