1 2 /* 3 Factorization code for BAIJ format. 4 */ 5 #include <../src/mat/impls/baij/seq/baij.h> 6 #include <petsc/private/kernels/blockinvert.h> 7 8 /* 9 Version for when blocks are 3 by 3 10 */ 11 #undef __FUNCT__ 12 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_inplace" 13 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_inplace(Mat C,Mat A,const MatFactorInfo *info) 14 { 15 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data; 16 IS isrow = b->row,isicol = b->icol; 17 PetscErrorCode ierr; 18 const PetscInt *r,*ic; 19 PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 20 PetscInt *ajtmpold,*ajtmp,nz,row,*ai=a->i,*aj=a->j; 21 PetscInt *diag_offset = b->diag,idx,*pj; 22 MatScalar *pv,*v,*rtmp,*pc,*w,*x; 23 MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4; 24 MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9; 25 MatScalar *ba = b->a,*aa = a->a; 26 PetscReal shift = info->shiftamount; 27 PetscBool allowzeropivot,zeropivotdetected; 28 29 PetscFunctionBegin; 30 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 31 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 32 ierr = PetscMalloc1(9*(n+1),&rtmp);CHKERRQ(ierr); 33 34 for (i=0; i<n; i++) { 35 nz = bi[i+1] - bi[i]; 36 ajtmp = bj + bi[i]; 37 for (j=0; j<nz; j++) { 38 x = rtmp + 9*ajtmp[j]; 39 x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = 0.0; 40 } 41 /* load in initial (unfactored row) */ 42 idx = r[i]; 43 nz = ai[idx+1] - ai[idx]; 44 ajtmpold = aj + ai[idx]; 45 v = aa + 9*ai[idx]; 46 for (j=0; j<nz; j++) { 47 x = rtmp + 9*ic[ajtmpold[j]]; 48 x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 49 x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8]; 50 v += 9; 51 } 52 row = *ajtmp++; 53 while (row < i) { 54 pc = rtmp + 9*row; 55 p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 56 p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8]; 57 if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || 58 p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0) { 59 pv = ba + 9*diag_offset[row]; 60 pj = bj + diag_offset[row] + 1; 61 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 62 x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 63 pc[0] = m1 = p1*x1 + p4*x2 + p7*x3; 64 pc[1] = m2 = p2*x1 + p5*x2 + p8*x3; 65 pc[2] = m3 = p3*x1 + p6*x2 + p9*x3; 66 67 pc[3] = m4 = p1*x4 + p4*x5 + p7*x6; 68 pc[4] = m5 = p2*x4 + p5*x5 + p8*x6; 69 pc[5] = m6 = p3*x4 + p6*x5 + p9*x6; 70 71 pc[6] = m7 = p1*x7 + p4*x8 + p7*x9; 72 pc[7] = m8 = p2*x7 + p5*x8 + p8*x9; 73 pc[8] = m9 = p3*x7 + p6*x8 + p9*x9; 74 nz = bi[row+1] - diag_offset[row] - 1; 75 pv += 9; 76 for (j=0; j<nz; j++) { 77 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 78 x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 79 x = rtmp + 9*pj[j]; 80 x[0] -= m1*x1 + m4*x2 + m7*x3; 81 x[1] -= m2*x1 + m5*x2 + m8*x3; 82 x[2] -= m3*x1 + m6*x2 + m9*x3; 83 84 x[3] -= m1*x4 + m4*x5 + m7*x6; 85 x[4] -= m2*x4 + m5*x5 + m8*x6; 86 x[5] -= m3*x4 + m6*x5 + m9*x6; 87 88 x[6] -= m1*x7 + m4*x8 + m7*x9; 89 x[7] -= m2*x7 + m5*x8 + m8*x9; 90 x[8] -= m3*x7 + m6*x8 + m9*x9; 91 pv += 9; 92 } 93 ierr = PetscLogFlops(54.0*nz+36.0);CHKERRQ(ierr); 94 } 95 row = *ajtmp++; 96 } 97 /* finished row so stick it into b->a */ 98 pv = ba + 9*bi[i]; 99 pj = bj + bi[i]; 100 nz = bi[i+1] - bi[i]; 101 for (j=0; j<nz; j++) { 102 x = rtmp + 9*pj[j]; 103 pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 104 pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8]; 105 pv += 9; 106 } 107 /* invert diagonal block */ 108 w = ba + 9*diag_offset[i]; 109 allowzeropivot = PetscNot(A->erroriffailure); 110 ierr = PetscKernel_A_gets_inverse_A_3(w,shift,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr); 111 if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 112 } 113 114 ierr = PetscFree(rtmp);CHKERRQ(ierr); 115 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 116 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 117 118 C->ops->solve = MatSolve_SeqBAIJ_3_inplace; 119 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_inplace; 120 C->assembled = PETSC_TRUE; 121 122 ierr = PetscLogFlops(1.333333333333*3*3*3*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 123 PetscFunctionReturn(0); 124 } 125 126 /* MatLUFactorNumeric_SeqBAIJ_3 - 127 copied from MatLUFactorNumeric_SeqBAIJ_N_inplace() and manually re-implemented 128 PetscKernel_A_gets_A_times_B() 129 PetscKernel_A_gets_A_minus_B_times_C() 130 PetscKernel_A_gets_inverse_A() 131 */ 132 #undef __FUNCT__ 133 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3" 134 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3(Mat B,Mat A,const MatFactorInfo *info) 135 { 136 Mat C =B; 137 Mat_SeqBAIJ *a =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data; 138 IS isrow = b->row,isicol = b->icol; 139 PetscErrorCode ierr; 140 const PetscInt *r,*ic; 141 PetscInt i,j,k,nz,nzL,row; 142 const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 143 const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2; 144 MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a; 145 PetscInt flg; 146 PetscReal shift = info->shiftamount; 147 PetscBool allowzeropivot,zeropivotdetected; 148 149 PetscFunctionBegin; 150 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 151 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 152 153 /* generate work space needed by the factorization */ 154 ierr = PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);CHKERRQ(ierr); 155 ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 156 157 for (i=0; i<n; i++) { 158 /* zero rtmp */ 159 /* L part */ 160 nz = bi[i+1] - bi[i]; 161 bjtmp = bj + bi[i]; 162 for (j=0; j<nz; j++) { 163 ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 164 } 165 166 /* U part */ 167 nz = bdiag[i] - bdiag[i+1]; 168 bjtmp = bj + bdiag[i+1]+1; 169 for (j=0; j<nz; j++) { 170 ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 171 } 172 173 /* load in initial (unfactored row) */ 174 nz = ai[r[i]+1] - ai[r[i]]; 175 ajtmp = aj + ai[r[i]]; 176 v = aa + bs2*ai[r[i]]; 177 for (j=0; j<nz; j++) { 178 ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 179 } 180 181 /* elimination */ 182 bjtmp = bj + bi[i]; 183 nzL = bi[i+1] - bi[i]; 184 for (k = 0; k < nzL; k++) { 185 row = bjtmp[k]; 186 pc = rtmp + bs2*row; 187 for (flg=0,j=0; j<bs2; j++) { 188 if (pc[j]!=0.0) { 189 flg = 1; 190 break; 191 } 192 } 193 if (flg) { 194 pv = b->a + bs2*bdiag[row]; 195 /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 196 ierr = PetscKernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr); 197 198 pj = b->j + bdiag[row+1] + 1; /* beginning of U(row,:) */ 199 pv = b->a + bs2*(bdiag[row+1]+1); 200 nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */ 201 for (j=0; j<nz; j++) { 202 /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 203 /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 204 v = rtmp + bs2*pj[j]; 205 ierr = PetscKernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr); 206 pv += bs2; 207 } 208 ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 209 } 210 } 211 212 /* finished row so stick it into b->a */ 213 /* L part */ 214 pv = b->a + bs2*bi[i]; 215 pj = b->j + bi[i]; 216 nz = bi[i+1] - bi[i]; 217 for (j=0; j<nz; j++) { 218 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 219 } 220 221 /* Mark diagonal and invert diagonal for simplier triangular solves */ 222 pv = b->a + bs2*bdiag[i]; 223 pj = b->j + bdiag[i]; 224 ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 225 /* ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ 226 allowzeropivot = PetscNot(A->erroriffailure); 227 ierr = PetscKernel_A_gets_inverse_A_3(pv,shift,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr); 228 if (zeropivotdetected) B->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 229 230 /* U part */ 231 pj = b->j + bdiag[i+1] + 1; 232 pv = b->a + bs2*(bdiag[i+1]+1); 233 nz = bdiag[i] - bdiag[i+1] - 1; 234 for (j=0; j<nz; j++) { 235 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 236 } 237 } 238 239 ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr); 240 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 241 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 242 243 C->ops->solve = MatSolve_SeqBAIJ_3; 244 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3; 245 C->assembled = PETSC_TRUE; 246 247 ierr = PetscLogFlops(1.333333333333*3*3*3*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 248 PetscFunctionReturn(0); 249 } 250 251 #undef __FUNCT__ 252 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace" 253 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info) 254 { 255 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data; 256 PetscErrorCode ierr; 257 PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 258 PetscInt *ajtmpold,*ajtmp,nz,row; 259 PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; 260 MatScalar *pv,*v,*rtmp,*pc,*w,*x; 261 MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4; 262 MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9; 263 MatScalar *ba = b->a,*aa = a->a; 264 PetscReal shift = info->shiftamount; 265 PetscBool allowzeropivot,zeropivotdetected; 266 267 PetscFunctionBegin; 268 ierr = PetscMalloc1(9*(n+1),&rtmp);CHKERRQ(ierr); 269 270 for (i=0; i<n; i++) { 271 nz = bi[i+1] - bi[i]; 272 ajtmp = bj + bi[i]; 273 for (j=0; j<nz; j++) { 274 x = rtmp+9*ajtmp[j]; 275 x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = 0.0; 276 } 277 /* load in initial (unfactored row) */ 278 nz = ai[i+1] - ai[i]; 279 ajtmpold = aj + ai[i]; 280 v = aa + 9*ai[i]; 281 for (j=0; j<nz; j++) { 282 x = rtmp+9*ajtmpold[j]; 283 x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 284 x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8]; 285 v += 9; 286 } 287 row = *ajtmp++; 288 while (row < i) { 289 pc = rtmp + 9*row; 290 p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 291 p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8]; 292 if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || 293 p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0) { 294 pv = ba + 9*diag_offset[row]; 295 pj = bj + diag_offset[row] + 1; 296 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 297 x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 298 pc[0] = m1 = p1*x1 + p4*x2 + p7*x3; 299 pc[1] = m2 = p2*x1 + p5*x2 + p8*x3; 300 pc[2] = m3 = p3*x1 + p6*x2 + p9*x3; 301 302 pc[3] = m4 = p1*x4 + p4*x5 + p7*x6; 303 pc[4] = m5 = p2*x4 + p5*x5 + p8*x6; 304 pc[5] = m6 = p3*x4 + p6*x5 + p9*x6; 305 306 pc[6] = m7 = p1*x7 + p4*x8 + p7*x9; 307 pc[7] = m8 = p2*x7 + p5*x8 + p8*x9; 308 pc[8] = m9 = p3*x7 + p6*x8 + p9*x9; 309 310 nz = bi[row+1] - diag_offset[row] - 1; 311 pv += 9; 312 for (j=0; j<nz; j++) { 313 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 314 x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 315 x = rtmp + 9*pj[j]; 316 x[0] -= m1*x1 + m4*x2 + m7*x3; 317 x[1] -= m2*x1 + m5*x2 + m8*x3; 318 x[2] -= m3*x1 + m6*x2 + m9*x3; 319 320 x[3] -= m1*x4 + m4*x5 + m7*x6; 321 x[4] -= m2*x4 + m5*x5 + m8*x6; 322 x[5] -= m3*x4 + m6*x5 + m9*x6; 323 324 x[6] -= m1*x7 + m4*x8 + m7*x9; 325 x[7] -= m2*x7 + m5*x8 + m8*x9; 326 x[8] -= m3*x7 + m6*x8 + m9*x9; 327 pv += 9; 328 } 329 ierr = PetscLogFlops(54.0*nz+36.0);CHKERRQ(ierr); 330 } 331 row = *ajtmp++; 332 } 333 /* finished row so stick it into b->a */ 334 pv = ba + 9*bi[i]; 335 pj = bj + bi[i]; 336 nz = bi[i+1] - bi[i]; 337 for (j=0; j<nz; j++) { 338 x = rtmp+9*pj[j]; 339 pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 340 pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8]; 341 pv += 9; 342 } 343 /* invert diagonal block */ 344 w = ba + 9*diag_offset[i]; 345 allowzeropivot = PetscNot(A->erroriffailure); 346 ierr = PetscKernel_A_gets_inverse_A_3(w,shift,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr); 347 if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 348 } 349 350 ierr = PetscFree(rtmp);CHKERRQ(ierr); 351 352 C->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering_inplace; 353 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering_inplace; 354 C->assembled = PETSC_TRUE; 355 356 ierr = PetscLogFlops(1.333333333333*3*3*3*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 357 PetscFunctionReturn(0); 358 } 359 360 /* 361 MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering - 362 copied from MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace() 363 */ 364 #undef __FUNCT__ 365 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering" 366 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info) 367 { 368 Mat C =B; 369 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data; 370 PetscErrorCode ierr; 371 PetscInt i,j,k,nz,nzL,row; 372 const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 373 const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2; 374 MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a; 375 PetscInt flg; 376 PetscReal shift = info->shiftamount; 377 PetscBool allowzeropivot,zeropivotdetected; 378 379 PetscFunctionBegin; 380 /* generate work space needed by the factorization */ 381 ierr = PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);CHKERRQ(ierr); 382 ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 383 384 for (i=0; i<n; i++) { 385 /* zero rtmp */ 386 /* L part */ 387 nz = bi[i+1] - bi[i]; 388 bjtmp = bj + bi[i]; 389 for (j=0; j<nz; j++) { 390 ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 391 } 392 393 /* U part */ 394 nz = bdiag[i] - bdiag[i+1]; 395 bjtmp = bj + bdiag[i+1] + 1; 396 for (j=0; j<nz; j++) { 397 ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 398 } 399 400 /* load in initial (unfactored row) */ 401 nz = ai[i+1] - ai[i]; 402 ajtmp = aj + ai[i]; 403 v = aa + bs2*ai[i]; 404 for (j=0; j<nz; j++) { 405 ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 406 } 407 408 /* elimination */ 409 bjtmp = bj + bi[i]; 410 nzL = bi[i+1] - bi[i]; 411 for (k=0; k<nzL; k++) { 412 row = bjtmp[k]; 413 pc = rtmp + bs2*row; 414 for (flg=0,j=0; j<bs2; j++) { 415 if (pc[j]!=0.0) { 416 flg = 1; 417 break; 418 } 419 } 420 if (flg) { 421 pv = b->a + bs2*bdiag[row]; 422 /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 423 ierr = PetscKernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr); 424 425 pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */ 426 pv = b->a + bs2*(bdiag[row+1]+1); 427 nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */ 428 for (j=0; j<nz; j++) { 429 /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 430 /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 431 v = rtmp + bs2*pj[j]; 432 ierr = PetscKernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr); 433 pv += bs2; 434 } 435 ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 436 } 437 } 438 439 /* finished row so stick it into b->a */ 440 /* L part */ 441 pv = b->a + bs2*bi[i]; 442 pj = b->j + bi[i]; 443 nz = bi[i+1] - bi[i]; 444 for (j=0; j<nz; j++) { 445 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 446 } 447 448 /* Mark diagonal and invert diagonal for simplier triangular solves */ 449 pv = b->a + bs2*bdiag[i]; 450 pj = b->j + bdiag[i]; 451 ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 452 /* ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ 453 allowzeropivot = PetscNot(A->erroriffailure); 454 ierr = PetscKernel_A_gets_inverse_A_3(pv,shift,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr); 455 if (zeropivotdetected) B->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 456 457 /* U part */ 458 pv = b->a + bs2*(bdiag[i+1]+1); 459 pj = b->j + bdiag[i+1]+1; 460 nz = bdiag[i] - bdiag[i+1] - 1; 461 for (j=0; j<nz; j++) { 462 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 463 } 464 } 465 ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr); 466 467 C->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering; 468 C->ops->forwardsolve = MatForwardSolve_SeqBAIJ_3_NaturalOrdering; 469 C->ops->backwardsolve = MatBackwardSolve_SeqBAIJ_3_NaturalOrdering; 470 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering; 471 C->assembled = PETSC_TRUE; 472 473 ierr = PetscLogFlops(1.333333333333*3*3*3*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 474 PetscFunctionReturn(0); 475 } 476 477