xref: /petsc/src/mat/impls/baij/seq/baijfact3.c (revision 6ba06ab77c48d2bbbcba7e829ffcd04868efd3d5)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
383287d42SBarry Smith /*
483287d42SBarry Smith     Factorization code for BAIJ format.
583287d42SBarry Smith */
67c4f633dSBarry Smith #include "../src/mat/impls/baij/seq/baij.h"
7c60f0209SBarry Smith #include "../src/mat/blockinvert.h"
883287d42SBarry Smith 
9db4efbfdSBarry Smith #undef __FUNCT__
104dd39f65SShri Abhyankar #define __FUNCT__ "MatSeqBAIJSetNumericFactorization"
11db4efbfdSBarry Smith /*
12db4efbfdSBarry Smith    This is used to set the numeric factorization for both LU and ILU symbolic factorization
13db4efbfdSBarry Smith */
144dd39f65SShri Abhyankar PetscErrorCode MatSeqBAIJSetNumericFactorization(Mat fact,PetscTruth natural)
15ae3d28f0SHong Zhang {
16ae3d28f0SHong Zhang   PetscFunctionBegin;
176929473cSShri Abhyankar   if(natural){
186929473cSShri Abhyankar     switch (fact->rmap->bs){
19048b5e81SShri Abhyankar     case 1:
20048b5e81SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
21048b5e81SShri Abhyankar       break;
226929473cSShri Abhyankar     case 2:
234dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering;
246929473cSShri Abhyankar       break;
256929473cSShri Abhyankar     case 3:
264dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering;
276929473cSShri Abhyankar       break;
286929473cSShri Abhyankar     case 4:
294dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering;
306929473cSShri Abhyankar       break;
316929473cSShri Abhyankar     case 5:
324dd39f65SShri Abhyankar        fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering;
336929473cSShri Abhyankar        break;
346929473cSShri Abhyankar     case 6:
354dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering;
366929473cSShri Abhyankar       break;
376929473cSShri Abhyankar     case 7:
384dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering;
396929473cSShri Abhyankar       break;
4029a97285SShri Abhyankar     case 15:
4129a97285SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_15_NaturalOrdering;
4229a97285SShri Abhyankar       break;
436929473cSShri Abhyankar     default:
444dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N;
456929473cSShri Abhyankar       break;
466929473cSShri Abhyankar     }
47*6ba06ab7SHong Zhang   } else {
48ae3d28f0SHong Zhang     switch (fact->rmap->bs){
49048b5e81SShri Abhyankar     case 1:
50048b5e81SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
51048b5e81SShri Abhyankar       break;
52ae3d28f0SHong Zhang     case 2:
534dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
54ae3d28f0SHong Zhang       break;
55ae3d28f0SHong Zhang     case 3:
564dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
57ae3d28f0SHong Zhang       break;
58ae3d28f0SHong Zhang     case 4:
594dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
60ae3d28f0SHong Zhang       break;
61ae3d28f0SHong Zhang     case 5:
624dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
63ae3d28f0SHong Zhang       break;
64ae3d28f0SHong Zhang     case 6:
654dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6;
66ae3d28f0SHong Zhang       break;
67ae3d28f0SHong Zhang     case 7:
684dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7;
69ae3d28f0SHong Zhang       break;
70ae3d28f0SHong Zhang     default:
714dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N;
72ae3d28f0SHong Zhang       break;
73ae3d28f0SHong Zhang     }
746929473cSShri Abhyankar   }
75ae3d28f0SHong Zhang   PetscFunctionReturn(0);
76ae3d28f0SHong Zhang }
77ae3d28f0SHong Zhang 
788b1456e3SHong Zhang PetscErrorCode MatSeqBAIJSetNumericFactorization_inplace(Mat inA,PetscTruth natural)
79db4efbfdSBarry Smith {
80db4efbfdSBarry Smith   PetscFunctionBegin;
81db4efbfdSBarry Smith   if (natural) {
82db4efbfdSBarry Smith     switch (inA->rmap->bs) {
83db4efbfdSBarry Smith     case 1:
8406e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1_inplace;
85db4efbfdSBarry Smith       break;
86db4efbfdSBarry Smith     case 2:
8706e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace;
88db4efbfdSBarry Smith       break;
89db4efbfdSBarry Smith     case 3:
9006e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace;
91db4efbfdSBarry Smith       break;
92db4efbfdSBarry Smith     case 4:
9365460251SBarry Smith #if defined(PETSC_USE_SCALAR_MAT_SINGLE)
94db4efbfdSBarry Smith       {
95db4efbfdSBarry Smith         PetscTruth  sse_enabled_local;
96db4efbfdSBarry Smith         PetscErrorCode ierr;
97db4efbfdSBarry Smith         ierr = PetscSSEIsEnabled(inA->comm,&sse_enabled_local,PETSC_NULL);CHKERRQ(ierr);
98db4efbfdSBarry Smith         if (sse_enabled_local) {
99db4efbfdSBarry Smith #  if defined(PETSC_HAVE_SSE)
100db4efbfdSBarry Smith           int i,*AJ=a->j,nz=a->nz,n=a->mbs;
101db4efbfdSBarry Smith           if (n==(unsigned short)n) {
102db4efbfdSBarry Smith             unsigned short *aj=(unsigned short *)AJ;
103db4efbfdSBarry Smith             for (i=0;i<nz;i++) {
104db4efbfdSBarry Smith               aj[i] = (unsigned short)AJ[i];
105db4efbfdSBarry Smith             }
106db4efbfdSBarry Smith             inA->ops->setunfactored   = MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj;
107db4efbfdSBarry Smith             inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE_usj;
108db4efbfdSBarry Smith             ierr = PetscInfo(inA,"Using special SSE, in-place natural ordering, ushort j index factor BS=4\n");CHKERRQ(ierr);
109db4efbfdSBarry Smith           } else {
110db4efbfdSBarry Smith         /* Scale the column indices for easier indexing in MatSolve. */
111db4efbfdSBarry Smith /*            for (i=0;i<nz;i++) { */
112db4efbfdSBarry Smith /*              AJ[i] = AJ[i]*4; */
113db4efbfdSBarry Smith /*            } */
114db4efbfdSBarry Smith             inA->ops->setunfactored   = MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE;
115db4efbfdSBarry Smith             inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE;
116db4efbfdSBarry Smith             ierr = PetscInfo(inA,"Using special SSE, in-place natural ordering, int j index factor BS=4\n");CHKERRQ(ierr);
117db4efbfdSBarry Smith           }
118db4efbfdSBarry Smith #  else
119db4efbfdSBarry Smith         /* This should never be reached.  If so, problem in PetscSSEIsEnabled. */
120e32f2f54SBarry Smith           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSE Hardware unavailable");
121db4efbfdSBarry Smith #  endif
122db4efbfdSBarry Smith         } else {
12306e38f1dSHong Zhang           inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_inplace;
124db4efbfdSBarry Smith         }
125db4efbfdSBarry Smith       }
126db4efbfdSBarry Smith #else
12706e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_inplace;
128db4efbfdSBarry Smith #endif
129db4efbfdSBarry Smith       break;
130db4efbfdSBarry Smith     case 5:
13106e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering_inplace;
132db4efbfdSBarry Smith       break;
133db4efbfdSBarry Smith     case 6:
13406e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering_inplace;
135db4efbfdSBarry Smith       break;
136db4efbfdSBarry Smith     case 7:
13706e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering_inplace;
138db4efbfdSBarry Smith       break;
139db4efbfdSBarry Smith     default:
14006e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N_inplace;
141db4efbfdSBarry Smith       break;
142db4efbfdSBarry Smith     }
143db4efbfdSBarry Smith   } else {
144db4efbfdSBarry Smith     switch (inA->rmap->bs) {
145db4efbfdSBarry Smith     case 1:
14606e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1_inplace;
147db4efbfdSBarry Smith       break;
148db4efbfdSBarry Smith     case 2:
14906e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_inplace;
150db4efbfdSBarry Smith       break;
151db4efbfdSBarry Smith     case 3:
15206e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_inplace;
153db4efbfdSBarry Smith       break;
154db4efbfdSBarry Smith     case 4:
15506e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_inplace;
156db4efbfdSBarry Smith       break;
157db4efbfdSBarry Smith     case 5:
15806e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_inplace;
159db4efbfdSBarry Smith       break;
160db4efbfdSBarry Smith     case 6:
16106e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_inplace;
162db4efbfdSBarry Smith       break;
163db4efbfdSBarry Smith     case 7:
16406e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_inplace;
165db4efbfdSBarry Smith       break;
166db4efbfdSBarry Smith     default:
16706e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N_inplace;
168db4efbfdSBarry Smith       break;
169db4efbfdSBarry Smith     }
170db4efbfdSBarry Smith   }
171db4efbfdSBarry Smith   PetscFunctionReturn(0);
172db4efbfdSBarry Smith }
173db4efbfdSBarry Smith 
17483287d42SBarry Smith /*
17583287d42SBarry Smith     The symbolic factorization code is identical to that for AIJ format,
17683287d42SBarry Smith   except for very small changes since this is now a SeqBAIJ datastructure.
17783287d42SBarry Smith   NOT good code reuse.
17883287d42SBarry Smith */
179c7272abeSHong Zhang #include "petscbt.h"
180c7272abeSHong Zhang #include "../src/mat/utils/freespace.h"
181c7272abeSHong Zhang 
1824a2ae208SSatish Balay #undef __FUNCT__
1834dd39f65SShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_SeqBAIJ"
1844dd39f65SShri Abhyankar PetscErrorCode MatLUFactorSymbolic_SeqBAIJ(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
18583287d42SBarry Smith {
18683287d42SBarry Smith   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data,*b;
187c7272abeSHong Zhang   PetscInt           n=a->mbs,bs = A->rmap->bs,bs2=a->bs2;
188c7272abeSHong Zhang   PetscTruth         row_identity,col_identity,both_identity;
18983287d42SBarry Smith   IS                 isicol;
1906849ba73SBarry Smith   PetscErrorCode     ierr;
191c7272abeSHong Zhang   const PetscInt     *r,*ic;
192c7272abeSHong Zhang   PetscInt           i,*ai=a->i,*aj=a->j;
193c7272abeSHong Zhang   PetscInt           *bi,*bj,*ajtmp;
194c7272abeSHong Zhang   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
195c7272abeSHong Zhang   PetscReal          f;
196c7272abeSHong Zhang   PetscInt           nlnk,*lnk,k,**bi_ptr;
197c7272abeSHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
198c7272abeSHong Zhang   PetscBT            lnkbt;
19983287d42SBarry Smith 
20083287d42SBarry Smith   PetscFunctionBegin;
201e32f2f54SBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"matrix must be square");
202*6ba06ab7SHong Zhang   if (bs>1){  /* check shifttype */
203*6ba06ab7SHong Zhang     if (info->shifttype == MAT_SHIFT_NONZERO || info->shifttype == MAT_SHIFT_POSITIVE_DEFINITE)
204*6ba06ab7SHong Zhang       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only MAT_SHIFT_NONE and MAT_SHIFT_INBLOCKS are supported for BAIJ matrix");
205*6ba06ab7SHong Zhang   }
206*6ba06ab7SHong Zhang 
20783287d42SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
20883287d42SBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
20983287d42SBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
21083287d42SBarry Smith 
211bc74c81fSJed Brown   /* get new row and diagonal pointers, must be allocated separately because they will be given to the Mat_SeqAIJ and freed separately */
212bc74c81fSJed Brown   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
213bc74c81fSJed Brown   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
214a6df321fSShri Abhyankar   bi[0] = bdiag[0] = 0;
215b2b2dd24SShri Abhyankar 
216b2b2dd24SShri Abhyankar   /* linked list for storing column indices of the active row */
217b2b2dd24SShri Abhyankar   nlnk = n + 1;
218b2b2dd24SShri Abhyankar   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
219b2b2dd24SShri Abhyankar 
220b2b2dd24SShri Abhyankar   ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr);
221b2b2dd24SShri Abhyankar 
222b2b2dd24SShri Abhyankar   /* initial FreeSpace size is f*(ai[n]+1) */
223b2b2dd24SShri Abhyankar   f = info->fill;
224b2b2dd24SShri Abhyankar   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
225b2b2dd24SShri Abhyankar   current_space = free_space;
226b2b2dd24SShri Abhyankar 
227b2b2dd24SShri Abhyankar   for (i=0; i<n; i++) {
228b2b2dd24SShri Abhyankar     /* copy previous fill into linked list */
229b2b2dd24SShri Abhyankar     nzi = 0;
230b2b2dd24SShri Abhyankar     nnz = ai[r[i]+1] - ai[r[i]];
231e32f2f54SBarry Smith     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);
232b2b2dd24SShri Abhyankar     ajtmp = aj + ai[r[i]];
233b2b2dd24SShri Abhyankar     ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
234b2b2dd24SShri Abhyankar     nzi += nlnk;
235b2b2dd24SShri Abhyankar 
236b2b2dd24SShri Abhyankar     /* add pivot rows into linked list */
237b2b2dd24SShri Abhyankar     row = lnk[n];
238b2b2dd24SShri Abhyankar     while (row < i) {
239b2b2dd24SShri Abhyankar       nzbd    = bdiag[row] + 1; /* num of entries in the row with column index <= row */
240b2b2dd24SShri Abhyankar       ajtmp   = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */
241b2b2dd24SShri Abhyankar       ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
242b2b2dd24SShri Abhyankar       nzi += nlnk;
243b2b2dd24SShri Abhyankar       row  = lnk[row];
244b2b2dd24SShri Abhyankar     }
245b2b2dd24SShri Abhyankar     bi[i+1] = bi[i] + nzi;
246b2b2dd24SShri Abhyankar     im[i]   = nzi;
247b2b2dd24SShri Abhyankar 
248b2b2dd24SShri Abhyankar     /* mark bdiag */
249b2b2dd24SShri Abhyankar     nzbd = 0;
250b2b2dd24SShri Abhyankar     nnz  = nzi;
251b2b2dd24SShri Abhyankar     k    = lnk[n];
252b2b2dd24SShri Abhyankar     while (nnz-- && k < i){
253b2b2dd24SShri Abhyankar       nzbd++;
254b2b2dd24SShri Abhyankar       k = lnk[k];
255b2b2dd24SShri Abhyankar     }
2562ce24eb6SHong Zhang     bdiag[i] = nzbd; /* note : bdaig[i] = nnzL as input for PetscFreeSpaceContiguous_LU() */
257b2b2dd24SShri Abhyankar 
258b2b2dd24SShri Abhyankar     /* if free space is not available, make more free space */
259b2b2dd24SShri Abhyankar     if (current_space->local_remaining<nzi) {
260b2b2dd24SShri Abhyankar       nnz = 2*(n - i)*nzi; /* estimated and max additional space needed */
261b2b2dd24SShri Abhyankar       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
262b2b2dd24SShri Abhyankar       reallocs++;
263b2b2dd24SShri Abhyankar     }
264b2b2dd24SShri Abhyankar 
265b2b2dd24SShri Abhyankar     /* copy data into free space, then initialize lnk */
266b2b2dd24SShri Abhyankar     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
267b2b2dd24SShri Abhyankar     bi_ptr[i] = current_space->array;
268b2b2dd24SShri Abhyankar     current_space->array           += nzi;
269b2b2dd24SShri Abhyankar     current_space->local_used      += nzi;
270b2b2dd24SShri Abhyankar     current_space->local_remaining -= nzi;
271b2b2dd24SShri Abhyankar   }
272b2b2dd24SShri Abhyankar #if defined(PETSC_USE_INFO)
273b2b2dd24SShri Abhyankar   if (ai[n] != 0) {
274b2b2dd24SShri Abhyankar     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
275b2b2dd24SShri Abhyankar     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
276b2b2dd24SShri Abhyankar     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
277b2b2dd24SShri Abhyankar     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
278b2b2dd24SShri Abhyankar     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
279b2b2dd24SShri Abhyankar   } else {
280b2b2dd24SShri Abhyankar     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
281b2b2dd24SShri Abhyankar   }
282b2b2dd24SShri Abhyankar #endif
283b2b2dd24SShri Abhyankar 
284b2b2dd24SShri Abhyankar   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
285b2b2dd24SShri Abhyankar   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
286b2b2dd24SShri Abhyankar 
287b2b2dd24SShri Abhyankar   /* destroy list of free space and other temporary array(s) */
288b2b2dd24SShri Abhyankar   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
2892ce24eb6SHong Zhang   ierr = PetscFreeSpaceContiguous_LU(&free_space,bj,n,bi,bdiag);CHKERRQ(ierr);
290b2b2dd24SShri Abhyankar   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
291b2b2dd24SShri Abhyankar   ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr);
292b2b2dd24SShri Abhyankar 
293b2b2dd24SShri Abhyankar   /* put together the new matrix */
294b2b2dd24SShri Abhyankar   ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
295b2b2dd24SShri Abhyankar   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
296b2b2dd24SShri Abhyankar   b    = (Mat_SeqBAIJ*)(B)->data;
297b2b2dd24SShri Abhyankar   b->free_a       = PETSC_TRUE;
298b2b2dd24SShri Abhyankar   b->free_ij      = PETSC_TRUE;
299b2b2dd24SShri Abhyankar   b->singlemalloc = PETSC_FALSE;
300b2b2dd24SShri Abhyankar   ierr          = PetscMalloc((bdiag[0]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
301b2b2dd24SShri Abhyankar   b->j          = bj;
302b2b2dd24SShri Abhyankar   b->i          = bi;
303b2b2dd24SShri Abhyankar   b->diag       = bdiag;
304b2b2dd24SShri Abhyankar   b->free_diag  = PETSC_TRUE;
305b2b2dd24SShri Abhyankar   b->ilen       = 0;
306b2b2dd24SShri Abhyankar   b->imax       = 0;
307b2b2dd24SShri Abhyankar   b->row        = isrow;
308b2b2dd24SShri Abhyankar   b->col        = iscol;
309b2b2dd24SShri Abhyankar   b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
310b2b2dd24SShri Abhyankar   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
311b2b2dd24SShri Abhyankar   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
312b2b2dd24SShri Abhyankar   b->icol       = isicol;
313b2b2dd24SShri Abhyankar   ierr = PetscMalloc((bs*n+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
314b2b2dd24SShri Abhyankar   ierr = PetscLogObjectMemory(B,(bdiag[0]+1)*(sizeof(PetscInt)+sizeof(PetscScalar)*bs2));CHKERRQ(ierr);
315b2b2dd24SShri Abhyankar 
316b2b2dd24SShri Abhyankar   b->maxnz = b->nz = bdiag[0]+1;
317d5f3da31SBarry Smith   B->factortype            =  MAT_FACTOR_LU;
318b2b2dd24SShri Abhyankar   B->info.factor_mallocs   = reallocs;
319b2b2dd24SShri Abhyankar   B->info.fill_ratio_given = f;
320b2b2dd24SShri Abhyankar 
321b2b2dd24SShri Abhyankar   if (ai[n] != 0) {
322b2b2dd24SShri Abhyankar     B->info.fill_ratio_needed = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]);
323b2b2dd24SShri Abhyankar   } else {
324b2b2dd24SShri Abhyankar     B->info.fill_ratio_needed = 0.0;
325b2b2dd24SShri Abhyankar   }
326b2b2dd24SShri Abhyankar 
327b2b2dd24SShri Abhyankar   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
328b2b2dd24SShri Abhyankar   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
329b2b2dd24SShri Abhyankar   both_identity = (PetscTruth) (row_identity && col_identity);
3304dd39f65SShri Abhyankar   ierr = MatSeqBAIJSetNumericFactorization(B,both_identity);CHKERRQ(ierr);
331b2b2dd24SShri Abhyankar   PetscFunctionReturn(0);
332b2b2dd24SShri Abhyankar  }
333b2b2dd24SShri Abhyankar 
334b2b2dd24SShri Abhyankar #undef __FUNCT__
33506e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorSymbolic_SeqBAIJ_inplace"
33606e38f1dSHong Zhang PetscErrorCode MatLUFactorSymbolic_SeqBAIJ_inplace(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
337faca2338SShri Abhyankar {
338faca2338SShri Abhyankar   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data,*b;
339faca2338SShri Abhyankar   PetscInt           n=a->mbs,bs = A->rmap->bs,bs2=a->bs2;
340faca2338SShri Abhyankar   PetscTruth         row_identity,col_identity,both_identity;
341faca2338SShri Abhyankar   IS                 isicol;
342faca2338SShri Abhyankar   PetscErrorCode     ierr;
343faca2338SShri Abhyankar   const PetscInt     *r,*ic;
344faca2338SShri Abhyankar   PetscInt           i,*ai=a->i,*aj=a->j;
345faca2338SShri Abhyankar   PetscInt           *bi,*bj,*ajtmp;
346faca2338SShri Abhyankar   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
347faca2338SShri Abhyankar   PetscReal          f;
348faca2338SShri Abhyankar   PetscInt           nlnk,*lnk,k,**bi_ptr;
349faca2338SShri Abhyankar   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
350faca2338SShri Abhyankar   PetscBT            lnkbt;
351faca2338SShri Abhyankar 
352faca2338SShri Abhyankar   PetscFunctionBegin;
353e32f2f54SBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"matrix must be square");
354faca2338SShri Abhyankar   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
355faca2338SShri Abhyankar   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
356faca2338SShri Abhyankar   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
357faca2338SShri Abhyankar 
358bc74c81fSJed Brown   /* get new row and diagonal pointers, must be allocated separately because they will be given to the Mat_SeqAIJ and freed separately */
359bc74c81fSJed Brown   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
360bc74c81fSJed Brown   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
361bc74c81fSJed Brown 
362a6df321fSShri Abhyankar   bi[0] = bdiag[0] = 0;
363c7272abeSHong Zhang 
364c7272abeSHong Zhang   /* linked list for storing column indices of the active row */
365c7272abeSHong Zhang   nlnk = n + 1;
366c7272abeSHong Zhang   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
367c7272abeSHong Zhang 
368c7272abeSHong Zhang   ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr);
369c7272abeSHong Zhang 
370c7272abeSHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
371c7272abeSHong Zhang   f = info->fill;
372c7272abeSHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
373c7272abeSHong Zhang   current_space = free_space;
37483287d42SBarry Smith 
37583287d42SBarry Smith   for (i=0; i<n; i++) {
376c7272abeSHong Zhang     /* copy previous fill into linked list */
37783287d42SBarry Smith     nzi = 0;
378c7272abeSHong Zhang     nnz = ai[r[i]+1] - ai[r[i]];
379e32f2f54SBarry Smith     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);
380c7272abeSHong Zhang     ajtmp = aj + ai[r[i]];
381c7272abeSHong Zhang     ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
382c7272abeSHong Zhang     nzi += nlnk;
383c7272abeSHong Zhang 
384c7272abeSHong Zhang     /* add pivot rows into linked list */
385c7272abeSHong Zhang     row = lnk[n];
386c7272abeSHong Zhang     while (row < i) {
387c7272abeSHong Zhang       nzbd    = bdiag[row] - bi[row] + 1; /* num of entries in the row with column index <= row */
388c7272abeSHong Zhang       ajtmp   = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */
389c7272abeSHong Zhang       ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
390c7272abeSHong Zhang       nzi += nlnk;
391c7272abeSHong Zhang       row  = lnk[row];
39283287d42SBarry Smith     }
393c7272abeSHong Zhang     bi[i+1] = bi[i] + nzi;
394c7272abeSHong Zhang     im[i]   = nzi;
395c7272abeSHong Zhang 
396c7272abeSHong Zhang     /* mark bdiag */
397c7272abeSHong Zhang     nzbd = 0;
398c7272abeSHong Zhang     nnz  = nzi;
399c7272abeSHong Zhang     k    = lnk[n];
400c7272abeSHong Zhang     while (nnz-- && k < i){
401c7272abeSHong Zhang       nzbd++;
402c7272abeSHong Zhang       k = lnk[k];
403c7272abeSHong Zhang     }
404c7272abeSHong Zhang     bdiag[i] = bi[i] + nzbd;
405c7272abeSHong Zhang 
406c7272abeSHong Zhang     /* if free space is not available, make more free space */
407c7272abeSHong Zhang     if (current_space->local_remaining<nzi) {
408c7272abeSHong Zhang       nnz = (n - i)*nzi; /* estimated and max additional space needed */
409c7272abeSHong Zhang       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
410c7272abeSHong Zhang       reallocs++;
41183287d42SBarry Smith     }
41283287d42SBarry Smith 
413c7272abeSHong Zhang     /* copy data into free space, then initialize lnk */
414c7272abeSHong Zhang     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
415c7272abeSHong Zhang     bi_ptr[i] = current_space->array;
416c7272abeSHong Zhang     current_space->array           += nzi;
417c7272abeSHong Zhang     current_space->local_used      += nzi;
418c7272abeSHong Zhang     current_space->local_remaining -= nzi;
419c7272abeSHong Zhang   }
4206cf91177SBarry Smith #if defined(PETSC_USE_INFO)
42183287d42SBarry Smith   if (ai[n] != 0) {
422c7272abeSHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
423ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
424ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
425ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
426ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
42783287d42SBarry Smith   } else {
428c7272abeSHong Zhang     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
42983287d42SBarry Smith   }
43063ba0a88SBarry Smith #endif
43183287d42SBarry Smith 
43283287d42SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
43383287d42SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
43483287d42SBarry Smith 
435c7272abeSHong Zhang   /* destroy list of free space and other temporary array(s) */
436c7272abeSHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
437c7272abeSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
438c7272abeSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
439c7272abeSHong Zhang   ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr);
44083287d42SBarry Smith 
44183287d42SBarry Smith   /* put together the new matrix */
442719d5645SBarry Smith   ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
443719d5645SBarry Smith   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
444719d5645SBarry Smith   b    = (Mat_SeqBAIJ*)(B)->data;
445e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
446e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
447c7272abeSHong Zhang   b->singlemalloc = PETSC_FALSE;
448c7272abeSHong Zhang   ierr          = PetscMalloc((bi[n]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
449c7272abeSHong Zhang   b->j          = bj;
450c7272abeSHong Zhang   b->i          = bi;
451c7272abeSHong Zhang   b->diag       = bdiag;
4527f53bb6cSHong Zhang   b->free_diag  = PETSC_TRUE;
45383287d42SBarry Smith   b->ilen       = 0;
45483287d42SBarry Smith   b->imax       = 0;
45583287d42SBarry Smith   b->row        = isrow;
45683287d42SBarry Smith   b->col        = iscol;
457d3d32019SBarry Smith   b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
45883287d42SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
45983287d42SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
46083287d42SBarry Smith   b->icol       = isicol;
46187828ca2SBarry Smith   ierr = PetscMalloc((bs*n+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
462c7272abeSHong Zhang   ierr = PetscLogObjectMemory(B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)*bs2));CHKERRQ(ierr);
46383287d42SBarry Smith 
464c7272abeSHong Zhang   b->maxnz = b->nz = bi[n] ;
465d5f3da31SBarry Smith   (B)->factortype            =  MAT_FACTOR_LU;
466719d5645SBarry Smith   (B)->info.factor_mallocs   = reallocs;
467719d5645SBarry Smith   (B)->info.fill_ratio_given = f;
468c7272abeSHong Zhang 
46983287d42SBarry Smith   if (ai[n] != 0) {
470c7272abeSHong Zhang     (B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
47183287d42SBarry Smith   } else {
472719d5645SBarry Smith     (B)->info.fill_ratio_needed = 0.0;
47383287d42SBarry Smith   }
474c7272abeSHong Zhang 
475db4efbfdSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
476db4efbfdSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
4777d18ce8fSMatthew Knepley   both_identity = (PetscTruth) (row_identity && col_identity);
4788b1456e3SHong Zhang   ierr = MatSeqBAIJSetNumericFactorization_inplace(B,both_identity);CHKERRQ(ierr);
47983287d42SBarry Smith   PetscFunctionReturn(0);
48083287d42SBarry Smith  }
481db4efbfdSBarry Smith 
482