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