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