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 EXTERN_C_BEGIN 255 #undef __FUNCT__ 256 #define __FUNCT__ "MatGetFactor_seqbaij_petsc" 257 PetscErrorCode MatGetFactor_seqbaij_petsc(Mat A,MatFactorType ftype,Mat *B) 258 { 259 PetscInt n = A->rmap.n; 260 PetscErrorCode ierr; 261 262 PetscFunctionBegin; 263 ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr); 264 ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr); 265 if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) { 266 ierr = MatSetType(*B,MATSEQBAIJ);CHKERRQ(ierr); 267 } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) { 268 ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr); 269 ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 270 (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqBAIJ; 271 (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqBAIJ; 272 } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported"); 273 PetscFunctionReturn(0); 274 } 275 EXTERN_C_END 276 277 /* ----------------------------------------------------------- */ 278 #undef __FUNCT__ 279 #define __FUNCT__ "MatLUFactor_SeqBAIJ" 280 PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,MatFactorInfo *info) 281 { 282 PetscErrorCode ierr; 283 Mat C; 284 285 PetscFunctionBegin; 286 ierr = MatGetFactor(A,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr); 287 ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr); 288 ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr); 289 ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 290 ierr = PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr); 291 PetscFunctionReturn(0); 292 } 293 294 #include "src/mat/impls/sbaij/seq/sbaij.h" 295 #undef __FUNCT__ 296 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N" 297 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat A,MatFactorInfo *info,Mat *B) 298 { 299 PetscErrorCode ierr; 300 Mat C = *B; 301 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 302 Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 303 IS ip=b->row; 304 PetscInt *rip,i,j,mbs=a->mbs,bs=A->rmap.bs,*bi=b->i,*bj=b->j,*bcol; 305 PetscInt *ai=a->i,*aj=a->j; 306 PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 307 MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; 308 PetscReal zeropivot,rs,shiftnz; 309 PetscReal shiftpd; 310 ChShift_Ctx sctx; 311 PetscInt newshift; 312 313 PetscFunctionBegin; 314 if (bs > 1) { 315 if (!a->sbaijMat){ 316 ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 317 } 318 ierr = (a->sbaijMat)->ops->choleskyfactornumeric(a->sbaijMat,info,B);CHKERRQ(ierr); 319 ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr); 320 a->sbaijMat = PETSC_NULL; 321 PetscFunctionReturn(0); 322 } 323 324 /* initialization */ 325 shiftnz = info->shiftnz; 326 shiftpd = info->shiftpd; 327 zeropivot = info->zeropivot; 328 329 ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 330 nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); 331 ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 332 jl = il + mbs; 333 rtmp = (MatScalar*)(jl + mbs); 334 335 sctx.shift_amount = 0; 336 sctx.nshift = 0; 337 do { 338 sctx.chshift = PETSC_FALSE; 339 for (i=0; i<mbs; i++) { 340 rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 341 } 342 343 for (k = 0; k<mbs; k++){ 344 bval = ba + bi[k]; 345 /* initialize k-th row by the perm[k]-th row of A */ 346 jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 347 for (j = jmin; j < jmax; j++){ 348 col = rip[aj[j]]; 349 if (col >= k){ /* only take upper triangular entry */ 350 rtmp[col] = aa[j]; 351 *bval++ = 0.0; /* for in-place factorization */ 352 } 353 } 354 355 /* shift the diagonal of the matrix */ 356 if (sctx.nshift) rtmp[k] += sctx.shift_amount; 357 358 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 359 dk = rtmp[k]; 360 i = jl[k]; /* first row to be added to k_th row */ 361 362 while (i < k){ 363 nexti = jl[i]; /* next row to be added to k_th row */ 364 365 /* compute multiplier, update diag(k) and U(i,k) */ 366 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 367 uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 368 dk += uikdi*ba[ili]; 369 ba[ili] = uikdi; /* -U(i,k) */ 370 371 /* add multiple of row i to k-th row */ 372 jmin = ili + 1; jmax = bi[i+1]; 373 if (jmin < jmax){ 374 for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 375 /* update il and jl for row i */ 376 il[i] = jmin; 377 j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 378 } 379 i = nexti; 380 } 381 382 /* shift the diagonals when zero pivot is detected */ 383 /* compute rs=sum of abs(off-diagonal) */ 384 rs = 0.0; 385 jmin = bi[k]+1; 386 nz = bi[k+1] - jmin; 387 if (nz){ 388 bcol = bj + jmin; 389 while (nz--){ 390 rs += PetscAbsScalar(rtmp[*bcol]); 391 bcol++; 392 } 393 } 394 395 sctx.rs = rs; 396 sctx.pv = dk; 397 ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr); 398 if (newshift == 1) break; 399 400 /* copy data into U(k,:) */ 401 ba[bi[k]] = 1.0/dk; /* U(k,k) */ 402 jmin = bi[k]+1; jmax = bi[k+1]; 403 if (jmin < jmax) { 404 for (j=jmin; j<jmax; j++){ 405 col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 406 } 407 /* add the k-th row into il and jl */ 408 il[k] = jmin; 409 i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 410 } 411 } 412 } while (sctx.chshift); 413 ierr = PetscFree(il);CHKERRQ(ierr); 414 415 ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 416 C->factor = MAT_FACTOR_CHOLESKY; 417 C->assembled = PETSC_TRUE; 418 C->preallocated = PETSC_TRUE; 419 ierr = PetscLogFlops(C->rmap.N);CHKERRQ(ierr); 420 if (sctx.nshift){ 421 if (shiftnz) { 422 ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 423 } else if (shiftpd) { 424 ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 425 } 426 } 427 PetscFunctionReturn(0); 428 } 429 430 #undef __FUNCT__ 431 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering" 432 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *fact) 433 { 434 Mat C = *fact; 435 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 436 Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 437 PetscErrorCode ierr; 438 PetscInt i,j,am=a->mbs; 439 PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 440 PetscInt k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz; 441 MatScalar *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval; 442 PetscReal zeropivot,rs,shiftnz; 443 PetscReal shiftpd; 444 ChShift_Ctx sctx; 445 PetscInt newshift; 446 447 PetscFunctionBegin; 448 /* initialization */ 449 shiftnz = info->shiftnz; 450 shiftpd = info->shiftpd; 451 zeropivot = info->zeropivot; 452 453 nz = (2*am+1)*sizeof(PetscInt)+am*sizeof(MatScalar); 454 ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 455 jl = il + am; 456 rtmp = (MatScalar*)(jl + am); 457 458 sctx.shift_amount = 0; 459 sctx.nshift = 0; 460 do { 461 sctx.chshift = PETSC_FALSE; 462 for (i=0; i<am; i++) { 463 rtmp[i] = 0.0; jl[i] = am; il[0] = 0; 464 } 465 466 for (k = 0; k<am; k++){ 467 /* initialize k-th row with elements nonzero in row perm(k) of A */ 468 nz = ai[k+1] - ai[k]; 469 acol = aj + ai[k]; 470 aval = aa + ai[k]; 471 bval = ba + bi[k]; 472 while (nz -- ){ 473 if (*acol < k) { /* skip lower triangular entries */ 474 acol++; aval++; 475 } else { 476 rtmp[*acol++] = *aval++; 477 *bval++ = 0.0; /* for in-place factorization */ 478 } 479 } 480 481 /* shift the diagonal of the matrix */ 482 if (sctx.nshift) rtmp[k] += sctx.shift_amount; 483 484 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 485 dk = rtmp[k]; 486 i = jl[k]; /* first row to be added to k_th row */ 487 488 while (i < k){ 489 nexti = jl[i]; /* next row to be added to k_th row */ 490 /* compute multiplier, update D(k) and U(i,k) */ 491 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 492 uikdi = - ba[ili]*ba[bi[i]]; 493 dk += uikdi*ba[ili]; 494 ba[ili] = uikdi; /* -U(i,k) */ 495 496 /* add multiple of row i to k-th row ... */ 497 jmin = ili + 1; 498 nz = bi[i+1] - jmin; 499 if (nz > 0){ 500 bcol = bj + jmin; 501 bval = ba + jmin; 502 while (nz --) rtmp[*bcol++] += uikdi*(*bval++); 503 /* update il and jl for i-th row */ 504 il[i] = jmin; 505 j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 506 } 507 i = nexti; 508 } 509 510 /* shift the diagonals when zero pivot is detected */ 511 /* compute rs=sum of abs(off-diagonal) */ 512 rs = 0.0; 513 jmin = bi[k]+1; 514 nz = bi[k+1] - jmin; 515 if (nz){ 516 bcol = bj + jmin; 517 while (nz--){ 518 rs += PetscAbsScalar(rtmp[*bcol]); 519 bcol++; 520 } 521 } 522 523 sctx.rs = rs; 524 sctx.pv = dk; 525 ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr); 526 if (newshift == 1) break; /* sctx.shift_amount is updated */ 527 528 /* copy data into U(k,:) */ 529 ba[bi[k]] = 1.0/dk; 530 jmin = bi[k]+1; 531 nz = bi[k+1] - jmin; 532 if (nz){ 533 bcol = bj + jmin; 534 bval = ba + jmin; 535 while (nz--){ 536 *bval++ = rtmp[*bcol]; 537 rtmp[*bcol++] = 0.0; 538 } 539 /* add k-th row into il and jl */ 540 il[k] = jmin; 541 i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 542 } 543 } 544 } while (sctx.chshift); 545 ierr = PetscFree(il);CHKERRQ(ierr); 546 547 C->factor = MAT_FACTOR_CHOLESKY; 548 C->assembled = PETSC_TRUE; 549 C->preallocated = PETSC_TRUE; 550 ierr = PetscLogFlops(C->rmap.N);CHKERRQ(ierr); 551 if (sctx.nshift){ 552 if (shiftnz) { 553 ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 554 } else if (shiftpd) { 555 ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 556 } 557 } 558 PetscFunctionReturn(0); 559 } 560 561 #include "petscbt.h" 562 #include "src/mat/utils/freespace.h" 563 #undef __FUNCT__ 564 #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ" 565 PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 566 { 567 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 568 Mat_SeqSBAIJ *b; 569 Mat B; 570 PetscErrorCode ierr; 571 PetscTruth perm_identity; 572 PetscInt reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap.bs,*ui; 573 PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow; 574 PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr; 575 PetscReal fill=info->fill,levels=info->levels; 576 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 577 PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 578 PetscBT lnkbt; 579 580 PetscFunctionBegin; 581 if (bs > 1){ 582 if (!a->sbaijMat){ 583 ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 584 } 585 (*fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */ 586 ierr = MatICCFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr); 587 PetscFunctionReturn(0); 588 } 589 590 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 591 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 592 593 /* special case that simply copies fill pattern */ 594 if (!levels && perm_identity) { 595 ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr); 596 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 597 for (i=0; i<am; i++) { 598 ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */ 599 } 600 B = *fact; 601 ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr); 602 603 604 b = (Mat_SeqSBAIJ*)B->data; 605 uj = b->j; 606 for (i=0; i<am; i++) { 607 aj = a->j + a->diag[i]; 608 for (j=0; j<ui[i]; j++){ 609 *uj++ = *aj++; 610 } 611 b->ilen[i] = ui[i]; 612 } 613 ierr = PetscFree(ui);CHKERRQ(ierr); 614 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 615 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 616 617 B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 618 B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 619 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 620 PetscFunctionReturn(0); 621 } 622 623 /* initialization */ 624 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 625 ui[0] = 0; 626 ierr = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr); 627 628 /* jl: linked list for storing indices of the pivot rows 629 il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 630 ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr); 631 il = jl + am; 632 uj_ptr = (PetscInt**)(il + am); 633 uj_lvl_ptr = (PetscInt**)(uj_ptr + am); 634 for (i=0; i<am; i++){ 635 jl[i] = am; il[i] = 0; 636 } 637 638 /* create and initialize a linked list for storing column indices of the active row k */ 639 nlnk = am + 1; 640 ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 641 642 /* initial FreeSpace size is fill*(ai[am]+1) */ 643 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 644 current_space = free_space; 645 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr); 646 current_space_lvl = free_space_lvl; 647 648 for (k=0; k<am; k++){ /* for each active row k */ 649 /* initialize lnk by the column indices of row rip[k] of A */ 650 nzk = 0; 651 ncols = ai[rip[k]+1] - ai[rip[k]]; 652 ncols_upper = 0; 653 cols = cols_lvl + am; 654 for (j=0; j<ncols; j++){ 655 i = rip[*(aj + ai[rip[k]] + j)]; 656 if (i >= k){ /* only take upper triangular entry */ 657 cols[ncols_upper] = i; 658 cols_lvl[ncols_upper] = -1; /* initialize level for nonzero entries */ 659 ncols_upper++; 660 } 661 } 662 ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 663 nzk += nlnk; 664 665 /* update lnk by computing fill-in for each pivot row to be merged in */ 666 prow = jl[k]; /* 1st pivot row */ 667 668 while (prow < k){ 669 nextprow = jl[prow]; 670 671 /* merge prow into k-th row */ 672 jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 673 jmax = ui[prow+1]; 674 ncols = jmax-jmin; 675 i = jmin - ui[prow]; 676 cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 677 for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j); 678 ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 679 nzk += nlnk; 680 681 /* update il and jl for prow */ 682 if (jmin < jmax){ 683 il[prow] = jmin; 684 j = *cols; jl[prow] = jl[j]; jl[j] = prow; 685 } 686 prow = nextprow; 687 } 688 689 /* if free space is not available, make more free space */ 690 if (current_space->local_remaining<nzk) { 691 i = am - k + 1; /* num of unfactored rows */ 692 i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 693 ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 694 ierr = PetscFreeSpaceGet(i,¤t_space_lvl);CHKERRQ(ierr); 695 reallocs++; 696 } 697 698 /* copy data into free_space and free_space_lvl, then initialize lnk */ 699 ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 700 701 /* add the k-th row into il and jl */ 702 if (nzk-1 > 0){ 703 i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 704 jl[k] = jl[i]; jl[i] = k; 705 il[k] = ui[k] + 1; 706 } 707 uj_ptr[k] = current_space->array; 708 uj_lvl_ptr[k] = current_space_lvl->array; 709 710 current_space->array += nzk; 711 current_space->local_used += nzk; 712 current_space->local_remaining -= nzk; 713 714 current_space_lvl->array += nzk; 715 current_space_lvl->local_used += nzk; 716 current_space_lvl->local_remaining -= nzk; 717 718 ui[k+1] = ui[k] + nzk; 719 } 720 721 #if defined(PETSC_USE_INFO) 722 if (ai[am] != 0) { 723 PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]); 724 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 725 ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 726 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 727 } else { 728 ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 729 } 730 #endif 731 732 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 733 ierr = PetscFree(jl);CHKERRQ(ierr); 734 ierr = PetscFree(cols_lvl);CHKERRQ(ierr); 735 736 /* destroy list of free space and other temporary array(s) */ 737 ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 738 ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 739 ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 740 ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 741 742 /* put together the new matrix in MATSEQSBAIJ format */ 743 B = *fact; 744 ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 745 746 b = (Mat_SeqSBAIJ*)B->data; 747 b->singlemalloc = PETSC_FALSE; 748 b->free_a = PETSC_TRUE; 749 b->free_ij = PETSC_TRUE; 750 ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 751 b->j = uj; 752 b->i = ui; 753 b->diag = 0; 754 b->ilen = 0; 755 b->imax = 0; 756 b->row = perm; 757 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 758 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 759 b->icol = perm; 760 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 761 ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 762 ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 763 b->maxnz = b->nz = ui[am]; 764 765 B->factor = MAT_FACTOR_CHOLESKY; 766 B->info.factor_mallocs = reallocs; 767 B->info.fill_ratio_given = fill; 768 if (ai[am] != 0) { 769 B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 770 } else { 771 B->info.fill_ratio_needed = 0.0; 772 } 773 if (perm_identity){ 774 B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 775 B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 776 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 777 } else { 778 (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 779 } 780 PetscFunctionReturn(0); 781 } 782 783 #undef __FUNCT__ 784 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ" 785 PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 786 { 787 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 788 Mat_SeqSBAIJ *b; 789 Mat B; 790 PetscErrorCode ierr; 791 PetscTruth perm_identity; 792 PetscReal fill = info->fill; 793 PetscInt *rip,i,mbs=a->mbs,bs=A->rmap.bs,*ai=a->i,*aj=a->j,reallocs=0,prow; 794 PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 795 PetscInt nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr; 796 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 797 PetscBT lnkbt; 798 799 PetscFunctionBegin; 800 if (bs > 1) { /* convert to seqsbaij */ 801 if (!a->sbaijMat){ 802 ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 803 } 804 (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */ 805 ierr = MatCholeskyFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr); 806 PetscFunctionReturn(0); 807 } 808 809 /* check whether perm is the identity mapping */ 810 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 811 if (!perm_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported"); 812 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 813 814 /* initialization */ 815 ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 816 ui[0] = 0; 817 818 /* jl: linked list for storing indices of the pivot rows 819 il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */ 820 ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr); 821 il = jl + mbs; 822 cols = il + mbs; 823 ui_ptr = (PetscInt**)(cols + mbs); 824 for (i=0; i<mbs; i++){ 825 jl[i] = mbs; il[i] = 0; 826 } 827 828 /* create and initialize a linked list for storing column indices of the active row k */ 829 nlnk = mbs + 1; 830 ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 831 832 /* initial FreeSpace size is fill*(ai[mbs]+1) */ 833 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr); 834 current_space = free_space; 835 836 for (k=0; k<mbs; k++){ /* for each active row k */ 837 /* initialize lnk by the column indices of row rip[k] of A */ 838 nzk = 0; 839 ncols = ai[rip[k]+1] - ai[rip[k]]; 840 ncols_upper = 0; 841 for (j=0; j<ncols; j++){ 842 i = rip[*(aj + ai[rip[k]] + j)]; 843 if (i >= k){ /* only take upper triangular entry */ 844 cols[ncols_upper] = i; 845 ncols_upper++; 846 } 847 } 848 ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 849 nzk += nlnk; 850 851 /* update lnk by computing fill-in for each pivot row to be merged in */ 852 prow = jl[k]; /* 1st pivot row */ 853 854 while (prow < k){ 855 nextprow = jl[prow]; 856 /* merge prow into k-th row */ 857 jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */ 858 jmax = ui[prow+1]; 859 ncols = jmax-jmin; 860 uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */ 861 ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 862 nzk += nlnk; 863 864 /* update il and jl for prow */ 865 if (jmin < jmax){ 866 il[prow] = jmin; 867 j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 868 } 869 prow = nextprow; 870 } 871 872 /* if free space is not available, make more free space */ 873 if (current_space->local_remaining<nzk) { 874 i = mbs - k + 1; /* num of unfactored rows */ 875 i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 876 ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 877 reallocs++; 878 } 879 880 /* copy data into free space, then initialize lnk */ 881 ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 882 883 /* add the k-th row into il and jl */ 884 if (nzk-1 > 0){ 885 i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */ 886 jl[k] = jl[i]; jl[i] = k; 887 il[k] = ui[k] + 1; 888 } 889 ui_ptr[k] = current_space->array; 890 current_space->array += nzk; 891 current_space->local_used += nzk; 892 current_space->local_remaining -= nzk; 893 894 ui[k+1] = ui[k] + nzk; 895 } 896 897 #if defined(PETSC_USE_INFO) 898 if (ai[mbs] != 0) { 899 PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 900 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 901 ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 902 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 903 } else { 904 ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 905 } 906 #endif 907 908 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 909 ierr = PetscFree(jl);CHKERRQ(ierr); 910 911 /* destroy list of free space and other temporary array(s) */ 912 ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 913 ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 914 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 915 916 /* put together the new matrix in MATSEQSBAIJ format */ 917 B = *fact; 918 ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 919 920 b = (Mat_SeqSBAIJ*)B->data; 921 b->singlemalloc = PETSC_FALSE; 922 b->free_a = PETSC_TRUE; 923 b->free_ij = PETSC_TRUE; 924 ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 925 b->j = uj; 926 b->i = ui; 927 b->diag = 0; 928 b->ilen = 0; 929 b->imax = 0; 930 b->row = perm; 931 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 932 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 933 b->icol = perm; 934 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 935 ierr = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 936 ierr = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 937 b->maxnz = b->nz = ui[mbs]; 938 939 B->factor = MAT_FACTOR_CHOLESKY; 940 B->info.factor_mallocs = reallocs; 941 B->info.fill_ratio_given = fill; 942 if (ai[mbs] != 0) { 943 B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 944 } else { 945 B->info.fill_ratio_needed = 0.0; 946 } 947 if (perm_identity){ 948 B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 949 B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 950 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 951 } else { 952 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 953 } 954 PetscFunctionReturn(0); 955 } 956