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