xref: /petsc/src/mat/impls/baij/seq/baijfact3.c (revision ece542f93d9ecceaf66cc7a4435fb274c8ff4ac6)
1be1d678aSKris Buschelman 
283287d42SBarry Smith /*
383287d42SBarry Smith     Factorization code for BAIJ format.
483287d42SBarry Smith */
5c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
606873bf2SBarry Smith #include <petsc-private/kernels/blockinvert.h>
783287d42SBarry Smith 
8db4efbfdSBarry Smith #undef __FUNCT__
94dd39f65SShri Abhyankar #define __FUNCT__ "MatSeqBAIJSetNumericFactorization"
10db4efbfdSBarry Smith /*
11db4efbfdSBarry Smith    This is used to set the numeric factorization for both LU and ILU symbolic factorization
12db4efbfdSBarry Smith */
13ace3abfcSBarry Smith PetscErrorCode MatSeqBAIJSetNumericFactorization(Mat fact,PetscBool natural)
14ae3d28f0SHong Zhang {
15ae3d28f0SHong Zhang   PetscFunctionBegin;
166929473cSShri Abhyankar   if (natural) {
176929473cSShri Abhyankar     switch (fact->rmap->bs) {
18048b5e81SShri Abhyankar     case 1:
19048b5e81SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
20048b5e81SShri Abhyankar       break;
216929473cSShri Abhyankar     case 2:
224dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering;
236929473cSShri Abhyankar       break;
246929473cSShri Abhyankar     case 3:
254dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering;
266929473cSShri Abhyankar       break;
276929473cSShri Abhyankar     case 4:
284dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering;
296929473cSShri Abhyankar       break;
306929473cSShri Abhyankar     case 5:
314dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering;
326929473cSShri Abhyankar       break;
336929473cSShri Abhyankar     case 6:
344dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering;
356929473cSShri Abhyankar       break;
366929473cSShri Abhyankar     case 7:
374dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering;
386929473cSShri Abhyankar       break;
3929a97285SShri Abhyankar     case 15:
4029a97285SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_15_NaturalOrdering;
4129a97285SShri Abhyankar       break;
426929473cSShri Abhyankar     default:
434dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N;
446929473cSShri Abhyankar       break;
456929473cSShri Abhyankar     }
466ba06ab7SHong Zhang   } else {
47ae3d28f0SHong Zhang     switch (fact->rmap->bs) {
48048b5e81SShri Abhyankar     case 1:
49048b5e81SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
50048b5e81SShri Abhyankar       break;
51ae3d28f0SHong Zhang     case 2:
524dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
53ae3d28f0SHong Zhang       break;
54ae3d28f0SHong Zhang     case 3:
554dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
56ae3d28f0SHong Zhang       break;
57ae3d28f0SHong Zhang     case 4:
584dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
59ae3d28f0SHong Zhang       break;
60ae3d28f0SHong Zhang     case 5:
614dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
62ae3d28f0SHong Zhang       break;
63ae3d28f0SHong Zhang     case 6:
644dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6;
65ae3d28f0SHong Zhang       break;
66ae3d28f0SHong Zhang     case 7:
674dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7;
68ae3d28f0SHong Zhang       break;
69ae3d28f0SHong Zhang     default:
704dd39f65SShri Abhyankar       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N;
71ae3d28f0SHong Zhang       break;
72ae3d28f0SHong Zhang     }
736929473cSShri Abhyankar   }
74ae3d28f0SHong Zhang   PetscFunctionReturn(0);
75ae3d28f0SHong Zhang }
76ae3d28f0SHong Zhang 
771008731dSBarry Smith #undef __FUNCT__
781008731dSBarry Smith #define __FUNCT__ "MatSeqBAIJSetNumericFactorization_inplace"
79ace3abfcSBarry Smith PetscErrorCode MatSeqBAIJSetNumericFactorization_inplace(Mat inA,PetscBool natural)
80db4efbfdSBarry Smith {
81db4efbfdSBarry Smith   PetscFunctionBegin;
82db4efbfdSBarry Smith   if (natural) {
83db4efbfdSBarry Smith     switch (inA->rmap->bs) {
84db4efbfdSBarry Smith     case 1:
8506e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1_inplace;
86db4efbfdSBarry Smith       break;
87db4efbfdSBarry Smith     case 2:
8806e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace;
89db4efbfdSBarry Smith       break;
90db4efbfdSBarry Smith     case 3:
9106e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace;
92db4efbfdSBarry Smith       break;
93db4efbfdSBarry Smith     case 4:
94ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
95db4efbfdSBarry Smith       {
96ace3abfcSBarry Smith         PetscBool      sse_enabled_local;
97db4efbfdSBarry Smith         PetscErrorCode ierr;
980298fd71SBarry Smith         ierr = PetscSSEIsEnabled(inA->comm,&sse_enabled_local,NULL);CHKERRQ(ierr);
99db4efbfdSBarry Smith         if (sse_enabled_local) {
100db4efbfdSBarry Smith #  if defined(PETSC_HAVE_SSE)
101db4efbfdSBarry Smith           int i,*AJ=a->j,nz=a->nz,n=a->mbs;
102db4efbfdSBarry Smith           if (n==(unsigned short)n) {
103db4efbfdSBarry Smith             unsigned short *aj=(unsigned short*)AJ;
10426fbe8dcSKarl Rupp             for (i=0; i<nz; i++) aj[i] = (unsigned short)AJ[i];
10526fbe8dcSKarl Rupp 
106db4efbfdSBarry Smith             inA->ops->setunfactored   = MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj;
107db4efbfdSBarry Smith             inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE_usj;
10826fbe8dcSKarl Rupp 
109db4efbfdSBarry Smith             ierr = PetscInfo(inA,"Using special SSE, in-place natural ordering, ushort j index factor BS=4\n");CHKERRQ(ierr);
110db4efbfdSBarry Smith           } else {
111db4efbfdSBarry Smith             /* Scale the column indices for easier indexing in MatSolve. */
112db4efbfdSBarry Smith /*            for (i=0;i<nz;i++) { */
113db4efbfdSBarry Smith /*              AJ[i] = AJ[i]*4; */
114db4efbfdSBarry Smith /*            } */
115db4efbfdSBarry Smith             inA->ops->setunfactored   = MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE;
116db4efbfdSBarry Smith             inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE;
11726fbe8dcSKarl Rupp 
118db4efbfdSBarry Smith             ierr = PetscInfo(inA,"Using special SSE, in-place natural ordering, int j index factor BS=4\n");CHKERRQ(ierr);
119db4efbfdSBarry Smith           }
120db4efbfdSBarry Smith #  else
121db4efbfdSBarry Smith           /* This should never be reached.  If so, problem in PetscSSEIsEnabled. */
122e32f2f54SBarry Smith           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSE Hardware unavailable");
123db4efbfdSBarry Smith #  endif
124db4efbfdSBarry Smith         } else {
12506e38f1dSHong Zhang           inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_inplace;
126db4efbfdSBarry Smith         }
127db4efbfdSBarry Smith       }
128db4efbfdSBarry Smith #else
12906e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_inplace;
130db4efbfdSBarry Smith #endif
131db4efbfdSBarry Smith       break;
132db4efbfdSBarry Smith     case 5:
13306e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering_inplace;
134db4efbfdSBarry Smith       break;
135db4efbfdSBarry Smith     case 6:
13606e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering_inplace;
137db4efbfdSBarry Smith       break;
138db4efbfdSBarry Smith     case 7:
13906e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering_inplace;
140db4efbfdSBarry Smith       break;
141db4efbfdSBarry Smith     default:
14206e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N_inplace;
143db4efbfdSBarry Smith       break;
144db4efbfdSBarry Smith     }
145db4efbfdSBarry Smith   } else {
146db4efbfdSBarry Smith     switch (inA->rmap->bs) {
147db4efbfdSBarry Smith     case 1:
14806e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1_inplace;
149db4efbfdSBarry Smith       break;
150db4efbfdSBarry Smith     case 2:
15106e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_inplace;
152db4efbfdSBarry Smith       break;
153db4efbfdSBarry Smith     case 3:
15406e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_inplace;
155db4efbfdSBarry Smith       break;
156db4efbfdSBarry Smith     case 4:
15706e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_inplace;
158db4efbfdSBarry Smith       break;
159db4efbfdSBarry Smith     case 5:
16006e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_inplace;
161db4efbfdSBarry Smith       break;
162db4efbfdSBarry Smith     case 6:
16306e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_inplace;
164db4efbfdSBarry Smith       break;
165db4efbfdSBarry Smith     case 7:
16606e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_inplace;
167db4efbfdSBarry Smith       break;
168db4efbfdSBarry Smith     default:
16906e38f1dSHong Zhang       inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N_inplace;
170db4efbfdSBarry Smith       break;
171db4efbfdSBarry Smith     }
172db4efbfdSBarry Smith   }
173db4efbfdSBarry Smith   PetscFunctionReturn(0);
174db4efbfdSBarry Smith }
175db4efbfdSBarry Smith 
17683287d42SBarry Smith /*
17783287d42SBarry Smith     The symbolic factorization code is identical to that for AIJ format,
17883287d42SBarry Smith   except for very small changes since this is now a SeqBAIJ datastructure.
17983287d42SBarry Smith   NOT good code reuse.
18083287d42SBarry Smith */
181c6db04a5SJed Brown #include <petscbt.h>
182c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
183c7272abeSHong Zhang 
1844a2ae208SSatish Balay #undef __FUNCT__
1854dd39f65SShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_SeqBAIJ"
1864dd39f65SShri Abhyankar PetscErrorCode MatLUFactorSymbolic_SeqBAIJ(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
18783287d42SBarry Smith {
18883287d42SBarry Smith   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data,*b;
189c7272abeSHong Zhang   PetscInt           n  =a->mbs,bs = A->rmap->bs,bs2=a->bs2;
190ace3abfcSBarry Smith   PetscBool          row_identity,col_identity,both_identity;
19183287d42SBarry Smith   IS                 isicol;
1926849ba73SBarry Smith   PetscErrorCode     ierr;
193c7272abeSHong Zhang   const PetscInt     *r,*ic;
194c7272abeSHong Zhang   PetscInt           i,*ai=a->i,*aj=a->j;
195c7272abeSHong Zhang   PetscInt           *bi,*bj,*ajtmp;
196c7272abeSHong Zhang   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
197c7272abeSHong Zhang   PetscReal          f;
198c7272abeSHong Zhang   PetscInt           nlnk,*lnk,k,**bi_ptr;
1990298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
200c7272abeSHong Zhang   PetscBT            lnkbt;
201*ece542f9SHong Zhang   PetscBool          missing;
20283287d42SBarry Smith 
20383287d42SBarry Smith   PetscFunctionBegin;
204e32f2f54SBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"matrix must be square");
205*ece542f9SHong Zhang   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
206*ece542f9SHong Zhang   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
207*ece542f9SHong Zhang 
2086ba06ab7SHong Zhang   if (bs>1) {  /* check shifttype */
209f23aa3ddSBarry 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");
2106ba06ab7SHong Zhang   }
2116ba06ab7SHong Zhang 
21283287d42SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
21383287d42SBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
21483287d42SBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
21583287d42SBarry Smith 
216bc74c81fSJed Brown   /* get new row and diagonal pointers, must be allocated separately because they will be given to the Mat_SeqAIJ and freed separately */
217785e854fSJed Brown   ierr  = PetscMalloc1((n+1),&bi);CHKERRQ(ierr);
218785e854fSJed Brown   ierr  = PetscMalloc1((n+1),&bdiag);CHKERRQ(ierr);
219a6df321fSShri Abhyankar   bi[0] = bdiag[0] = 0;
220b2b2dd24SShri Abhyankar 
221b2b2dd24SShri Abhyankar   /* linked list for storing column indices of the active row */
222b2b2dd24SShri Abhyankar   nlnk = n + 1;
223b2b2dd24SShri Abhyankar   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
224b2b2dd24SShri Abhyankar 
225dcca6d9dSJed Brown   ierr = PetscMalloc2(n+1,&bi_ptr,n+1,&im);CHKERRQ(ierr);
226b2b2dd24SShri Abhyankar 
227b2b2dd24SShri Abhyankar   /* initial FreeSpace size is f*(ai[n]+1) */
228b2b2dd24SShri Abhyankar   f    = info->fill;
229b2b2dd24SShri Abhyankar   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
23026fbe8dcSKarl Rupp 
231b2b2dd24SShri Abhyankar   current_space = free_space;
232b2b2dd24SShri Abhyankar 
233b2b2dd24SShri Abhyankar   for (i=0; i<n; i++) {
234b2b2dd24SShri Abhyankar     /* copy previous fill into linked list */
235b2b2dd24SShri Abhyankar     nzi = 0;
236b2b2dd24SShri Abhyankar     nnz = ai[r[i]+1] - ai[r[i]];
237b2b2dd24SShri Abhyankar     ajtmp = aj + ai[r[i]];
238b2b2dd24SShri Abhyankar     ierr  = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
239b2b2dd24SShri Abhyankar     nzi  += nlnk;
240b2b2dd24SShri Abhyankar 
241b2b2dd24SShri Abhyankar     /* add pivot rows into linked list */
242b2b2dd24SShri Abhyankar     row = lnk[n];
243b2b2dd24SShri Abhyankar     while (row < i) {
244b2b2dd24SShri Abhyankar       nzbd  = bdiag[row] + 1;   /* num of entries in the row with column index <= row */
245b2b2dd24SShri Abhyankar       ajtmp = bi_ptr[row] + nzbd;   /* points to the entry next to the diagonal */
246b2b2dd24SShri Abhyankar       ierr  = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
247b2b2dd24SShri Abhyankar       nzi  += nlnk;
248b2b2dd24SShri Abhyankar       row   = lnk[row];
249b2b2dd24SShri Abhyankar     }
250b2b2dd24SShri Abhyankar     bi[i+1] = bi[i] + nzi;
251b2b2dd24SShri Abhyankar     im[i]   = nzi;
252b2b2dd24SShri Abhyankar 
253b2b2dd24SShri Abhyankar     /* mark bdiag */
254b2b2dd24SShri Abhyankar     nzbd = 0;
255b2b2dd24SShri Abhyankar     nnz  = nzi;
256b2b2dd24SShri Abhyankar     k    = lnk[n];
257b2b2dd24SShri Abhyankar     while (nnz-- && k < i) {
258b2b2dd24SShri Abhyankar       nzbd++;
259b2b2dd24SShri Abhyankar       k = lnk[k];
260b2b2dd24SShri Abhyankar     }
2612ce24eb6SHong Zhang     bdiag[i] = nzbd; /* note : bdaig[i] = nnzL as input for PetscFreeSpaceContiguous_LU() */
262b2b2dd24SShri Abhyankar 
263b2b2dd24SShri Abhyankar     /* if free space is not available, make more free space */
264b2b2dd24SShri Abhyankar     if (current_space->local_remaining<nzi) {
265b2b2dd24SShri Abhyankar       nnz  = 2*(n - i)*nzi; /* estimated and max additional space needed */
266b2b2dd24SShri Abhyankar       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
267b2b2dd24SShri Abhyankar       reallocs++;
268b2b2dd24SShri Abhyankar     }
269b2b2dd24SShri Abhyankar 
270b2b2dd24SShri Abhyankar     /* copy data into free space, then initialize lnk */
271b2b2dd24SShri Abhyankar     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
27226fbe8dcSKarl Rupp 
273b2b2dd24SShri Abhyankar     bi_ptr[i]                       = current_space->array;
274b2b2dd24SShri Abhyankar     current_space->array           += nzi;
275b2b2dd24SShri Abhyankar     current_space->local_used      += nzi;
276b2b2dd24SShri Abhyankar     current_space->local_remaining -= nzi;
277b2b2dd24SShri Abhyankar   }
278b2b2dd24SShri Abhyankar 
279b2b2dd24SShri Abhyankar   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
280b2b2dd24SShri Abhyankar   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
281b2b2dd24SShri Abhyankar 
2829263d837SHong Zhang   /* copy free_space into bj and free free_space; set bi, bj, bdiag in new datastructure; */
283785e854fSJed Brown   ierr = PetscMalloc1((bi[n]+1),&bj);CHKERRQ(ierr);
2842ce24eb6SHong Zhang   ierr = PetscFreeSpaceContiguous_LU(&free_space,bj,n,bi,bdiag);CHKERRQ(ierr);
285b2b2dd24SShri Abhyankar   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
286b2b2dd24SShri Abhyankar   ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr);
287b2b2dd24SShri Abhyankar 
288b2b2dd24SShri Abhyankar   /* put together the new matrix */
2890298fd71SBarry Smith   ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
2903bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)B,(PetscObject)isicol);CHKERRQ(ierr);
291b2b2dd24SShri Abhyankar   b    = (Mat_SeqBAIJ*)(B)->data;
29226fbe8dcSKarl Rupp 
293b2b2dd24SShri Abhyankar   b->free_a       = PETSC_TRUE;
294b2b2dd24SShri Abhyankar   b->free_ij      = PETSC_TRUE;
295b2b2dd24SShri Abhyankar   b->singlemalloc = PETSC_FALSE;
29626fbe8dcSKarl Rupp 
297b2b2dd24SShri Abhyankar   ierr             = PetscMalloc((bdiag[0]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
298b2b2dd24SShri Abhyankar   b->j             = bj;
299b2b2dd24SShri Abhyankar   b->i             = bi;
300b2b2dd24SShri Abhyankar   b->diag          = bdiag;
301b2b2dd24SShri Abhyankar   b->free_diag     = PETSC_TRUE;
302b2b2dd24SShri Abhyankar   b->ilen          = 0;
303b2b2dd24SShri Abhyankar   b->imax          = 0;
304b2b2dd24SShri Abhyankar   b->row           = isrow;
305b2b2dd24SShri Abhyankar   b->col           = iscol;
306b2b2dd24SShri Abhyankar   b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
30726fbe8dcSKarl Rupp 
308b2b2dd24SShri Abhyankar   ierr    = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
309b2b2dd24SShri Abhyankar   ierr    = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
310b2b2dd24SShri Abhyankar   b->icol = isicol;
311785e854fSJed Brown   ierr    = PetscMalloc1((bs*n+bs),&b->solve_work);CHKERRQ(ierr);
3123bb1ff40SBarry Smith   ierr    = PetscLogObjectMemory((PetscObject)B,(bdiag[0]+1)*(sizeof(PetscInt)+sizeof(PetscScalar)*bs2));CHKERRQ(ierr);
313b2b2dd24SShri Abhyankar 
314b2b2dd24SShri Abhyankar   b->maxnz = b->nz = bdiag[0]+1;
31526fbe8dcSKarl Rupp 
316d5f3da31SBarry Smith   B->factortype            =  MAT_FACTOR_LU;
317b2b2dd24SShri Abhyankar   B->info.factor_mallocs   = reallocs;
318b2b2dd24SShri Abhyankar   B->info.fill_ratio_given = f;
319b2b2dd24SShri Abhyankar 
320b2b2dd24SShri Abhyankar   if (ai[n] != 0) {
321b2b2dd24SShri Abhyankar     B->info.fill_ratio_needed = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]);
322b2b2dd24SShri Abhyankar   } else {
323b2b2dd24SShri Abhyankar     B->info.fill_ratio_needed = 0.0;
324b2b2dd24SShri Abhyankar   }
3259263d837SHong Zhang #if defined(PETSC_USE_INFO)
3269263d837SHong Zhang   if (ai[n] != 0) {
3279263d837SHong Zhang     PetscReal af = B->info.fill_ratio_needed;
32857622a8eSBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);CHKERRQ(ierr);
32957622a8eSBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
33057622a8eSBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g);\n",(double)af);CHKERRQ(ierr);
3319263d837SHong Zhang     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
3329263d837SHong Zhang   } else {
3339263d837SHong Zhang     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
3349263d837SHong Zhang   }
3359263d837SHong Zhang #endif
336b2b2dd24SShri Abhyankar 
337b2b2dd24SShri Abhyankar   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
338b2b2dd24SShri Abhyankar   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
33926fbe8dcSKarl Rupp 
340ace3abfcSBarry Smith   both_identity = (PetscBool) (row_identity && col_identity);
34126fbe8dcSKarl Rupp 
3424dd39f65SShri Abhyankar   ierr = MatSeqBAIJSetNumericFactorization(B,both_identity);CHKERRQ(ierr);
343b2b2dd24SShri Abhyankar   PetscFunctionReturn(0);
344b2b2dd24SShri Abhyankar }
345b2b2dd24SShri Abhyankar 
346b2b2dd24SShri Abhyankar #undef __FUNCT__
34706e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorSymbolic_SeqBAIJ_inplace"
34806e38f1dSHong Zhang PetscErrorCode MatLUFactorSymbolic_SeqBAIJ_inplace(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
349faca2338SShri Abhyankar {
350faca2338SShri Abhyankar   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data,*b;
351faca2338SShri Abhyankar   PetscInt           n  =a->mbs,bs = A->rmap->bs,bs2=a->bs2;
352ace3abfcSBarry Smith   PetscBool          row_identity,col_identity,both_identity;
353faca2338SShri Abhyankar   IS                 isicol;
354faca2338SShri Abhyankar   PetscErrorCode     ierr;
355faca2338SShri Abhyankar   const PetscInt     *r,*ic;
356faca2338SShri Abhyankar   PetscInt           i,*ai=a->i,*aj=a->j;
357faca2338SShri Abhyankar   PetscInt           *bi,*bj,*ajtmp;
358faca2338SShri Abhyankar   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
359faca2338SShri Abhyankar   PetscReal          f;
360faca2338SShri Abhyankar   PetscInt           nlnk,*lnk,k,**bi_ptr;
3610298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
362faca2338SShri Abhyankar   PetscBT            lnkbt;
363*ece542f9SHong Zhang   PetscBool          missing;
364faca2338SShri Abhyankar 
365faca2338SShri Abhyankar   PetscFunctionBegin;
366e32f2f54SBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"matrix must be square");
367*ece542f9SHong Zhang   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
368*ece542f9SHong Zhang   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
369*ece542f9SHong Zhang 
370faca2338SShri Abhyankar   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
371faca2338SShri Abhyankar   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
372faca2338SShri Abhyankar   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
373faca2338SShri Abhyankar 
374bc74c81fSJed Brown   /* get new row and diagonal pointers, must be allocated separately because they will be given to the Mat_SeqAIJ and freed separately */
375785e854fSJed Brown   ierr = PetscMalloc1((n+1),&bi);CHKERRQ(ierr);
376785e854fSJed Brown   ierr = PetscMalloc1((n+1),&bdiag);CHKERRQ(ierr);
377bc74c81fSJed Brown 
378a6df321fSShri Abhyankar   bi[0] = bdiag[0] = 0;
379c7272abeSHong Zhang 
380c7272abeSHong Zhang   /* linked list for storing column indices of the active row */
381c7272abeSHong Zhang   nlnk = n + 1;
382c7272abeSHong Zhang   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
383c7272abeSHong Zhang 
384dcca6d9dSJed Brown   ierr = PetscMalloc2(n+1,&bi_ptr,n+1,&im);CHKERRQ(ierr);
385c7272abeSHong Zhang 
386c7272abeSHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
387c7272abeSHong Zhang   f             = info->fill;
388c7272abeSHong Zhang   ierr          = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
389c7272abeSHong Zhang   current_space = free_space;
39083287d42SBarry Smith 
39183287d42SBarry Smith   for (i=0; i<n; i++) {
392c7272abeSHong Zhang     /* copy previous fill into linked list */
39383287d42SBarry Smith     nzi = 0;
394c7272abeSHong Zhang     nnz = ai[r[i]+1] - ai[r[i]];
395c7272abeSHong Zhang     ajtmp = aj + ai[r[i]];
396c7272abeSHong Zhang     ierr  = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
397c7272abeSHong Zhang     nzi  += nlnk;
398c7272abeSHong Zhang 
399c7272abeSHong Zhang     /* add pivot rows into linked list */
400c7272abeSHong Zhang     row = lnk[n];
401c7272abeSHong Zhang     while (row < i) {
402c7272abeSHong Zhang       nzbd  = bdiag[row] - bi[row] + 1;   /* num of entries in the row with column index <= row */
403c7272abeSHong Zhang       ajtmp = bi_ptr[row] + nzbd;   /* points to the entry next to the diagonal */
404c7272abeSHong Zhang       ierr  = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
405c7272abeSHong Zhang       nzi  += nlnk;
406c7272abeSHong Zhang       row   = lnk[row];
40783287d42SBarry Smith     }
408c7272abeSHong Zhang     bi[i+1] = bi[i] + nzi;
409c7272abeSHong Zhang     im[i]   = nzi;
410c7272abeSHong Zhang 
411c7272abeSHong Zhang     /* mark bdiag */
412c7272abeSHong Zhang     nzbd = 0;
413c7272abeSHong Zhang     nnz  = nzi;
414c7272abeSHong Zhang     k    = lnk[n];
415c7272abeSHong Zhang     while (nnz-- && k < i) {
416c7272abeSHong Zhang       nzbd++;
417c7272abeSHong Zhang       k = lnk[k];
418c7272abeSHong Zhang     }
419c7272abeSHong Zhang     bdiag[i] = bi[i] + nzbd;
420c7272abeSHong Zhang 
421c7272abeSHong Zhang     /* if free space is not available, make more free space */
422c7272abeSHong Zhang     if (current_space->local_remaining<nzi) {
423c7272abeSHong Zhang       nnz  = (n - i)*nzi; /* estimated and max additional space needed */
424c7272abeSHong Zhang       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
425c7272abeSHong Zhang       reallocs++;
42683287d42SBarry Smith     }
42783287d42SBarry Smith 
428c7272abeSHong Zhang     /* copy data into free space, then initialize lnk */
429c7272abeSHong Zhang     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
43026fbe8dcSKarl Rupp 
431c7272abeSHong Zhang     bi_ptr[i]                       = current_space->array;
432c7272abeSHong Zhang     current_space->array           += nzi;
433c7272abeSHong Zhang     current_space->local_used      += nzi;
434c7272abeSHong Zhang     current_space->local_remaining -= nzi;
435c7272abeSHong Zhang   }
4366cf91177SBarry Smith #if defined(PETSC_USE_INFO)
43783287d42SBarry Smith   if (ai[n] != 0) {
438c7272abeSHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
43957622a8eSBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);CHKERRQ(ierr);
44057622a8eSBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
44157622a8eSBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g);\n",(double)af);CHKERRQ(ierr);
442ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
44383287d42SBarry Smith   } else {
444c7272abeSHong Zhang     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
44583287d42SBarry Smith   }
44663ba0a88SBarry Smith #endif
44783287d42SBarry Smith 
44883287d42SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
44983287d42SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
45083287d42SBarry Smith 
451c7272abeSHong Zhang   /* destroy list of free space and other temporary array(s) */
452785e854fSJed Brown   ierr = PetscMalloc1((bi[n]+1),&bj);CHKERRQ(ierr);
453c7272abeSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
454c7272abeSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
455c7272abeSHong Zhang   ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr);
45683287d42SBarry Smith 
45783287d42SBarry Smith   /* put together the new matrix */
4580298fd71SBarry Smith   ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
4593bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)B,(PetscObject)isicol);CHKERRQ(ierr);
460719d5645SBarry Smith   b    = (Mat_SeqBAIJ*)(B)->data;
46126fbe8dcSKarl Rupp 
462e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
463e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
464c7272abeSHong Zhang   b->singlemalloc = PETSC_FALSE;
46526fbe8dcSKarl Rupp 
466c7272abeSHong Zhang   ierr             = PetscMalloc((bi[n]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
467c7272abeSHong Zhang   b->j             = bj;
468c7272abeSHong Zhang   b->i             = bi;
469c7272abeSHong Zhang   b->diag          = bdiag;
4707f53bb6cSHong Zhang   b->free_diag     = PETSC_TRUE;
47183287d42SBarry Smith   b->ilen          = 0;
47283287d42SBarry Smith   b->imax          = 0;
47383287d42SBarry Smith   b->row           = isrow;
47483287d42SBarry Smith   b->col           = iscol;
475d3d32019SBarry Smith   b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
47626fbe8dcSKarl Rupp 
47783287d42SBarry Smith   ierr    = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
47883287d42SBarry Smith   ierr    = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
47983287d42SBarry Smith   b->icol = isicol;
48026fbe8dcSKarl Rupp 
481785e854fSJed Brown   ierr = PetscMalloc1((bs*n+bs),&b->solve_work);CHKERRQ(ierr);
4823bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)*bs2));CHKERRQ(ierr);
48383287d42SBarry Smith 
484c7272abeSHong Zhang   b->maxnz = b->nz = bi[n];
48526fbe8dcSKarl Rupp 
486d5f3da31SBarry Smith   (B)->factortype            =  MAT_FACTOR_LU;
487719d5645SBarry Smith   (B)->info.factor_mallocs   = reallocs;
488719d5645SBarry Smith   (B)->info.fill_ratio_given = f;
489c7272abeSHong Zhang 
49083287d42SBarry Smith   if (ai[n] != 0) {
491c7272abeSHong Zhang     (B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
49283287d42SBarry Smith   } else {
493719d5645SBarry Smith     (B)->info.fill_ratio_needed = 0.0;
49483287d42SBarry Smith   }
495c7272abeSHong Zhang 
496db4efbfdSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
497db4efbfdSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
49826fbe8dcSKarl Rupp 
499ace3abfcSBarry Smith   both_identity = (PetscBool) (row_identity && col_identity);
50026fbe8dcSKarl Rupp 
5018b1456e3SHong Zhang   ierr = MatSeqBAIJSetNumericFactorization_inplace(B,both_identity);CHKERRQ(ierr);
50283287d42SBarry Smith   PetscFunctionReturn(0);
50383287d42SBarry Smith }
504db4efbfdSBarry Smith 
505