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