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