1 2 #include <../src/mat/impls/sbaij/seq/sbaij.h> 3 #include <petsc/private/kernels/blockinvert.h> 4 5 /* Version for when blocks are 3 by 3 */ 6 #undef __FUNCT__ 7 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_3" 8 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat C,Mat A,const MatFactorInfo *info) 9 { 10 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data; 11 IS perm = b->row; 12 PetscErrorCode ierr; 13 const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j; 14 PetscInt *a2anew,i,j,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 PetscReal shift = info->shiftamount; 18 19 PetscFunctionBegin; 20 /* initialization */ 21 ierr = PetscCalloc1(9*mbs,&rtmp);CHKERRQ(ierr); 22 ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr); 23 for (i=0; i<mbs; i++) { 24 jl[i] = mbs; il[0] = 0; 25 } 26 ierr = PetscMalloc2(9,&dk,9,&uik);CHKERRQ(ierr); 27 ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr); 28 29 /* check permutation */ 30 if (!a->permute) { 31 ai = a->i; aj = a->j; aa = a->a; 32 } else { 33 ai = a->inew; aj = a->jnew; 34 ierr = PetscMalloc1(9*ai[mbs],&aa);CHKERRQ(ierr); 35 ierr = PetscMemcpy(aa,a->a,9*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); 36 ierr = PetscMalloc1(ai[mbs],&a2anew);CHKERRQ(ierr); 37 ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr); 38 39 for (i=0; i<mbs; i++) { 40 jmin = ai[i]; jmax = ai[i+1]; 41 for (j=jmin; j<jmax; j++) { 42 while (a2anew[j] != j) { 43 k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k; 44 for (k1=0; k1<9; k1++) { 45 dk[k1] = aa[k*9+k1]; 46 aa[k*9+k1] = aa[j*9+k1]; 47 aa[j*9+k1] = dk[k1]; 48 } 49 } 50 /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */ 51 if (i > aj[j]) { 52 /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */ 53 ap = aa + j*9; /* ptr to the beginning of j-th block of aa */ 54 for (k=0; k<9; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */ 55 for (k=0; k<3; k++) { /* j-th block of aa <- dk^T */ 56 for (k1=0; k1<3; k1++) *ap++ = dk[k + 3*k1]; 57 } 58 } 59 } 60 } 61 ierr = PetscFree(a2anew);CHKERRQ(ierr); 62 } 63 64 /* for each row k */ 65 for (k = 0; k<mbs; k++) { 66 67 /*initialize k-th row with elements nonzero in row perm(k) of A */ 68 jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1]; 69 if (jmin < jmax) { 70 ap = aa + jmin*9; 71 for (j = jmin; j < jmax; j++) { 72 vj = perm_ptr[aj[j]]; /* block col. index */ 73 rtmp_ptr = rtmp + vj*9; 74 for (i=0; i<9; i++) *rtmp_ptr++ = *ap++; 75 } 76 } 77 78 /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 79 ierr = PetscMemcpy(dk,rtmp+k*9,9*sizeof(MatScalar));CHKERRQ(ierr); 80 i = jl[k]; /* first row to be added to k_th row */ 81 82 while (i < mbs) { 83 nexti = jl[i]; /* next row to be added to k_th row */ 84 85 /* compute multiplier */ 86 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 87 88 /* uik = -inv(Di)*U_bar(i,k) */ 89 diag = ba + i*9; 90 u = ba + ili*9; 91 92 uik[0] = -(diag[0]*u[0] + diag[3]*u[1] + diag[6]*u[2]); 93 uik[1] = -(diag[1]*u[0] + diag[4]*u[1] + diag[7]*u[2]); 94 uik[2] = -(diag[2]*u[0] + diag[5]*u[1] + diag[8]*u[2]); 95 96 uik[3] = -(diag[0]*u[3] + diag[3]*u[4] + diag[6]*u[5]); 97 uik[4] = -(diag[1]*u[3] + diag[4]*u[4] + diag[7]*u[5]); 98 uik[5] = -(diag[2]*u[3] + diag[5]*u[4] + diag[8]*u[5]); 99 100 uik[6] = -(diag[0]*u[6] + diag[3]*u[7] + diag[6]*u[8]); 101 uik[7] = -(diag[1]*u[6] + diag[4]*u[7] + diag[7]*u[8]); 102 uik[8] = -(diag[2]*u[6] + diag[5]*u[7] + diag[8]*u[8]); 103 104 /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 105 dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2]; 106 dk[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2]; 107 dk[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2]; 108 109 dk[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5]; 110 dk[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5]; 111 dk[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5]; 112 113 dk[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8]; 114 dk[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8]; 115 dk[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8]; 116 117 ierr = PetscLogFlops(27.0*4.0);CHKERRQ(ierr); 118 119 /* update -U(i,k) */ 120 ierr = PetscMemcpy(ba+ili*9,uik,9*sizeof(MatScalar));CHKERRQ(ierr); 121 122 /* add multiple of row i to k-th row ... */ 123 jmin = ili + 1; jmax = bi[i+1]; 124 if (jmin < jmax) { 125 for (j=jmin; j<jmax; j++) { 126 /* rtmp += -U(i,k)^T * U_bar(i,j) */ 127 rtmp_ptr = rtmp + bj[j]*9; 128 u = ba + j*9; 129 rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2]; 130 rtmp_ptr[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2]; 131 rtmp_ptr[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2]; 132 133 rtmp_ptr[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5]; 134 rtmp_ptr[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5]; 135 rtmp_ptr[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5]; 136 137 rtmp_ptr[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8]; 138 rtmp_ptr[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8]; 139 rtmp_ptr[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8]; 140 } 141 ierr = PetscLogFlops(2.0*27.0*(jmax-jmin));CHKERRQ(ierr); 142 143 /* ... add i to row list for next nonzero entry */ 144 il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 145 j = bj[jmin]; 146 jl[i] = jl[j]; jl[j] = i; /* update jl */ 147 } 148 i = nexti; 149 } 150 151 /* save nonzero entries in k-th row of U ... */ 152 153 /* invert diagonal block */ 154 diag = ba+k*9; 155 ierr = PetscMemcpy(diag,dk,9*sizeof(MatScalar));CHKERRQ(ierr); 156 ierr = PetscKernel_A_gets_inverse_A_3(diag,shift);CHKERRQ(ierr); 157 158 jmin = bi[k]; jmax = bi[k+1]; 159 if (jmin < jmax) { 160 for (j=jmin; j<jmax; j++) { 161 vj = bj[j]; /* block col. index of U */ 162 u = ba + j*9; 163 rtmp_ptr = rtmp + vj*9; 164 for (k1=0; k1<9; k1++) { 165 *u++ = *rtmp_ptr; 166 *rtmp_ptr++ = 0.0; 167 } 168 } 169 170 /* ... add k to row list for first nonzero entry in k-th row */ 171 il[k] = jmin; 172 i = bj[jmin]; 173 jl[k] = jl[i]; jl[i] = k; 174 } 175 } 176 177 ierr = PetscFree(rtmp);CHKERRQ(ierr); 178 ierr = PetscFree2(il,jl);CHKERRQ(ierr); 179 ierr = PetscFree2(dk,uik);CHKERRQ(ierr); 180 if (a->permute) { 181 ierr = PetscFree(aa);CHKERRQ(ierr); 182 } 183 184 ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr); 185 186 C->ops->solve = MatSolve_SeqSBAIJ_3_inplace; 187 C->ops->solvetranspose = MatSolve_SeqSBAIJ_3_inplace; 188 C->assembled = PETSC_TRUE; 189 C->preallocated = PETSC_TRUE; 190 191 ierr = PetscLogFlops(1.3333*27*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 192 PetscFunctionReturn(0); 193 } 194