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/inline/ilu.h" 8 9 /* ------------------------------------------------------------*/ 10 /* 11 Version for when blocks are 2 by 2 12 */ 13 #undef __FUNCT__ 14 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2" 15 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat A,MatFactorInfo *info,Mat *B) 16 { 17 Mat C = *B; 18 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 19 IS isrow = b->row,isicol = b->icol; 20 PetscErrorCode ierr; 21 PetscInt *r,*ic,i,j,n = a->mbs,*bi = b->i,*bj = b->j; 22 PetscInt *ajtmpold,*ajtmp,nz,row; 23 PetscInt *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj; 24 MatScalar *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4; 25 MatScalar p1,p2,p3,p4; 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(4*(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+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 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 + 4*ai[idx]; 45 for (j=0; j<nz; j++) { 46 x = rtmp+4*ic[ajtmpold[j]]; 47 x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 48 v += 4; 49 } 50 row = *ajtmp++; 51 while (row < i) { 52 pc = rtmp + 4*row; 53 p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 54 if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) { 55 pv = ba + 4*diag_offset[row]; 56 pj = bj + diag_offset[row] + 1; 57 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 58 pc[0] = m1 = p1*x1 + p3*x2; 59 pc[1] = m2 = p2*x1 + p4*x2; 60 pc[2] = m3 = p1*x3 + p3*x4; 61 pc[3] = m4 = p2*x3 + p4*x4; 62 nz = bi[row+1] - diag_offset[row] - 1; 63 pv += 4; 64 for (j=0; j<nz; j++) { 65 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 66 x = rtmp + 4*pj[j]; 67 x[0] -= m1*x1 + m3*x2; 68 x[1] -= m2*x1 + m4*x2; 69 x[2] -= m1*x3 + m3*x4; 70 x[3] -= m2*x3 + m4*x4; 71 pv += 4; 72 } 73 ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); 74 } 75 row = *ajtmp++; 76 } 77 /* finished row so stick it into b->a */ 78 pv = ba + 4*bi[i]; 79 pj = bj + bi[i]; 80 nz = bi[i+1] - bi[i]; 81 for (j=0; j<nz; j++) { 82 x = rtmp+4*pj[j]; 83 pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 84 pv += 4; 85 } 86 /* invert diagonal block */ 87 w = ba + 4*diag_offset[i]; 88 ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr); 89 } 90 91 ierr = PetscFree(rtmp);CHKERRQ(ierr); 92 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 93 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 94 C->factor = MAT_FACTOR_LU; 95 C->assembled = PETSC_TRUE; 96 ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 97 PetscFunctionReturn(0); 98 } 99 /* 100 Version for when blocks are 2 by 2 Using natural ordering 101 */ 102 #undef __FUNCT__ 103 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering" 104 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B) 105 { 106 Mat C = *B; 107 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 108 PetscErrorCode ierr; 109 PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 110 PetscInt *ajtmpold,*ajtmp,nz,row; 111 PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; 112 MatScalar *pv,*v,*rtmp,*pc,*w,*x; 113 MatScalar p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4; 114 MatScalar *ba = b->a,*aa = a->a; 115 PetscReal shift = info->shiftinblocks; 116 117 PetscFunctionBegin; 118 ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 119 120 for (i=0; i<n; i++) { 121 nz = bi[i+1] - bi[i]; 122 ajtmp = bj + bi[i]; 123 for (j=0; j<nz; j++) { 124 x = rtmp+4*ajtmp[j]; 125 x[0] = x[1] = x[2] = x[3] = 0.0; 126 } 127 /* load in initial (unfactored row) */ 128 nz = ai[i+1] - ai[i]; 129 ajtmpold = aj + ai[i]; 130 v = aa + 4*ai[i]; 131 for (j=0; j<nz; j++) { 132 x = rtmp+4*ajtmpold[j]; 133 x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 134 v += 4; 135 } 136 row = *ajtmp++; 137 while (row < i) { 138 pc = rtmp + 4*row; 139 p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 140 if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) { 141 pv = ba + 4*diag_offset[row]; 142 pj = bj + diag_offset[row] + 1; 143 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 144 pc[0] = m1 = p1*x1 + p3*x2; 145 pc[1] = m2 = p2*x1 + p4*x2; 146 pc[2] = m3 = p1*x3 + p3*x4; 147 pc[3] = m4 = p2*x3 + p4*x4; 148 nz = bi[row+1] - diag_offset[row] - 1; 149 pv += 4; 150 for (j=0; j<nz; j++) { 151 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 152 x = rtmp + 4*pj[j]; 153 x[0] -= m1*x1 + m3*x2; 154 x[1] -= m2*x1 + m4*x2; 155 x[2] -= m1*x3 + m3*x4; 156 x[3] -= m2*x3 + m4*x4; 157 pv += 4; 158 } 159 ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); 160 } 161 row = *ajtmp++; 162 } 163 /* finished row so stick it into b->a */ 164 pv = ba + 4*bi[i]; 165 pj = bj + bi[i]; 166 nz = bi[i+1] - bi[i]; 167 for (j=0; j<nz; j++) { 168 x = rtmp+4*pj[j]; 169 pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 170 pv += 4; 171 } 172 /* invert diagonal block */ 173 w = ba + 4*diag_offset[i]; 174 ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr); 175 } 176 177 ierr = PetscFree(rtmp);CHKERRQ(ierr); 178 C->factor = MAT_FACTOR_LU; 179 C->assembled = PETSC_TRUE; 180 ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 181 PetscFunctionReturn(0); 182 } 183 184 /* ----------------------------------------------------------- */ 185 /* 186 Version for when blocks are 1 by 1. 187 */ 188 #undef __FUNCT__ 189 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1" 190 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat A,MatFactorInfo *info,Mat *B) 191 { 192 Mat C = *B; 193 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 194 IS isrow = b->row,isicol = b->icol; 195 PetscErrorCode ierr; 196 PetscInt *r,*ic,i,j,n = a->mbs,*bi = b->i,*bj = b->j; 197 PetscInt *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j; 198 PetscInt *diag_offset = b->diag,diag,*pj; 199 MatScalar *pv,*v,*rtmp,multiplier,*pc; 200 MatScalar *ba = b->a,*aa = a->a; 201 202 PetscFunctionBegin; 203 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 204 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 205 ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 206 207 for (i=0; i<n; i++) { 208 nz = bi[i+1] - bi[i]; 209 ajtmp = bj + bi[i]; 210 for (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0; 211 212 /* load in initial (unfactored row) */ 213 nz = ai[r[i]+1] - ai[r[i]]; 214 ajtmpold = aj + ai[r[i]]; 215 v = aa + ai[r[i]]; 216 for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] = v[j]; 217 218 row = *ajtmp++; 219 while (row < i) { 220 pc = rtmp + row; 221 if (*pc != 0.0) { 222 pv = ba + diag_offset[row]; 223 pj = bj + diag_offset[row] + 1; 224 multiplier = *pc * *pv++; 225 *pc = multiplier; 226 nz = bi[row+1] - diag_offset[row] - 1; 227 for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j]; 228 ierr = PetscLogFlops(1+2*nz);CHKERRQ(ierr); 229 } 230 row = *ajtmp++; 231 } 232 /* finished row so stick it into b->a */ 233 pv = ba + bi[i]; 234 pj = bj + bi[i]; 235 nz = bi[i+1] - bi[i]; 236 for (j=0; j<nz; j++) {pv[j] = rtmp[pj[j]];} 237 diag = diag_offset[i] - bi[i]; 238 /* check pivot entry for current row */ 239 if (pv[diag] == 0.0) { 240 SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot"); 241 } 242 pv[diag] = 1.0/pv[diag]; 243 } 244 245 ierr = PetscFree(rtmp);CHKERRQ(ierr); 246 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 247 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 248 C->factor = MAT_FACTOR_LU; 249 C->assembled = PETSC_TRUE; 250 ierr = PetscLogFlops(C->cmap.n);CHKERRQ(ierr); 251 PetscFunctionReturn(0); 252 } 253 254 #undef __FUNCT__ 255 #define __FUNCT__ "MatGetFactor_seqbaij_petsc" 256 PetscErrorCode MatGetFactor_seqbaij_petsc(Mat A,MatFactorType ftype,Mat *B) 257 { 258 PetscInt n = A->rmap.n; 259 PetscErrorCode ierr; 260 261 PetscFunctionBegin; 262 ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr); 263 ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr); 264 if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) { 265 ierr = MatSetType(*B,MATSEQBAIJ);CHKERRQ(ierr); 266 } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) { 267 ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr); 268 ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 269 (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqBAIJ; 270 (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqBAIJ; 271 } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported"); 272 PetscFunctionReturn(0); 273 } 274 275 276 /* ----------------------------------------------------------- */ 277 #undef __FUNCT__ 278 #define __FUNCT__ "MatLUFactor_SeqBAIJ" 279 PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,MatFactorInfo *info) 280 { 281 PetscErrorCode ierr; 282 Mat C; 283 284 PetscFunctionBegin; 285 ierr = MatGetFactor(A,"petsc",MAT_FACTOR_LU,&C);CHKERRQ(ierr); 286 ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr); 287 ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr); 288 ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 289 ierr = PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr); 290 PetscFunctionReturn(0); 291 } 292 293 #include "src/mat/impls/sbaij/seq/sbaij.h" 294 #undef __FUNCT__ 295 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N" 296 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat A,MatFactorInfo *info,Mat *B) 297 { 298 PetscErrorCode ierr; 299 Mat C = *B; 300 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 301 Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 302 IS ip=b->row; 303 PetscInt *rip,i,j,mbs=a->mbs,bs=A->rmap.bs,*bi=b->i,*bj=b->j,*bcol; 304 PetscInt *ai=a->i,*aj=a->j; 305 PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 306 MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; 307 PetscReal zeropivot,rs,shiftnz; 308 PetscReal shiftpd; 309 ChShift_Ctx sctx; 310 PetscInt newshift; 311 312 PetscFunctionBegin; 313 if (bs > 1) { 314 if (!a->sbaijMat){ 315 ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 316 } 317 ierr = (a->sbaijMat)->ops->choleskyfactornumeric(a->sbaijMat,info,B);CHKERRQ(ierr); 318 ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr); 319 a->sbaijMat = PETSC_NULL; 320 PetscFunctionReturn(0); 321 } 322 323 /* initialization */ 324 shiftnz = info->shiftnz; 325 shiftpd = info->shiftpd; 326 zeropivot = info->zeropivot; 327 328 ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 329 nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); 330 ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 331 jl = il + mbs; 332 rtmp = (MatScalar*)(jl + mbs); 333 334 sctx.shift_amount = 0; 335 sctx.nshift = 0; 336 do { 337 sctx.chshift = PETSC_FALSE; 338 for (i=0; i<mbs; i++) { 339 rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 340 } 341 342 for (k = 0; k<mbs; k++){ 343 bval = ba + bi[k]; 344 /* initialize k-th row by the perm[k]-th row of A */ 345 jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 346 for (j = jmin; j < jmax; j++){ 347 col = rip[aj[j]]; 348 if (col >= k){ /* only take upper triangular entry */ 349 rtmp[col] = aa[j]; 350 *bval++ = 0.0; /* for in-place factorization */ 351 } 352 } 353 354 /* shift the diagonal of the matrix */ 355 if (sctx.nshift) rtmp[k] += sctx.shift_amount; 356 357 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 358 dk = rtmp[k]; 359 i = jl[k]; /* first row to be added to k_th row */ 360 361 while (i < k){ 362 nexti = jl[i]; /* next row to be added to k_th row */ 363 364 /* compute multiplier, update diag(k) and U(i,k) */ 365 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 366 uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 367 dk += uikdi*ba[ili]; 368 ba[ili] = uikdi; /* -U(i,k) */ 369 370 /* add multiple of row i to k-th row */ 371 jmin = ili + 1; jmax = bi[i+1]; 372 if (jmin < jmax){ 373 for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 374 /* update il and jl for row i */ 375 il[i] = jmin; 376 j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 377 } 378 i = nexti; 379 } 380 381 /* shift the diagonals when zero pivot is detected */ 382 /* compute rs=sum of abs(off-diagonal) */ 383 rs = 0.0; 384 jmin = bi[k]+1; 385 nz = bi[k+1] - jmin; 386 if (nz){ 387 bcol = bj + jmin; 388 while (nz--){ 389 rs += PetscAbsScalar(rtmp[*bcol]); 390 bcol++; 391 } 392 } 393 394 sctx.rs = rs; 395 sctx.pv = dk; 396 ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr); 397 if (newshift == 1) break; 398 399 /* copy data into U(k,:) */ 400 ba[bi[k]] = 1.0/dk; /* U(k,k) */ 401 jmin = bi[k]+1; jmax = bi[k+1]; 402 if (jmin < jmax) { 403 for (j=jmin; j<jmax; j++){ 404 col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 405 } 406 /* add the k-th row into il and jl */ 407 il[k] = jmin; 408 i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 409 } 410 } 411 } while (sctx.chshift); 412 ierr = PetscFree(il);CHKERRQ(ierr); 413 414 ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 415 C->factor = MAT_FACTOR_CHOLESKY; 416 C->assembled = PETSC_TRUE; 417 C->preallocated = PETSC_TRUE; 418 ierr = PetscLogFlops(C->rmap.N);CHKERRQ(ierr); 419 if (sctx.nshift){ 420 if (shiftnz) { 421 ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 422 } else if (shiftpd) { 423 ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 424 } 425 } 426 PetscFunctionReturn(0); 427 } 428 429 #undef __FUNCT__ 430 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering" 431 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *fact) 432 { 433 Mat C = *fact; 434 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 435 Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 436 PetscErrorCode ierr; 437 PetscInt i,j,am=a->mbs; 438 PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 439 PetscInt k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz; 440 MatScalar *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval; 441 PetscReal zeropivot,rs,shiftnz; 442 PetscReal shiftpd; 443 ChShift_Ctx sctx; 444 PetscInt newshift; 445 446 PetscFunctionBegin; 447 /* initialization */ 448 shiftnz = info->shiftnz; 449 shiftpd = info->shiftpd; 450 zeropivot = info->zeropivot; 451 452 nz = (2*am+1)*sizeof(PetscInt)+am*sizeof(MatScalar); 453 ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 454 jl = il + am; 455 rtmp = (MatScalar*)(jl + am); 456 457 sctx.shift_amount = 0; 458 sctx.nshift = 0; 459 do { 460 sctx.chshift = PETSC_FALSE; 461 for (i=0; i<am; i++) { 462 rtmp[i] = 0.0; jl[i] = am; il[0] = 0; 463 } 464 465 for (k = 0; k<am; k++){ 466 /* initialize k-th row with elements nonzero in row perm(k) of A */ 467 nz = ai[k+1] - ai[k]; 468 acol = aj + ai[k]; 469 aval = aa + ai[k]; 470 bval = ba + bi[k]; 471 while (nz -- ){ 472 if (*acol < k) { /* skip lower triangular entries */ 473 acol++; aval++; 474 } else { 475 rtmp[*acol++] = *aval++; 476 *bval++ = 0.0; /* for in-place factorization */ 477 } 478 } 479 480 /* shift the diagonal of the matrix */ 481 if (sctx.nshift) rtmp[k] += sctx.shift_amount; 482 483 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 484 dk = rtmp[k]; 485 i = jl[k]; /* first row to be added to k_th row */ 486 487 while (i < k){ 488 nexti = jl[i]; /* next row to be added to k_th row */ 489 /* compute multiplier, update D(k) and U(i,k) */ 490 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 491 uikdi = - ba[ili]*ba[bi[i]]; 492 dk += uikdi*ba[ili]; 493 ba[ili] = uikdi; /* -U(i,k) */ 494 495 /* add multiple of row i to k-th row ... */ 496 jmin = ili + 1; 497 nz = bi[i+1] - jmin; 498 if (nz > 0){ 499 bcol = bj + jmin; 500 bval = ba + jmin; 501 while (nz --) rtmp[*bcol++] += uikdi*(*bval++); 502 /* update il and jl for i-th row */ 503 il[i] = jmin; 504 j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 505 } 506 i = nexti; 507 } 508 509 /* shift the diagonals when zero pivot is detected */ 510 /* compute rs=sum of abs(off-diagonal) */ 511 rs = 0.0; 512 jmin = bi[k]+1; 513 nz = bi[k+1] - jmin; 514 if (nz){ 515 bcol = bj + jmin; 516 while (nz--){ 517 rs += PetscAbsScalar(rtmp[*bcol]); 518 bcol++; 519 } 520 } 521 522 sctx.rs = rs; 523 sctx.pv = dk; 524 ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr); 525 if (newshift == 1) break; /* sctx.shift_amount is updated */ 526 527 /* copy data into U(k,:) */ 528 ba[bi[k]] = 1.0/dk; 529 jmin = bi[k]+1; 530 nz = bi[k+1] - jmin; 531 if (nz){ 532 bcol = bj + jmin; 533 bval = ba + jmin; 534 while (nz--){ 535 *bval++ = rtmp[*bcol]; 536 rtmp[*bcol++] = 0.0; 537 } 538 /* add k-th row into il and jl */ 539 il[k] = jmin; 540 i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 541 } 542 } 543 } while (sctx.chshift); 544 ierr = PetscFree(il);CHKERRQ(ierr); 545 546 C->factor = MAT_FACTOR_CHOLESKY; 547 C->assembled = PETSC_TRUE; 548 C->preallocated = PETSC_TRUE; 549 ierr = PetscLogFlops(C->rmap.N);CHKERRQ(ierr); 550 if (sctx.nshift){ 551 if (shiftnz) { 552 ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 553 } else if (shiftpd) { 554 ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 555 } 556 } 557 PetscFunctionReturn(0); 558 } 559 560 #include "petscbt.h" 561 #include "src/mat/utils/freespace.h" 562 #undef __FUNCT__ 563 #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ" 564 PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 565 { 566 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 567 Mat_SeqSBAIJ *b; 568 Mat B; 569 PetscErrorCode ierr; 570 PetscTruth perm_identity; 571 PetscInt reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap.bs,*ui; 572 PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow; 573 PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr; 574 PetscReal fill=info->fill,levels=info->levels; 575 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 576 PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 577 PetscBT lnkbt; 578 579 PetscFunctionBegin; 580 if (bs > 1){ 581 if (!a->sbaijMat){ 582 ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 583 } 584 (*fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */ 585 ierr = MatICCFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr); 586 PetscFunctionReturn(0); 587 } 588 589 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 590 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 591 592 /* special case that simply copies fill pattern */ 593 if (!levels && perm_identity) { 594 ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr); 595 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 596 for (i=0; i<am; i++) { 597 ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */ 598 } 599 B = *fact; 600 ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr); 601 602 603 b = (Mat_SeqSBAIJ*)B->data; 604 uj = b->j; 605 for (i=0; i<am; i++) { 606 aj = a->j + a->diag[i]; 607 for (j=0; j<ui[i]; j++){ 608 *uj++ = *aj++; 609 } 610 b->ilen[i] = ui[i]; 611 } 612 ierr = PetscFree(ui);CHKERRQ(ierr); 613 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 614 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 615 616 B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 617 B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 618 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 619 PetscFunctionReturn(0); 620 } 621 622 /* initialization */ 623 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 624 ui[0] = 0; 625 ierr = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr); 626 627 /* jl: linked list for storing indices of the pivot rows 628 il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 629 ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr); 630 il = jl + am; 631 uj_ptr = (PetscInt**)(il + am); 632 uj_lvl_ptr = (PetscInt**)(uj_ptr + am); 633 for (i=0; i<am; i++){ 634 jl[i] = am; il[i] = 0; 635 } 636 637 /* create and initialize a linked list for storing column indices of the active row k */ 638 nlnk = am + 1; 639 ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 640 641 /* initial FreeSpace size is fill*(ai[am]+1) */ 642 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 643 current_space = free_space; 644 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr); 645 current_space_lvl = free_space_lvl; 646 647 for (k=0; k<am; k++){ /* for each active row k */ 648 /* initialize lnk by the column indices of row rip[k] of A */ 649 nzk = 0; 650 ncols = ai[rip[k]+1] - ai[rip[k]]; 651 ncols_upper = 0; 652 cols = cols_lvl + am; 653 for (j=0; j<ncols; j++){ 654 i = rip[*(aj + ai[rip[k]] + j)]; 655 if (i >= k){ /* only take upper triangular entry */ 656 cols[ncols_upper] = i; 657 cols_lvl[ncols_upper] = -1; /* initialize level for nonzero entries */ 658 ncols_upper++; 659 } 660 } 661 ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 662 nzk += nlnk; 663 664 /* update lnk by computing fill-in for each pivot row to be merged in */ 665 prow = jl[k]; /* 1st pivot row */ 666 667 while (prow < k){ 668 nextprow = jl[prow]; 669 670 /* merge prow into k-th row */ 671 jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 672 jmax = ui[prow+1]; 673 ncols = jmax-jmin; 674 i = jmin - ui[prow]; 675 cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 676 for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j); 677 ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 678 nzk += nlnk; 679 680 /* update il and jl for prow */ 681 if (jmin < jmax){ 682 il[prow] = jmin; 683 j = *cols; jl[prow] = jl[j]; jl[j] = prow; 684 } 685 prow = nextprow; 686 } 687 688 /* if free space is not available, make more free space */ 689 if (current_space->local_remaining<nzk) { 690 i = am - k + 1; /* num of unfactored rows */ 691 i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 692 ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 693 ierr = PetscFreeSpaceGet(i,¤t_space_lvl);CHKERRQ(ierr); 694 reallocs++; 695 } 696 697 /* copy data into free_space and free_space_lvl, then initialize lnk */ 698 ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 699 700 /* add the k-th row into il and jl */ 701 if (nzk-1 > 0){ 702 i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 703 jl[k] = jl[i]; jl[i] = k; 704 il[k] = ui[k] + 1; 705 } 706 uj_ptr[k] = current_space->array; 707 uj_lvl_ptr[k] = current_space_lvl->array; 708 709 current_space->array += nzk; 710 current_space->local_used += nzk; 711 current_space->local_remaining -= nzk; 712 713 current_space_lvl->array += nzk; 714 current_space_lvl->local_used += nzk; 715 current_space_lvl->local_remaining -= nzk; 716 717 ui[k+1] = ui[k] + nzk; 718 } 719 720 #if defined(PETSC_USE_INFO) 721 if (ai[am] != 0) { 722 PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]); 723 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 724 ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 725 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 726 } else { 727 ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 728 } 729 #endif 730 731 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 732 ierr = PetscFree(jl);CHKERRQ(ierr); 733 ierr = PetscFree(cols_lvl);CHKERRQ(ierr); 734 735 /* destroy list of free space and other temporary array(s) */ 736 ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 737 ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 738 ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 739 ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 740 741 /* put together the new matrix in MATSEQSBAIJ format */ 742 B = *fact; 743 ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 744 745 b = (Mat_SeqSBAIJ*)B->data; 746 b->singlemalloc = PETSC_FALSE; 747 b->free_a = PETSC_TRUE; 748 b->free_ij = PETSC_TRUE; 749 ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 750 b->j = uj; 751 b->i = ui; 752 b->diag = 0; 753 b->ilen = 0; 754 b->imax = 0; 755 b->row = perm; 756 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 757 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 758 b->icol = perm; 759 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 760 ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 761 ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 762 b->maxnz = b->nz = ui[am]; 763 764 B->factor = MAT_FACTOR_CHOLESKY; 765 B->info.factor_mallocs = reallocs; 766 B->info.fill_ratio_given = fill; 767 if (ai[am] != 0) { 768 B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 769 } else { 770 B->info.fill_ratio_needed = 0.0; 771 } 772 if (perm_identity){ 773 B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 774 B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 775 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 776 } else { 777 (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 778 } 779 PetscFunctionReturn(0); 780 } 781 782 #undef __FUNCT__ 783 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ" 784 PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 785 { 786 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 787 Mat_SeqSBAIJ *b; 788 Mat B; 789 PetscErrorCode ierr; 790 PetscTruth perm_identity; 791 PetscReal fill = info->fill; 792 PetscInt *rip,i,mbs=a->mbs,bs=A->rmap.bs,*ai=a->i,*aj=a->j,reallocs=0,prow; 793 PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 794 PetscInt nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr; 795 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 796 PetscBT lnkbt; 797 798 PetscFunctionBegin; 799 if (bs > 1) { /* convert to seqsbaij */ 800 if (!a->sbaijMat){ 801 ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 802 } 803 (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */ 804 ierr = MatCholeskyFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr); 805 PetscFunctionReturn(0); 806 } 807 808 /* check whether perm is the identity mapping */ 809 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 810 if (!perm_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported"); 811 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 812 813 /* initialization */ 814 ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 815 ui[0] = 0; 816 817 /* jl: linked list for storing indices of the pivot rows 818 il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */ 819 ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr); 820 il = jl + mbs; 821 cols = il + mbs; 822 ui_ptr = (PetscInt**)(cols + mbs); 823 for (i=0; i<mbs; i++){ 824 jl[i] = mbs; il[i] = 0; 825 } 826 827 /* create and initialize a linked list for storing column indices of the active row k */ 828 nlnk = mbs + 1; 829 ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 830 831 /* initial FreeSpace size is fill*(ai[mbs]+1) */ 832 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr); 833 current_space = free_space; 834 835 for (k=0; k<mbs; k++){ /* for each active row k */ 836 /* initialize lnk by the column indices of row rip[k] of A */ 837 nzk = 0; 838 ncols = ai[rip[k]+1] - ai[rip[k]]; 839 ncols_upper = 0; 840 for (j=0; j<ncols; j++){ 841 i = rip[*(aj + ai[rip[k]] + j)]; 842 if (i >= k){ /* only take upper triangular entry */ 843 cols[ncols_upper] = i; 844 ncols_upper++; 845 } 846 } 847 ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 848 nzk += nlnk; 849 850 /* update lnk by computing fill-in for each pivot row to be merged in */ 851 prow = jl[k]; /* 1st pivot row */ 852 853 while (prow < k){ 854 nextprow = jl[prow]; 855 /* merge prow into k-th row */ 856 jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */ 857 jmax = ui[prow+1]; 858 ncols = jmax-jmin; 859 uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */ 860 ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 861 nzk += nlnk; 862 863 /* update il and jl for prow */ 864 if (jmin < jmax){ 865 il[prow] = jmin; 866 j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 867 } 868 prow = nextprow; 869 } 870 871 /* if free space is not available, make more free space */ 872 if (current_space->local_remaining<nzk) { 873 i = mbs - k + 1; /* num of unfactored rows */ 874 i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 875 ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 876 reallocs++; 877 } 878 879 /* copy data into free space, then initialize lnk */ 880 ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 881 882 /* add the k-th row into il and jl */ 883 if (nzk-1 > 0){ 884 i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */ 885 jl[k] = jl[i]; jl[i] = k; 886 il[k] = ui[k] + 1; 887 } 888 ui_ptr[k] = current_space->array; 889 current_space->array += nzk; 890 current_space->local_used += nzk; 891 current_space->local_remaining -= nzk; 892 893 ui[k+1] = ui[k] + nzk; 894 } 895 896 #if defined(PETSC_USE_INFO) 897 if (ai[mbs] != 0) { 898 PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 899 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 900 ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 901 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 902 } else { 903 ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 904 } 905 #endif 906 907 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 908 ierr = PetscFree(jl);CHKERRQ(ierr); 909 910 /* destroy list of free space and other temporary array(s) */ 911 ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 912 ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 913 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 914 915 /* put together the new matrix in MATSEQSBAIJ format */ 916 B = *fact; 917 ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 918 919 b = (Mat_SeqSBAIJ*)B->data; 920 b->singlemalloc = PETSC_FALSE; 921 b->free_a = PETSC_TRUE; 922 b->free_ij = PETSC_TRUE; 923 ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 924 b->j = uj; 925 b->i = ui; 926 b->diag = 0; 927 b->ilen = 0; 928 b->imax = 0; 929 b->row = perm; 930 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 931 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 932 b->icol = perm; 933 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 934 ierr = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 935 ierr = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 936 b->maxnz = b->nz = ui[mbs]; 937 938 B->factor = MAT_FACTOR_CHOLESKY; 939 B->info.factor_mallocs = reallocs; 940 B->info.fill_ratio_given = fill; 941 if (ai[mbs] != 0) { 942 B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 943 } else { 944 B->info.fill_ratio_needed = 0.0; 945 } 946 if (perm_identity){ 947 B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 948 B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 949 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 950 } else { 951 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 952 } 953 PetscFunctionReturn(0); 954 } 955