1 #define PETSCMAT_DLL 2 3 #include "../src/mat/impls/sbaij/seq/sbaij.h" 4 #include "../src/mat/blockinvert.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->shiftamount; 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 = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr); 27 for (i=0; i<mbs; i++) { 28 jl[i] = mbs; il[0] = 0; 29 } 30 ierr = PetscMalloc2(9,MatScalar,&dk,9,MatScalar,&uik);CHKERRQ(ierr); 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 ierr = PetscLogFlops(27.0*4.0);CHKERRQ(ierr); 87 88 /* update -U(i,k) */ 89 ierr = PetscMemcpy(ba+ili*9,uik,9*sizeof(MatScalar));CHKERRQ(ierr); 90 91 /* add multiple of row i to k-th row ... */ 92 jmin = ili + 1; jmax = bi[i+1]; 93 if (jmin < jmax){ 94 for (j=jmin; j<jmax; j++) { 95 /* rtmp += -U(i,k)^T * U_bar(i,j) */ 96 rtmp_ptr = rtmp + bj[j]*9; 97 u = ba + j*9; 98 rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2]; 99 rtmp_ptr[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2]; 100 rtmp_ptr[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2]; 101 102 rtmp_ptr[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5]; 103 rtmp_ptr[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5]; 104 rtmp_ptr[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5]; 105 106 rtmp_ptr[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8]; 107 rtmp_ptr[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8]; 108 rtmp_ptr[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8]; 109 } 110 ierr = PetscLogFlops(2.0*27.0*(jmax-jmin));CHKERRQ(ierr); 111 112 /* ... add i to row list for next nonzero entry */ 113 il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 114 j = bj[jmin]; 115 jl[i] = jl[j]; jl[j] = i; /* update jl */ 116 } 117 i = nexti; 118 } 119 120 /* save nonzero entries in k-th row of U ... */ 121 122 /* invert diagonal block */ 123 diag = ba+k*9; 124 ierr = PetscMemcpy(diag,dk,9*sizeof(MatScalar));CHKERRQ(ierr); 125 ierr = Kernel_A_gets_inverse_A_3(diag,shift);CHKERRQ(ierr); 126 127 jmin = bi[k]; jmax = bi[k+1]; 128 if (jmin < jmax) { 129 for (j=jmin; j<jmax; j++){ 130 vj = bj[j]; /* block col. index of U */ 131 u = ba + j*9; 132 rtmp_ptr = rtmp + vj*9; 133 for (k1=0; k1<9; k1++){ 134 *u++ = *rtmp_ptr; 135 *rtmp_ptr++ = 0.0; 136 } 137 } 138 139 /* ... add k to row list for first nonzero entry in k-th row */ 140 il[k] = jmin; 141 i = bj[jmin]; 142 jl[k] = jl[i]; jl[i] = k; 143 } 144 } 145 146 ierr = PetscFree(rtmp);CHKERRQ(ierr); 147 ierr = PetscFree2(il,jl);CHKERRQ(ierr); 148 ierr = PetscFree2(dk,uik);CHKERRQ(ierr); 149 150 C->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering_inplace; 151 C->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering_inplace; 152 C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace; 153 C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace; 154 155 C->assembled = PETSC_TRUE; 156 C->preallocated = PETSC_TRUE; 157 ierr = PetscLogFlops(1.3333*27*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 158 PetscFunctionReturn(0); 159 } 160