xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact.c (revision 2a6744eb01855f5aa328eb8fdf4b0d01e72ad151)
1 #define PETSCMAT_DLL
2 
3 #include "src/mat/impls/baij/seq/baij.h"
4 #include "src/mat/impls/sbaij/seq/sbaij.h"
5 #include "src/inline/ilu.h"
6 #include "include/petscis.h"
7 
8 #if !defined(PETSC_USE_COMPLEX)
9 /*
10   input:
11    F -- numeric factor
12   output:
13    nneg, nzero, npos: matrix inertia
14 */
15 
16 #undef __FUNCT__
17 #define __FUNCT__ "MatGetInertia_SeqSBAIJ"
18 PetscErrorCode MatGetInertia_SeqSBAIJ(Mat F,PetscInt *nneig,PetscInt *nzero,PetscInt *npos)
19 {
20   Mat_SeqSBAIJ *fact_ptr = (Mat_SeqSBAIJ*)F->data;
21   PetscScalar  *dd = fact_ptr->a;
22   PetscInt     mbs=fact_ptr->mbs,bs=F->rmap.bs,i,nneig_tmp,npos_tmp,*fi = fact_ptr->i;
23 
24   PetscFunctionBegin;
25   if (bs != 1) SETERRQ1(PETSC_ERR_SUP,"No support for bs: %D >1 yet",bs);
26   nneig_tmp = 0; npos_tmp = 0;
27   for (i=0; i<mbs; i++){
28     if (PetscRealPart(dd[*fi]) > 0.0){
29       npos_tmp++;
30     } else if (PetscRealPart(dd[*fi]) < 0.0){
31       nneig_tmp++;
32     }
33     fi++;
34   }
35   if (nneig) *nneig = nneig_tmp;
36   if (npos)  *npos  = npos_tmp;
37   if (nzero) *nzero = mbs - nneig_tmp - npos_tmp;
38 
39   PetscFunctionReturn(0);
40 }
41 #endif /* !defined(PETSC_USE_COMPLEX) */
42 
43 /*
44   Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP.
45   Use Modified Sparse Row (MSR) storage for u and ju. See page 85, "Iterative Methods ..." by Saad.
46 */
47 #undef __FUNCT__
48 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ_MSR"
49 PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(Mat A,IS perm,MatFactorInfo *info,Mat *B)
50 {
51   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b;
52   PetscErrorCode ierr;
53   PetscInt       *rip,i,mbs = a->mbs,*ai,*aj;
54   PetscInt       *jutmp,bs = A->rmap.bs,bs2=a->bs2;
55   PetscInt       m,reallocs = 0,prow;
56   PetscInt       *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd;
57   PetscReal      f = info->fill;
58   PetscTruth     perm_identity;
59 
60   PetscFunctionBegin;
61   /* check whether perm is the identity mapping */
62   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
63   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
64 
65   if (perm_identity){ /* without permutation */
66     a->permute = PETSC_FALSE;
67     ai = a->i; aj = a->j;
68   } else {            /* non-trivial permutation */
69     a->permute = PETSC_TRUE;
70     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
71     ai = a->inew; aj = a->jnew;
72   }
73 
74   /* initialization */
75   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&iu);CHKERRQ(ierr);
76   umax  = (PetscInt)(f*ai[mbs] + 1); umax += mbs + 1;
77   ierr  = PetscMalloc(umax*sizeof(PetscInt),&ju);CHKERRQ(ierr);
78   iu[0] = mbs+1;
79   juidx = mbs + 1; /* index for ju */
80   ierr  = PetscMalloc(2*mbs*sizeof(PetscInt),&jl);CHKERRQ(ierr); /* linked list for pivot row */
81   q     = jl + mbs;   /* linked list for col index */
82   for (i=0; i<mbs; i++){
83     jl[i] = mbs;
84     q[i] = 0;
85   }
86 
87   /* for each row k */
88   for (k=0; k<mbs; k++){
89     for (i=0; i<mbs; i++) q[i] = 0;  /* to be removed! */
90     nzk  = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */
91     q[k] = mbs;
92     /* initialize nonzero structure of k-th row to row rip[k] of A */
93     jmin = ai[rip[k]] +1; /* exclude diag[k] */
94     jmax = ai[rip[k]+1];
95     for (j=jmin; j<jmax; j++){
96       vj = rip[aj[j]]; /* col. value */
97       if(vj > k){
98         qm = k;
99         do {
100           m  = qm; qm = q[m];
101         } while(qm < vj);
102         if (qm == vj) {
103           SETERRQ(PETSC_ERR_PLIB,"Duplicate entry in A\n");
104         }
105         nzk++;
106         q[m]  = vj;
107         q[vj] = qm;
108       } /* if(vj > k) */
109     } /* for (j=jmin; j<jmax; j++) */
110 
111     /* modify nonzero structure of k-th row by computing fill-in
112        for each row i to be merged in */
113     prow = k;
114     prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */
115 
116     while (prow < k){
117       /* merge row prow into k-th row */
118       jmin = iu[prow] + 1; jmax = iu[prow+1];
119       qm = k;
120       for (j=jmin; j<jmax; j++){
121         vj = ju[j];
122         do {
123           m = qm; qm = q[m];
124         } while (qm < vj);
125         if (qm != vj){
126          nzk++; q[m] = vj; q[vj] = qm; qm = vj;
127         }
128       }
129       prow = jl[prow]; /* next pivot row */
130     }
131 
132     /* add k to row list for first nonzero element in k-th row */
133     if (nzk > 0){
134       i = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */
135       jl[k] = jl[i]; jl[i] = k;
136     }
137     iu[k+1] = iu[k] + nzk;
138 
139     /* allocate more space to ju if needed */
140     if (iu[k+1] > umax) {
141       /* estimate how much additional space we will need */
142       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
143       /* just double the memory each time */
144       maxadd = umax;
145       if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2;
146       umax += maxadd;
147 
148       /* allocate a longer ju */
149       ierr = PetscMalloc(umax*sizeof(PetscInt),&jutmp);CHKERRQ(ierr);
150       ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(PetscInt));CHKERRQ(ierr);
151       ierr = PetscFree(ju);CHKERRQ(ierr);
152       ju   = jutmp;
153       reallocs++; /* count how many times we realloc */
154     }
155 
156     /* save nonzero structure of k-th row in ju */
157     i=k;
158     while (nzk --) {
159       i           = q[i];
160       ju[juidx++] = i;
161     }
162   }
163 
164 #if defined(PETSC_USE_INFO)
165   if (ai[mbs] != 0) {
166     PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
167     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
168     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
169     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
170     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
171   } else {
172     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
173   }
174 #endif
175 
176   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
177   ierr = PetscFree(jl);CHKERRQ(ierr);
178 
179   /* put together the new matrix */
180   ierr = MatCreate(A->comm,B);CHKERRQ(ierr);
181   ierr = MatSetSizes(*B,bs*mbs,bs*mbs,bs*mbs,bs*mbs);CHKERRQ(ierr);
182   ierr = MatSetType(*B,A->type_name);CHKERRQ(ierr);
183   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
184 
185   /* ierr = PetscLogObjectParent(*B,iperm);CHKERRQ(ierr); */
186   b = (Mat_SeqSBAIJ*)(*B)->data;
187   b->singlemalloc = PETSC_FALSE;
188   b->free_a       = PETSC_TRUE;
189   b->free_ij       = PETSC_TRUE;
190   ierr = PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
191   b->j    = ju;
192   b->i    = iu;
193   b->diag = 0;
194   b->ilen = 0;
195   b->imax = 0;
196   b->row  = perm;
197   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
198   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
199   b->icol = perm;
200   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
201   ierr    = PetscMalloc((bs*mbs+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
202   /* In b structure:  Free imax, ilen, old a, old j.
203      Allocate idnew, solve_work, new a, new j */
204   ierr = PetscLogObjectMemory(*B,(iu[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
205   b->maxnz = b->nz = iu[mbs];
206 
207   (*B)->factor                 = FACTOR_CHOLESKY;
208   (*B)->info.factor_mallocs    = reallocs;
209   (*B)->info.fill_ratio_given  = f;
210   if (ai[mbs] != 0) {
211     (*B)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
212   } else {
213     (*B)->info.fill_ratio_needed = 0.0;
214   }
215 
216   if (perm_identity){
217     switch (bs) {
218       case 1:
219         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering;
220         (*B)->ops->solve                 = MatSolve_SeqSBAIJ_1_NaturalOrdering;
221         ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=1\n");CHKERRQ(ierr);
222         break;
223       case 2:
224         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
225         (*B)->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
226         ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=2\n");CHKERRQ(ierr);
227         break;
228       case 3:
229         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
230         (*B)->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
231         ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=3\n");CHKERRQ(ierr);
232         break;
233       case 4:
234         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
235         (*B)->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
236         ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=4\n");CHKERRQ(ierr);
237         break;
238       case 5:
239         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
240         (*B)->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
241         ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=5\n");CHKERRQ(ierr);
242         break;
243       case 6:
244         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
245         (*B)->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
246         ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=6\n");CHKERRQ(ierr);
247         break;
248       case 7:
249         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
250         (*B)->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
251         ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=7\n");CHKERRQ(ierr);
252       break;
253       default:
254         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering;
255         (*B)->ops->solve           = MatSolve_SeqSBAIJ_N_NaturalOrdering;
256         ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS>7\n");CHKERRQ(ierr);
257       break;
258     }
259   }
260   PetscFunctionReturn(0);
261 }
262 /*
263     Symbolic U^T*D*U factorization for SBAIJ format.
264 */
265 #include "petscbt.h"
266 #include "src/mat/utils/freespace.h"
267 #undef __FUNCT__
268 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ"
269 PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
270 {
271   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data;
272   Mat_SeqSBAIJ       *b;
273   Mat                B;
274   PetscErrorCode     ierr;
275   PetscTruth         perm_identity;
276   PetscReal          fill = info->fill;
277   PetscInt           *rip,i,mbs=a->mbs,bs=A->rmap.bs,*ai,*aj,reallocs=0,prow;
278   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
279   PetscInt           nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr;
280   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
281   PetscBT            lnkbt;
282 
283   PetscFunctionBegin;
284   /*
285    This code originally uses Modified Sparse Row (MSR) storage
286    (see page 85, "Iterative Methods ..." by Saad) for the output matrix B - bad choise!
287    Then it is rewritten so the factor B takes seqsbaij format. However the associated
288    MatCholeskyFactorNumeric_() have not been modified for the cases of bs>1 or !perm_identity,
289    thus the original code in MSR format is still used for these cases.
290    The code below should replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR() whenever
291    MatCholeskyFactorNumeric_() is modified for using sbaij symbolic factor.
292   */
293   if (bs > 1){
294     ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(A,perm,info,fact);CHKERRQ(ierr);
295     PetscFunctionReturn(0);
296   }
297 
298   /* check whether perm is the identity mapping */
299   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
300 
301   if (perm_identity){
302     a->permute = PETSC_FALSE;
303     ai = a->i; aj = a->j;
304   } else {
305     a->permute = PETSC_TRUE;
306     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
307     ai = a->inew; aj = a->jnew;
308   }
309   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
310 
311   /* initialization */
312   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
313   ui[0] = 0;
314 
315   /* jl: linked list for storing indices of the pivot rows
316      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
317   ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
318   il     = jl + mbs;
319   cols   = il + mbs;
320   ui_ptr = (PetscInt**)(cols + mbs);
321 
322   for (i=0; i<mbs; i++){
323     jl[i] = mbs; il[i] = 0;
324   }
325 
326   /* create and initialize a linked list for storing column indices of the active row k */
327   nlnk = mbs + 1;
328   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
329 
330   /* initial FreeSpace size is fill*(ai[mbs]+1) */
331   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
332   current_space = free_space;
333 
334   for (k=0; k<mbs; k++){  /* for each active row k */
335     /* initialize lnk by the column indices of row rip[k] of A */
336     nzk   = 0;
337     ncols = ai[rip[k]+1] - ai[rip[k]];
338     for (j=0; j<ncols; j++){
339       i = *(aj + ai[rip[k]] + j);
340       cols[j] = rip[i];
341     }
342     ierr = PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
343     nzk += nlnk;
344 
345     /* update lnk by computing fill-in for each pivot row to be merged in */
346     prow = jl[k]; /* 1st pivot row */
347 
348     while (prow < k){
349       nextprow = jl[prow];
350       /* merge prow into k-th row */
351       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
352       jmax = ui[prow+1];
353       ncols = jmax-jmin;
354       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
355       ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
356       nzk += nlnk;
357 
358       /* update il and jl for prow */
359       if (jmin < jmax){
360         il[prow] = jmin;
361         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
362       }
363       prow = nextprow;
364     }
365 
366     /* if free space is not available, make more free space */
367     if (current_space->local_remaining<nzk) {
368       i = mbs - k + 1; /* num of unfactored rows */
369       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
370       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
371       reallocs++;
372     }
373 
374     /* copy data into free space, then initialize lnk */
375     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
376 
377     /* add the k-th row into il and jl */
378     if (nzk-1 > 0){
379       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
380       jl[k] = jl[i]; jl[i] = k;
381       il[k] = ui[k] + 1;
382     }
383     ui_ptr[k] = current_space->array;
384     current_space->array           += nzk;
385     current_space->local_used      += nzk;
386     current_space->local_remaining -= nzk;
387 
388     ui[k+1] = ui[k] + nzk;
389   }
390 
391 #if defined(PETSC_USE_INFO)
392   if (ai[mbs] != 0) {
393     PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
394     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
395     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
396     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
397   } else {
398     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
399   }
400 #endif
401 
402   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
403   ierr = PetscFree(jl);CHKERRQ(ierr);
404 
405   /* destroy list of free space and other temporary array(s) */
406   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
407   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
408   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
409 
410   /* put together the new matrix in MATSEQSBAIJ format */
411   ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr);
412   ierr = MatSetSizes(*fact,mbs,mbs,mbs,mbs);CHKERRQ(ierr);
413   B    = *fact;
414   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
415   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
416 
417   b = (Mat_SeqSBAIJ*)B->data;
418   b->singlemalloc = PETSC_FALSE;
419   b->free_a       = PETSC_TRUE;
420   b->free_ij      = PETSC_TRUE;
421   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
422   b->j    = uj;
423   b->i    = ui;
424   b->diag = 0;
425   b->ilen = 0;
426   b->imax = 0;
427   b->row  = perm;
428   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
429   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
430   b->icol = perm;
431   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
432   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
433   ierr    = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
434   b->maxnz = b->nz = ui[mbs];
435 
436   B->factor                 = FACTOR_CHOLESKY;
437   B->info.factor_mallocs    = reallocs;
438   B->info.fill_ratio_given  = fill;
439   if (ai[mbs] != 0) {
440     B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
441   } else {
442     B->info.fill_ratio_needed = 0.0;
443   }
444 
445   if (perm_identity){
446     switch (bs) {
447       case 1:
448         B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering;
449         B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
450         ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=1\n");CHKERRQ(ierr);
451         break;
452       case 2:
453         B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
454         B->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
455         ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=2\n");CHKERRQ(ierr);
456         break;
457       case 3:
458         B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
459         B->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
460         ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=3\n");CHKERRQ(ierr);
461         break;
462       case 4:
463         B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
464         B->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
465         ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=4\n");CHKERRQ(ierr);
466         break;
467       case 5:
468         B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
469         B->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
470         ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=5\n");CHKERRQ(ierr);
471         break;
472       case 6:
473         B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
474         B->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
475         ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=6\n");CHKERRQ(ierr);
476         break;
477       case 7:
478         B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
479         B->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
480         ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=7\n");CHKERRQ(ierr);
481       break;
482       default:
483         B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering;
484         B->ops->solve           = MatSolve_SeqSBAIJ_N_NaturalOrdering;
485         ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS>7\n");CHKERRQ(ierr);
486       break;
487     }
488   }
489   PetscFunctionReturn(0);
490 }
491 #undef __FUNCT__
492 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N"
493 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat A,MatFactorInfo *info,Mat *B)
494 {
495   Mat            C = *B;
496   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
497   IS             perm = b->row;
498   PetscErrorCode ierr;
499   PetscInt       *perm_ptr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
500   PetscInt       *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
501   PetscInt       bs=A->rmap.bs,bs2 = a->bs2;
502   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
503   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
504   MatScalar      *work;
505   PetscInt       *pivots;
506 
507   PetscFunctionBegin;
508   /* initialization */
509   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
510   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
511   ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr);
512   jl   = il + mbs;
513   for (i=0; i<mbs; i++) {
514     jl[i] = mbs; il[0] = 0;
515   }
516   ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr);
517   uik  = dk + bs2;
518   work = uik + bs2;
519   ierr = PetscMalloc(bs*sizeof(PetscInt),&pivots);CHKERRQ(ierr);
520 
521   ierr  = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
522 
523   /* check permutation */
524   if (!a->permute){
525     ai = a->i; aj = a->j; aa = a->a;
526   } else {
527     ai   = a->inew; aj = a->jnew;
528     ierr = PetscMalloc(bs2*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
529     ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
530     ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr);
531     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
532 
533     for (i=0; i<mbs; i++){
534       jmin = ai[i]; jmax = ai[i+1];
535       for (j=jmin; j<jmax; j++){
536         while (a2anew[j] != j){
537           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
538           for (k1=0; k1<bs2; k1++){
539             dk[k1]       = aa[k*bs2+k1];
540             aa[k*bs2+k1] = aa[j*bs2+k1];
541             aa[j*bs2+k1] = dk[k1];
542           }
543         }
544         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
545         if (i > aj[j]){
546           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
547           ap = aa + j*bs2;                     /* ptr to the beginning of j-th block of aa */
548           for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
549           for (k=0; k<bs; k++){               /* j-th block of aa <- dk^T */
550             for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1];
551           }
552         }
553       }
554     }
555     ierr = PetscFree(a2anew);CHKERRQ(ierr);
556   }
557 
558   /* for each row k */
559   for (k = 0; k<mbs; k++){
560 
561     /*initialize k-th row with elements nonzero in row perm(k) of A */
562     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
563 
564     ap = aa + jmin*bs2;
565     for (j = jmin; j < jmax; j++){
566       vj = perm_ptr[aj[j]];         /* block col. index */
567       rtmp_ptr = rtmp + vj*bs2;
568       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
569     }
570 
571     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
572     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
573     i = jl[k]; /* first row to be added to k_th row  */
574 
575     while (i < k){
576       nexti = jl[i]; /* next row to be added to k_th row */
577 
578       /* compute multiplier */
579       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
580 
581       /* uik = -inv(Di)*U_bar(i,k) */
582       diag = ba + i*bs2;
583       u    = ba + ili*bs2;
584       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
585       Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
586 
587       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
588       Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
589 
590       /* update -U(i,k) */
591       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
592 
593       /* add multiple of row i to k-th row ... */
594       jmin = ili + 1; jmax = bi[i+1];
595       if (jmin < jmax){
596         for (j=jmin; j<jmax; j++) {
597           /* rtmp += -U(i,k)^T * U_bar(i,j) */
598           rtmp_ptr = rtmp + bj[j]*bs2;
599           u = ba + j*bs2;
600           Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
601         }
602 
603         /* ... add i to row list for next nonzero entry */
604         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
605         j     = bj[jmin];
606         jl[i] = jl[j]; jl[j] = i; /* update jl */
607       }
608       i = nexti;
609     }
610 
611     /* save nonzero entries in k-th row of U ... */
612 
613     /* invert diagonal block */
614     diag = ba+k*bs2;
615     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
616     ierr = Kernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr);
617 
618     jmin = bi[k]; jmax = bi[k+1];
619     if (jmin < jmax) {
620       for (j=jmin; j<jmax; j++){
621          vj = bj[j];           /* block col. index of U */
622          u   = ba + j*bs2;
623          rtmp_ptr = rtmp + vj*bs2;
624          for (k1=0; k1<bs2; k1++){
625            *u++        = *rtmp_ptr;
626            *rtmp_ptr++ = 0.0;
627          }
628       }
629 
630       /* ... add k to row list for first nonzero entry in k-th row */
631       il[k] = jmin;
632       i     = bj[jmin];
633       jl[k] = jl[i]; jl[i] = k;
634     }
635   }
636 
637   ierr = PetscFree(rtmp);CHKERRQ(ierr);
638   ierr = PetscFree(il);CHKERRQ(ierr);
639   ierr = PetscFree(dk);CHKERRQ(ierr);
640   ierr = PetscFree(pivots);CHKERRQ(ierr);
641   if (a->permute){
642     ierr = PetscFree(aa);CHKERRQ(ierr);
643   }
644 
645   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
646   C->factor       = FACTOR_CHOLESKY;
647   C->assembled    = PETSC_TRUE;
648   C->preallocated = PETSC_TRUE;
649   ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
650   PetscFunctionReturn(0);
651 }
652 
653 #undef __FUNCT__
654 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering"
655 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B)
656 {
657   Mat            C = *B;
658   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
659   PetscErrorCode ierr;
660   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
661   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
662   PetscInt       bs=A->rmap.bs,bs2 = a->bs2;
663   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
664   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
665   MatScalar      *work;
666   PetscInt       *pivots;
667 
668   PetscFunctionBegin;
669   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
670   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
671   ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr);
672   jl   = il + mbs;
673   for (i=0; i<mbs; i++) {
674     jl[i] = mbs; il[0] = 0;
675   }
676   ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr);
677   uik  = dk + bs2;
678   work = uik + bs2;
679   ierr = PetscMalloc(bs*sizeof(PetscInt),&pivots);CHKERRQ(ierr);
680 
681   ai = a->i; aj = a->j; aa = a->a;
682 
683   /* for each row k */
684   for (k = 0; k<mbs; k++){
685 
686     /*initialize k-th row with elements nonzero in row k of A */
687     jmin = ai[k]; jmax = ai[k+1];
688     ap = aa + jmin*bs2;
689     for (j = jmin; j < jmax; j++){
690       vj = aj[j];         /* block col. index */
691       rtmp_ptr = rtmp + vj*bs2;
692       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
693     }
694 
695     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
696     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
697     i = jl[k]; /* first row to be added to k_th row  */
698 
699     while (i < k){
700       nexti = jl[i]; /* next row to be added to k_th row */
701 
702       /* compute multiplier */
703       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
704 
705       /* uik = -inv(Di)*U_bar(i,k) */
706       diag = ba + i*bs2;
707       u    = ba + ili*bs2;
708       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
709       Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
710 
711       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
712       Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
713 
714       /* update -U(i,k) */
715       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
716 
717       /* add multiple of row i to k-th row ... */
718       jmin = ili + 1; jmax = bi[i+1];
719       if (jmin < jmax){
720         for (j=jmin; j<jmax; j++) {
721           /* rtmp += -U(i,k)^T * U_bar(i,j) */
722           rtmp_ptr = rtmp + bj[j]*bs2;
723           u = ba + j*bs2;
724           Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
725         }
726 
727         /* ... add i to row list for next nonzero entry */
728         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
729         j     = bj[jmin];
730         jl[i] = jl[j]; jl[j] = i; /* update jl */
731       }
732       i = nexti;
733     }
734 
735     /* save nonzero entries in k-th row of U ... */
736 
737     /* invert diagonal block */
738     diag = ba+k*bs2;
739     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
740     ierr = Kernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr);
741 
742     jmin = bi[k]; jmax = bi[k+1];
743     if (jmin < jmax) {
744       for (j=jmin; j<jmax; j++){
745          vj = bj[j];           /* block col. index of U */
746          u   = ba + j*bs2;
747          rtmp_ptr = rtmp + vj*bs2;
748          for (k1=0; k1<bs2; k1++){
749            *u++        = *rtmp_ptr;
750            *rtmp_ptr++ = 0.0;
751          }
752       }
753 
754       /* ... add k to row list for first nonzero entry in k-th row */
755       il[k] = jmin;
756       i     = bj[jmin];
757       jl[k] = jl[i]; jl[i] = k;
758     }
759   }
760 
761   ierr = PetscFree(rtmp);CHKERRQ(ierr);
762   ierr = PetscFree(il);CHKERRQ(ierr);
763   ierr = PetscFree(dk);CHKERRQ(ierr);
764   ierr = PetscFree(pivots);CHKERRQ(ierr);
765 
766   C->factor    = FACTOR_CHOLESKY;
767   C->assembled = PETSC_TRUE;
768   C->preallocated = PETSC_TRUE;
769   ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
770   PetscFunctionReturn(0);
771 }
772 
773 /*
774     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
775     Version for blocks 2 by 2.
776 */
777 #undef __FUNCT__
778 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2"
779 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat A,MatFactorInfo *info,Mat *B)
780 {
781   Mat            C = *B;
782   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
783   IS             perm = b->row;
784   PetscErrorCode ierr;
785   PetscInt       *perm_ptr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
786   PetscInt       *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
787   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
788   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
789 
790   PetscFunctionBegin;
791   /* initialization */
792   /* il and jl record the first nonzero element in each row of the accessing
793      window U(0:k, k:mbs-1).
794      jl:    list of rows to be added to uneliminated rows
795             i>= k: jl(i) is the first row to be added to row i
796             i<  k: jl(i) is the row following row i in some list of rows
797             jl(i) = mbs indicates the end of a list
798      il(i): points to the first nonzero element in columns k,...,mbs-1 of
799             row i of U */
800   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
801   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
802   ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr);
803   jl   = il + mbs;
804   for (i=0; i<mbs; i++) {
805     jl[i] = mbs; il[0] = 0;
806   }
807   ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr);
808   uik  = dk + 4;
809   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
810 
811   /* check permutation */
812   if (!a->permute){
813     ai = a->i; aj = a->j; aa = a->a;
814   } else {
815     ai   = a->inew; aj = a->jnew;
816     ierr = PetscMalloc(4*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
817     ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
818     ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr);
819     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
820 
821     for (i=0; i<mbs; i++){
822       jmin = ai[i]; jmax = ai[i+1];
823       for (j=jmin; j<jmax; j++){
824         while (a2anew[j] != j){
825           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
826           for (k1=0; k1<4; k1++){
827             dk[k1]       = aa[k*4+k1];
828             aa[k*4+k1] = aa[j*4+k1];
829             aa[j*4+k1] = dk[k1];
830           }
831         }
832         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
833         if (i > aj[j]){
834           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
835           ap = aa + j*4;     /* ptr to the beginning of the block */
836           dk[1] = ap[1];     /* swap ap[1] and ap[2] */
837           ap[1] = ap[2];
838           ap[2] = dk[1];
839         }
840       }
841     }
842     ierr = PetscFree(a2anew);CHKERRQ(ierr);
843   }
844 
845   /* for each row k */
846   for (k = 0; k<mbs; k++){
847 
848     /*initialize k-th row with elements nonzero in row perm(k) of A */
849     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
850     ap = aa + jmin*4;
851     for (j = jmin; j < jmax; j++){
852       vj = perm_ptr[aj[j]];         /* block col. index */
853       rtmp_ptr = rtmp + vj*4;
854       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
855     }
856 
857     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
858     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
859     i = jl[k]; /* first row to be added to k_th row  */
860 
861     while (i < k){
862       nexti = jl[i]; /* next row to be added to k_th row */
863 
864       /* compute multiplier */
865       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
866 
867       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
868       diag = ba + i*4;
869       u    = ba + ili*4;
870       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
871       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
872       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
873       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
874 
875       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
876       dk[0] += uik[0]*u[0] + uik[1]*u[1];
877       dk[1] += uik[2]*u[0] + uik[3]*u[1];
878       dk[2] += uik[0]*u[2] + uik[1]*u[3];
879       dk[3] += uik[2]*u[2] + uik[3]*u[3];
880 
881       /* update -U(i,k): ba[ili] = uik */
882       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
883 
884       /* add multiple of row i to k-th row ... */
885       jmin = ili + 1; jmax = bi[i+1];
886       if (jmin < jmax){
887         for (j=jmin; j<jmax; j++) {
888           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
889           rtmp_ptr = rtmp + bj[j]*4;
890           u = ba + j*4;
891           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
892           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
893           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
894           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
895         }
896 
897         /* ... add i to row list for next nonzero entry */
898         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
899         j     = bj[jmin];
900         jl[i] = jl[j]; jl[j] = i; /* update jl */
901       }
902       i = nexti;
903     }
904 
905     /* save nonzero entries in k-th row of U ... */
906 
907     /* invert diagonal block */
908     diag = ba+k*4;
909     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
910     ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr);
911 
912     jmin = bi[k]; jmax = bi[k+1];
913     if (jmin < jmax) {
914       for (j=jmin; j<jmax; j++){
915          vj = bj[j];           /* block col. index of U */
916          u   = ba + j*4;
917          rtmp_ptr = rtmp + vj*4;
918          for (k1=0; k1<4; k1++){
919            *u++        = *rtmp_ptr;
920            *rtmp_ptr++ = 0.0;
921          }
922       }
923 
924       /* ... add k to row list for first nonzero entry in k-th row */
925       il[k] = jmin;
926       i     = bj[jmin];
927       jl[k] = jl[i]; jl[i] = k;
928     }
929   }
930 
931   ierr = PetscFree(rtmp);CHKERRQ(ierr);
932   ierr = PetscFree(il);CHKERRQ(ierr);
933   ierr = PetscFree(dk);CHKERRQ(ierr);
934   if (a->permute) {
935     ierr = PetscFree(aa);CHKERRQ(ierr);
936   }
937   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
938   C->factor    = FACTOR_CHOLESKY;
939   C->assembled = PETSC_TRUE;
940   C->preallocated = PETSC_TRUE;
941   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
942   PetscFunctionReturn(0);
943 }
944 
945 /*
946       Version for when blocks are 2 by 2 Using natural ordering
947 */
948 #undef __FUNCT__
949 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering"
950 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B)
951 {
952   Mat            C = *B;
953   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
954   PetscErrorCode ierr;
955   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
956   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
957   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
958   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
959 
960   PetscFunctionBegin;
961   /* initialization */
962   /* il and jl record the first nonzero element in each row of the accessing
963      window U(0:k, k:mbs-1).
964      jl:    list of rows to be added to uneliminated rows
965             i>= k: jl(i) is the first row to be added to row i
966             i<  k: jl(i) is the row following row i in some list of rows
967             jl(i) = mbs indicates the end of a list
968      il(i): points to the first nonzero element in columns k,...,mbs-1 of
969             row i of U */
970   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
971   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
972   ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr);
973   jl   = il + mbs;
974   for (i=0; i<mbs; i++) {
975     jl[i] = mbs; il[0] = 0;
976   }
977   ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr);
978   uik  = dk + 4;
979 
980   ai = a->i; aj = a->j; aa = a->a;
981 
982   /* for each row k */
983   for (k = 0; k<mbs; k++){
984 
985     /*initialize k-th row with elements nonzero in row k of A */
986     jmin = ai[k]; jmax = ai[k+1];
987     ap = aa + jmin*4;
988     for (j = jmin; j < jmax; j++){
989       vj = aj[j];         /* block col. index */
990       rtmp_ptr = rtmp + vj*4;
991       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
992     }
993 
994     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
995     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
996     i = jl[k]; /* first row to be added to k_th row  */
997 
998     while (i < k){
999       nexti = jl[i]; /* next row to be added to k_th row */
1000 
1001       /* compute multiplier */
1002       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1003 
1004       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
1005       diag = ba + i*4;
1006       u    = ba + ili*4;
1007       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
1008       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
1009       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
1010       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
1011 
1012       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
1013       dk[0] += uik[0]*u[0] + uik[1]*u[1];
1014       dk[1] += uik[2]*u[0] + uik[3]*u[1];
1015       dk[2] += uik[0]*u[2] + uik[1]*u[3];
1016       dk[3] += uik[2]*u[2] + uik[3]*u[3];
1017 
1018       /* update -U(i,k): ba[ili] = uik */
1019       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
1020 
1021       /* add multiple of row i to k-th row ... */
1022       jmin = ili + 1; jmax = bi[i+1];
1023       if (jmin < jmax){
1024         for (j=jmin; j<jmax; j++) {
1025           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
1026           rtmp_ptr = rtmp + bj[j]*4;
1027           u = ba + j*4;
1028           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
1029           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
1030           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
1031           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
1032         }
1033 
1034         /* ... add i to row list for next nonzero entry */
1035         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
1036         j     = bj[jmin];
1037         jl[i] = jl[j]; jl[j] = i; /* update jl */
1038       }
1039       i = nexti;
1040     }
1041 
1042     /* save nonzero entries in k-th row of U ... */
1043 
1044     /* invert diagonal block */
1045     diag = ba+k*4;
1046     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
1047     ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr);
1048 
1049     jmin = bi[k]; jmax = bi[k+1];
1050     if (jmin < jmax) {
1051       for (j=jmin; j<jmax; j++){
1052          vj = bj[j];           /* block col. index of U */
1053          u   = ba + j*4;
1054          rtmp_ptr = rtmp + vj*4;
1055          for (k1=0; k1<4; k1++){
1056            *u++        = *rtmp_ptr;
1057            *rtmp_ptr++ = 0.0;
1058          }
1059       }
1060 
1061       /* ... add k to row list for first nonzero entry in k-th row */
1062       il[k] = jmin;
1063       i     = bj[jmin];
1064       jl[k] = jl[i]; jl[i] = k;
1065     }
1066   }
1067 
1068   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1069   ierr = PetscFree(il);CHKERRQ(ierr);
1070   ierr = PetscFree(dk);CHKERRQ(ierr);
1071 
1072   C->factor    = FACTOR_CHOLESKY;
1073   C->assembled = PETSC_TRUE;
1074   C->preallocated = PETSC_TRUE;
1075   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
1076   PetscFunctionReturn(0);
1077 }
1078 
1079 /*
1080     Numeric U^T*D*U factorization for SBAIJ format.
1081     Version for blocks are 1 by 1.
1082 */
1083 #undef __FUNCT__
1084 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1"
1085 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat A,MatFactorInfo *info,Mat *B)
1086 {
1087   Mat            C = *B;
1088   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data;
1089   IS             ip=b->row;
1090   PetscErrorCode ierr;
1091   PetscInt       *rip,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j,*bcol;
1092   PetscInt       *ai,*aj,*a2anew;
1093   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1094   MatScalar      *rtmp,*ba=b->a,*bval,*aa,dk,uikdi;
1095   PetscReal      zeropivot,rs,shiftnz;
1096   PetscReal      shiftpd;
1097   ChShift_Ctx    sctx;
1098   PetscInt       newshift;
1099 
1100   PetscFunctionBegin;
1101   /* initialization */
1102   shiftnz   = info->shiftnz;
1103   shiftpd   = info->shiftpd;
1104   zeropivot = info->zeropivot;
1105 
1106   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
1107   if (!a->permute){
1108     ai = a->i; aj = a->j; aa = a->a;
1109   } else {
1110     ai = a->inew; aj = a->jnew;
1111     nz = ai[mbs];
1112     ierr = PetscMalloc(nz*sizeof(MatScalar),&aa);CHKERRQ(ierr);
1113     a2anew = a->a2anew;
1114     bval   = a->a;
1115     for (j=0; j<nz; j++){
1116       aa[a2anew[j]] = *(bval++);
1117     }
1118   }
1119 
1120   /* initialization */
1121   /* il and jl record the first nonzero element in each row of the accessing
1122      window U(0:k, k:mbs-1).
1123      jl:    list of rows to be added to uneliminated rows
1124             i>= k: jl(i) is the first row to be added to row i
1125             i<  k: jl(i) is the row following row i in some list of rows
1126             jl(i) = mbs indicates the end of a list
1127      il(i): points to the first nonzero element in columns k,...,mbs-1 of
1128             row i of U */
1129   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
1130   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
1131   jl   = il + mbs;
1132   rtmp = (MatScalar*)(jl + mbs);
1133 
1134   sctx.shift_amount = 0;
1135   sctx.nshift       = 0;
1136   do {
1137     sctx.chshift = PETSC_FALSE;
1138     for (i=0; i<mbs; i++) {
1139       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1140     }
1141 
1142     for (k = 0; k<mbs; k++){
1143       /*initialize k-th row by the perm[k]-th row of A */
1144       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
1145       bval = ba + bi[k];
1146       for (j = jmin; j < jmax; j++){
1147         col = rip[aj[j]];
1148         rtmp[col] = aa[j];
1149         *bval++  = 0.0; /* for in-place factorization */
1150       }
1151 
1152       /* shift the diagonal of the matrix */
1153       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
1154 
1155       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1156       dk = rtmp[k];
1157       i = jl[k]; /* first row to be added to k_th row  */
1158 
1159       while (i < k){
1160         nexti = jl[i]; /* next row to be added to k_th row */
1161 
1162         /* compute multiplier, update diag(k) and U(i,k) */
1163         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1164         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
1165         dk += uikdi*ba[ili];
1166         ba[ili] = uikdi; /* -U(i,k) */
1167 
1168         /* add multiple of row i to k-th row */
1169         jmin = ili + 1; jmax = bi[i+1];
1170         if (jmin < jmax){
1171           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1172           /* update il and jl for row i */
1173           il[i] = jmin;
1174           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1175         }
1176         i = nexti;
1177       }
1178 
1179       /* shift the diagonals when zero pivot is detected */
1180       /* compute rs=sum of abs(off-diagonal) */
1181       rs   = 0.0;
1182       jmin = bi[k]+1;
1183       nz   = bi[k+1] - jmin;
1184       if (nz){
1185         bcol = bj + jmin;
1186         while (nz--){
1187           rs += PetscAbsScalar(rtmp[*bcol]);
1188           bcol++;
1189         }
1190       }
1191 
1192       sctx.rs = rs;
1193       sctx.pv = dk;
1194       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
1195       if (newshift == 1) break;    /* sctx.shift_amount is updated */
1196 
1197       /* copy data into U(k,:) */
1198       ba[bi[k]] = 1.0/dk; /* U(k,k) */
1199       jmin = bi[k]+1; jmax = bi[k+1];
1200       if (jmin < jmax) {
1201         for (j=jmin; j<jmax; j++){
1202           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
1203         }
1204         /* add the k-th row into il and jl */
1205         il[k] = jmin;
1206         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
1207       }
1208     }
1209   } while (sctx.chshift);
1210   ierr = PetscFree(il);CHKERRQ(ierr);
1211   if (a->permute){ierr = PetscFree(aa);CHKERRQ(ierr);}
1212 
1213   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
1214   C->factor       = FACTOR_CHOLESKY;
1215   C->assembled    = PETSC_TRUE;
1216   C->preallocated = PETSC_TRUE;
1217   ierr = PetscLogFlops(C->rmap.N);CHKERRQ(ierr);
1218     if (sctx.nshift){
1219     if (shiftnz) {
1220       ierr = PetscInfo2(0,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1221     } else if (shiftpd) {
1222       ierr = PetscInfo2(0,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1223     }
1224   }
1225   PetscFunctionReturn(0);
1226 }
1227 
1228 /*
1229   Version for when blocks are 1 by 1 Using natural ordering
1230 */
1231 #undef __FUNCT__
1232 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering"
1233 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B)
1234 {
1235   Mat            C = *B;
1236   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data;
1237   PetscErrorCode ierr;
1238   PetscInt       i,j,mbs = a->mbs;
1239   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
1240   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
1241   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
1242   PetscReal      zeropivot,rs,shiftnz;
1243   PetscReal      shiftpd;
1244   ChShift_Ctx    sctx;
1245   PetscInt       newshift;
1246 
1247   PetscFunctionBegin;
1248   /* initialization */
1249   shiftnz   = info->shiftnz;
1250   shiftpd   = info->shiftpd;
1251   zeropivot = info->zeropivot;
1252 
1253   /* il and jl record the first nonzero element in each row of the accessing
1254      window U(0:k, k:mbs-1).
1255      jl:    list of rows to be added to uneliminated rows
1256             i>= k: jl(i) is the first row to be added to row i
1257             i<  k: jl(i) is the row following row i in some list of rows
1258             jl(i) = mbs indicates the end of a list
1259      il(i): points to the first nonzero element in U(i,k:mbs-1)
1260   */
1261   ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
1262   ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr);
1263   jl   = il + mbs;
1264 
1265   sctx.shift_amount = 0;
1266   sctx.nshift       = 0;
1267   do {
1268     sctx.chshift = PETSC_FALSE;
1269     for (i=0; i<mbs; i++) {
1270       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1271     }
1272 
1273     for (k = 0; k<mbs; k++){
1274       /*initialize k-th row with elements nonzero in row perm(k) of A */
1275       nz   = ai[k+1] - ai[k];
1276       acol = aj + ai[k];
1277       aval = aa + ai[k];
1278       bval = ba + bi[k];
1279       while (nz -- ){
1280         rtmp[*acol++] = *aval++;
1281         *bval++       = 0.0; /* for in-place factorization */
1282       }
1283 
1284       /* shift the diagonal of the matrix */
1285       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
1286 
1287       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1288       dk = rtmp[k];
1289       i  = jl[k]; /* first row to be added to k_th row  */
1290 
1291       while (i < k){
1292         nexti = jl[i]; /* next row to be added to k_th row */
1293         /* compute multiplier, update D(k) and U(i,k) */
1294         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1295         uikdi = - ba[ili]*ba[bi[i]];
1296         dk   += uikdi*ba[ili];
1297         ba[ili] = uikdi; /* -U(i,k) */
1298 
1299         /* add multiple of row i to k-th row ... */
1300         jmin = ili + 1;
1301         nz   = bi[i+1] - jmin;
1302         if (nz > 0){
1303           bcol = bj + jmin;
1304           bval = ba + jmin;
1305           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
1306           /* update il and jl for i-th row */
1307           il[i] = jmin;
1308           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1309         }
1310         i = nexti;
1311       }
1312 
1313       /* shift the diagonals when zero pivot is detected */
1314       /* compute rs=sum of abs(off-diagonal) */
1315       rs   = 0.0;
1316       jmin = bi[k]+1;
1317       nz   = bi[k+1] - jmin;
1318       if (nz){
1319         bcol = bj + jmin;
1320         while (nz--){
1321           rs += PetscAbsScalar(rtmp[*bcol]);
1322           bcol++;
1323         }
1324       }
1325 
1326       sctx.rs = rs;
1327       sctx.pv = dk;
1328       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
1329       if (newshift == 1) break;    /* sctx.shift_amount is updated */
1330 
1331       /* copy data into U(k,:) */
1332       ba[bi[k]] = 1.0/dk;
1333       jmin      = bi[k]+1;
1334       nz        = bi[k+1] - jmin;
1335       if (nz){
1336         bcol = bj + jmin;
1337         bval = ba + jmin;
1338         while (nz--){
1339           *bval++       = rtmp[*bcol];
1340           rtmp[*bcol++] = 0.0;
1341         }
1342         /* add k-th row into il and jl */
1343         il[k] = jmin;
1344         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
1345       }
1346     } /* end of for (k = 0; k<mbs; k++) */
1347   } while (sctx.chshift);
1348   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1349   ierr = PetscFree(il);CHKERRQ(ierr);
1350 
1351   C->factor       = FACTOR_CHOLESKY;
1352   C->assembled    = PETSC_TRUE;
1353   C->preallocated = PETSC_TRUE;
1354   ierr = PetscLogFlops(C->rmap.N);CHKERRQ(ierr);
1355   if (sctx.nshift){
1356     if (shiftnz) {
1357       ierr = PetscInfo2(0,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1358     } else if (shiftpd) {
1359       ierr = PetscInfo2(0,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1360     }
1361   }
1362   PetscFunctionReturn(0);
1363 }
1364 
1365 #undef __FUNCT__
1366 #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ"
1367 PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,MatFactorInfo *info)
1368 {
1369   PetscErrorCode ierr;
1370   Mat            C;
1371 
1372   PetscFunctionBegin;
1373   ierr = MatCholeskyFactorSymbolic(A,perm,info,&C);CHKERRQ(ierr);
1374   ierr = MatCholeskyFactorNumeric(A,info,&C);CHKERRQ(ierr);
1375   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
1376   PetscFunctionReturn(0);
1377 }
1378 
1379 
1380