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