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