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