xref: /petsc/src/mat/impls/baij/seq/baijfact3.c (revision e32f2f54e699d0aa6e733466c00da7e34666fe5e)
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){
196929473cSShri Abhyankar     case 2:
204dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering;
216929473cSShri Abhyankar       break;
226929473cSShri Abhyankar     case 3:
234dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering;
246929473cSShri Abhyankar       break;
256929473cSShri Abhyankar     case 4:
264dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering;
276929473cSShri Abhyankar       break;
286929473cSShri Abhyankar     case 5:
294dd39f65SShri Abhyankar        fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering;
306929473cSShri Abhyankar        break;
316929473cSShri Abhyankar     case 6:
324dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering;
336929473cSShri Abhyankar       break;
346929473cSShri Abhyankar     case 7:
354dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering;
366929473cSShri Abhyankar       break;
3729a97285SShri Abhyankar     case 15:
3829a97285SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_15_NaturalOrdering;
3929a97285SShri Abhyankar       break;
406929473cSShri Abhyankar     default:
414dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N;
426929473cSShri Abhyankar       break;
436929473cSShri Abhyankar     }
446929473cSShri Abhyankar   }
456929473cSShri Abhyankar   else{
46ae3d28f0SHong Zhang     switch (fact->rmap->bs){
47ae3d28f0SHong Zhang     case 2:
484dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
49ae3d28f0SHong Zhang       break;
50ae3d28f0SHong Zhang     case 3:
514dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
52ae3d28f0SHong Zhang       break;
53ae3d28f0SHong Zhang     case 4:
544dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
55ae3d28f0SHong Zhang       break;
56ae3d28f0SHong Zhang     case 5:
574dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
58ae3d28f0SHong Zhang       break;
59ae3d28f0SHong Zhang     case 6:
604dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6;
61ae3d28f0SHong Zhang       break;
62ae3d28f0SHong Zhang     case 7:
634dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7;
64ae3d28f0SHong Zhang       break;
65ae3d28f0SHong Zhang     default:
664dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N;
67ae3d28f0SHong Zhang       break;
68ae3d28f0SHong Zhang     }
696929473cSShri Abhyankar   }
70ae3d28f0SHong Zhang   PetscFunctionReturn(0);
71ae3d28f0SHong Zhang }
72ae3d28f0SHong Zhang 
738b1456e3SHong Zhang PetscErrorCode MatSeqBAIJSetNumericFactorization_inplace(Mat inA,PetscTruth natural)
74db4efbfdSBarry Smith {
75db4efbfdSBarry Smith   PetscFunctionBegin;
76db4efbfdSBarry Smith   if (natural) {
77db4efbfdSBarry Smith     switch (inA->rmap->bs) {
78db4efbfdSBarry Smith     case 1:
7906e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1_inplace;
80db4efbfdSBarry Smith       break;
81db4efbfdSBarry Smith     case 2:
8206e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace;
83db4efbfdSBarry Smith       break;
84db4efbfdSBarry Smith     case 3:
8506e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace;
86db4efbfdSBarry Smith       break;
87db4efbfdSBarry Smith     case 4:
8865460251SBarry Smith #if defined(PETSC_USE_SCALAR_MAT_SINGLE)
89db4efbfdSBarry Smith       {
90db4efbfdSBarry Smith         PetscTruth  sse_enabled_local;
91db4efbfdSBarry Smith         PetscErrorCode ierr;
92db4efbfdSBarry Smith         ierr = PetscSSEIsEnabled(inA->comm,&sse_enabled_local,PETSC_NULL);CHKERRQ(ierr);
93db4efbfdSBarry Smith         if (sse_enabled_local) {
94db4efbfdSBarry Smith #  if defined(PETSC_HAVE_SSE)
95db4efbfdSBarry Smith           int i,*AJ=a->j,nz=a->nz,n=a->mbs;
96db4efbfdSBarry Smith           if (n==(unsigned short)n) {
97db4efbfdSBarry Smith             unsigned short *aj=(unsigned short *)AJ;
98db4efbfdSBarry Smith             for (i=0;i<nz;i++) {
99db4efbfdSBarry Smith               aj[i] = (unsigned short)AJ[i];
100db4efbfdSBarry Smith             }
101db4efbfdSBarry Smith             inA->ops->setunfactored   = MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj;
102db4efbfdSBarry Smith             inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE_usj;
103db4efbfdSBarry Smith             ierr = PetscInfo(inA,"Using special SSE, in-place natural ordering, ushort j index factor BS=4\n");CHKERRQ(ierr);
104db4efbfdSBarry Smith           } else {
105db4efbfdSBarry Smith         /* Scale the column indices for easier indexing in MatSolve. */
106db4efbfdSBarry Smith /*            for (i=0;i<nz;i++) { */
107db4efbfdSBarry Smith /*              AJ[i] = AJ[i]*4; */
108db4efbfdSBarry Smith /*            } */
109db4efbfdSBarry Smith             inA->ops->setunfactored   = MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE;
110db4efbfdSBarry Smith             inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE;
111db4efbfdSBarry Smith             ierr = PetscInfo(inA,"Using special SSE, in-place natural ordering, int j index factor BS=4\n");CHKERRQ(ierr);
112db4efbfdSBarry Smith           }
113db4efbfdSBarry Smith #  else
114db4efbfdSBarry Smith         /* This should never be reached.  If so, problem in PetscSSEIsEnabled. */
115*e32f2f54SBarry Smith           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSE Hardware unavailable");
116db4efbfdSBarry Smith #  endif
117db4efbfdSBarry Smith         } else {
11806e38f1dSHong Zhang           inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_inplace;
119db4efbfdSBarry Smith         }
120db4efbfdSBarry Smith       }
121db4efbfdSBarry Smith #else
12206e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_inplace;
123db4efbfdSBarry Smith #endif
124db4efbfdSBarry Smith       break;
125db4efbfdSBarry Smith     case 5:
12606e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering_inplace;
127db4efbfdSBarry Smith       break;
128db4efbfdSBarry Smith     case 6:
12906e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering_inplace;
130db4efbfdSBarry Smith       break;
131db4efbfdSBarry Smith     case 7:
13206e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering_inplace;
133db4efbfdSBarry Smith       break;
134db4efbfdSBarry Smith     default:
13506e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N_inplace;
136db4efbfdSBarry Smith       break;
137db4efbfdSBarry Smith     }
138db4efbfdSBarry Smith   } else {
139db4efbfdSBarry Smith     switch (inA->rmap->bs) {
140db4efbfdSBarry Smith     case 1:
14106e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1_inplace;
142db4efbfdSBarry Smith       break;
143db4efbfdSBarry Smith     case 2:
14406e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_inplace;
145db4efbfdSBarry Smith       break;
146db4efbfdSBarry Smith     case 3:
14706e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_inplace;
148db4efbfdSBarry Smith       break;
149db4efbfdSBarry Smith     case 4:
15006e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_inplace;
151db4efbfdSBarry Smith       break;
152db4efbfdSBarry Smith     case 5:
15306e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_inplace;
154db4efbfdSBarry Smith       break;
155db4efbfdSBarry Smith     case 6:
15606e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_inplace;
157db4efbfdSBarry Smith       break;
158db4efbfdSBarry Smith     case 7:
15906e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_inplace;
160db4efbfdSBarry Smith       break;
161db4efbfdSBarry Smith     default:
16206e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N_inplace;
163db4efbfdSBarry Smith       break;
164db4efbfdSBarry Smith     }
165db4efbfdSBarry Smith   }
166db4efbfdSBarry Smith   PetscFunctionReturn(0);
167db4efbfdSBarry Smith }
168db4efbfdSBarry Smith 
16983287d42SBarry Smith /*
17083287d42SBarry Smith     The symbolic factorization code is identical to that for AIJ format,
17183287d42SBarry Smith   except for very small changes since this is now a SeqBAIJ datastructure.
17283287d42SBarry Smith   NOT good code reuse.
17383287d42SBarry Smith */
174c7272abeSHong Zhang #include "petscbt.h"
175c7272abeSHong Zhang #include "../src/mat/utils/freespace.h"
176c7272abeSHong Zhang 
1774a2ae208SSatish Balay #undef __FUNCT__
1784dd39f65SShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_SeqBAIJ"
1794dd39f65SShri Abhyankar PetscErrorCode MatLUFactorSymbolic_SeqBAIJ(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
18083287d42SBarry Smith {
18183287d42SBarry Smith   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data,*b;
182c7272abeSHong Zhang   PetscInt           n=a->mbs,bs = A->rmap->bs,bs2=a->bs2;
183c7272abeSHong Zhang   PetscTruth         row_identity,col_identity,both_identity;
18483287d42SBarry Smith   IS                 isicol;
1856849ba73SBarry Smith   PetscErrorCode     ierr;
186c7272abeSHong Zhang   const PetscInt     *r,*ic;
187c7272abeSHong Zhang   PetscInt           i,*ai=a->i,*aj=a->j;
188c7272abeSHong Zhang   PetscInt           *bi,*bj,*ajtmp;
189c7272abeSHong Zhang   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
190c7272abeSHong Zhang   PetscReal          f;
191c7272abeSHong Zhang   PetscInt           nlnk,*lnk,k,**bi_ptr;
192c7272abeSHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
193c7272abeSHong Zhang   PetscBT            lnkbt;
19483287d42SBarry Smith 
19583287d42SBarry Smith   PetscFunctionBegin;
196*e32f2f54SBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"matrix must be square");
19783287d42SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
19883287d42SBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
19983287d42SBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
20083287d42SBarry Smith 
201bc74c81fSJed Brown   /* get new row and diagonal pointers, must be allocated separately because they will be given to the Mat_SeqAIJ and freed separately */
202bc74c81fSJed Brown   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
203bc74c81fSJed Brown   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
204a6df321fSShri Abhyankar   bi[0] = bdiag[0] = 0;
205b2b2dd24SShri Abhyankar 
206b2b2dd24SShri Abhyankar   /* linked list for storing column indices of the active row */
207b2b2dd24SShri Abhyankar   nlnk = n + 1;
208b2b2dd24SShri Abhyankar   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
209b2b2dd24SShri Abhyankar 
210b2b2dd24SShri Abhyankar   ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr);
211b2b2dd24SShri Abhyankar 
212b2b2dd24SShri Abhyankar   /* initial FreeSpace size is f*(ai[n]+1) */
213b2b2dd24SShri Abhyankar   f = info->fill;
214b2b2dd24SShri Abhyankar   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
215b2b2dd24SShri Abhyankar   current_space = free_space;
216b2b2dd24SShri Abhyankar 
217b2b2dd24SShri Abhyankar   for (i=0; i<n; i++) {
218b2b2dd24SShri Abhyankar     /* copy previous fill into linked list */
219b2b2dd24SShri Abhyankar     nzi = 0;
220b2b2dd24SShri Abhyankar     nnz = ai[r[i]+1] - ai[r[i]];
221*e32f2f54SBarry 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);
222b2b2dd24SShri Abhyankar     ajtmp = aj + ai[r[i]];
223b2b2dd24SShri Abhyankar     ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
224b2b2dd24SShri Abhyankar     nzi += nlnk;
225b2b2dd24SShri Abhyankar 
226b2b2dd24SShri Abhyankar     /* add pivot rows into linked list */
227b2b2dd24SShri Abhyankar     row = lnk[n];
228b2b2dd24SShri Abhyankar     while (row < i) {
229b2b2dd24SShri Abhyankar       nzbd    = bdiag[row] + 1; /* num of entries in the row with column index <= row */
230b2b2dd24SShri Abhyankar       ajtmp   = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */
231b2b2dd24SShri Abhyankar       ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
232b2b2dd24SShri Abhyankar       nzi += nlnk;
233b2b2dd24SShri Abhyankar       row  = lnk[row];
234b2b2dd24SShri Abhyankar     }
235b2b2dd24SShri Abhyankar     bi[i+1] = bi[i] + nzi;
236b2b2dd24SShri Abhyankar     im[i]   = nzi;
237b2b2dd24SShri Abhyankar 
238b2b2dd24SShri Abhyankar     /* mark bdiag */
239b2b2dd24SShri Abhyankar     nzbd = 0;
240b2b2dd24SShri Abhyankar     nnz  = nzi;
241b2b2dd24SShri Abhyankar     k    = lnk[n];
242b2b2dd24SShri Abhyankar     while (nnz-- && k < i){
243b2b2dd24SShri Abhyankar       nzbd++;
244b2b2dd24SShri Abhyankar       k = lnk[k];
245b2b2dd24SShri Abhyankar     }
2462ce24eb6SHong Zhang     bdiag[i] = nzbd; /* note : bdaig[i] = nnzL as input for PetscFreeSpaceContiguous_LU() */
247b2b2dd24SShri Abhyankar 
248b2b2dd24SShri Abhyankar     /* if free space is not available, make more free space */
249b2b2dd24SShri Abhyankar     if (current_space->local_remaining<nzi) {
250b2b2dd24SShri Abhyankar       nnz = 2*(n - i)*nzi; /* estimated and max additional space needed */
251b2b2dd24SShri Abhyankar       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
252b2b2dd24SShri Abhyankar       reallocs++;
253b2b2dd24SShri Abhyankar     }
254b2b2dd24SShri Abhyankar 
255b2b2dd24SShri Abhyankar     /* copy data into free space, then initialize lnk */
256b2b2dd24SShri Abhyankar     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
257b2b2dd24SShri Abhyankar     bi_ptr[i] = current_space->array;
258b2b2dd24SShri Abhyankar     current_space->array           += nzi;
259b2b2dd24SShri Abhyankar     current_space->local_used      += nzi;
260b2b2dd24SShri Abhyankar     current_space->local_remaining -= nzi;
261b2b2dd24SShri Abhyankar   }
262b2b2dd24SShri Abhyankar #if defined(PETSC_USE_INFO)
263b2b2dd24SShri Abhyankar   if (ai[n] != 0) {
264b2b2dd24SShri Abhyankar     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
265b2b2dd24SShri Abhyankar     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
266b2b2dd24SShri Abhyankar     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
267b2b2dd24SShri Abhyankar     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
268b2b2dd24SShri Abhyankar     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
269b2b2dd24SShri Abhyankar   } else {
270b2b2dd24SShri Abhyankar     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
271b2b2dd24SShri Abhyankar   }
272b2b2dd24SShri Abhyankar #endif
273b2b2dd24SShri Abhyankar 
274b2b2dd24SShri Abhyankar   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
275b2b2dd24SShri Abhyankar   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
276b2b2dd24SShri Abhyankar 
277b2b2dd24SShri Abhyankar   /* destroy list of free space and other temporary array(s) */
278b2b2dd24SShri Abhyankar   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
2792ce24eb6SHong Zhang   ierr = PetscFreeSpaceContiguous_LU(&free_space,bj,n,bi,bdiag);CHKERRQ(ierr);
280b2b2dd24SShri Abhyankar   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
281b2b2dd24SShri Abhyankar   ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr);
282b2b2dd24SShri Abhyankar 
283b2b2dd24SShri Abhyankar   /* put together the new matrix */
284b2b2dd24SShri Abhyankar   ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
285b2b2dd24SShri Abhyankar   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
286b2b2dd24SShri Abhyankar   b    = (Mat_SeqBAIJ*)(B)->data;
287b2b2dd24SShri Abhyankar   b->free_a       = PETSC_TRUE;
288b2b2dd24SShri Abhyankar   b->free_ij      = PETSC_TRUE;
289b2b2dd24SShri Abhyankar   b->singlemalloc = PETSC_FALSE;
290b2b2dd24SShri Abhyankar   ierr          = PetscMalloc((bdiag[0]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
291b2b2dd24SShri Abhyankar   b->j          = bj;
292b2b2dd24SShri Abhyankar   b->i          = bi;
293b2b2dd24SShri Abhyankar   b->diag       = bdiag;
294b2b2dd24SShri Abhyankar   b->free_diag  = PETSC_TRUE;
295b2b2dd24SShri Abhyankar   b->ilen       = 0;
296b2b2dd24SShri Abhyankar   b->imax       = 0;
297b2b2dd24SShri Abhyankar   b->row        = isrow;
298b2b2dd24SShri Abhyankar   b->col        = iscol;
299b2b2dd24SShri Abhyankar   b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
300b2b2dd24SShri Abhyankar   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
301b2b2dd24SShri Abhyankar   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
302b2b2dd24SShri Abhyankar   b->icol       = isicol;
303b2b2dd24SShri Abhyankar   ierr = PetscMalloc((bs*n+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
304b2b2dd24SShri Abhyankar   ierr = PetscLogObjectMemory(B,(bdiag[0]+1)*(sizeof(PetscInt)+sizeof(PetscScalar)*bs2));CHKERRQ(ierr);
305b2b2dd24SShri Abhyankar 
306b2b2dd24SShri Abhyankar   b->maxnz = b->nz = bdiag[0]+1;
307d5f3da31SBarry Smith   B->factortype            =  MAT_FACTOR_LU;
308b2b2dd24SShri Abhyankar   B->info.factor_mallocs   = reallocs;
309b2b2dd24SShri Abhyankar   B->info.fill_ratio_given = f;
310b2b2dd24SShri Abhyankar 
311b2b2dd24SShri Abhyankar   if (ai[n] != 0) {
312b2b2dd24SShri Abhyankar     B->info.fill_ratio_needed = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]);
313b2b2dd24SShri Abhyankar   } else {
314b2b2dd24SShri Abhyankar     B->info.fill_ratio_needed = 0.0;
315b2b2dd24SShri Abhyankar   }
316b2b2dd24SShri Abhyankar 
317b2b2dd24SShri Abhyankar   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
318b2b2dd24SShri Abhyankar   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
319b2b2dd24SShri Abhyankar   both_identity = (PetscTruth) (row_identity && col_identity);
3204dd39f65SShri Abhyankar   ierr = MatSeqBAIJSetNumericFactorization(B,both_identity);CHKERRQ(ierr);
321b2b2dd24SShri Abhyankar   PetscFunctionReturn(0);
322b2b2dd24SShri Abhyankar  }
323b2b2dd24SShri Abhyankar 
324b2b2dd24SShri Abhyankar #undef __FUNCT__
32506e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorSymbolic_SeqBAIJ_inplace"
32606e38f1dSHong Zhang PetscErrorCode MatLUFactorSymbolic_SeqBAIJ_inplace(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
327faca2338SShri Abhyankar {
328faca2338SShri Abhyankar   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data,*b;
329faca2338SShri Abhyankar   PetscInt           n=a->mbs,bs = A->rmap->bs,bs2=a->bs2;
330faca2338SShri Abhyankar   PetscTruth         row_identity,col_identity,both_identity;
331faca2338SShri Abhyankar   IS                 isicol;
332faca2338SShri Abhyankar   PetscErrorCode     ierr;
333faca2338SShri Abhyankar   const PetscInt     *r,*ic;
334faca2338SShri Abhyankar   PetscInt           i,*ai=a->i,*aj=a->j;
335faca2338SShri Abhyankar   PetscInt           *bi,*bj,*ajtmp;
336faca2338SShri Abhyankar   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
337faca2338SShri Abhyankar   PetscReal          f;
338faca2338SShri Abhyankar   PetscInt           nlnk,*lnk,k,**bi_ptr;
339faca2338SShri Abhyankar   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
340faca2338SShri Abhyankar   PetscBT            lnkbt;
341faca2338SShri Abhyankar 
342faca2338SShri Abhyankar   PetscFunctionBegin;
343*e32f2f54SBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"matrix must be square");
344faca2338SShri Abhyankar   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
345faca2338SShri Abhyankar   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
346faca2338SShri Abhyankar   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
347faca2338SShri Abhyankar 
348bc74c81fSJed Brown   /* get new row and diagonal pointers, must be allocated separately because they will be given to the Mat_SeqAIJ and freed separately */
349bc74c81fSJed Brown   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
350bc74c81fSJed Brown   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
351bc74c81fSJed Brown 
352a6df321fSShri Abhyankar   bi[0] = bdiag[0] = 0;
353c7272abeSHong Zhang 
354c7272abeSHong Zhang   /* linked list for storing column indices of the active row */
355c7272abeSHong Zhang   nlnk = n + 1;
356c7272abeSHong Zhang   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
357c7272abeSHong Zhang 
358c7272abeSHong Zhang   ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr);
359c7272abeSHong Zhang 
360c7272abeSHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
361c7272abeSHong Zhang   f = info->fill;
362c7272abeSHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
363c7272abeSHong Zhang   current_space = free_space;
36483287d42SBarry Smith 
36583287d42SBarry Smith   for (i=0; i<n; i++) {
366c7272abeSHong Zhang     /* copy previous fill into linked list */
36783287d42SBarry Smith     nzi = 0;
368c7272abeSHong Zhang     nnz = ai[r[i]+1] - ai[r[i]];
369*e32f2f54SBarry 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);
370c7272abeSHong Zhang     ajtmp = aj + ai[r[i]];
371c7272abeSHong Zhang     ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
372c7272abeSHong Zhang     nzi += nlnk;
373c7272abeSHong Zhang 
374c7272abeSHong Zhang     /* add pivot rows into linked list */
375c7272abeSHong Zhang     row = lnk[n];
376c7272abeSHong Zhang     while (row < i) {
377c7272abeSHong Zhang       nzbd    = bdiag[row] - bi[row] + 1; /* num of entries in the row with column index <= row */
378c7272abeSHong Zhang       ajtmp   = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */
379c7272abeSHong Zhang       ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
380c7272abeSHong Zhang       nzi += nlnk;
381c7272abeSHong Zhang       row  = lnk[row];
38283287d42SBarry Smith     }
383c7272abeSHong Zhang     bi[i+1] = bi[i] + nzi;
384c7272abeSHong Zhang     im[i]   = nzi;
385c7272abeSHong Zhang 
386c7272abeSHong Zhang     /* mark bdiag */
387c7272abeSHong Zhang     nzbd = 0;
388c7272abeSHong Zhang     nnz  = nzi;
389c7272abeSHong Zhang     k    = lnk[n];
390c7272abeSHong Zhang     while (nnz-- && k < i){
391c7272abeSHong Zhang       nzbd++;
392c7272abeSHong Zhang       k = lnk[k];
393c7272abeSHong Zhang     }
394c7272abeSHong Zhang     bdiag[i] = bi[i] + nzbd;
395c7272abeSHong Zhang 
396c7272abeSHong Zhang     /* if free space is not available, make more free space */
397c7272abeSHong Zhang     if (current_space->local_remaining<nzi) {
398c7272abeSHong Zhang       nnz = (n - i)*nzi; /* estimated and max additional space needed */
399c7272abeSHong Zhang       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
400c7272abeSHong Zhang       reallocs++;
40183287d42SBarry Smith     }
40283287d42SBarry Smith 
403c7272abeSHong Zhang     /* copy data into free space, then initialize lnk */
404c7272abeSHong Zhang     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
405c7272abeSHong Zhang     bi_ptr[i] = current_space->array;
406c7272abeSHong Zhang     current_space->array           += nzi;
407c7272abeSHong Zhang     current_space->local_used      += nzi;
408c7272abeSHong Zhang     current_space->local_remaining -= nzi;
409c7272abeSHong Zhang   }
4106cf91177SBarry Smith #if defined(PETSC_USE_INFO)
41183287d42SBarry Smith   if (ai[n] != 0) {
412c7272abeSHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
413ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
414ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
415ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
416ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
41783287d42SBarry Smith   } else {
418c7272abeSHong Zhang     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
41983287d42SBarry Smith   }
42063ba0a88SBarry Smith #endif
42183287d42SBarry Smith 
42283287d42SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
42383287d42SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
42483287d42SBarry Smith 
425c7272abeSHong Zhang   /* destroy list of free space and other temporary array(s) */
426c7272abeSHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
427c7272abeSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
428c7272abeSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
429c7272abeSHong Zhang   ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr);
43083287d42SBarry Smith 
43183287d42SBarry Smith   /* put together the new matrix */
432719d5645SBarry Smith   ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
433719d5645SBarry Smith   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
434719d5645SBarry Smith   b    = (Mat_SeqBAIJ*)(B)->data;
435e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
436e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
437c7272abeSHong Zhang   b->singlemalloc = PETSC_FALSE;
438c7272abeSHong Zhang   ierr          = PetscMalloc((bi[n]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
439c7272abeSHong Zhang   b->j          = bj;
440c7272abeSHong Zhang   b->i          = bi;
441c7272abeSHong Zhang   b->diag       = bdiag;
4427f53bb6cSHong Zhang   b->free_diag  = PETSC_TRUE;
44383287d42SBarry Smith   b->ilen       = 0;
44483287d42SBarry Smith   b->imax       = 0;
44583287d42SBarry Smith   b->row        = isrow;
44683287d42SBarry Smith   b->col        = iscol;
447d3d32019SBarry Smith   b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
44883287d42SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
44983287d42SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
45083287d42SBarry Smith   b->icol       = isicol;
45187828ca2SBarry Smith   ierr = PetscMalloc((bs*n+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
452c7272abeSHong Zhang   ierr = PetscLogObjectMemory(B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)*bs2));CHKERRQ(ierr);
45383287d42SBarry Smith 
454c7272abeSHong Zhang   b->maxnz = b->nz = bi[n] ;
455d5f3da31SBarry Smith   (B)->factortype            =  MAT_FACTOR_LU;
456719d5645SBarry Smith   (B)->info.factor_mallocs   = reallocs;
457719d5645SBarry Smith   (B)->info.fill_ratio_given = f;
458c7272abeSHong Zhang 
45983287d42SBarry Smith   if (ai[n] != 0) {
460c7272abeSHong Zhang     (B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
46183287d42SBarry Smith   } else {
462719d5645SBarry Smith     (B)->info.fill_ratio_needed = 0.0;
46383287d42SBarry Smith   }
464c7272abeSHong Zhang 
465db4efbfdSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
466db4efbfdSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
4677d18ce8fSMatthew Knepley   both_identity = (PetscTruth) (row_identity && col_identity);
4688b1456e3SHong Zhang   ierr = MatSeqBAIJSetNumericFactorization_inplace(B,both_identity);CHKERRQ(ierr);
46983287d42SBarry Smith   PetscFunctionReturn(0);
47083287d42SBarry Smith  }
471db4efbfdSBarry Smith 
472