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