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