1 #define PETSCMAT_DLL 2 3 /* 4 Factorization code for BAIJ format. 5 */ 6 #include "../src/mat/impls/baij/seq/baij.h" 7 #include "../src/mat/blockinvert.h" 8 9 #undef __FUNCT__ 10 #define __FUNCT__ "MatSeqBAIJSetNumericFactorization" 11 /* 12 This is used to set the numeric factorization for both LU and ILU symbolic factorization 13 */ 14 PetscErrorCode MatSeqBAIJSetNumericFactorization(Mat fact,PetscTruth natural) 15 { 16 PetscFunctionBegin; 17 if(natural){ 18 switch (fact->rmap->bs){ 19 case 2: 20 fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering; 21 break; 22 case 3: 23 fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering; 24 break; 25 case 4: 26 fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering; 27 break; 28 case 5: 29 fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering; 30 break; 31 case 6: 32 fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering; 33 break; 34 case 7: 35 fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering; 36 break; 37 default: 38 fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N; 39 break; 40 } 41 } 42 else{ 43 switch (fact->rmap->bs){ 44 case 2: 45 fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 46 break; 47 case 3: 48 fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 49 break; 50 case 4: 51 fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 52 break; 53 case 5: 54 fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 55 break; 56 case 6: 57 fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6; 58 break; 59 case 7: 60 fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7; 61 break; 62 case 15: 63 fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_15; 64 break; 65 default: 66 fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N; 67 break; 68 } 69 } 70 PetscFunctionReturn(0); 71 } 72 73 PetscErrorCode MatSeqBAIJSetNumericFactorization_inplace(Mat inA,PetscTruth natural) 74 { 75 PetscFunctionBegin; 76 if (natural) { 77 switch (inA->rmap->bs) { 78 case 1: 79 inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1_inplace; 80 break; 81 case 2: 82 inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace; 83 break; 84 case 3: 85 inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace; 86 break; 87 case 4: 88 #if defined(PETSC_USE_SCALAR_MAT_SINGLE) 89 { 90 PetscTruth sse_enabled_local; 91 PetscErrorCode ierr; 92 ierr = PetscSSEIsEnabled(inA->comm,&sse_enabled_local,PETSC_NULL);CHKERRQ(ierr); 93 if (sse_enabled_local) { 94 # if defined(PETSC_HAVE_SSE) 95 int i,*AJ=a->j,nz=a->nz,n=a->mbs; 96 if (n==(unsigned short)n) { 97 unsigned short *aj=(unsigned short *)AJ; 98 for (i=0;i<nz;i++) { 99 aj[i] = (unsigned short)AJ[i]; 100 } 101 inA->ops->setunfactored = MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj; 102 inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE_usj; 103 ierr = PetscInfo(inA,"Using special SSE, in-place natural ordering, ushort j index factor BS=4\n");CHKERRQ(ierr); 104 } else { 105 /* Scale the column indices for easier indexing in MatSolve. */ 106 /* for (i=0;i<nz;i++) { */ 107 /* AJ[i] = AJ[i]*4; */ 108 /* } */ 109 inA->ops->setunfactored = MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE; 110 inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE; 111 ierr = PetscInfo(inA,"Using special SSE, in-place natural ordering, int j index factor BS=4\n");CHKERRQ(ierr); 112 } 113 # else 114 /* This should never be reached. If so, problem in PetscSSEIsEnabled. */ 115 SETERRQ(PETSC_ERR_SUP,"SSE Hardware unavailable"); 116 # endif 117 } else { 118 inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_inplace; 119 } 120 } 121 #else 122 inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_inplace; 123 #endif 124 break; 125 case 5: 126 inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering_inplace; 127 break; 128 case 6: 129 inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering_inplace; 130 break; 131 case 7: 132 inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering_inplace; 133 break; 134 default: 135 inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N_inplace; 136 break; 137 } 138 } else { 139 switch (inA->rmap->bs) { 140 case 1: 141 inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1_inplace; 142 break; 143 case 2: 144 inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_inplace; 145 break; 146 case 3: 147 inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_inplace; 148 break; 149 case 4: 150 inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_inplace; 151 break; 152 case 5: 153 inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_inplace; 154 break; 155 case 6: 156 inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_inplace; 157 break; 158 case 7: 159 inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_inplace; 160 break; 161 default: 162 inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N_inplace; 163 break; 164 } 165 } 166 PetscFunctionReturn(0); 167 } 168 169 /* 170 The symbolic factorization code is identical to that for AIJ format, 171 except for very small changes since this is now a SeqBAIJ datastructure. 172 NOT good code reuse. 173 */ 174 #include "petscbt.h" 175 #include "../src/mat/utils/freespace.h" 176 177 #undef __FUNCT__ 178 #define __FUNCT__ "MatLUFactorSymbolic_SeqBAIJ" 179 PetscErrorCode MatLUFactorSymbolic_SeqBAIJ(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 180 { 181 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b; 182 PetscInt n=a->mbs,bs = A->rmap->bs,bs2=a->bs2; 183 PetscTruth row_identity,col_identity,both_identity; 184 IS isicol; 185 PetscErrorCode ierr; 186 const PetscInt *r,*ic; 187 PetscInt i,*ai=a->i,*aj=a->j; 188 PetscInt *bi,*bj,*ajtmp; 189 PetscInt *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im; 190 PetscReal f; 191 PetscInt nlnk,*lnk,k,**bi_ptr; 192 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 193 PetscBT lnkbt; 194 PetscTruth newdatastruct=PETSC_FALSE; 195 196 PetscFunctionBegin; 197 ierr = PetscOptionsGetTruth(PETSC_NULL,"-lu_old",&newdatastruct,PETSC_NULL);CHKERRQ(ierr); 198 if(newdatastruct){ 199 ierr = MatLUFactorSymbolic_SeqBAIJ_inplace(B,A,isrow,iscol,info);CHKERRQ(ierr); 200 PetscFunctionReturn(0); 201 } 202 if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square"); 203 ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 204 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 205 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 206 207 /* get new row and diagonal pointers, must be allocated separately because they will be given to the Mat_SeqAIJ and freed separately */ 208 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 209 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 210 bi[0] = bdiag[0] = 0; 211 212 /* linked list for storing column indices of the active row */ 213 nlnk = n + 1; 214 ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 215 216 ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr); 217 218 /* initial FreeSpace size is f*(ai[n]+1) */ 219 f = info->fill; 220 ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr); 221 current_space = free_space; 222 223 for (i=0; i<n; i++) { 224 /* copy previous fill into linked list */ 225 nzi = 0; 226 nnz = ai[r[i]+1] - ai[r[i]]; 227 if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i); 228 ajtmp = aj + ai[r[i]]; 229 ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 230 nzi += nlnk; 231 232 /* add pivot rows into linked list */ 233 row = lnk[n]; 234 while (row < i) { 235 nzbd = bdiag[row] + 1; /* num of entries in the row with column index <= row */ 236 ajtmp = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */ 237 ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr); 238 nzi += nlnk; 239 row = lnk[row]; 240 } 241 bi[i+1] = bi[i] + nzi; 242 im[i] = nzi; 243 244 /* mark bdiag */ 245 nzbd = 0; 246 nnz = nzi; 247 k = lnk[n]; 248 while (nnz-- && k < i){ 249 nzbd++; 250 k = lnk[k]; 251 } 252 bdiag[i] = nzbd; /* note : bdaig[i] = nnzL as input for PetscFreeSpaceContiguous_LU() */ 253 254 /* if free space is not available, make more free space */ 255 if (current_space->local_remaining<nzi) { 256 nnz = 2*(n - i)*nzi; /* estimated and max additional space needed */ 257 ierr = PetscFreeSpaceGet(nnz,¤t_space);CHKERRQ(ierr); 258 reallocs++; 259 } 260 261 /* copy data into free space, then initialize lnk */ 262 ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 263 bi_ptr[i] = current_space->array; 264 current_space->array += nzi; 265 current_space->local_used += nzi; 266 current_space->local_remaining -= nzi; 267 } 268 #if defined(PETSC_USE_INFO) 269 if (ai[n] != 0) { 270 PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 271 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr); 272 ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 273 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr); 274 ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 275 } else { 276 ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr); 277 } 278 #endif 279 280 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 281 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 282 283 /* destroy list of free space and other temporary array(s) */ 284 ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 285 ierr = PetscFreeSpaceContiguous_LU(&free_space,bj,n,bi,bdiag);CHKERRQ(ierr); 286 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 287 ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr); 288 289 /* put together the new matrix */ 290 ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 291 ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr); 292 b = (Mat_SeqBAIJ*)(B)->data; 293 b->free_a = PETSC_TRUE; 294 b->free_ij = PETSC_TRUE; 295 b->singlemalloc = PETSC_FALSE; 296 ierr = PetscMalloc((bdiag[0]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr); 297 b->j = bj; 298 b->i = bi; 299 b->diag = bdiag; 300 b->free_diag = PETSC_TRUE; 301 b->ilen = 0; 302 b->imax = 0; 303 b->row = isrow; 304 b->col = iscol; 305 b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 306 ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 307 ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 308 b->icol = isicol; 309 ierr = PetscMalloc((bs*n+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 310 ierr = PetscLogObjectMemory(B,(bdiag[0]+1)*(sizeof(PetscInt)+sizeof(PetscScalar)*bs2));CHKERRQ(ierr); 311 312 b->maxnz = b->nz = bdiag[0]+1; 313 B->factor = MAT_FACTOR_LU; 314 B->info.factor_mallocs = reallocs; 315 B->info.fill_ratio_given = f; 316 317 if (ai[n] != 0) { 318 B->info.fill_ratio_needed = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]); 319 } else { 320 B->info.fill_ratio_needed = 0.0; 321 } 322 323 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 324 ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 325 both_identity = (PetscTruth) (row_identity && col_identity); 326 ierr = MatSeqBAIJSetNumericFactorization(B,both_identity);CHKERRQ(ierr); 327 PetscFunctionReturn(0); 328 } 329 330 #undef __FUNCT__ 331 #define __FUNCT__ "MatLUFactorSymbolic_SeqBAIJ_inplace" 332 PetscErrorCode MatLUFactorSymbolic_SeqBAIJ_inplace(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 333 { 334 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b; 335 PetscInt n=a->mbs,bs = A->rmap->bs,bs2=a->bs2; 336 PetscTruth row_identity,col_identity,both_identity; 337 IS isicol; 338 PetscErrorCode ierr; 339 const PetscInt *r,*ic; 340 PetscInt i,*ai=a->i,*aj=a->j; 341 PetscInt *bi,*bj,*ajtmp; 342 PetscInt *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im; 343 PetscReal f; 344 PetscInt nlnk,*lnk,k,**bi_ptr; 345 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 346 PetscBT lnkbt; 347 348 PetscFunctionBegin; 349 if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square"); 350 ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 351 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 352 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 353 354 /* get new row and diagonal pointers, must be allocated separately because they will be given to the Mat_SeqAIJ and freed separately */ 355 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 356 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 357 358 bi[0] = bdiag[0] = 0; 359 360 /* linked list for storing column indices of the active row */ 361 nlnk = n + 1; 362 ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 363 364 ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr); 365 366 /* initial FreeSpace size is f*(ai[n]+1) */ 367 f = info->fill; 368 ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr); 369 current_space = free_space; 370 371 for (i=0; i<n; i++) { 372 /* copy previous fill into linked list */ 373 nzi = 0; 374 nnz = ai[r[i]+1] - ai[r[i]]; 375 if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i); 376 ajtmp = aj + ai[r[i]]; 377 ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 378 nzi += nlnk; 379 380 /* add pivot rows into linked list */ 381 row = lnk[n]; 382 while (row < i) { 383 nzbd = bdiag[row] - bi[row] + 1; /* num of entries in the row with column index <= row */ 384 ajtmp = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */ 385 ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr); 386 nzi += nlnk; 387 row = lnk[row]; 388 } 389 bi[i+1] = bi[i] + nzi; 390 im[i] = nzi; 391 392 /* mark bdiag */ 393 nzbd = 0; 394 nnz = nzi; 395 k = lnk[n]; 396 while (nnz-- && k < i){ 397 nzbd++; 398 k = lnk[k]; 399 } 400 bdiag[i] = bi[i] + nzbd; 401 402 /* if free space is not available, make more free space */ 403 if (current_space->local_remaining<nzi) { 404 nnz = (n - i)*nzi; /* estimated and max additional space needed */ 405 ierr = PetscFreeSpaceGet(nnz,¤t_space);CHKERRQ(ierr); 406 reallocs++; 407 } 408 409 /* copy data into free space, then initialize lnk */ 410 ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 411 bi_ptr[i] = current_space->array; 412 current_space->array += nzi; 413 current_space->local_used += nzi; 414 current_space->local_remaining -= nzi; 415 } 416 #if defined(PETSC_USE_INFO) 417 if (ai[n] != 0) { 418 PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 419 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr); 420 ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 421 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr); 422 ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 423 } else { 424 ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr); 425 } 426 #endif 427 428 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 429 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 430 431 /* destroy list of free space and other temporary array(s) */ 432 ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 433 ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 434 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 435 ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr); 436 437 /* put together the new matrix */ 438 ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 439 ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr); 440 b = (Mat_SeqBAIJ*)(B)->data; 441 b->free_a = PETSC_TRUE; 442 b->free_ij = PETSC_TRUE; 443 b->singlemalloc = PETSC_FALSE; 444 ierr = PetscMalloc((bi[n]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr); 445 b->j = bj; 446 b->i = bi; 447 b->diag = bdiag; 448 b->free_diag = PETSC_TRUE; 449 b->ilen = 0; 450 b->imax = 0; 451 b->row = isrow; 452 b->col = iscol; 453 b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 454 ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 455 ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 456 b->icol = isicol; 457 ierr = PetscMalloc((bs*n+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 458 ierr = PetscLogObjectMemory(B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)*bs2));CHKERRQ(ierr); 459 460 b->maxnz = b->nz = bi[n] ; 461 (B)->factor = MAT_FACTOR_LU; 462 (B)->info.factor_mallocs = reallocs; 463 (B)->info.fill_ratio_given = f; 464 465 if (ai[n] != 0) { 466 (B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]); 467 } else { 468 (B)->info.fill_ratio_needed = 0.0; 469 } 470 471 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 472 ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 473 both_identity = (PetscTruth) (row_identity && col_identity); 474 ierr = MatSeqBAIJSetNumericFactorization_inplace(B,both_identity);CHKERRQ(ierr); 475 PetscFunctionReturn(0); 476 } 477 478