xref: /petsc/src/mat/impls/baij/seq/baijfact3.c (revision 06e38f1d1dda9d815f70edc0143bb77eead08adb)
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__
10ae3d28f0SHong Zhang #define __FUNCT__ "MatSeqBAIJSetNumericFactorization_newdatastruct"
11db4efbfdSBarry Smith /*
12db4efbfdSBarry Smith    This is used to set the numeric factorization for both LU and ILU symbolic factorization
13db4efbfdSBarry Smith */
14ae3d28f0SHong Zhang PetscErrorCode MatSeqBAIJSetNumericFactorization_newdatastruct(Mat fact,PetscTruth natural)
15ae3d28f0SHong Zhang {
16ae3d28f0SHong Zhang   PetscFunctionBegin;
176929473cSShri Abhyankar   if(natural){
186929473cSShri Abhyankar     switch (fact->rmap->bs){
196929473cSShri Abhyankar     case 2:
206929473cSShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_newdatastruct;
216929473cSShri Abhyankar       break;
226929473cSShri Abhyankar     case 3:
236929473cSShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_newdatastruct;
246929473cSShri Abhyankar       break;
256929473cSShri Abhyankar     case 4:
266929473cSShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_newdatastruct;
276929473cSShri Abhyankar       break;
286929473cSShri Abhyankar     case 5:
296929473cSShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering_newdatastruct;
306929473cSShri Abhyankar       break;
316929473cSShri Abhyankar     case 6:
326929473cSShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering_newdatastruct;
336929473cSShri Abhyankar       break;
346929473cSShri Abhyankar     case 7:
356929473cSShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering_newdatastruct;
366929473cSShri Abhyankar       break;
376929473cSShri Abhyankar     default:
386929473cSShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N_newdatastruct;
396929473cSShri Abhyankar       break;
406929473cSShri Abhyankar     }
416929473cSShri Abhyankar   }
426929473cSShri Abhyankar   else{
43ae3d28f0SHong Zhang     switch (fact->rmap->bs){
44ae3d28f0SHong Zhang     case 2:
45ae3d28f0SHong Zhang       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_newdatastruct;
46ae3d28f0SHong Zhang       break;
47ae3d28f0SHong Zhang     case 3:
48ae3d28f0SHong Zhang       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_newdatastruct;
49ae3d28f0SHong Zhang       break;
50ae3d28f0SHong Zhang     case 4:
51ae3d28f0SHong Zhang       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_newdatastruct;
52ae3d28f0SHong Zhang       break;
53ae3d28f0SHong Zhang     case 5:
54ae3d28f0SHong Zhang       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_newdatastruct;
55ae3d28f0SHong Zhang       break;
56ae3d28f0SHong Zhang     case 6:
57ae3d28f0SHong Zhang       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_newdatastruct;
58ae3d28f0SHong Zhang       break;
59ae3d28f0SHong Zhang     case 7:
60ae3d28f0SHong Zhang       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_newdatastruct;
61ae3d28f0SHong Zhang       break;
62ae3d28f0SHong Zhang     default:
63ae3d28f0SHong Zhang       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N_newdatastruct;
64ae3d28f0SHong Zhang       break;
65ae3d28f0SHong Zhang     }
666929473cSShri Abhyankar   }
67ae3d28f0SHong Zhang   PetscFunctionReturn(0);
68ae3d28f0SHong Zhang }
69ae3d28f0SHong Zhang 
70db4efbfdSBarry Smith PetscErrorCode MatSeqBAIJSetNumericFactorization(Mat inA,PetscTruth natural)
71db4efbfdSBarry Smith {
72db4efbfdSBarry Smith   PetscFunctionBegin;
73db4efbfdSBarry Smith   if (natural) {
74db4efbfdSBarry Smith     switch (inA->rmap->bs) {
75db4efbfdSBarry Smith     case 1:
76*06e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1_inplace;
77db4efbfdSBarry Smith       break;
78db4efbfdSBarry Smith     case 2:
79*06e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace;
80db4efbfdSBarry Smith       break;
81db4efbfdSBarry Smith     case 3:
82*06e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace;
83db4efbfdSBarry Smith       break;
84db4efbfdSBarry Smith     case 4:
8565460251SBarry Smith #if defined(PETSC_USE_SCALAR_MAT_SINGLE)
86db4efbfdSBarry Smith       {
87db4efbfdSBarry Smith         PetscTruth  sse_enabled_local;
88db4efbfdSBarry Smith         PetscErrorCode ierr;
89db4efbfdSBarry Smith         ierr = PetscSSEIsEnabled(inA->comm,&sse_enabled_local,PETSC_NULL);CHKERRQ(ierr);
90db4efbfdSBarry Smith         if (sse_enabled_local) {
91db4efbfdSBarry Smith #  if defined(PETSC_HAVE_SSE)
92db4efbfdSBarry Smith           int i,*AJ=a->j,nz=a->nz,n=a->mbs;
93db4efbfdSBarry Smith           if (n==(unsigned short)n) {
94db4efbfdSBarry Smith             unsigned short *aj=(unsigned short *)AJ;
95db4efbfdSBarry Smith             for (i=0;i<nz;i++) {
96db4efbfdSBarry Smith               aj[i] = (unsigned short)AJ[i];
97db4efbfdSBarry Smith             }
98db4efbfdSBarry Smith             inA->ops->setunfactored   = MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj;
99db4efbfdSBarry Smith             inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE_usj;
100db4efbfdSBarry Smith             ierr = PetscInfo(inA,"Using special SSE, in-place natural ordering, ushort j index factor BS=4\n");CHKERRQ(ierr);
101db4efbfdSBarry Smith           } else {
102db4efbfdSBarry Smith         /* Scale the column indices for easier indexing in MatSolve. */
103db4efbfdSBarry Smith /*            for (i=0;i<nz;i++) { */
104db4efbfdSBarry Smith /*              AJ[i] = AJ[i]*4; */
105db4efbfdSBarry Smith /*            } */
106db4efbfdSBarry Smith             inA->ops->setunfactored   = MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE;
107db4efbfdSBarry Smith             inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE;
108db4efbfdSBarry Smith             ierr = PetscInfo(inA,"Using special SSE, in-place natural ordering, int j index factor BS=4\n");CHKERRQ(ierr);
109db4efbfdSBarry Smith           }
110db4efbfdSBarry Smith #  else
111db4efbfdSBarry Smith         /* This should never be reached.  If so, problem in PetscSSEIsEnabled. */
112db4efbfdSBarry Smith           SETERRQ(PETSC_ERR_SUP,"SSE Hardware unavailable");
113db4efbfdSBarry Smith #  endif
114db4efbfdSBarry Smith         } else {
115*06e38f1dSHong Zhang           inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_inplace;
116db4efbfdSBarry Smith         }
117db4efbfdSBarry Smith       }
118db4efbfdSBarry Smith #else
119*06e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_inplace;
120db4efbfdSBarry Smith #endif
121db4efbfdSBarry Smith       break;
122db4efbfdSBarry Smith     case 5:
123*06e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering_inplace;
124db4efbfdSBarry Smith       break;
125db4efbfdSBarry Smith     case 6:
126*06e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering_inplace;
127db4efbfdSBarry Smith       break;
128db4efbfdSBarry Smith     case 7:
129*06e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering_inplace;
130db4efbfdSBarry Smith       break;
131db4efbfdSBarry Smith     default:
132*06e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N_inplace;
133db4efbfdSBarry Smith       break;
134db4efbfdSBarry Smith     }
135db4efbfdSBarry Smith   } else {
136db4efbfdSBarry Smith     switch (inA->rmap->bs) {
137db4efbfdSBarry Smith     case 1:
138*06e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1_inplace;
139db4efbfdSBarry Smith       break;
140db4efbfdSBarry Smith     case 2:
141*06e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_inplace;
142db4efbfdSBarry Smith       break;
143db4efbfdSBarry Smith     case 3:
144*06e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_inplace;
145db4efbfdSBarry Smith       break;
146db4efbfdSBarry Smith     case 4:
147*06e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_inplace;
148db4efbfdSBarry Smith       break;
149db4efbfdSBarry Smith     case 5:
150*06e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_inplace;
151db4efbfdSBarry Smith       break;
152db4efbfdSBarry Smith     case 6:
153*06e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_inplace;
154db4efbfdSBarry Smith       break;
155db4efbfdSBarry Smith     case 7:
156*06e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_inplace;
157db4efbfdSBarry Smith       break;
158db4efbfdSBarry Smith     default:
159*06e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N_inplace;
160db4efbfdSBarry Smith       break;
161db4efbfdSBarry Smith     }
162db4efbfdSBarry Smith   }
163db4efbfdSBarry Smith   PetscFunctionReturn(0);
164db4efbfdSBarry Smith }
165db4efbfdSBarry Smith 
16683287d42SBarry Smith /*
16783287d42SBarry Smith     The symbolic factorization code is identical to that for AIJ format,
16883287d42SBarry Smith   except for very small changes since this is now a SeqBAIJ datastructure.
16983287d42SBarry Smith   NOT good code reuse.
17083287d42SBarry Smith */
171c7272abeSHong Zhang #include "petscbt.h"
172c7272abeSHong Zhang #include "../src/mat/utils/freespace.h"
173c7272abeSHong Zhang 
1744a2ae208SSatish Balay #undef __FUNCT__
175faca2338SShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_SeqBAIJ_newdatastruct"
176faca2338SShri Abhyankar PetscErrorCode MatLUFactorSymbolic_SeqBAIJ_newdatastruct(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
17783287d42SBarry Smith {
17883287d42SBarry Smith   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data,*b;
179c7272abeSHong Zhang   PetscInt           n=a->mbs,bs = A->rmap->bs,bs2=a->bs2;
180c7272abeSHong Zhang   PetscTruth         row_identity,col_identity,both_identity;
18183287d42SBarry Smith   IS                 isicol;
1826849ba73SBarry Smith   PetscErrorCode     ierr;
183c7272abeSHong Zhang   const PetscInt     *r,*ic;
184c7272abeSHong Zhang   PetscInt           i,*ai=a->i,*aj=a->j;
185c7272abeSHong Zhang   PetscInt           *bi,*bj,*ajtmp;
186c7272abeSHong Zhang   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
187c7272abeSHong Zhang   PetscReal          f;
188c7272abeSHong Zhang   PetscInt           nlnk,*lnk,k,**bi_ptr;
189c7272abeSHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
190c7272abeSHong Zhang   PetscBT            lnkbt;
191*06e38f1dSHong Zhang   PetscTruth         newdatastruct=PETSC_FALSE;
19283287d42SBarry Smith 
19383287d42SBarry Smith   PetscFunctionBegin;
194*06e38f1dSHong Zhang   ierr = PetscOptionsGetTruth(PETSC_NULL,"-lu_old",&newdatastruct,PETSC_NULL);CHKERRQ(ierr);
195*06e38f1dSHong Zhang   if(newdatastruct){
196*06e38f1dSHong Zhang     ierr = MatLUFactorSymbolic_SeqBAIJ_inplace(B,A,isrow,iscol,info);CHKERRQ(ierr);
197*06e38f1dSHong Zhang     PetscFunctionReturn(0);
198*06e38f1dSHong Zhang   }
199d0f46423SBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
20083287d42SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
20183287d42SBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
20283287d42SBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
20383287d42SBarry Smith 
204bc74c81fSJed Brown   /* get new row and diagonal pointers, must be allocated separately because they will be given to the Mat_SeqAIJ and freed separately */
205bc74c81fSJed Brown   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
206bc74c81fSJed Brown   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
207a6df321fSShri Abhyankar   bi[0] = bdiag[0] = 0;
208b2b2dd24SShri Abhyankar 
209b2b2dd24SShri Abhyankar   /* linked list for storing column indices of the active row */
210b2b2dd24SShri Abhyankar   nlnk = n + 1;
211b2b2dd24SShri Abhyankar   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
212b2b2dd24SShri Abhyankar 
213b2b2dd24SShri Abhyankar   ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr);
214b2b2dd24SShri Abhyankar 
215b2b2dd24SShri Abhyankar   /* initial FreeSpace size is f*(ai[n]+1) */
216b2b2dd24SShri Abhyankar   f = info->fill;
217b2b2dd24SShri Abhyankar   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
218b2b2dd24SShri Abhyankar   current_space = free_space;
219b2b2dd24SShri Abhyankar 
220b2b2dd24SShri Abhyankar   for (i=0; i<n; i++) {
221b2b2dd24SShri Abhyankar     /* copy previous fill into linked list */
222b2b2dd24SShri Abhyankar     nzi = 0;
223b2b2dd24SShri Abhyankar     nnz = ai[r[i]+1] - ai[r[i]];
224b2b2dd24SShri Abhyankar     if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
225b2b2dd24SShri Abhyankar     ajtmp = aj + ai[r[i]];
226b2b2dd24SShri Abhyankar     ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
227b2b2dd24SShri Abhyankar     nzi += nlnk;
228b2b2dd24SShri Abhyankar 
229b2b2dd24SShri Abhyankar     /* add pivot rows into linked list */
230b2b2dd24SShri Abhyankar     row = lnk[n];
231b2b2dd24SShri Abhyankar     while (row < i) {
232b2b2dd24SShri Abhyankar       nzbd    = bdiag[row] + 1; /* num of entries in the row with column index <= row */
233b2b2dd24SShri Abhyankar       ajtmp   = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */
234b2b2dd24SShri Abhyankar       ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
235b2b2dd24SShri Abhyankar       nzi += nlnk;
236b2b2dd24SShri Abhyankar       row  = lnk[row];
237b2b2dd24SShri Abhyankar     }
238b2b2dd24SShri Abhyankar     bi[i+1] = bi[i] + nzi;
239b2b2dd24SShri Abhyankar     im[i]   = nzi;
240b2b2dd24SShri Abhyankar 
241b2b2dd24SShri Abhyankar     /* mark bdiag */
242b2b2dd24SShri Abhyankar     nzbd = 0;
243b2b2dd24SShri Abhyankar     nnz  = nzi;
244b2b2dd24SShri Abhyankar     k    = lnk[n];
245b2b2dd24SShri Abhyankar     while (nnz-- && k < i){
246b2b2dd24SShri Abhyankar       nzbd++;
247b2b2dd24SShri Abhyankar       k = lnk[k];
248b2b2dd24SShri Abhyankar     }
249b2b2dd24SShri Abhyankar     bdiag[i] = nzbd; /* note : bdaig[i] = nnzL as input for PetscFreeSpaceContiguous_LU_v2() */
250b2b2dd24SShri Abhyankar 
251b2b2dd24SShri Abhyankar     /* if free space is not available, make more free space */
252b2b2dd24SShri Abhyankar     if (current_space->local_remaining<nzi) {
253b2b2dd24SShri Abhyankar       nnz = 2*(n - i)*nzi; /* estimated and max additional space needed */
254b2b2dd24SShri Abhyankar       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
255b2b2dd24SShri Abhyankar       reallocs++;
256b2b2dd24SShri Abhyankar     }
257b2b2dd24SShri Abhyankar 
258b2b2dd24SShri Abhyankar     /* copy data into free space, then initialize lnk */
259b2b2dd24SShri Abhyankar     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
260b2b2dd24SShri Abhyankar     bi_ptr[i] = current_space->array;
261b2b2dd24SShri Abhyankar     current_space->array           += nzi;
262b2b2dd24SShri Abhyankar     current_space->local_used      += nzi;
263b2b2dd24SShri Abhyankar     current_space->local_remaining -= nzi;
264b2b2dd24SShri Abhyankar   }
265b2b2dd24SShri Abhyankar #if defined(PETSC_USE_INFO)
266b2b2dd24SShri Abhyankar   if (ai[n] != 0) {
267b2b2dd24SShri Abhyankar     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
268b2b2dd24SShri Abhyankar     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
269b2b2dd24SShri Abhyankar     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
270b2b2dd24SShri Abhyankar     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
271b2b2dd24SShri Abhyankar     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
272b2b2dd24SShri Abhyankar   } else {
273b2b2dd24SShri Abhyankar     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
274b2b2dd24SShri Abhyankar   }
275b2b2dd24SShri Abhyankar #endif
276b2b2dd24SShri Abhyankar 
277b2b2dd24SShri Abhyankar   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
278b2b2dd24SShri Abhyankar   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
279b2b2dd24SShri Abhyankar 
280b2b2dd24SShri Abhyankar   /* destroy list of free space and other temporary array(s) */
281b2b2dd24SShri Abhyankar   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
282b2b2dd24SShri Abhyankar   ierr = PetscFreeSpaceContiguous_LU_v2(&free_space,bj,n,bi,bdiag);CHKERRQ(ierr);
283b2b2dd24SShri Abhyankar   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
284b2b2dd24SShri Abhyankar   ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr);
285b2b2dd24SShri Abhyankar 
286b2b2dd24SShri Abhyankar   /* put together the new matrix */
287b2b2dd24SShri Abhyankar   ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
288b2b2dd24SShri Abhyankar   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
289b2b2dd24SShri Abhyankar   b    = (Mat_SeqBAIJ*)(B)->data;
290b2b2dd24SShri Abhyankar   b->free_a       = PETSC_TRUE;
291b2b2dd24SShri Abhyankar   b->free_ij      = PETSC_TRUE;
292b2b2dd24SShri Abhyankar   b->singlemalloc = PETSC_FALSE;
293b2b2dd24SShri Abhyankar   ierr          = PetscMalloc((bdiag[0]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
294b2b2dd24SShri Abhyankar   b->j          = bj;
295b2b2dd24SShri Abhyankar   b->i          = bi;
296b2b2dd24SShri Abhyankar   b->diag       = bdiag;
297b2b2dd24SShri Abhyankar   b->free_diag  = PETSC_TRUE;
298b2b2dd24SShri Abhyankar   b->ilen       = 0;
299b2b2dd24SShri Abhyankar   b->imax       = 0;
300b2b2dd24SShri Abhyankar   b->row        = isrow;
301b2b2dd24SShri Abhyankar   b->col        = iscol;
302b2b2dd24SShri Abhyankar   b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
303b2b2dd24SShri Abhyankar   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
304b2b2dd24SShri Abhyankar   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
305b2b2dd24SShri Abhyankar   b->icol       = isicol;
306b2b2dd24SShri Abhyankar   ierr = PetscMalloc((bs*n+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
307b2b2dd24SShri Abhyankar   ierr = PetscLogObjectMemory(B,(bdiag[0]+1)*(sizeof(PetscInt)+sizeof(PetscScalar)*bs2));CHKERRQ(ierr);
308b2b2dd24SShri Abhyankar 
309b2b2dd24SShri Abhyankar   b->maxnz = b->nz = bdiag[0]+1;
310b2b2dd24SShri Abhyankar   B->factor                =  MAT_FACTOR_LU;
311b2b2dd24SShri Abhyankar   B->info.factor_mallocs   = reallocs;
312b2b2dd24SShri Abhyankar   B->info.fill_ratio_given = f;
313b2b2dd24SShri Abhyankar 
314b2b2dd24SShri Abhyankar   if (ai[n] != 0) {
315b2b2dd24SShri Abhyankar     B->info.fill_ratio_needed = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]);
316b2b2dd24SShri Abhyankar   } else {
317b2b2dd24SShri Abhyankar     B->info.fill_ratio_needed = 0.0;
318b2b2dd24SShri Abhyankar   }
319b2b2dd24SShri Abhyankar 
320b2b2dd24SShri Abhyankar   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
321b2b2dd24SShri Abhyankar   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
322b2b2dd24SShri Abhyankar   both_identity = (PetscTruth) (row_identity && col_identity);
3236ccb8380SShri Abhyankar   ierr = MatSeqBAIJSetNumericFactorization_newdatastruct(B,both_identity);CHKERRQ(ierr);
324b2b2dd24SShri Abhyankar   PetscFunctionReturn(0);
325b2b2dd24SShri Abhyankar  }
326b2b2dd24SShri Abhyankar 
327b2b2dd24SShri Abhyankar #undef __FUNCT__
328*06e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorSymbolic_SeqBAIJ_inplace"
329*06e38f1dSHong Zhang PetscErrorCode MatLUFactorSymbolic_SeqBAIJ_inplace(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
330faca2338SShri Abhyankar {
331faca2338SShri Abhyankar   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data,*b;
332faca2338SShri Abhyankar   PetscInt           n=a->mbs,bs = A->rmap->bs,bs2=a->bs2;
333faca2338SShri Abhyankar   PetscTruth         row_identity,col_identity,both_identity;
334faca2338SShri Abhyankar   IS                 isicol;
335faca2338SShri Abhyankar   PetscErrorCode     ierr;
336faca2338SShri Abhyankar   const PetscInt     *r,*ic;
337faca2338SShri Abhyankar   PetscInt           i,*ai=a->i,*aj=a->j;
338faca2338SShri Abhyankar   PetscInt           *bi,*bj,*ajtmp;
339faca2338SShri Abhyankar   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
340faca2338SShri Abhyankar   PetscReal          f;
341faca2338SShri Abhyankar   PetscInt           nlnk,*lnk,k,**bi_ptr;
342faca2338SShri Abhyankar   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
343faca2338SShri Abhyankar   PetscBT            lnkbt;
344faca2338SShri Abhyankar 
345faca2338SShri Abhyankar   PetscFunctionBegin;
346faca2338SShri Abhyankar   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
347faca2338SShri Abhyankar   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
348faca2338SShri Abhyankar   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
349faca2338SShri Abhyankar   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
350faca2338SShri Abhyankar 
351bc74c81fSJed Brown   /* get new row and diagonal pointers, must be allocated separately because they will be given to the Mat_SeqAIJ and freed separately */
352bc74c81fSJed Brown   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
353bc74c81fSJed Brown   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
354bc74c81fSJed Brown 
355a6df321fSShri Abhyankar   bi[0] = bdiag[0] = 0;
356c7272abeSHong Zhang 
357c7272abeSHong Zhang   /* linked list for storing column indices of the active row */
358c7272abeSHong Zhang   nlnk = n + 1;
359c7272abeSHong Zhang   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
360c7272abeSHong Zhang 
361c7272abeSHong Zhang   ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr);
362c7272abeSHong Zhang 
363c7272abeSHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
364c7272abeSHong Zhang   f = info->fill;
365c7272abeSHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
366c7272abeSHong Zhang   current_space = free_space;
36783287d42SBarry Smith 
36883287d42SBarry Smith   for (i=0; i<n; i++) {
369c7272abeSHong Zhang     /* copy previous fill into linked list */
37083287d42SBarry Smith     nzi = 0;
371c7272abeSHong Zhang     nnz = ai[r[i]+1] - ai[r[i]];
372c7272abeSHong Zhang     if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
373c7272abeSHong Zhang     ajtmp = aj + ai[r[i]];
374c7272abeSHong Zhang     ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
375c7272abeSHong Zhang     nzi += nlnk;
376c7272abeSHong Zhang 
377c7272abeSHong Zhang     /* add pivot rows into linked list */
378c7272abeSHong Zhang     row = lnk[n];
379c7272abeSHong Zhang     while (row < i) {
380c7272abeSHong Zhang       nzbd    = bdiag[row] - bi[row] + 1; /* num of entries in the row with column index <= row */
381c7272abeSHong Zhang       ajtmp   = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */
382c7272abeSHong Zhang       ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
383c7272abeSHong Zhang       nzi += nlnk;
384c7272abeSHong Zhang       row  = lnk[row];
38583287d42SBarry Smith     }
386c7272abeSHong Zhang     bi[i+1] = bi[i] + nzi;
387c7272abeSHong Zhang     im[i]   = nzi;
388c7272abeSHong Zhang 
389c7272abeSHong Zhang     /* mark bdiag */
390c7272abeSHong Zhang     nzbd = 0;
391c7272abeSHong Zhang     nnz  = nzi;
392c7272abeSHong Zhang     k    = lnk[n];
393c7272abeSHong Zhang     while (nnz-- && k < i){
394c7272abeSHong Zhang       nzbd++;
395c7272abeSHong Zhang       k = lnk[k];
396c7272abeSHong Zhang     }
397c7272abeSHong Zhang     bdiag[i] = bi[i] + nzbd;
398c7272abeSHong Zhang 
399c7272abeSHong Zhang     /* if free space is not available, make more free space */
400c7272abeSHong Zhang     if (current_space->local_remaining<nzi) {
401c7272abeSHong Zhang       nnz = (n - i)*nzi; /* estimated and max additional space needed */
402c7272abeSHong Zhang       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
403c7272abeSHong Zhang       reallocs++;
40483287d42SBarry Smith     }
40583287d42SBarry Smith 
406c7272abeSHong Zhang     /* copy data into free space, then initialize lnk */
407c7272abeSHong Zhang     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
408c7272abeSHong Zhang     bi_ptr[i] = current_space->array;
409c7272abeSHong Zhang     current_space->array           += nzi;
410c7272abeSHong Zhang     current_space->local_used      += nzi;
411c7272abeSHong Zhang     current_space->local_remaining -= nzi;
412c7272abeSHong Zhang   }
4136cf91177SBarry Smith #if defined(PETSC_USE_INFO)
41483287d42SBarry Smith   if (ai[n] != 0) {
415c7272abeSHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
416ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
417ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
418ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
419ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
42083287d42SBarry Smith   } else {
421c7272abeSHong Zhang     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
42283287d42SBarry Smith   }
42363ba0a88SBarry Smith #endif
42483287d42SBarry Smith 
42583287d42SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
42683287d42SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
42783287d42SBarry Smith 
428c7272abeSHong Zhang   /* destroy list of free space and other temporary array(s) */
429c7272abeSHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
430c7272abeSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
431c7272abeSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
432c7272abeSHong Zhang   ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr);
43383287d42SBarry Smith 
43483287d42SBarry Smith   /* put together the new matrix */
435719d5645SBarry Smith   ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
436719d5645SBarry Smith   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
437719d5645SBarry Smith   b    = (Mat_SeqBAIJ*)(B)->data;
438e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
439e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
440c7272abeSHong Zhang   b->singlemalloc = PETSC_FALSE;
441c7272abeSHong Zhang   ierr          = PetscMalloc((bi[n]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
442c7272abeSHong Zhang   b->j          = bj;
443c7272abeSHong Zhang   b->i          = bi;
444c7272abeSHong Zhang   b->diag       = bdiag;
4457f53bb6cSHong Zhang   b->free_diag  = PETSC_TRUE;
44683287d42SBarry Smith   b->ilen       = 0;
44783287d42SBarry Smith   b->imax       = 0;
44883287d42SBarry Smith   b->row        = isrow;
44983287d42SBarry Smith   b->col        = iscol;
450d3d32019SBarry Smith   b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
45183287d42SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
45283287d42SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
45383287d42SBarry Smith   b->icol       = isicol;
45487828ca2SBarry Smith   ierr = PetscMalloc((bs*n+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
455c7272abeSHong Zhang   ierr = PetscLogObjectMemory(B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)*bs2));CHKERRQ(ierr);
45683287d42SBarry Smith 
457c7272abeSHong Zhang   b->maxnz = b->nz = bi[n] ;
458c7272abeSHong Zhang   (B)->factor                =  MAT_FACTOR_LU;
459719d5645SBarry Smith   (B)->info.factor_mallocs   = reallocs;
460719d5645SBarry Smith   (B)->info.fill_ratio_given = f;
461c7272abeSHong Zhang 
46283287d42SBarry Smith   if (ai[n] != 0) {
463c7272abeSHong Zhang     (B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
46483287d42SBarry Smith   } else {
465719d5645SBarry Smith     (B)->info.fill_ratio_needed = 0.0;
46683287d42SBarry Smith   }
467c7272abeSHong Zhang 
468db4efbfdSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
469db4efbfdSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
4707d18ce8fSMatthew Knepley   both_identity = (PetscTruth) (row_identity && col_identity);
47141df41f0SMatthew Knepley   ierr = MatSeqBAIJSetNumericFactorization(B,both_identity);CHKERRQ(ierr);
47283287d42SBarry Smith   PetscFunctionReturn(0);
47383287d42SBarry Smith  }
474db4efbfdSBarry Smith 
475