1 #define PETSCMAT_DLL 2 3 /* 4 Factorization code for BAIJ format. 5 */ 6 #include "../src/mat/impls/baij/seq/baij.h" 7 #include "../src/mat/blockinvert.h" 8 9 /* 10 Version for when blocks are 3 by 3 11 */ 12 #undef __FUNCT__ 13 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_inplace" 14 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_inplace(Mat C,Mat A,const MatFactorInfo *info) 15 { 16 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 17 IS isrow = b->row,isicol = b->icol; 18 PetscErrorCode ierr; 19 const PetscInt *r,*ic; 20 PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 21 PetscInt *ajtmpold,*ajtmp,nz,row,*ai=a->i,*aj=a->j; 22 PetscInt *diag_offset = b->diag,idx,*pj; 23 MatScalar *pv,*v,*rtmp,*pc,*w,*x; 24 MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4; 25 MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9; 26 MatScalar *ba = b->a,*aa = a->a; 27 PetscReal shift = info->shiftamount; 28 29 PetscFunctionBegin; 30 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 31 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 32 ierr = PetscMalloc(9*(n+1)*sizeof(MatScalar),&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 ierr = Kernel_A_gets_inverse_A_3(w,shift);CHKERRQ(ierr); 110 } 111 112 ierr = PetscFree(rtmp);CHKERRQ(ierr); 113 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 114 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 115 C->ops->solve = MatSolve_SeqBAIJ_3_inplace; 116 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_inplace; 117 C->assembled = PETSC_TRUE; 118 ierr = PetscLogFlops(1.333333333333*3*3*3*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 119 PetscFunctionReturn(0); 120 } 121 122 /* MatLUFactorNumeric_SeqBAIJ_3 - 123 copied from MatLUFactorNumeric_SeqBAIJ_N_inplace() and manually re-implemented 124 Kernel_A_gets_A_times_B() 125 Kernel_A_gets_A_minus_B_times_C() 126 Kernel_A_gets_inverse_A() 127 */ 128 #undef __FUNCT__ 129 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3" 130 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3(Mat B,Mat A,const MatFactorInfo *info) 131 { 132 Mat C=B; 133 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data; 134 IS isrow = b->row,isicol = b->icol; 135 PetscErrorCode ierr; 136 const PetscInt *r,*ic,*ics; 137 PetscInt i,j,k,nz,nzL,row; 138 const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 139 const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2; 140 MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a; 141 PetscInt flg; 142 PetscReal shift = info->shiftamount; 143 144 PetscFunctionBegin; 145 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 146 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 147 148 /* generate work space needed by the factorization */ 149 ierr = PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);CHKERRQ(ierr); 150 ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 151 ics = ic; 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++) { if (pc[j]!=0.0) { flg = 1; break; }} 184 if (flg) { 185 pv = b->a + bs2*bdiag[row]; 186 /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 187 ierr = Kernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr); 188 189 pj = b->j + bdiag[row+1] + 1; /* beginning of U(row,:) */ 190 pv = b->a + bs2*(bdiag[row+1]+1); 191 nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */ 192 for (j=0; j<nz; j++) { 193 /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 194 /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 195 v = rtmp + bs2*pj[j]; 196 ierr = Kernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr); 197 pv += bs2; 198 } 199 ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 200 } 201 } 202 203 /* finished row so stick it into b->a */ 204 /* L part */ 205 pv = b->a + bs2*bi[i] ; 206 pj = b->j + bi[i] ; 207 nz = bi[i+1] - bi[i]; 208 for (j=0; j<nz; j++) { 209 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 210 } 211 212 /* Mark diagonal and invert diagonal for simplier triangular solves */ 213 pv = b->a + bs2*bdiag[i]; 214 pj = b->j + bdiag[i]; 215 ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 216 /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ 217 ierr = Kernel_A_gets_inverse_A_3(pv,shift);CHKERRQ(ierr); 218 219 /* U part */ 220 pj = b->j + bdiag[i+1] + 1; 221 pv = b->a + bs2*(bdiag[i+1]+1); 222 nz = bdiag[i] - bdiag[i+1] - 1; 223 for (j=0; j<nz; j++){ 224 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 225 } 226 } 227 228 ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr); 229 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 230 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 231 C->ops->solve = MatSolve_SeqBAIJ_3; 232 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3; 233 234 C->assembled = PETSC_TRUE; 235 ierr = PetscLogFlops(1.333333333333*3*3*3*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 236 PetscFunctionReturn(0); 237 } 238 239 #undef __FUNCT__ 240 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace" 241 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info) 242 { 243 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 244 PetscErrorCode ierr; 245 PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 246 PetscInt *ajtmpold,*ajtmp,nz,row; 247 PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; 248 MatScalar *pv,*v,*rtmp,*pc,*w,*x; 249 MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4; 250 MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9; 251 MatScalar *ba = b->a,*aa = a->a; 252 PetscReal shift = info->shiftamount; 253 254 PetscFunctionBegin; 255 ierr = PetscMalloc(9*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 256 257 for (i=0; i<n; i++) { 258 nz = bi[i+1] - bi[i]; 259 ajtmp = bj + bi[i]; 260 for (j=0; j<nz; j++) { 261 x = rtmp+9*ajtmp[j]; 262 x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = 0.0; 263 } 264 /* load in initial (unfactored row) */ 265 nz = ai[i+1] - ai[i]; 266 ajtmpold = aj + ai[i]; 267 v = aa + 9*ai[i]; 268 for (j=0; j<nz; j++) { 269 x = rtmp+9*ajtmpold[j]; 270 x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 271 x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8]; 272 v += 9; 273 } 274 row = *ajtmp++; 275 while (row < i) { 276 pc = rtmp + 9*row; 277 p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 278 p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8]; 279 if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || 280 p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0) { 281 pv = ba + 9*diag_offset[row]; 282 pj = bj + diag_offset[row] + 1; 283 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 284 x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 285 pc[0] = m1 = p1*x1 + p4*x2 + p7*x3; 286 pc[1] = m2 = p2*x1 + p5*x2 + p8*x3; 287 pc[2] = m3 = p3*x1 + p6*x2 + p9*x3; 288 289 pc[3] = m4 = p1*x4 + p4*x5 + p7*x6; 290 pc[4] = m5 = p2*x4 + p5*x5 + p8*x6; 291 pc[5] = m6 = p3*x4 + p6*x5 + p9*x6; 292 293 pc[6] = m7 = p1*x7 + p4*x8 + p7*x9; 294 pc[7] = m8 = p2*x7 + p5*x8 + p8*x9; 295 pc[8] = m9 = p3*x7 + p6*x8 + p9*x9; 296 297 nz = bi[row+1] - diag_offset[row] - 1; 298 pv += 9; 299 for (j=0; j<nz; j++) { 300 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 301 x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 302 x = rtmp + 9*pj[j]; 303 x[0] -= m1*x1 + m4*x2 + m7*x3; 304 x[1] -= m2*x1 + m5*x2 + m8*x3; 305 x[2] -= m3*x1 + m6*x2 + m9*x3; 306 307 x[3] -= m1*x4 + m4*x5 + m7*x6; 308 x[4] -= m2*x4 + m5*x5 + m8*x6; 309 x[5] -= m3*x4 + m6*x5 + m9*x6; 310 311 x[6] -= m1*x7 + m4*x8 + m7*x9; 312 x[7] -= m2*x7 + m5*x8 + m8*x9; 313 x[8] -= m3*x7 + m6*x8 + m9*x9; 314 pv += 9; 315 } 316 ierr = PetscLogFlops(54.0*nz+36.0);CHKERRQ(ierr); 317 } 318 row = *ajtmp++; 319 } 320 /* finished row so stick it into b->a */ 321 pv = ba + 9*bi[i]; 322 pj = bj + bi[i]; 323 nz = bi[i+1] - bi[i]; 324 for (j=0; j<nz; j++) { 325 x = rtmp+9*pj[j]; 326 pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 327 pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8]; 328 pv += 9; 329 } 330 /* invert diagonal block */ 331 w = ba + 9*diag_offset[i]; 332 ierr = Kernel_A_gets_inverse_A_3(w,shift);CHKERRQ(ierr); 333 } 334 335 ierr = PetscFree(rtmp);CHKERRQ(ierr); 336 C->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering_inplace; 337 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering_inplace; 338 C->assembled = PETSC_TRUE; 339 ierr = PetscLogFlops(1.333333333333*3*3*3*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 340 PetscFunctionReturn(0); 341 } 342 343 /* 344 MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering - 345 copied from MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace() 346 */ 347 #undef __FUNCT__ 348 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering" 349 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info) 350 { 351 Mat C=B; 352 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data; 353 PetscErrorCode ierr; 354 PetscInt i,j,k,nz,nzL,row; 355 const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 356 const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2; 357 MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a; 358 PetscInt flg; 359 PetscReal shift = info->shiftamount; 360 361 PetscFunctionBegin; 362 /* generate work space needed by the factorization */ 363 ierr = PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);CHKERRQ(ierr); 364 ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 365 366 for (i=0; i<n; i++){ 367 /* zero rtmp */ 368 /* L part */ 369 nz = bi[i+1] - bi[i]; 370 bjtmp = bj + bi[i]; 371 for (j=0; j<nz; j++){ 372 ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 373 } 374 375 /* U part */ 376 nz = bdiag[i] - bdiag[i+1]; 377 bjtmp = bj + bdiag[i+1] + 1; 378 for (j=0; j<nz; j++){ 379 ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 380 } 381 382 /* load in initial (unfactored row) */ 383 nz = ai[i+1] - ai[i]; 384 ajtmp = aj + ai[i]; 385 v = aa + bs2*ai[i]; 386 for (j=0; j<nz; j++) { 387 ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 388 } 389 390 /* elimination */ 391 bjtmp = bj + bi[i]; 392 nzL = bi[i+1] - bi[i]; 393 for(k=0;k<nzL;k++){ 394 row = bjtmp[k]; 395 pc = rtmp + bs2*row; 396 for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }} 397 if (flg) { 398 pv = b->a + bs2*bdiag[row]; 399 /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 400 ierr = Kernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr); 401 402 pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */ 403 pv = b->a + bs2*(bdiag[row+1]+1); 404 nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */ 405 for (j=0; j<nz; j++) { 406 /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 407 /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 408 v = rtmp + bs2*pj[j]; 409 ierr = Kernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr); 410 pv += bs2; 411 } 412 ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 413 } 414 } 415 416 /* finished row so stick it into b->a */ 417 /* L part */ 418 pv = b->a + bs2*bi[i] ; 419 pj = b->j + bi[i] ; 420 nz = bi[i+1] - bi[i]; 421 for (j=0; j<nz; j++) { 422 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 423 } 424 425 /* Mark diagonal and invert diagonal for simplier triangular solves */ 426 pv = b->a + bs2*bdiag[i]; 427 pj = b->j + bdiag[i]; 428 ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 429 /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ 430 ierr = Kernel_A_gets_inverse_A_3(pv,shift);CHKERRQ(ierr); 431 432 /* U part */ 433 pv = b->a + bs2*(bdiag[i+1]+1); 434 pj = b->j + bdiag[i+1]+1; 435 nz = bdiag[i] - bdiag[i+1] - 1; 436 for (j=0; j<nz; j++){ 437 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 438 } 439 } 440 ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr); 441 C->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering; 442 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering; 443 C->assembled = PETSC_TRUE; 444 ierr = PetscLogFlops(1.333333333333*3*3*3*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 445 PetscFunctionReturn(0); 446 } 447 448