xref: /petsc/src/mat/impls/baij/seq/baijfact3.c (revision 96e086a24b9cdc04c368e2851f46791dcaa837d7)
1be1d678aSKris Buschelman 
283287d42SBarry Smith /*
383287d42SBarry Smith     Factorization code for BAIJ format.
483287d42SBarry Smith */
5c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
6af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h>
783287d42SBarry Smith 
8db4efbfdSBarry Smith /*
9db4efbfdSBarry Smith    This is used to set the numeric factorization for both LU and ILU symbolic factorization
10db4efbfdSBarry Smith */
11ace3abfcSBarry Smith PetscErrorCode MatSeqBAIJSetNumericFactorization(Mat fact,PetscBool natural)
12ae3d28f0SHong Zhang {
13ae3d28f0SHong Zhang   PetscFunctionBegin;
146929473cSShri Abhyankar   if (natural) {
156929473cSShri Abhyankar     switch (fact->rmap->bs) {
16048b5e81SShri Abhyankar     case 1:
17048b5e81SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
18048b5e81SShri Abhyankar       break;
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;
37*96e086a2SDaniel Kokron     case 9:
38*96e086a2SDaniel Kokron       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_9_NaturalOrdering;
39*96e086a2SDaniel Kokron       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     }
476ba06ab7SHong 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 
78ace3abfcSBarry Smith PetscErrorCode MatSeqBAIJSetNumericFactorization_inplace(Mat inA,PetscBool 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:
93ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
94db4efbfdSBarry Smith       {
95ace3abfcSBarry Smith         PetscBool      sse_enabled_local;
96db4efbfdSBarry Smith         PetscErrorCode ierr;
970298fd71SBarry Smith         ierr = PetscSSEIsEnabled(inA->comm,&sse_enabled_local,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;
10326fbe8dcSKarl Rupp             for (i=0; i<nz; i++) aj[i] = (unsigned short)AJ[i];
10426fbe8dcSKarl Rupp 
105db4efbfdSBarry Smith             inA->ops->setunfactored   = MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj;
106db4efbfdSBarry Smith             inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE_usj;
10726fbe8dcSKarl Rupp 
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;
11626fbe8dcSKarl Rupp 
117db4efbfdSBarry Smith             ierr = PetscInfo(inA,"Using special SSE, in-place natural ordering, int j index factor BS=4\n");CHKERRQ(ierr);
118db4efbfdSBarry Smith           }
119db4efbfdSBarry Smith #  else
120db4efbfdSBarry Smith           /* This should never be reached.  If so, problem in PetscSSEIsEnabled. */
121e32f2f54SBarry Smith           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSE Hardware unavailable");
122db4efbfdSBarry Smith #  endif
123db4efbfdSBarry Smith         } else {
12406e38f1dSHong Zhang           inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_inplace;
125db4efbfdSBarry Smith         }
126db4efbfdSBarry Smith       }
127db4efbfdSBarry Smith #else
12806e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_inplace;
129db4efbfdSBarry Smith #endif
130db4efbfdSBarry Smith       break;
131db4efbfdSBarry Smith     case 5:
13206e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering_inplace;
133db4efbfdSBarry Smith       break;
134db4efbfdSBarry Smith     case 6:
13506e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering_inplace;
136db4efbfdSBarry Smith       break;
137db4efbfdSBarry Smith     case 7:
13806e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering_inplace;
139db4efbfdSBarry Smith       break;
140db4efbfdSBarry Smith     default:
14106e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N_inplace;
142db4efbfdSBarry Smith       break;
143db4efbfdSBarry Smith     }
144db4efbfdSBarry Smith   } else {
145db4efbfdSBarry Smith     switch (inA->rmap->bs) {
146db4efbfdSBarry Smith     case 1:
14706e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1_inplace;
148db4efbfdSBarry Smith       break;
149db4efbfdSBarry Smith     case 2:
15006e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_inplace;
151db4efbfdSBarry Smith       break;
152db4efbfdSBarry Smith     case 3:
15306e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_inplace;
154db4efbfdSBarry Smith       break;
155db4efbfdSBarry Smith     case 4:
15606e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_inplace;
157db4efbfdSBarry Smith       break;
158db4efbfdSBarry Smith     case 5:
15906e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_inplace;
160db4efbfdSBarry Smith       break;
161db4efbfdSBarry Smith     case 6:
16206e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_inplace;
163db4efbfdSBarry Smith       break;
164db4efbfdSBarry Smith     case 7:
16506e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_inplace;
166db4efbfdSBarry Smith       break;
167db4efbfdSBarry Smith     default:
16806e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N_inplace;
169db4efbfdSBarry Smith       break;
170db4efbfdSBarry Smith     }
171db4efbfdSBarry Smith   }
172db4efbfdSBarry Smith   PetscFunctionReturn(0);
173db4efbfdSBarry Smith }
174db4efbfdSBarry Smith 
17583287d42SBarry Smith /*
17683287d42SBarry Smith     The symbolic factorization code is identical to that for AIJ format,
17783287d42SBarry Smith   except for very small changes since this is now a SeqBAIJ datastructure.
17883287d42SBarry Smith   NOT good code reuse.
17983287d42SBarry Smith */
180c6db04a5SJed Brown #include <petscbt.h>
181c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
182c7272abeSHong Zhang 
1834dd39f65SShri Abhyankar PetscErrorCode MatLUFactorSymbolic_SeqBAIJ(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
18483287d42SBarry Smith {
18583287d42SBarry Smith   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data,*b;
186c7272abeSHong Zhang   PetscInt           n  =a->mbs,bs = A->rmap->bs,bs2=a->bs2;
187ace3abfcSBarry Smith   PetscBool          row_identity,col_identity,both_identity;
18883287d42SBarry Smith   IS                 isicol;
1896849ba73SBarry Smith   PetscErrorCode     ierr;
190c7272abeSHong Zhang   const PetscInt     *r,*ic;
191c7272abeSHong Zhang   PetscInt           i,*ai=a->i,*aj=a->j;
192c7272abeSHong Zhang   PetscInt           *bi,*bj,*ajtmp;
193c7272abeSHong Zhang   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
194c7272abeSHong Zhang   PetscReal          f;
195c7272abeSHong Zhang   PetscInt           nlnk,*lnk,k,**bi_ptr;
1960298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
197c7272abeSHong Zhang   PetscBT            lnkbt;
198ece542f9SHong Zhang   PetscBool          missing;
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");
202ece542f9SHong Zhang   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
203ece542f9SHong Zhang   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
204ece542f9SHong Zhang 
2056ba06ab7SHong Zhang   if (bs>1) {  /* check shifttype */
206f23aa3ddSBarry Smith     if (info->shifttype == MAT_SHIFT_NONZERO || info->shifttype == MAT_SHIFT_POSITIVE_DEFINITE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only MAT_SHIFT_NONE and MAT_SHIFT_INBLOCKS are supported for BAIJ matrix");
2076ba06ab7SHong Zhang   }
2086ba06ab7SHong Zhang 
20983287d42SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
21083287d42SBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
21183287d42SBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
21283287d42SBarry Smith 
213bc74c81fSJed Brown   /* get new row and diagonal pointers, must be allocated separately because they will be given to the Mat_SeqAIJ and freed separately */
214854ce69bSBarry Smith   ierr  = PetscMalloc1(n+1,&bi);CHKERRQ(ierr);
215854ce69bSBarry Smith   ierr  = PetscMalloc1(n+1,&bdiag);CHKERRQ(ierr);
216a6df321fSShri Abhyankar   bi[0] = bdiag[0] = 0;
217b2b2dd24SShri Abhyankar 
218b2b2dd24SShri Abhyankar   /* linked list for storing column indices of the active row */
219b2b2dd24SShri Abhyankar   nlnk = n + 1;
220b2b2dd24SShri Abhyankar   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
221b2b2dd24SShri Abhyankar 
222dcca6d9dSJed Brown   ierr = PetscMalloc2(n+1,&bi_ptr,n+1,&im);CHKERRQ(ierr);
223b2b2dd24SShri Abhyankar 
224b2b2dd24SShri Abhyankar   /* initial FreeSpace size is f*(ai[n]+1) */
225b2b2dd24SShri Abhyankar   f    = info->fill;
226f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space);CHKERRQ(ierr);
22726fbe8dcSKarl Rupp 
228b2b2dd24SShri Abhyankar   current_space = free_space;
229b2b2dd24SShri Abhyankar 
230b2b2dd24SShri Abhyankar   for (i=0; i<n; i++) {
231b2b2dd24SShri Abhyankar     /* copy previous fill into linked list */
232b2b2dd24SShri Abhyankar     nzi = 0;
233b2b2dd24SShri Abhyankar     nnz = ai[r[i]+1] - ai[r[i]];
234b2b2dd24SShri Abhyankar     ajtmp = aj + ai[r[i]];
235b2b2dd24SShri Abhyankar     ierr  = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
236b2b2dd24SShri Abhyankar     nzi  += nlnk;
237b2b2dd24SShri Abhyankar 
238b2b2dd24SShri Abhyankar     /* add pivot rows into linked list */
239b2b2dd24SShri Abhyankar     row = lnk[n];
240b2b2dd24SShri Abhyankar     while (row < i) {
241b2b2dd24SShri Abhyankar       nzbd  = bdiag[row] + 1;   /* num of entries in the row with column index <= row */
242b2b2dd24SShri Abhyankar       ajtmp = bi_ptr[row] + nzbd;   /* points to the entry next to the diagonal */
243b2b2dd24SShri Abhyankar       ierr  = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
244b2b2dd24SShri Abhyankar       nzi  += nlnk;
245b2b2dd24SShri Abhyankar       row   = lnk[row];
246b2b2dd24SShri Abhyankar     }
247b2b2dd24SShri Abhyankar     bi[i+1] = bi[i] + nzi;
248b2b2dd24SShri Abhyankar     im[i]   = nzi;
249b2b2dd24SShri Abhyankar 
250b2b2dd24SShri Abhyankar     /* mark bdiag */
251b2b2dd24SShri Abhyankar     nzbd = 0;
252b2b2dd24SShri Abhyankar     nnz  = nzi;
253b2b2dd24SShri Abhyankar     k    = lnk[n];
254b2b2dd24SShri Abhyankar     while (nnz-- && k < i) {
255b2b2dd24SShri Abhyankar       nzbd++;
256b2b2dd24SShri Abhyankar       k = lnk[k];
257b2b2dd24SShri Abhyankar     }
2582ce24eb6SHong Zhang     bdiag[i] = nzbd; /* note : bdaig[i] = nnzL as input for PetscFreeSpaceContiguous_LU() */
259b2b2dd24SShri Abhyankar 
260b2b2dd24SShri Abhyankar     /* if free space is not available, make more free space */
261b2b2dd24SShri Abhyankar     if (current_space->local_remaining<nzi) {
262f91af8c7SBarry Smith       nnz  = PetscIntMultTruncate(2,PetscIntMultTruncate(n - i,nzi)); /* estimated and max additional space needed */
263b2b2dd24SShri Abhyankar       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
264b2b2dd24SShri Abhyankar       reallocs++;
265b2b2dd24SShri Abhyankar     }
266b2b2dd24SShri Abhyankar 
267b2b2dd24SShri Abhyankar     /* copy data into free space, then initialize lnk */
268b2b2dd24SShri Abhyankar     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
26926fbe8dcSKarl Rupp 
270b2b2dd24SShri Abhyankar     bi_ptr[i]                       = current_space->array;
271b2b2dd24SShri Abhyankar     current_space->array           += nzi;
272b2b2dd24SShri Abhyankar     current_space->local_used      += nzi;
273b2b2dd24SShri Abhyankar     current_space->local_remaining -= nzi;
274b2b2dd24SShri Abhyankar   }
275b2b2dd24SShri Abhyankar 
276b2b2dd24SShri Abhyankar   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
277b2b2dd24SShri Abhyankar   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
278b2b2dd24SShri Abhyankar 
2799263d837SHong Zhang   /* copy free_space into bj and free free_space; set bi, bj, bdiag in new datastructure; */
280854ce69bSBarry Smith   ierr = PetscMalloc1(bi[n]+1,&bj);CHKERRQ(ierr);
2812ce24eb6SHong Zhang   ierr = PetscFreeSpaceContiguous_LU(&free_space,bj,n,bi,bdiag);CHKERRQ(ierr);
282b2b2dd24SShri Abhyankar   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
283b2b2dd24SShri Abhyankar   ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr);
284b2b2dd24SShri Abhyankar 
285b2b2dd24SShri Abhyankar   /* put together the new matrix */
286367daffbSBarry Smith   ierr = MatSeqBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
2873bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)B,(PetscObject)isicol);CHKERRQ(ierr);
288b2b2dd24SShri Abhyankar   b    = (Mat_SeqBAIJ*)(B)->data;
28926fbe8dcSKarl Rupp 
290b2b2dd24SShri Abhyankar   b->free_a       = PETSC_TRUE;
291b2b2dd24SShri Abhyankar   b->free_ij      = PETSC_TRUE;
292b2b2dd24SShri Abhyankar   b->singlemalloc = PETSC_FALSE;
29326fbe8dcSKarl Rupp 
294854ce69bSBarry Smith   ierr             = PetscMalloc1((bdiag[0]+1)*bs2,&b->a);CHKERRQ(ierr);
295b2b2dd24SShri Abhyankar   b->j             = bj;
296b2b2dd24SShri Abhyankar   b->i             = bi;
297b2b2dd24SShri Abhyankar   b->diag          = bdiag;
298b2b2dd24SShri Abhyankar   b->free_diag     = PETSC_TRUE;
299b2b2dd24SShri Abhyankar   b->ilen          = 0;
300b2b2dd24SShri Abhyankar   b->imax          = 0;
301b2b2dd24SShri Abhyankar   b->row           = isrow;
302b2b2dd24SShri Abhyankar   b->col           = iscol;
303b2b2dd24SShri Abhyankar   b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
30426fbe8dcSKarl Rupp 
305b2b2dd24SShri Abhyankar   ierr    = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
306b2b2dd24SShri Abhyankar   ierr    = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
307b2b2dd24SShri Abhyankar   b->icol = isicol;
308854ce69bSBarry Smith   ierr    = PetscMalloc1(bs*n+bs,&b->solve_work);CHKERRQ(ierr);
3093bb1ff40SBarry Smith   ierr    = PetscLogObjectMemory((PetscObject)B,(bdiag[0]+1)*(sizeof(PetscInt)+sizeof(PetscScalar)*bs2));CHKERRQ(ierr);
310b2b2dd24SShri Abhyankar 
311b2b2dd24SShri Abhyankar   b->maxnz = b->nz = bdiag[0]+1;
31226fbe8dcSKarl Rupp 
313d5f3da31SBarry Smith   B->factortype            =  MAT_FACTOR_LU;
314b2b2dd24SShri Abhyankar   B->info.factor_mallocs   = reallocs;
315b2b2dd24SShri Abhyankar   B->info.fill_ratio_given = f;
316b2b2dd24SShri Abhyankar 
317b2b2dd24SShri Abhyankar   if (ai[n] != 0) {
318b2b2dd24SShri Abhyankar     B->info.fill_ratio_needed = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]);
319b2b2dd24SShri Abhyankar   } else {
320b2b2dd24SShri Abhyankar     B->info.fill_ratio_needed = 0.0;
321b2b2dd24SShri Abhyankar   }
3229263d837SHong Zhang #if defined(PETSC_USE_INFO)
3239263d837SHong Zhang   if (ai[n] != 0) {
3249263d837SHong Zhang     PetscReal af = B->info.fill_ratio_needed;
32557622a8eSBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);CHKERRQ(ierr);
32657622a8eSBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
32757622a8eSBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g);\n",(double)af);CHKERRQ(ierr);
3289263d837SHong Zhang     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
3299263d837SHong Zhang   } else {
3309263d837SHong Zhang     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
3319263d837SHong Zhang   }
3329263d837SHong Zhang #endif
333b2b2dd24SShri Abhyankar 
334b2b2dd24SShri Abhyankar   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
335b2b2dd24SShri Abhyankar   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
33626fbe8dcSKarl Rupp 
337ace3abfcSBarry Smith   both_identity = (PetscBool) (row_identity && col_identity);
33826fbe8dcSKarl Rupp 
3394dd39f65SShri Abhyankar   ierr = MatSeqBAIJSetNumericFactorization(B,both_identity);CHKERRQ(ierr);
340b2b2dd24SShri Abhyankar   PetscFunctionReturn(0);
341b2b2dd24SShri Abhyankar }
342b2b2dd24SShri Abhyankar 
34306e38f1dSHong Zhang PetscErrorCode MatLUFactorSymbolic_SeqBAIJ_inplace(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
344faca2338SShri Abhyankar {
345faca2338SShri Abhyankar   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data,*b;
346faca2338SShri Abhyankar   PetscInt           n  =a->mbs,bs = A->rmap->bs,bs2=a->bs2;
347ace3abfcSBarry Smith   PetscBool          row_identity,col_identity,both_identity;
348faca2338SShri Abhyankar   IS                 isicol;
349faca2338SShri Abhyankar   PetscErrorCode     ierr;
350faca2338SShri Abhyankar   const PetscInt     *r,*ic;
351faca2338SShri Abhyankar   PetscInt           i,*ai=a->i,*aj=a->j;
352faca2338SShri Abhyankar   PetscInt           *bi,*bj,*ajtmp;
353faca2338SShri Abhyankar   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
354faca2338SShri Abhyankar   PetscReal          f;
355faca2338SShri Abhyankar   PetscInt           nlnk,*lnk,k,**bi_ptr;
3560298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
357faca2338SShri Abhyankar   PetscBT            lnkbt;
358ece542f9SHong Zhang   PetscBool          missing;
359faca2338SShri Abhyankar 
360faca2338SShri Abhyankar   PetscFunctionBegin;
361e32f2f54SBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"matrix must be square");
362ece542f9SHong Zhang   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
363ece542f9SHong Zhang   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
364ece542f9SHong Zhang 
365faca2338SShri Abhyankar   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
366faca2338SShri Abhyankar   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
367faca2338SShri Abhyankar   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
368faca2338SShri Abhyankar 
369bc74c81fSJed Brown   /* get new row and diagonal pointers, must be allocated separately because they will be given to the Mat_SeqAIJ and freed separately */
370854ce69bSBarry Smith   ierr = PetscMalloc1(n+1,&bi);CHKERRQ(ierr);
371854ce69bSBarry Smith   ierr = PetscMalloc1(n+1,&bdiag);CHKERRQ(ierr);
372bc74c81fSJed Brown 
373a6df321fSShri Abhyankar   bi[0] = bdiag[0] = 0;
374c7272abeSHong Zhang 
375c7272abeSHong Zhang   /* linked list for storing column indices of the active row */
376c7272abeSHong Zhang   nlnk = n + 1;
377c7272abeSHong Zhang   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
378c7272abeSHong Zhang 
379dcca6d9dSJed Brown   ierr = PetscMalloc2(n+1,&bi_ptr,n+1,&im);CHKERRQ(ierr);
380c7272abeSHong Zhang 
381c7272abeSHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
382c7272abeSHong Zhang   f             = info->fill;
383f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space);CHKERRQ(ierr);
384c7272abeSHong Zhang   current_space = free_space;
38583287d42SBarry Smith 
38683287d42SBarry Smith   for (i=0; i<n; i++) {
387c7272abeSHong Zhang     /* copy previous fill into linked list */
38883287d42SBarry Smith     nzi = 0;
389c7272abeSHong Zhang     nnz = ai[r[i]+1] - ai[r[i]];
390c7272abeSHong Zhang     ajtmp = aj + ai[r[i]];
391c7272abeSHong Zhang     ierr  = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
392c7272abeSHong Zhang     nzi  += nlnk;
393c7272abeSHong Zhang 
394c7272abeSHong Zhang     /* add pivot rows into linked list */
395c7272abeSHong Zhang     row = lnk[n];
396c7272abeSHong Zhang     while (row < i) {
397c7272abeSHong Zhang       nzbd  = bdiag[row] - bi[row] + 1;   /* num of entries in the row with column index <= row */
398c7272abeSHong Zhang       ajtmp = bi_ptr[row] + nzbd;   /* points to the entry next to the diagonal */
399c7272abeSHong Zhang       ierr  = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
400c7272abeSHong Zhang       nzi  += nlnk;
401c7272abeSHong Zhang       row   = lnk[row];
40283287d42SBarry Smith     }
403c7272abeSHong Zhang     bi[i+1] = bi[i] + nzi;
404c7272abeSHong Zhang     im[i]   = nzi;
405c7272abeSHong Zhang 
406c7272abeSHong Zhang     /* mark bdiag */
407c7272abeSHong Zhang     nzbd = 0;
408c7272abeSHong Zhang     nnz  = nzi;
409c7272abeSHong Zhang     k    = lnk[n];
410c7272abeSHong Zhang     while (nnz-- && k < i) {
411c7272abeSHong Zhang       nzbd++;
412c7272abeSHong Zhang       k = lnk[k];
413c7272abeSHong Zhang     }
414c7272abeSHong Zhang     bdiag[i] = bi[i] + nzbd;
415c7272abeSHong Zhang 
416c7272abeSHong Zhang     /* if free space is not available, make more free space */
417c7272abeSHong Zhang     if (current_space->local_remaining<nzi) {
418f91af8c7SBarry Smith       nnz  = PetscIntMultTruncate(n - i,nzi); /* estimated and max additional space needed */
419c7272abeSHong Zhang       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
420c7272abeSHong Zhang       reallocs++;
42183287d42SBarry Smith     }
42283287d42SBarry Smith 
423c7272abeSHong Zhang     /* copy data into free space, then initialize lnk */
424c7272abeSHong Zhang     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
42526fbe8dcSKarl Rupp 
426c7272abeSHong Zhang     bi_ptr[i]                       = current_space->array;
427c7272abeSHong Zhang     current_space->array           += nzi;
428c7272abeSHong Zhang     current_space->local_used      += nzi;
429c7272abeSHong Zhang     current_space->local_remaining -= nzi;
430c7272abeSHong Zhang   }
4316cf91177SBarry Smith #if defined(PETSC_USE_INFO)
43283287d42SBarry Smith   if (ai[n] != 0) {
433c7272abeSHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
43457622a8eSBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);CHKERRQ(ierr);
43557622a8eSBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
43657622a8eSBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g);\n",(double)af);CHKERRQ(ierr);
437ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
43883287d42SBarry Smith   } else {
439c7272abeSHong Zhang     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
44083287d42SBarry Smith   }
44163ba0a88SBarry Smith #endif
44283287d42SBarry Smith 
44383287d42SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
44483287d42SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
44583287d42SBarry Smith 
446c7272abeSHong Zhang   /* destroy list of free space and other temporary array(s) */
447854ce69bSBarry Smith   ierr = PetscMalloc1(bi[n]+1,&bj);CHKERRQ(ierr);
448c7272abeSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
449c7272abeSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
450c7272abeSHong Zhang   ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr);
45183287d42SBarry Smith 
45283287d42SBarry Smith   /* put together the new matrix */
453367daffbSBarry Smith   ierr = MatSeqBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
4543bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)B,(PetscObject)isicol);CHKERRQ(ierr);
455719d5645SBarry Smith   b    = (Mat_SeqBAIJ*)(B)->data;
45626fbe8dcSKarl Rupp 
457e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
458e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
459c7272abeSHong Zhang   b->singlemalloc = PETSC_FALSE;
46026fbe8dcSKarl Rupp 
461854ce69bSBarry Smith   ierr             = PetscMalloc1((bi[n]+1)*bs2,&b->a);CHKERRQ(ierr);
462c7272abeSHong Zhang   b->j             = bj;
463c7272abeSHong Zhang   b->i             = bi;
464c7272abeSHong Zhang   b->diag          = bdiag;
4657f53bb6cSHong Zhang   b->free_diag     = PETSC_TRUE;
46683287d42SBarry Smith   b->ilen          = 0;
46783287d42SBarry Smith   b->imax          = 0;
46883287d42SBarry Smith   b->row           = isrow;
46983287d42SBarry Smith   b->col           = iscol;
470d3d32019SBarry Smith   b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
47126fbe8dcSKarl Rupp 
47283287d42SBarry Smith   ierr    = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
47383287d42SBarry Smith   ierr    = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
47483287d42SBarry Smith   b->icol = isicol;
47526fbe8dcSKarl Rupp 
476854ce69bSBarry Smith   ierr = PetscMalloc1(bs*n+bs,&b->solve_work);CHKERRQ(ierr);
4773bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)*bs2));CHKERRQ(ierr);
47883287d42SBarry Smith 
479c7272abeSHong Zhang   b->maxnz = b->nz = bi[n];
48026fbe8dcSKarl Rupp 
481d5f3da31SBarry Smith   (B)->factortype            =  MAT_FACTOR_LU;
482719d5645SBarry Smith   (B)->info.factor_mallocs   = reallocs;
483719d5645SBarry Smith   (B)->info.fill_ratio_given = f;
484c7272abeSHong Zhang 
48583287d42SBarry Smith   if (ai[n] != 0) {
486c7272abeSHong Zhang     (B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
48783287d42SBarry Smith   } else {
488719d5645SBarry Smith     (B)->info.fill_ratio_needed = 0.0;
48983287d42SBarry Smith   }
490c7272abeSHong Zhang 
491db4efbfdSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
492db4efbfdSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
49326fbe8dcSKarl Rupp 
494ace3abfcSBarry Smith   both_identity = (PetscBool) (row_identity && col_identity);
49526fbe8dcSKarl Rupp 
4968b1456e3SHong Zhang   ierr = MatSeqBAIJSetNumericFactorization_inplace(B,both_identity);CHKERRQ(ierr);
49783287d42SBarry Smith   PetscFunctionReturn(0);
49883287d42SBarry Smith }
499db4efbfdSBarry Smith 
500