1 #include "sbaij.h" 2 #include "src/inline/ilu.h" 3 4 /* 5 Version for when blocks are 3 by 3 Using natural ordering 6 */ 7 #undef __FUNC__ 8 #define __FUNC__ "MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering" 9 int MatCholeskyFactorNumeric_SeqSBAIJ_3_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(9*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 22 ierr = PetscMemzero(rtmp,9*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(18*sizeof(MatScalar),&dk);CHKERRQ(ierr); 29 uik = dk + 9; 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*9; 40 for (j = jmin; j < jmax; j++){ 41 vj = aj[j]; /* block col. index */ 42 rtmp_ptr = rtmp + vj*9; 43 for (i=0; i<9; 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*9,9*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*9; 59 u = ba + ili*9; 60 61 uik[0] = -(diag[0]*u[0] + diag[3]*u[1] + diag[6]*u[2]); 62 uik[1] = -(diag[1]*u[0] + diag[4]*u[1] + diag[7]*u[2]); 63 uik[2] = -(diag[2]*u[0] + diag[5]*u[1] + diag[8]*u[2]); 64 65 uik[3] = -(diag[0]*u[3] + diag[3]*u[4] + diag[6]*u[5]); 66 uik[4] = -(diag[1]*u[3] + diag[4]*u[4] + diag[7]*u[5]); 67 uik[5] = -(diag[2]*u[3] + diag[5]*u[4] + diag[8]*u[5]); 68 69 uik[6] = -(diag[0]*u[6] + diag[3]*u[7] + diag[6]*u[8]); 70 uik[7] = -(diag[1]*u[6] + diag[4]*u[7] + diag[7]*u[8]); 71 uik[8] = -(diag[2]*u[6] + diag[5]*u[7] + diag[8]*u[8]); 72 73 /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 74 dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2]; 75 dk[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2]; 76 dk[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2]; 77 78 dk[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5]; 79 dk[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5]; 80 dk[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5]; 81 82 dk[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8]; 83 dk[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8]; 84 dk[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8]; 85 86 /* update -U(i,k) */ 87 ierr = PetscMemcpy(ba+ili*9,uik,9*sizeof(MatScalar));CHKERRQ(ierr); 88 89 /* add multiple of row i to k-th row ... */ 90 jmin = ili + 1; jmax = bi[i+1]; 91 if (jmin < jmax){ 92 for (j=jmin; j<jmax; j++) { 93 /* rtmp += -U(i,k)^T * U_bar(i,j) */ 94 rtmp_ptr = rtmp + bj[j]*9; 95 u = ba + j*9; 96 rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2]; 97 rtmp_ptr[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2]; 98 rtmp_ptr[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2]; 99 100 rtmp_ptr[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5]; 101 rtmp_ptr[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5]; 102 rtmp_ptr[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5]; 103 104 rtmp_ptr[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8]; 105 rtmp_ptr[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8]; 106 rtmp_ptr[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8]; 107 } 108 109 /* ... add i to row list for next nonzero entry */ 110 il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 111 j = bj[jmin]; 112 jl[i] = jl[j]; jl[j] = i; /* update jl */ 113 } 114 i = nexti; 115 } 116 117 /* save nonzero entries in k-th row of U ... */ 118 119 /* invert diagonal block */ 120 diag = ba+k*9; 121 ierr = PetscMemcpy(diag,dk,9*sizeof(MatScalar));CHKERRQ(ierr); 122 ierr = Kernel_A_gets_inverse_A_3(diag);CHKERRQ(ierr); 123 124 jmin = bi[k]; jmax = bi[k+1]; 125 if (jmin < jmax) { 126 for (j=jmin; j<jmax; j++){ 127 vj = bj[j]; /* block col. index of U */ 128 u = ba + j*9; 129 rtmp_ptr = rtmp + vj*9; 130 for (k1=0; k1<9; k1++){ 131 *u++ = *rtmp_ptr; 132 *rtmp_ptr++ = 0.0; 133 } 134 } 135 136 /* ... add k to row list for first nonzero entry in k-th row */ 137 il[k] = jmin; 138 i = bj[jmin]; 139 jl[k] = jl[i]; jl[i] = k; 140 } 141 } 142 143 ierr = PetscFree(rtmp);CHKERRQ(ierr); 144 ierr = PetscFree(il);CHKERRQ(ierr); 145 ierr = PetscFree(dk);CHKERRQ(ierr); 146 147 C->factor = FACTOR_CHOLESKY; 148 C->assembled = PETSC_TRUE; 149 C->preallocated = PETSC_TRUE; 150 PetscLogFlops(1.3333*27*b->mbs); /* from inverting diagonal blocks */ 151 PetscFunctionReturn(0); 152 } 153