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