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