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