xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact.c (revision 5f8bbcca2696d61632b174eeae8a5433a128d797)
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(PetscRealIntMultTruncate(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    = PetscIntMultTruncate(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 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format");
428   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
429 
430   /* initialization */
431   ierr  = PetscMalloc1(mbs+1,&ui);CHKERRQ(ierr);
432   ui[0] = 0;
433 
434   /* jl: linked list for storing indices of the pivot rows
435      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
436   ierr = PetscMalloc4(mbs,&ui_ptr,mbs,&il,mbs,&jl,mbs,&cols);CHKERRQ(ierr);
437   for (i=0; i<mbs; i++) {
438     jl[i] = mbs; il[i] = 0;
439   }
440 
441   /* create and initialize a linked list for storing column indices of the active row k */
442   nlnk = mbs + 1;
443   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
444 
445   /* initial FreeSpace size is fill*(ai[mbs]+1) */
446   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,ai[mbs]+1),&free_space);CHKERRQ(ierr);
447   current_space = free_space;
448 
449   for (k=0; k<mbs; k++) {  /* for each active row k */
450     /* initialize lnk by the column indices of row rip[k] of A */
451     nzk   = 0;
452     ncols = ai[rip[k]+1] - ai[rip[k]];
453     for (j=0; j<ncols; j++) {
454       i       = *(aj + ai[rip[k]] + j);
455       cols[j] = rip[i];
456     }
457     ierr = PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
458     nzk += nlnk;
459 
460     /* update lnk by computing fill-in for each pivot row to be merged in */
461     prow = jl[k]; /* 1st pivot row */
462 
463     while (prow < k) {
464       nextprow = jl[prow];
465       /* merge prow into k-th row */
466       jmin   = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
467       jmax   = ui[prow+1];
468       ncols  = jmax-jmin;
469       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
470       ierr   = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
471       nzk   += nlnk;
472 
473       /* update il and jl for prow */
474       if (jmin < jmax) {
475         il[prow] = jmin;
476 
477         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
478       }
479       prow = nextprow;
480     }
481 
482     /* if free space is not available, make more free space */
483     if (current_space->local_remaining<nzk) {
484       i    = mbs - k + 1; /* num of unfactored rows */
485       i    = PetscMin(PetscIntMultTruncate(i,nzk), PetscIntMultTruncate(i,i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
486       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
487       reallocs++;
488     }
489 
490     /* copy data into free space, then initialize lnk */
491     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
492 
493     /* add the k-th row into il and jl */
494     if (nzk-1 > 0) {
495       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
496       jl[k] = jl[i]; jl[i] = k;
497       il[k] = ui[k] + 1;
498     }
499     ui_ptr[k] = current_space->array;
500 
501     current_space->array           += nzk;
502     current_space->local_used      += nzk;
503     current_space->local_remaining -= nzk;
504 
505     ui[k+1] = ui[k] + nzk;
506   }
507 
508   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
509   ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr);
510 
511   /* destroy list of free space and other temporary array(s) */
512   ierr = PetscMalloc1(ui[mbs]+1,&uj);CHKERRQ(ierr);
513   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
514   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
515 
516   /* put together the new matrix in MATSEQSBAIJ format */
517   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
518 
519   b               = (Mat_SeqSBAIJ*)fact->data;
520   b->singlemalloc = PETSC_FALSE;
521   b->free_a       = PETSC_TRUE;
522   b->free_ij      = PETSC_TRUE;
523 
524   ierr = PetscMalloc1(ui[mbs]+1,&b->a);CHKERRQ(ierr);
525 
526   b->j    = uj;
527   b->i    = ui;
528   b->diag = 0;
529   b->ilen = 0;
530   b->imax = 0;
531   b->row  = perm;
532 
533   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
534 
535   ierr     = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
536   b->icol  = perm;
537   ierr     = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
538   ierr     = PetscMalloc1(mbs+1,&b->solve_work);CHKERRQ(ierr);
539   ierr     = PetscLogObjectMemory((PetscObject)fact,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
540   b->maxnz = b->nz = ui[mbs];
541 
542   fact->info.factor_mallocs   = reallocs;
543   fact->info.fill_ratio_given = fill;
544   if (ai[mbs] != 0) {
545     fact->info.fill_ratio_needed = ((PetscReal)ui[mbs])/ai[mbs];
546   } else {
547     fact->info.fill_ratio_needed = 0.0;
548   }
549 #if defined(PETSC_USE_INFO)
550   if (ai[mbs] != 0) {
551     PetscReal af = fact->info.fill_ratio_needed;
552     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);CHKERRQ(ierr);
553     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
554     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);CHKERRQ(ierr);
555   } else {
556     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
557   }
558 #endif
559   ierr = MatSeqSBAIJSetNumericFactorization_inplace(fact,perm_identity);CHKERRQ(ierr);
560   PetscFunctionReturn(0);
561 }
562 
563 #undef __FUNCT__
564 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N"
565 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat C,Mat A,const MatFactorInfo *info)
566 {
567   Mat_SeqSBAIJ   *a   = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
568   IS             perm = b->row;
569   PetscErrorCode ierr;
570   const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j;
571   PetscInt       i,j;
572   PetscInt       *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
573   PetscInt       bs  =A->rmap->bs,bs2 = a->bs2;
574   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
575   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
576   MatScalar      *work;
577   PetscInt       *pivots;
578   PetscBool      allowzeropivot,zeropivotdetected;
579 
580   PetscFunctionBegin;
581   /* initialization */
582   ierr = PetscCalloc1(bs2*mbs,&rtmp);CHKERRQ(ierr);
583   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
584   allowzeropivot = PetscNot(A->erroriffailure);
585 
586   for (i=0; i<mbs; i++) {
587     jl[i] = mbs; il[0] = 0;
588   }
589   ierr = PetscMalloc3(bs2,&dk,bs2,&uik,bs,&work);CHKERRQ(ierr);
590   ierr = PetscMalloc1(bs,&pivots);CHKERRQ(ierr);
591 
592   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
593 
594   /* check permutation */
595   if (!a->permute) {
596     ai = a->i; aj = a->j; aa = a->a;
597   } else {
598     ai   = a->inew; aj = a->jnew;
599     ierr = PetscMalloc1(bs2*ai[mbs],&aa);CHKERRQ(ierr);
600     ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
601     ierr = PetscMalloc1(ai[mbs],&a2anew);CHKERRQ(ierr);
602     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
603 
604     for (i=0; i<mbs; i++) {
605       jmin = ai[i]; jmax = ai[i+1];
606       for (j=jmin; j<jmax; j++) {
607         while (a2anew[j] != j) {
608           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
609           for (k1=0; k1<bs2; k1++) {
610             dk[k1]       = aa[k*bs2+k1];
611             aa[k*bs2+k1] = aa[j*bs2+k1];
612             aa[j*bs2+k1] = dk[k1];
613           }
614         }
615         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
616         if (i > aj[j]) {
617           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
618           ap = aa + j*bs2;                     /* ptr to the beginning of j-th block of aa */
619           for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
620           for (k=0; k<bs; k++) {               /* j-th block of aa <- dk^T */
621             for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1];
622           }
623         }
624       }
625     }
626     ierr = PetscFree(a2anew);CHKERRQ(ierr);
627   }
628 
629   /* for each row k */
630   for (k = 0; k<mbs; k++) {
631 
632     /*initialize k-th row with elements nonzero in row perm(k) of A */
633     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
634 
635     ap = aa + jmin*bs2;
636     for (j = jmin; j < jmax; j++) {
637       vj       = perm_ptr[aj[j]];   /* block col. index */
638       rtmp_ptr = rtmp + vj*bs2;
639       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
640     }
641 
642     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
643     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
644     i    = jl[k]; /* first row to be added to k_th row  */
645 
646     while (i < k) {
647       nexti = jl[i]; /* next row to be added to k_th row */
648 
649       /* compute multiplier */
650       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
651 
652       /* uik = -inv(Di)*U_bar(i,k) */
653       diag = ba + i*bs2;
654       u    = ba + ili*bs2;
655       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
656       PetscKernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
657 
658       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
659       PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
660       ierr = PetscLogFlops(4.0*bs*bs2);CHKERRQ(ierr);
661 
662       /* update -U(i,k) */
663       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
664 
665       /* add multiple of row i to k-th row ... */
666       jmin = ili + 1; jmax = bi[i+1];
667       if (jmin < jmax) {
668         for (j=jmin; j<jmax; j++) {
669           /* rtmp += -U(i,k)^T * U_bar(i,j) */
670           rtmp_ptr = rtmp + bj[j]*bs2;
671           u        = ba + j*bs2;
672           PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
673         }
674         ierr = PetscLogFlops(2.0*bs*bs2*(jmax-jmin));CHKERRQ(ierr);
675 
676         /* ... add i to row list for next nonzero entry */
677         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
678         j     = bj[jmin];
679         jl[i] = jl[j]; jl[j] = i; /* update jl */
680       }
681       i = nexti;
682     }
683 
684     /* save nonzero entries in k-th row of U ... */
685 
686     /* invert diagonal block */
687     diag = ba+k*bs2;
688     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
689 
690     ierr = PetscKernel_A_gets_inverse_A(bs,diag,pivots,work,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr);
691     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
692 
693     jmin = bi[k]; jmax = bi[k+1];
694     if (jmin < jmax) {
695       for (j=jmin; j<jmax; j++) {
696         vj       = bj[j];      /* block col. index of U */
697         u        = ba + j*bs2;
698         rtmp_ptr = rtmp + vj*bs2;
699         for (k1=0; k1<bs2; k1++) {
700           *u++        = *rtmp_ptr;
701           *rtmp_ptr++ = 0.0;
702         }
703       }
704 
705       /* ... add k to row list for first nonzero entry in k-th row */
706       il[k] = jmin;
707       i     = bj[jmin];
708       jl[k] = jl[i]; jl[i] = k;
709     }
710   }
711 
712   ierr = PetscFree(rtmp);CHKERRQ(ierr);
713   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
714   ierr = PetscFree3(dk,uik,work);CHKERRQ(ierr);
715   ierr = PetscFree(pivots);CHKERRQ(ierr);
716   if (a->permute) {
717     ierr = PetscFree(aa);CHKERRQ(ierr);
718   }
719 
720   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
721 
722   C->ops->solve          = MatSolve_SeqSBAIJ_N_inplace;
723   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_inplace;
724   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_inplace;
725   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_inplace;
726 
727   C->assembled    = PETSC_TRUE;
728   C->preallocated = PETSC_TRUE;
729 
730   ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
731   PetscFunctionReturn(0);
732 }
733 
734 #undef __FUNCT__
735 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering"
736 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
737 {
738   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
739   PetscErrorCode ierr;
740   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
741   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
742   PetscInt       bs  =A->rmap->bs,bs2 = a->bs2;
743   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
744   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
745   MatScalar      *work;
746   PetscInt       *pivots;
747   PetscBool      allowzeropivot,zeropivotdetected;
748 
749   PetscFunctionBegin;
750   ierr = PetscCalloc1(bs2*mbs,&rtmp);CHKERRQ(ierr);
751   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
752   for (i=0; i<mbs; i++) {
753     jl[i] = mbs; il[0] = 0;
754   }
755   ierr = PetscMalloc3(bs2,&dk,bs2,&uik,bs,&work);CHKERRQ(ierr);
756   ierr = PetscMalloc1(bs,&pivots);CHKERRQ(ierr);
757   allowzeropivot = PetscNot(A->erroriffailure);
758 
759   ai = a->i; aj = a->j; aa = a->a;
760 
761   /* for each row k */
762   for (k = 0; k<mbs; k++) {
763 
764     /*initialize k-th row with elements nonzero in row k of A */
765     jmin = ai[k]; jmax = ai[k+1];
766     ap   = aa + jmin*bs2;
767     for (j = jmin; j < jmax; j++) {
768       vj       = aj[j];   /* block col. index */
769       rtmp_ptr = rtmp + vj*bs2;
770       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
771     }
772 
773     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
774     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
775     i    = jl[k]; /* first row to be added to k_th row  */
776 
777     while (i < k) {
778       nexti = jl[i]; /* next row to be added to k_th row */
779 
780       /* compute multiplier */
781       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
782 
783       /* uik = -inv(Di)*U_bar(i,k) */
784       diag = ba + i*bs2;
785       u    = ba + ili*bs2;
786       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
787       PetscKernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
788 
789       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
790       PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
791       ierr = PetscLogFlops(2.0*bs*bs2);CHKERRQ(ierr);
792 
793       /* update -U(i,k) */
794       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
795 
796       /* add multiple of row i to k-th row ... */
797       jmin = ili + 1; jmax = bi[i+1];
798       if (jmin < jmax) {
799         for (j=jmin; j<jmax; j++) {
800           /* rtmp += -U(i,k)^T * U_bar(i,j) */
801           rtmp_ptr = rtmp + bj[j]*bs2;
802           u        = ba + j*bs2;
803           PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
804         }
805         ierr = PetscLogFlops(2.0*bs*bs2*(jmax-jmin));CHKERRQ(ierr);
806 
807         /* ... add i to row list for next nonzero entry */
808         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
809         j     = bj[jmin];
810         jl[i] = jl[j]; jl[j] = i; /* update jl */
811       }
812       i = nexti;
813     }
814 
815     /* save nonzero entries in k-th row of U ... */
816 
817     /* invert diagonal block */
818     diag = ba+k*bs2;
819     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
820 
821     ierr = PetscKernel_A_gets_inverse_A(bs,diag,pivots,work,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr);
822     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
823 
824     jmin = bi[k]; jmax = bi[k+1];
825     if (jmin < jmax) {
826       for (j=jmin; j<jmax; j++) {
827         vj       = bj[j];      /* block col. index of U */
828         u        = ba + j*bs2;
829         rtmp_ptr = rtmp + vj*bs2;
830         for (k1=0; k1<bs2; k1++) {
831           *u++        = *rtmp_ptr;
832           *rtmp_ptr++ = 0.0;
833         }
834       }
835 
836       /* ... add k to row list for first nonzero entry in k-th row */
837       il[k] = jmin;
838       i     = bj[jmin];
839       jl[k] = jl[i]; jl[i] = k;
840     }
841   }
842 
843   ierr = PetscFree(rtmp);CHKERRQ(ierr);
844   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
845   ierr = PetscFree3(dk,uik,work);CHKERRQ(ierr);
846   ierr = PetscFree(pivots);CHKERRQ(ierr);
847 
848   C->ops->solve          = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
849   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
850   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
851   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
852   C->assembled           = PETSC_TRUE;
853   C->preallocated        = PETSC_TRUE;
854 
855   ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
856   PetscFunctionReturn(0);
857 }
858 
859 /*
860     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
861     Version for blocks 2 by 2.
862 */
863 #undef __FUNCT__
864 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2"
865 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat C,Mat A,const MatFactorInfo *info)
866 {
867   Mat_SeqSBAIJ   *a   = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
868   IS             perm = b->row;
869   PetscErrorCode ierr;
870   const PetscInt *ai,*aj,*perm_ptr;
871   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
872   PetscInt       *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
873   MatScalar      *ba = b->a,*aa,*ap;
874   MatScalar      *u,*diag,*rtmp,*rtmp_ptr,dk[4],uik[4];
875   PetscReal      shift = info->shiftamount;
876   PetscBool      allowzeropivot,zeropivotdetected;
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     allowzeropivot = PetscNot(A->erroriffailure);
998     ierr = PetscKernel_A_gets_inverse_A_2(diag,shift,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr);
999     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1000 
1001     jmin = bi[k]; jmax = bi[k+1];
1002     if (jmin < jmax) {
1003       for (j=jmin; j<jmax; j++) {
1004         vj       = bj[j];      /* block col. index of U */
1005         u        = ba + j*4;
1006         rtmp_ptr = rtmp + vj*4;
1007         for (k1=0; k1<4; k1++) {
1008           *u++        = *rtmp_ptr;
1009           *rtmp_ptr++ = 0.0;
1010         }
1011       }
1012 
1013       /* ... add k to row list for first nonzero entry in k-th row */
1014       il[k] = jmin;
1015       i     = bj[jmin];
1016       jl[k] = jl[i]; jl[i] = k;
1017     }
1018   }
1019 
1020   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1021   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
1022   if (a->permute) {
1023     ierr = PetscFree(aa);CHKERRQ(ierr);
1024   }
1025   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
1026 
1027   C->ops->solve          = MatSolve_SeqSBAIJ_2_inplace;
1028   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_inplace;
1029   C->assembled           = PETSC_TRUE;
1030   C->preallocated        = PETSC_TRUE;
1031 
1032   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
1033   PetscFunctionReturn(0);
1034 }
1035 
1036 /*
1037       Version for when blocks are 2 by 2 Using natural ordering
1038 */
1039 #undef __FUNCT__
1040 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering"
1041 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
1042 {
1043   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
1044   PetscErrorCode ierr;
1045   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
1046   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
1047   MatScalar      *ba = b->a,*aa,*ap,dk[8],uik[8];
1048   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
1049   PetscReal      shift = info->shiftamount;
1050   PetscBool      allowzeropivot,zeropivotdetected;
1051 
1052   PetscFunctionBegin;
1053   /* initialization */
1054   /* il and jl record the first nonzero element in each row of the accessing
1055      window U(0:k, k:mbs-1).
1056      jl:    list of rows to be added to uneliminated rows
1057             i>= k: jl(i) is the first row to be added to row i
1058             i<  k: jl(i) is the row following row i in some list of rows
1059             jl(i) = mbs indicates the end of a list
1060      il(i): points to the first nonzero element in columns k,...,mbs-1 of
1061             row i of U */
1062   ierr = PetscCalloc1(4*mbs,&rtmp);CHKERRQ(ierr);
1063   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
1064   for (i=0; i<mbs; i++) {
1065     jl[i] = mbs; il[0] = 0;
1066   }
1067   ai = a->i; aj = a->j; aa = a->a;
1068 
1069   /* for each row k */
1070   for (k = 0; k<mbs; k++) {
1071 
1072     /*initialize k-th row with elements nonzero in row k of A */
1073     jmin = ai[k]; jmax = ai[k+1];
1074     ap   = aa + jmin*4;
1075     for (j = jmin; j < jmax; j++) {
1076       vj       = aj[j];   /* block col. index */
1077       rtmp_ptr = rtmp + vj*4;
1078       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
1079     }
1080 
1081     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
1082     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
1083     i    = jl[k]; /* first row to be added to k_th row  */
1084 
1085     while (i < k) {
1086       nexti = jl[i]; /* next row to be added to k_th row */
1087 
1088       /* compute multiplier */
1089       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1090 
1091       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
1092       diag   = ba + i*4;
1093       u      = ba + ili*4;
1094       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
1095       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
1096       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
1097       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
1098 
1099       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
1100       dk[0] += uik[0]*u[0] + uik[1]*u[1];
1101       dk[1] += uik[2]*u[0] + uik[3]*u[1];
1102       dk[2] += uik[0]*u[2] + uik[1]*u[3];
1103       dk[3] += uik[2]*u[2] + uik[3]*u[3];
1104 
1105       ierr = PetscLogFlops(16.0*2.0);CHKERRQ(ierr);
1106 
1107       /* update -U(i,k): ba[ili] = uik */
1108       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
1109 
1110       /* add multiple of row i to k-th row ... */
1111       jmin = ili + 1; jmax = bi[i+1];
1112       if (jmin < jmax) {
1113         for (j=jmin; j<jmax; j++) {
1114           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
1115           rtmp_ptr     = rtmp + bj[j]*4;
1116           u            = ba + j*4;
1117           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
1118           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
1119           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
1120           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
1121         }
1122         ierr = PetscLogFlops(16.0*(jmax-jmin));CHKERRQ(ierr);
1123 
1124         /* ... add i to row list for next nonzero entry */
1125         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
1126         j     = bj[jmin];
1127         jl[i] = jl[j]; jl[j] = i; /* update jl */
1128       }
1129       i = nexti;
1130     }
1131 
1132     /* save nonzero entries in k-th row of U ... */
1133 
1134     /* invert diagonal block */
1135     diag = ba+k*4;
1136     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
1137     allowzeropivot = PetscNot(A->erroriffailure);
1138     ierr = PetscKernel_A_gets_inverse_A_2(diag,shift,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr);
1139     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1140 
1141     jmin = bi[k]; jmax = bi[k+1];
1142     if (jmin < jmax) {
1143       for (j=jmin; j<jmax; j++) {
1144         vj       = bj[j];      /* block col. index of U */
1145         u        = ba + j*4;
1146         rtmp_ptr = rtmp + vj*4;
1147         for (k1=0; k1<4; k1++) {
1148           *u++        = *rtmp_ptr;
1149           *rtmp_ptr++ = 0.0;
1150         }
1151       }
1152 
1153       /* ... add k to row list for first nonzero entry in k-th row */
1154       il[k] = jmin;
1155       i     = bj[jmin];
1156       jl[k] = jl[i]; jl[i] = k;
1157     }
1158   }
1159 
1160   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1161   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
1162 
1163   C->ops->solve          = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
1164   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
1165   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
1166   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
1167   C->assembled           = PETSC_TRUE;
1168   C->preallocated        = PETSC_TRUE;
1169 
1170   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
1171   PetscFunctionReturn(0);
1172 }
1173 
1174 /*
1175     Numeric U^T*D*U factorization for SBAIJ format.
1176     Version for blocks are 1 by 1.
1177 */
1178 #undef __FUNCT__
1179 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace"
1180 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace(Mat C,Mat A,const MatFactorInfo *info)
1181 {
1182   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ*)C->data;
1183   IS             ip=b->row;
1184   PetscErrorCode ierr;
1185   const PetscInt *ai,*aj,*rip;
1186   PetscInt       *a2anew,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j,*bcol;
1187   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1188   MatScalar      *rtmp,*ba=b->a,*bval,*aa,dk,uikdi;
1189   PetscReal      rs;
1190   FactorShiftCtx sctx;
1191 
1192   PetscFunctionBegin;
1193   /* MatPivotSetUp(): initialize shift context sctx */
1194   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
1195 
1196   ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr);
1197   if (!a->permute) {
1198     ai = a->i; aj = a->j; aa = a->a;
1199   } else {
1200     ai     = a->inew; aj = a->jnew;
1201     nz     = ai[mbs];
1202     ierr   = PetscMalloc1(nz,&aa);CHKERRQ(ierr);
1203     a2anew = a->a2anew;
1204     bval   = a->a;
1205     for (j=0; j<nz; j++) {
1206       aa[a2anew[j]] = *(bval++);
1207     }
1208   }
1209 
1210   /* initialization */
1211   /* il and jl record the first nonzero element in each row of the accessing
1212      window U(0:k, k:mbs-1).
1213      jl:    list of rows to be added to uneliminated rows
1214             i>= k: jl(i) is the first row to be added to row i
1215             i<  k: jl(i) is the row following row i in some list of rows
1216             jl(i) = mbs indicates the end of a list
1217      il(i): points to the first nonzero element in columns k,...,mbs-1 of
1218             row i of U */
1219   ierr = PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&jl);CHKERRQ(ierr);
1220 
1221   do {
1222     sctx.newshift = PETSC_FALSE;
1223     for (i=0; i<mbs; i++) {
1224       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1225     }
1226 
1227     for (k = 0; k<mbs; k++) {
1228       /*initialize k-th row by the perm[k]-th row of A */
1229       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
1230       bval = ba + bi[k];
1231       for (j = jmin; j < jmax; j++) {
1232         col       = rip[aj[j]];
1233         rtmp[col] = aa[j];
1234         *bval++   = 0.0; /* for in-place factorization */
1235       }
1236 
1237       /* shift the diagonal of the matrix */
1238       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
1239 
1240       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1241       dk = rtmp[k];
1242       i  = jl[k]; /* first row to be added to k_th row  */
1243 
1244       while (i < k) {
1245         nexti = jl[i]; /* next row to be added to k_th row */
1246 
1247         /* compute multiplier, update diag(k) and U(i,k) */
1248         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
1249         uikdi   = -ba[ili]*ba[bi[i]]; /* diagonal(k) */
1250         dk     += uikdi*ba[ili];
1251         ba[ili] = uikdi; /* -U(i,k) */
1252 
1253         /* add multiple of row i to k-th row */
1254         jmin = ili + 1; jmax = bi[i+1];
1255         if (jmin < jmax) {
1256           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1257           ierr = PetscLogFlops(2.0*(jmax-jmin));CHKERRQ(ierr);
1258 
1259           /* update il and jl for row i */
1260           il[i] = jmin;
1261           j     = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1262         }
1263         i = nexti;
1264       }
1265 
1266       /* shift the diagonals when zero pivot is detected */
1267       /* compute rs=sum of abs(off-diagonal) */
1268       rs   = 0.0;
1269       jmin = bi[k]+1;
1270       nz   = bi[k+1] - jmin;
1271       if (nz) {
1272         bcol = bj + jmin;
1273         while (nz--) {
1274           rs += PetscAbsScalar(rtmp[*bcol]);
1275           bcol++;
1276         }
1277       }
1278 
1279       sctx.rs = rs;
1280       sctx.pv = dk;
1281       ierr    = MatPivotCheck(C,A,info,&sctx,k);CHKERRQ(ierr);
1282       if (sctx.newshift) break;    /* sctx.shift_amount is updated */
1283       dk = sctx.pv;
1284 
1285       /* copy data into U(k,:) */
1286       ba[bi[k]] = 1.0/dk; /* U(k,k) */
1287       jmin      = bi[k]+1; jmax = bi[k+1];
1288       if (jmin < jmax) {
1289         for (j=jmin; j<jmax; j++) {
1290           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
1291         }
1292         /* add the k-th row into il and jl */
1293         il[k] = jmin;
1294         i     = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
1295       }
1296     }
1297   } while (sctx.newshift);
1298   ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr);
1299   if (a->permute) {ierr = PetscFree(aa);CHKERRQ(ierr);}
1300 
1301   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
1302 
1303   C->ops->solve          = MatSolve_SeqSBAIJ_1_inplace;
1304   C->ops->solves         = MatSolves_SeqSBAIJ_1_inplace;
1305   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace;
1306   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_inplace;
1307   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_inplace;
1308   C->assembled           = PETSC_TRUE;
1309   C->preallocated        = PETSC_TRUE;
1310 
1311   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
1312   if (sctx.nshift) {
1313     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
1314       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
1315     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
1316       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
1317     }
1318   }
1319   PetscFunctionReturn(0);
1320 }
1321 
1322 /*
1323   Version for when blocks are 1 by 1 Using natural ordering under new datastructure
1324   Modified from MatCholeskyFactorNumeric_SeqAIJ()
1325 */
1326 #undef __FUNCT__
1327 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering"
1328 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
1329 {
1330   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data;
1331   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)B->data;
1332   PetscErrorCode ierr;
1333   PetscInt       i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bdiag=b->diag,*bjtmp;
1334   PetscInt       *ai=a->i,*aj=a->j,*ajtmp;
1335   PetscInt       k,jmin,jmax,*c2r,*il,col,nexti,ili,nz;
1336   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1337   FactorShiftCtx sctx;
1338   PetscReal      rs;
1339   MatScalar      d,*v;
1340 
1341   PetscFunctionBegin;
1342   ierr = PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&c2r);CHKERRQ(ierr);
1343 
1344   /* MatPivotSetUp(): initialize shift context sctx */
1345   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
1346 
1347   if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */
1348     sctx.shift_top = info->zeropivot;
1349 
1350     ierr = PetscMemzero(rtmp,mbs*sizeof(MatScalar));CHKERRQ(ierr);
1351 
1352     for (i=0; i<mbs; i++) {
1353       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
1354       d        = (aa)[a->diag[i]];
1355       rtmp[i] += -PetscRealPart(d);  /* diagonal entry */
1356       ajtmp    = aj + ai[i] + 1;     /* exclude diagonal */
1357       v        = aa + ai[i] + 1;
1358       nz       = ai[i+1] - ai[i] - 1;
1359       for (j=0; j<nz; j++) {
1360         rtmp[i]        += PetscAbsScalar(v[j]);
1361         rtmp[ajtmp[j]] += PetscAbsScalar(v[j]);
1362       }
1363       if (PetscRealPart(rtmp[i]) > sctx.shift_top) sctx.shift_top = PetscRealPart(rtmp[i]);
1364     }
1365     sctx.shift_top *= 1.1;
1366     sctx.nshift_max = 5;
1367     sctx.shift_lo   = 0.;
1368     sctx.shift_hi   = 1.;
1369   }
1370 
1371   /* allocate working arrays
1372      c2r: linked list, keep track of pivot rows for a given column. c2r[col]: head of the list for a given col
1373      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
1374   */
1375   do {
1376     sctx.newshift = PETSC_FALSE;
1377 
1378     for (i=0; i<mbs; i++) c2r[i] = mbs;
1379     if (mbs) il[0] = 0;
1380 
1381     for (k = 0; k<mbs; k++) {
1382       /* zero rtmp */
1383       nz    = bi[k+1] - bi[k];
1384       bjtmp = bj + bi[k];
1385       for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
1386 
1387       /* load in initial unfactored row */
1388       bval = ba + bi[k];
1389       jmin = ai[k]; jmax = ai[k+1];
1390       for (j = jmin; j < jmax; j++) {
1391         col       = aj[j];
1392         rtmp[col] = aa[j];
1393         *bval++   = 0.0; /* for in-place factorization */
1394       }
1395       /* shift the diagonal of the matrix: ZeropivotApply() */
1396       rtmp[k] += sctx.shift_amount;  /* shift the diagonal of the matrix */
1397 
1398       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1399       dk = rtmp[k];
1400       i  = c2r[k]; /* first row to be added to k_th row  */
1401 
1402       while (i < k) {
1403         nexti = c2r[i]; /* next row to be added to k_th row */
1404 
1405         /* compute multiplier, update diag(k) and U(i,k) */
1406         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
1407         uikdi   = -ba[ili]*ba[bdiag[i]]; /* diagonal(k) */
1408         dk     += uikdi*ba[ili]; /* update diag[k] */
1409         ba[ili] = uikdi; /* -U(i,k) */
1410 
1411         /* add multiple of row i to k-th row */
1412         jmin = ili + 1; jmax = bi[i+1];
1413         if (jmin < jmax) {
1414           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1415           /* update il and c2r for row i */
1416           il[i] = jmin;
1417           j     = bj[jmin]; c2r[i] = c2r[j]; c2r[j] = i;
1418         }
1419         i = nexti;
1420       }
1421 
1422       /* copy data into U(k,:) */
1423       rs   = 0.0;
1424       jmin = bi[k]; jmax = bi[k+1]-1;
1425       if (jmin < jmax) {
1426         for (j=jmin; j<jmax; j++) {
1427           col = bj[j]; ba[j] = rtmp[col]; rs += PetscAbsScalar(ba[j]);
1428         }
1429         /* add the k-th row into il and c2r */
1430         il[k] = jmin;
1431         i     = bj[jmin]; c2r[k] = c2r[i]; c2r[i] = k;
1432       }
1433 
1434       sctx.rs = rs;
1435       sctx.pv = dk;
1436       ierr    = MatPivotCheck(B,A,info,&sctx,k);CHKERRQ(ierr);
1437       if (sctx.newshift) break;
1438       dk = sctx.pv;
1439 
1440       ba[bdiag[k]] = 1.0/dk; /* U(k,k) */
1441     }
1442   } while (sctx.newshift);
1443 
1444   ierr = PetscFree3(rtmp,il,c2r);CHKERRQ(ierr);
1445 
1446   B->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1447   B->ops->solves         = MatSolves_SeqSBAIJ_1;
1448   B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1449   B->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering;
1450   B->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering;
1451 
1452   B->assembled    = PETSC_TRUE;
1453   B->preallocated = PETSC_TRUE;
1454 
1455   ierr = PetscLogFlops(B->rmap->n);CHKERRQ(ierr);
1456 
1457   /* MatPivotView() */
1458   if (sctx.nshift) {
1459     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
1460       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);
1461     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
1462       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
1463     } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) {
1464       ierr = PetscInfo2(A,"number of shift_inblocks applied %D, each shift_amount %g\n",sctx.nshift,(double)info->shiftamount);CHKERRQ(ierr);
1465     }
1466   }
1467   PetscFunctionReturn(0);
1468 }
1469 
1470 #undef __FUNCT__
1471 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace"
1472 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
1473 {
1474   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ*)C->data;
1475   PetscErrorCode ierr;
1476   PetscInt       i,j,mbs = a->mbs;
1477   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
1478   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
1479   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
1480   PetscReal      rs;
1481   FactorShiftCtx sctx;
1482 
1483   PetscFunctionBegin;
1484   /* MatPivotSetUp(): initialize shift context sctx */
1485   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
1486 
1487   /* initialization */
1488   /* il and jl record the first nonzero element in each row of the accessing
1489      window U(0:k, k:mbs-1).
1490      jl:    list of rows to be added to uneliminated rows
1491             i>= k: jl(i) is the first row to be added to row i
1492             i<  k: jl(i) is the row following row i in some list of rows
1493             jl(i) = mbs indicates the end of a list
1494      il(i): points to the first nonzero element in U(i,k:mbs-1)
1495   */
1496   ierr = PetscMalloc1(mbs,&rtmp);CHKERRQ(ierr);
1497   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
1498 
1499   do {
1500     sctx.newshift = PETSC_FALSE;
1501     for (i=0; i<mbs; i++) {
1502       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1503     }
1504 
1505     for (k = 0; k<mbs; k++) {
1506       /*initialize k-th row with elements nonzero in row perm(k) of A */
1507       nz   = ai[k+1] - ai[k];
1508       acol = aj + ai[k];
1509       aval = aa + ai[k];
1510       bval = ba + bi[k];
1511       while (nz--) {
1512         rtmp[*acol++] = *aval++;
1513         *bval++       = 0.0; /* for in-place factorization */
1514       }
1515 
1516       /* shift the diagonal of the matrix */
1517       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
1518 
1519       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1520       dk = rtmp[k];
1521       i  = jl[k]; /* first row to be added to k_th row  */
1522 
1523       while (i < k) {
1524         nexti = jl[i]; /* next row to be added to k_th row */
1525         /* compute multiplier, update D(k) and U(i,k) */
1526         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
1527         uikdi   = -ba[ili]*ba[bi[i]];
1528         dk     += uikdi*ba[ili];
1529         ba[ili] = uikdi; /* -U(i,k) */
1530 
1531         /* add multiple of row i to k-th row ... */
1532         jmin = ili + 1;
1533         nz   = bi[i+1] - jmin;
1534         if (nz > 0) {
1535           bcol = bj + jmin;
1536           bval = ba + jmin;
1537           ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr);
1538           while (nz--) rtmp[*bcol++] += uikdi*(*bval++);
1539 
1540           /* update il and jl for i-th row */
1541           il[i] = jmin;
1542           j     = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1543         }
1544         i = nexti;
1545       }
1546 
1547       /* shift the diagonals when zero pivot is detected */
1548       /* compute rs=sum of abs(off-diagonal) */
1549       rs   = 0.0;
1550       jmin = bi[k]+1;
1551       nz   = bi[k+1] - jmin;
1552       if (nz) {
1553         bcol = bj + jmin;
1554         while (nz--) {
1555           rs += PetscAbsScalar(rtmp[*bcol]);
1556           bcol++;
1557         }
1558       }
1559 
1560       sctx.rs = rs;
1561       sctx.pv = dk;
1562       ierr    = MatPivotCheck(C,A,info,&sctx,k);CHKERRQ(ierr);
1563       if (sctx.newshift) break;    /* sctx.shift_amount is updated */
1564       dk = sctx.pv;
1565 
1566       /* copy data into U(k,:) */
1567       ba[bi[k]] = 1.0/dk;
1568       jmin      = bi[k]+1;
1569       nz        = bi[k+1] - jmin;
1570       if (nz) {
1571         bcol = bj + jmin;
1572         bval = ba + jmin;
1573         while (nz--) {
1574           *bval++       = rtmp[*bcol];
1575           rtmp[*bcol++] = 0.0;
1576         }
1577         /* add k-th row into il and jl */
1578         il[k] = jmin;
1579         i     = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
1580       }
1581     } /* end of for (k = 0; k<mbs; k++) */
1582   } while (sctx.newshift);
1583   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1584   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
1585 
1586   C->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1587   C->ops->solves         = MatSolves_SeqSBAIJ_1_inplace;
1588   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1589   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1590   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1591 
1592   C->assembled    = PETSC_TRUE;
1593   C->preallocated = PETSC_TRUE;
1594 
1595   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
1596   if (sctx.nshift) {
1597     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
1598       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
1599     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
1600       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
1601     }
1602   }
1603   PetscFunctionReturn(0);
1604 }
1605 
1606 #undef __FUNCT__
1607 #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ"
1608 PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,const MatFactorInfo *info)
1609 {
1610   PetscErrorCode ierr;
1611   Mat            C;
1612 
1613   PetscFunctionBegin;
1614   ierr = MatGetFactor(A,"petsc",MAT_FACTOR_CHOLESKY,&C);CHKERRQ(ierr);
1615   ierr = MatCholeskyFactorSymbolic(C,A,perm,info);CHKERRQ(ierr);
1616   ierr = MatCholeskyFactorNumeric(C,A,info);CHKERRQ(ierr);
1617 
1618   A->ops->solve          = C->ops->solve;
1619   A->ops->solvetranspose = C->ops->solvetranspose;
1620 
1621   ierr = MatHeaderMerge(A,&C);CHKERRQ(ierr);
1622   PetscFunctionReturn(0);
1623 }
1624 
1625 
1626