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