1 2 #include <../src/mat/impls/sbaij/seq/sbaij.h> 3 #include <petsc/private/kernels/blockinvert.h> 4 5 /* 6 Version for when blocks are 3 by 3 Using natural ordering 7 */ 8 #undef __FUNCT__ 9 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering" 10 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info) 11 { 12 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data; 13 PetscErrorCode ierr; 14 PetscInt i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 15 PetscInt *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 PetscReal shift = info->shiftamount; 19 20 PetscFunctionBegin; 21 /* initialization */ 22 ierr = PetscCalloc1(9*mbs,&rtmp);CHKERRQ(ierr); 23 ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr); 24 for (i=0; i<mbs; i++) { 25 jl[i] = mbs; il[0] = 0; 26 } 27 ierr = PetscMalloc2(9,&dk,9,&uik);CHKERRQ(ierr); 28 ai = a->i; aj = a->j; aa = a->a; 29 30 /* for each row k */ 31 for (k = 0; k<mbs; k++) { 32 33 /*initialize k-th row with elements nonzero in row k of A */ 34 jmin = ai[k]; jmax = ai[k+1]; 35 if (jmin < jmax) { 36 ap = aa + jmin*9; 37 for (j = jmin; j < jmax; j++) { 38 vj = aj[j]; /* block col. index */ 39 rtmp_ptr = rtmp + vj*9; 40 for (i=0; i<9; i++) *rtmp_ptr++ = *ap++; 41 } 42 } 43 44 /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 45 ierr = PetscMemcpy(dk,rtmp+k*9,9*sizeof(MatScalar));CHKERRQ(ierr); 46 i = jl[k]; /* first row to be added to k_th row */ 47 48 while (i < mbs) { 49 nexti = jl[i]; /* next row to be added to k_th row */ 50 51 /* compute multiplier */ 52 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 53 54 /* uik = -inv(Di)*U_bar(i,k) */ 55 diag = ba + i*9; 56 u = ba + ili*9; 57 58 uik[0] = -(diag[0]*u[0] + diag[3]*u[1] + diag[6]*u[2]); 59 uik[1] = -(diag[1]*u[0] + diag[4]*u[1] + diag[7]*u[2]); 60 uik[2] = -(diag[2]*u[0] + diag[5]*u[1] + diag[8]*u[2]); 61 62 uik[3] = -(diag[0]*u[3] + diag[3]*u[4] + diag[6]*u[5]); 63 uik[4] = -(diag[1]*u[3] + diag[4]*u[4] + diag[7]*u[5]); 64 uik[5] = -(diag[2]*u[3] + diag[5]*u[4] + diag[8]*u[5]); 65 66 uik[6] = -(diag[0]*u[6] + diag[3]*u[7] + diag[6]*u[8]); 67 uik[7] = -(diag[1]*u[6] + diag[4]*u[7] + diag[7]*u[8]); 68 uik[8] = -(diag[2]*u[6] + diag[5]*u[7] + diag[8]*u[8]); 69 70 /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 71 dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2]; 72 dk[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2]; 73 dk[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2]; 74 75 dk[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5]; 76 dk[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5]; 77 dk[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5]; 78 79 dk[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8]; 80 dk[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8]; 81 dk[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8]; 82 83 ierr = PetscLogFlops(27.0*4.0);CHKERRQ(ierr); 84 85 /* update -U(i,k) */ 86 ierr = PetscMemcpy(ba+ili*9,uik,9*sizeof(MatScalar));CHKERRQ(ierr); 87 88 /* add multiple of row i to k-th row ... */ 89 jmin = ili + 1; jmax = bi[i+1]; 90 if (jmin < jmax) { 91 for (j=jmin; j<jmax; j++) { 92 /* rtmp += -U(i,k)^T * U_bar(i,j) */ 93 rtmp_ptr = rtmp + bj[j]*9; 94 u = ba + j*9; 95 rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2]; 96 rtmp_ptr[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2]; 97 rtmp_ptr[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2]; 98 99 rtmp_ptr[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5]; 100 rtmp_ptr[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5]; 101 rtmp_ptr[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5]; 102 103 rtmp_ptr[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8]; 104 rtmp_ptr[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8]; 105 rtmp_ptr[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8]; 106 } 107 ierr = PetscLogFlops(2.0*27.0*(jmax-jmin));CHKERRQ(ierr); 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 = PetscKernel_A_gets_inverse_A_3(diag,shift);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 = PetscFree2(il,jl);CHKERRQ(ierr); 145 ierr = PetscFree2(dk,uik);CHKERRQ(ierr); 146 147 C->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering_inplace; 148 C->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering_inplace; 149 C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace; 150 C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace; 151 152 C->assembled = PETSC_TRUE; 153 C->preallocated = PETSC_TRUE; 154 155 ierr = PetscLogFlops(1.3333*27*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 156 PetscFunctionReturn(0); 157 } 158