xref: /petsc/src/mat/impls/baij/seq/baijfact3.c (revision ee0bf03dc93d4c7c8a89c171f526ef710b24ee2e)
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, must be allocated separately because they will be given to the Mat_SeqAIJ and freed separately */
143   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
144   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
145   bi[0] = bdiag[0] = 0;
146 
147   /* linked list for storing column indices of the active row */
148   nlnk = n + 1;
149   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
150 
151   ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr);
152 
153   /* initial FreeSpace size is f*(ai[n]+1) */
154   f = info->fill;
155   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
156   current_space = free_space;
157 
158   for (i=0; i<n; i++) {
159     /* copy previous fill into linked list */
160     nzi = 0;
161     nnz = ai[r[i]+1] - ai[r[i]];
162     if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
163     ajtmp = aj + ai[r[i]];
164     ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
165     nzi += nlnk;
166 
167     /* add pivot rows into linked list */
168     row = lnk[n];
169     while (row < i) {
170       nzbd    = bdiag[row] + 1; /* num of entries in the row with column index <= row */
171       ajtmp   = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */
172       ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
173       nzi += nlnk;
174       row  = lnk[row];
175     }
176     bi[i+1] = bi[i] + nzi;
177     im[i]   = nzi;
178 
179     /* mark bdiag */
180     nzbd = 0;
181     nnz  = nzi;
182     k    = lnk[n];
183     while (nnz-- && k < i){
184       nzbd++;
185       k = lnk[k];
186     }
187     bdiag[i] = nzbd; /* note : bdaig[i] = nnzL as input for PetscFreeSpaceContiguous_LU_v2() */
188 
189     /* if free space is not available, make more free space */
190     if (current_space->local_remaining<nzi) {
191       nnz = 2*(n - i)*nzi; /* estimated and max additional space needed */
192       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
193       reallocs++;
194     }
195 
196     /* copy data into free space, then initialize lnk */
197     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
198     bi_ptr[i] = current_space->array;
199     current_space->array           += nzi;
200     current_space->local_used      += nzi;
201     current_space->local_remaining -= nzi;
202   }
203 #if defined(PETSC_USE_INFO)
204   if (ai[n] != 0) {
205     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
206     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
207     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
208     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
209     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
210   } else {
211     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
212   }
213 #endif
214 
215   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
216   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
217 
218   /* destroy list of free space and other temporary array(s) */
219   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
220   ierr = PetscFreeSpaceContiguous_LU_v2(&free_space,bj,n,bi,bdiag);CHKERRQ(ierr);
221   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
222   ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr);
223 
224   /* put together the new matrix */
225   ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
226   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
227   b    = (Mat_SeqBAIJ*)(B)->data;
228   b->free_a       = PETSC_TRUE;
229   b->free_ij      = PETSC_TRUE;
230   b->singlemalloc = PETSC_FALSE;
231   ierr          = PetscMalloc((bdiag[0]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
232   b->j          = bj;
233   b->i          = bi;
234   b->diag       = bdiag;
235   b->free_diag  = PETSC_TRUE;
236   b->ilen       = 0;
237   b->imax       = 0;
238   b->row        = isrow;
239   b->col        = iscol;
240   b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
241   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
242   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
243   b->icol       = isicol;
244   ierr = PetscMalloc((bs*n+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
245   ierr = PetscLogObjectMemory(B,(bdiag[0]+1)*(sizeof(PetscInt)+sizeof(PetscScalar)*bs2));CHKERRQ(ierr);
246 
247   b->maxnz = b->nz = bdiag[0]+1;
248   B->factor                =  MAT_FACTOR_LU;
249   B->info.factor_mallocs   = reallocs;
250   B->info.fill_ratio_given = f;
251 
252   if (ai[n] != 0) {
253     B->info.fill_ratio_needed = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]);
254   } else {
255     B->info.fill_ratio_needed = 0.0;
256   }
257 
258   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
259   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
260   both_identity = (PetscTruth) (row_identity && col_identity);
261   if(both_identity){
262      switch(bs){
263      case 2:
264         B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_newdatastruct;
265         B->ops->solve           = MatSolve_SeqBAIJ_2_NaturalOrdering_newdatastruct_v2;
266         break;
267      case 3:
268 	B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_newdatastruct;
269         B->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering_newdatastruct_v2;
270         break;
271      case 4:
272         B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_newdatastruct;
273         B->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering_newdatastruct_v2;
274         break;
275      case 5:
276         B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering_newdatastruct;
277         B->ops->solve = MatSolve_SeqBAIJ_5_NaturalOrdering_newdatastruct_v2;
278         break;
279      case 6:
280         B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering_newdatastruct;
281         B->ops->solve = MatSolve_SeqBAIJ_6_NaturalOrdering_newdatastruct_v2;
282         break;
283      case 7:
284         B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering_newdatastruct;
285         B->ops->solve = MatSolve_SeqBAIJ_7_NaturalOrdering_newdatastruct_v2;
286         break;
287      default:
288         B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N_newdatastruct;
289         B->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering_newdatastruct_v2;
290       }
291    }
292    else {
293      switch(bs){
294      case 2:
295         B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_newdatastruct;
296         B->ops->solve           = MatSolve_SeqBAIJ_2_newdatastruct_v2;
297         break;
298      case 3:
299 	B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_newdatastruct;
300         B->ops->solve = MatSolve_SeqBAIJ_3_newdatastruct_v2;
301         break;
302      case 4:
303         B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_newdatastruct;
304         B->ops->solve = MatSolve_SeqBAIJ_4_newdatastruct_v2;
305         break;
306      case 5:
307         B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_newdatastruct;
308         B->ops->solve = MatSolve_SeqBAIJ_5_newdatastruct_v2;
309         break;
310      case 6:
311         B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_newdatastruct;
312         B->ops->solve = MatSolve_SeqBAIJ_6_newdatastruct_v2;
313         break;
314      case 7:
315         B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_newdatastruct;
316         B->ops->solve = MatSolve_SeqBAIJ_7_newdatastruct_v2;
317         break;
318      default:
319         B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N_newdatastruct;
320         B->ops->solve = MatSolve_SeqBAIJ_N_newdatastruct_v2;
321       }
322    }
323   /*  ierr = MatSeqBAIJSetNumericFactorization(B,both_identity);CHKERRQ(ierr); */
324   PetscFunctionReturn(0);
325  }
326 
327 #undef __FUNCT__
328 #define __FUNCT__ "MatLUFactorSymbolic_SeqBAIJ"
329 PetscErrorCode MatLUFactorSymbolic_SeqBAIJ(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   PetscTruth         newdatastruct=PETSC_FALSE;
345 
346   PetscFunctionBegin;
347   ierr = PetscOptionsGetTruth(PETSC_NULL,"-lu_new",&newdatastruct,PETSC_NULL);CHKERRQ(ierr);
348   if(newdatastruct){
349     ierr = MatLUFactorSymbolic_SeqBAIJ_newdatastruct(B,A,isrow,iscol,info);CHKERRQ(ierr);
350     PetscFunctionReturn(0);
351   }
352 
353   if (A->rmap->N != A->cmap->N) SETERRQ(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_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,&current_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)->factor                =  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 = (PetscTruth) (row_identity && col_identity);
478   ierr = MatSeqBAIJSetNumericFactorization(B,both_identity);CHKERRQ(ierr);
479   PetscFunctionReturn(0);
480  }
481 
482