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