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