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