xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact.c (revision c07bf0741b2e8ac3bb29c3d024f5883d65444abc)
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 extern PetscErrorCode MatSeqSBAIJSetNumericFactorization(Mat,PetscTruth);
44 
45 /*
46   Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP.
47   Use Modified Sparse Row (MSR) storage for u and ju. See page 85, "Iterative Methods ..." by Saad.
48 */
49 #undef __FUNCT__
50 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ_MSR"
51 PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(Mat F,Mat A,IS perm,const MatFactorInfo *info)
52 {
53   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b;
54   PetscErrorCode ierr;
55   const PetscInt *rip,*ai,*aj;
56   PetscInt       i,mbs = a->mbs,*jutmp,bs = A->rmap->bs,bs2=a->bs2;
57   PetscInt       m,reallocs = 0,prow;
58   PetscInt       *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd;
59   PetscReal      f = info->fill;
60   PetscTruth     perm_identity;
61 
62   PetscFunctionBegin;
63   /* check whether perm is the identity mapping */
64   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
65   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
66 
67   if (perm_identity){ /* without permutation */
68     a->permute = PETSC_FALSE;
69     ai = a->i; aj = a->j;
70   } else {            /* non-trivial permutation */
71     a->permute = PETSC_TRUE;
72     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
73     ai = a->inew; aj = a->jnew;
74   }
75 
76   /* initialization */
77   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&iu);CHKERRQ(ierr);
78   umax  = (PetscInt)(f*ai[mbs] + 1); umax += mbs + 1;
79   ierr  = PetscMalloc(umax*sizeof(PetscInt),&ju);CHKERRQ(ierr);
80   iu[0] = mbs+1;
81   juidx = mbs + 1; /* index for ju */
82   /* jl linked list for pivot row -- linked list for col index */
83   ierr  = PetscMalloc2(mbs,PetscInt,&jl,mbs,PetscInt,&q);CHKERRQ(ierr);
84   for (i=0; i<mbs; i++){
85     jl[i] = mbs;
86     q[i] = 0;
87   }
88 
89   /* for each row k */
90   for (k=0; k<mbs; k++){
91     for (i=0; i<mbs; i++) q[i] = 0;  /* to be removed! */
92     nzk  = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */
93     q[k] = mbs;
94     /* initialize nonzero structure of k-th row to row rip[k] of A */
95     jmin = ai[rip[k]] +1; /* exclude diag[k] */
96     jmax = ai[rip[k]+1];
97     for (j=jmin; j<jmax; j++){
98       vj = rip[aj[j]]; /* col. value */
99       if(vj > k){
100         qm = k;
101         do {
102           m  = qm; qm = q[m];
103         } while(qm < vj);
104         if (qm == vj) {
105           SETERRQ(PETSC_ERR_PLIB,"Duplicate entry in A\n");
106         }
107         nzk++;
108         q[m]  = vj;
109         q[vj] = qm;
110       } /* if(vj > k) */
111     } /* for (j=jmin; j<jmax; j++) */
112 
113     /* modify nonzero structure of k-th row by computing fill-in
114        for each row i to be merged in */
115     prow = k;
116     prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */
117 
118     while (prow < k){
119       /* merge row prow into k-th row */
120       jmin = iu[prow] + 1; jmax = iu[prow+1];
121       qm = k;
122       for (j=jmin; j<jmax; j++){
123         vj = ju[j];
124         do {
125           m = qm; qm = q[m];
126         } while (qm < vj);
127         if (qm != vj){
128          nzk++; q[m] = vj; q[vj] = qm; qm = vj;
129         }
130       }
131       prow = jl[prow]; /* next pivot row */
132     }
133 
134     /* add k to row list for first nonzero element in k-th row */
135     if (nzk > 0){
136       i = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */
137       jl[k] = jl[i]; jl[i] = k;
138     }
139     iu[k+1] = iu[k] + nzk;
140 
141     /* allocate more space to ju if needed */
142     if (iu[k+1] > umax) {
143       /* estimate how much additional space we will need */
144       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
145       /* just double the memory each time */
146       maxadd = umax;
147       if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2;
148       umax += maxadd;
149 
150       /* allocate a longer ju */
151       ierr = PetscMalloc(umax*sizeof(PetscInt),&jutmp);CHKERRQ(ierr);
152       ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(PetscInt));CHKERRQ(ierr);
153       ierr = PetscFree(ju);CHKERRQ(ierr);
154       ju   = jutmp;
155       reallocs++; /* count how many times we realloc */
156     }
157 
158     /* save nonzero structure of k-th row in ju */
159     i=k;
160     while (nzk --) {
161       i           = q[i];
162       ju[juidx++] = i;
163     }
164   }
165 
166 #if defined(PETSC_USE_INFO)
167   if (ai[mbs] != 0) {
168     PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
169     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
170     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
171     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
172     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
173   } else {
174     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
175   }
176 #endif
177 
178   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
179   ierr = PetscFree2(jl,q);CHKERRQ(ierr);
180 
181   /* put together the new matrix */
182   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(F,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
183 
184   /* ierr = PetscLogObjectParent(B,iperm);CHKERRQ(ierr); */
185   b = (Mat_SeqSBAIJ*)(F)->data;
186   b->singlemalloc = PETSC_FALSE;
187   b->free_a       = PETSC_TRUE;
188   b->free_ij       = PETSC_TRUE;
189   ierr = PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
190   b->j    = ju;
191   b->i    = iu;
192   b->diag = 0;
193   b->ilen = 0;
194   b->imax = 0;
195   b->row  = perm;
196   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
197   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
198   b->icol = perm;
199   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
200   ierr    = PetscMalloc((bs*mbs+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
201   /* In b structure:  Free imax, ilen, old a, old j.
202      Allocate idnew, solve_work, new a, new j */
203   ierr = PetscLogObjectMemory(F,(iu[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
204   b->maxnz = b->nz = iu[mbs];
205 
206   (F)->info.factor_mallocs    = reallocs;
207   (F)->info.fill_ratio_given  = f;
208   if (ai[mbs] != 0) {
209     (F)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
210   } else {
211     (F)->info.fill_ratio_needed = 0.0;
212   }
213   ierr = MatSeqSBAIJSetNumericFactorization(F,perm_identity);CHKERRQ(ierr);
214   PetscFunctionReturn(0);
215 }
216 /*
217     Symbolic U^T*D*U factorization for SBAIJ format.
218 */
219 #include "petscbt.h"
220 #include "../src/mat/utils/freespace.h"
221 #undef __FUNCT__
222 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ"
223 PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
224 {
225   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data;
226   Mat_SeqSBAIJ       *b;
227   PetscErrorCode     ierr;
228   PetscTruth         perm_identity,missing;
229   PetscReal          fill = info->fill;
230   const PetscInt     *rip,*ai,*aj;
231   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,reallocs=0,prow,d;
232   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
233   PetscInt           nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr;
234   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
235   PetscBT            lnkbt;
236 
237   PetscFunctionBegin;
238   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
239   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
240 
241   /*
242    This code originally uses Modified Sparse Row (MSR) storage
243    (see page 85, "Iterative Methods ..." by Saad) for the output matrix B - bad choise!
244    Then it is rewritten so the factor B takes seqsbaij format. However the associated
245    MatCholeskyFactorNumeric_() have not been modified for the cases of bs>1 or !perm_identity,
246    thus the original code in MSR format is still used for these cases.
247    The code below should replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR() whenever
248    MatCholeskyFactorNumeric_() is modified for using sbaij symbolic factor.
249   */
250   if (bs > 1){
251     ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(fact,A,perm,info);CHKERRQ(ierr);
252     PetscFunctionReturn(0);
253   }
254 
255   /* check whether perm is the identity mapping */
256   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
257 
258   if (perm_identity){
259     a->permute = PETSC_FALSE;
260     ai = a->i; aj = a->j;
261   } else {
262     SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format");
263     /* There are bugs for reordeing. Needs further work.
264        MatReordering for sbaij cannot be efficient. User should use aij formt! */
265     a->permute = PETSC_TRUE;
266     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
267     ai = a->inew; aj = a->jnew;
268   }
269   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
270 
271   /* initialization */
272   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
273   ui[0] = 0;
274 
275   /* jl: linked list for storing indices of the pivot rows
276      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
277   ierr = PetscMalloc4(mbs,PetscInt*,&ui_ptr,mbs,PetscInt,&il,mbs,PetscInt,&jl,mbs,PetscInt,&cols);CHKERRQ(ierr);
278   for (i=0; i<mbs; i++){
279     jl[i] = mbs; il[i] = 0;
280   }
281 
282   /* create and initialize a linked list for storing column indices of the active row k */
283   nlnk = mbs + 1;
284   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
285 
286   /* initial FreeSpace size is fill*(ai[mbs]+1) */
287   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
288   current_space = free_space;
289 
290   for (k=0; k<mbs; k++){  /* for each active row k */
291     /* initialize lnk by the column indices of row rip[k] of A */
292     nzk   = 0;
293     ncols = ai[rip[k]+1] - ai[rip[k]];
294     for (j=0; j<ncols; j++){
295       i = *(aj + ai[rip[k]] + j);
296       cols[j] = rip[i];
297     }
298     ierr = PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
299     nzk += nlnk;
300 
301     /* update lnk by computing fill-in for each pivot row to be merged in */
302     prow = jl[k]; /* 1st pivot row */
303 
304     while (prow < k){
305       nextprow = jl[prow];
306       /* merge prow into k-th row */
307       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
308       jmax = ui[prow+1];
309       ncols = jmax-jmin;
310       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
311       ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
312       nzk += nlnk;
313 
314       /* update il and jl for prow */
315       if (jmin < jmax){
316         il[prow] = jmin;
317         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
318       }
319       prow = nextprow;
320     }
321 
322     /* if free space is not available, make more free space */
323     if (current_space->local_remaining<nzk) {
324       i = mbs - k + 1; /* num of unfactored rows */
325       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
326       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
327       reallocs++;
328     }
329 
330     /* copy data into free space, then initialize lnk */
331     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
332 
333     /* add the k-th row into il and jl */
334     if (nzk-1 > 0){
335       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
336       jl[k] = jl[i]; jl[i] = k;
337       il[k] = ui[k] + 1;
338     }
339     ui_ptr[k] = current_space->array;
340     current_space->array           += nzk;
341     current_space->local_used      += nzk;
342     current_space->local_remaining -= nzk;
343 
344     ui[k+1] = ui[k] + nzk;
345   }
346 
347 #if defined(PETSC_USE_INFO)
348   if (ai[mbs] != 0) {
349     PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
350     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
351     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
352     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
353   } else {
354     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
355   }
356 #endif
357 
358   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
359   ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr);
360 
361   /* destroy list of free space and other temporary array(s) */
362   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
363   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
364   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
365 
366   /* put together the new matrix in MATSEQSBAIJ format */
367   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(fact,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
368 
369   b = (Mat_SeqSBAIJ*)(fact)->data;
370   b->singlemalloc = PETSC_FALSE;
371   b->free_a       = PETSC_TRUE;
372   b->free_ij      = PETSC_TRUE;
373   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
374   b->j    = uj;
375   b->i    = ui;
376   b->diag = 0;
377   b->ilen = 0;
378   b->imax = 0;
379   b->row  = perm;
380   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
381   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
382   b->icol = perm;
383   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
384   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
385   ierr    = PetscLogObjectMemory(fact,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
386   b->maxnz = b->nz = ui[mbs];
387 
388   (fact)->info.factor_mallocs    = reallocs;
389   (fact)->info.fill_ratio_given  = fill;
390   if (ai[mbs] != 0) {
391     (fact)->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
392   } else {
393     (fact)->info.fill_ratio_needed = 0.0;
394   }
395   ierr = MatSeqSBAIJSetNumericFactorization(fact,perm_identity);CHKERRQ(ierr);
396   PetscFunctionReturn(0);
397 }
398 #undef __FUNCT__
399 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N"
400 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat C,Mat A,const MatFactorInfo *info)
401 {
402   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
403   IS             perm = b->row;
404   PetscErrorCode ierr;
405   const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j;
406   PetscInt       i,j;
407   PetscInt       *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
408   PetscInt       bs=A->rmap->bs,bs2 = a->bs2,bslog = 0;
409   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
410   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
411   MatScalar      *work;
412   PetscInt       *pivots;
413 
414   PetscFunctionBegin;
415   /* initialization */
416   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
417   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
418   ierr = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
419   for (i=0; i<mbs; i++) {
420     jl[i] = mbs; il[0] = 0;
421   }
422   ierr = PetscMalloc3(bs2,MatScalar,&dk,bs2,MatScalar,&uik,bs,MatScalar,&work);CHKERRQ(ierr);
423   ierr = PetscMalloc(bs*sizeof(PetscInt),&pivots);CHKERRQ(ierr);
424 
425   ierr  = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
426 
427   /* check permutation */
428   if (!a->permute){
429     ai = a->i; aj = a->j; aa = a->a;
430   } else {
431     ai   = a->inew; aj = a->jnew;
432     ierr = PetscMalloc(bs2*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
433     ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
434     ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr);
435     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
436 
437     /* flops in while loop */
438     bslog = 2*bs*bs2;
439 
440     for (i=0; i<mbs; i++){
441       jmin = ai[i]; jmax = ai[i+1];
442       for (j=jmin; j<jmax; j++){
443         while (a2anew[j] != j){
444           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
445           for (k1=0; k1<bs2; k1++){
446             dk[k1]       = aa[k*bs2+k1];
447             aa[k*bs2+k1] = aa[j*bs2+k1];
448             aa[j*bs2+k1] = dk[k1];
449           }
450         }
451         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
452         if (i > aj[j]){
453           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
454           ap = aa + j*bs2;                     /* ptr to the beginning of j-th block of aa */
455           for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
456           for (k=0; k<bs; k++){               /* j-th block of aa <- dk^T */
457             for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1];
458           }
459         }
460       }
461     }
462     ierr = PetscFree(a2anew);CHKERRQ(ierr);
463   }
464 
465   /* for each row k */
466   for (k = 0; k<mbs; k++){
467 
468     /*initialize k-th row with elements nonzero in row perm(k) of A */
469     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
470 
471     ap = aa + jmin*bs2;
472     for (j = jmin; j < jmax; j++){
473       vj = perm_ptr[aj[j]];         /* block col. index */
474       rtmp_ptr = rtmp + vj*bs2;
475       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
476     }
477 
478     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
479     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
480     i = jl[k]; /* first row to be added to k_th row  */
481 
482     while (i < k){
483       nexti = jl[i]; /* next row to be added to k_th row */
484 
485       /* compute multiplier */
486       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
487 
488       /* uik = -inv(Di)*U_bar(i,k) */
489       diag = ba + i*bs2;
490       u    = ba + ili*bs2;
491       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
492       Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
493 
494       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
495       Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
496       ierr = PetscLogFlops(bslog*2.0);CHKERRQ(ierr);
497 
498       /* update -U(i,k) */
499       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
500 
501       /* add multiple of row i to k-th row ... */
502       jmin = ili + 1; jmax = bi[i+1];
503       if (jmin < jmax){
504         for (j=jmin; j<jmax; j++) {
505           /* rtmp += -U(i,k)^T * U_bar(i,j) */
506           rtmp_ptr = rtmp + bj[j]*bs2;
507           u = ba + j*bs2;
508           Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
509         }
510         ierr = PetscLogFlops(bslog*(jmax-jmin));CHKERRQ(ierr);
511 
512         /* ... add i to row list for next nonzero entry */
513         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
514         j     = bj[jmin];
515         jl[i] = jl[j]; jl[j] = i; /* update jl */
516       }
517       i = nexti;
518     }
519 
520     /* save nonzero entries in k-th row of U ... */
521 
522     /* invert diagonal block */
523     diag = ba+k*bs2;
524     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
525     ierr = Kernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr);
526 
527     jmin = bi[k]; jmax = bi[k+1];
528     if (jmin < jmax) {
529       for (j=jmin; j<jmax; j++){
530          vj = bj[j];           /* block col. index of U */
531          u   = ba + j*bs2;
532          rtmp_ptr = rtmp + vj*bs2;
533          for (k1=0; k1<bs2; k1++){
534            *u++        = *rtmp_ptr;
535            *rtmp_ptr++ = 0.0;
536          }
537       }
538 
539       /* ... add k to row list for first nonzero entry in k-th row */
540       il[k] = jmin;
541       i     = bj[jmin];
542       jl[k] = jl[i]; jl[i] = k;
543     }
544   }
545 
546   ierr = PetscFree(rtmp);CHKERRQ(ierr);
547   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
548   ierr = PetscFree3(dk,uik,work);CHKERRQ(ierr);
549   ierr = PetscFree(pivots);CHKERRQ(ierr);
550   if (a->permute){
551     ierr = PetscFree(aa);CHKERRQ(ierr);
552   }
553 
554   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
555   C->ops->solve          = MatSolve_SeqSBAIJ_N_inplace;
556   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_inplace;
557   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_inplace;
558   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_inplace;
559 
560   C->assembled    = PETSC_TRUE;
561   C->preallocated = PETSC_TRUE;
562   ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
563   PetscFunctionReturn(0);
564 }
565 
566 #undef __FUNCT__
567 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering"
568 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
569 {
570   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
571   PetscErrorCode ierr;
572   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
573   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
574   PetscInt       bs=A->rmap->bs,bs2 = a->bs2,bslog;
575   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
576   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
577   MatScalar      *work;
578   PetscInt       *pivots;
579 
580   PetscFunctionBegin;
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   ai = a->i; aj = a->j; aa = a->a;
591 
592   /* flops in while loop */
593   bslog = 2*bs*bs2;
594 
595   /* for each row k */
596   for (k = 0; k<mbs; k++){
597 
598     /*initialize k-th row with elements nonzero in row k of A */
599     jmin = ai[k]; jmax = ai[k+1];
600     ap = aa + jmin*bs2;
601     for (j = jmin; j < jmax; j++){
602       vj = aj[j];         /* block col. index */
603       rtmp_ptr = rtmp + vj*bs2;
604       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
605     }
606 
607     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
608     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
609     i = jl[k]; /* first row to be added to k_th row  */
610 
611     while (i < k){
612       nexti = jl[i]; /* next row to be added to k_th row */
613 
614       /* compute multiplier */
615       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
616 
617       /* uik = -inv(Di)*U_bar(i,k) */
618       diag = ba + i*bs2;
619       u    = ba + ili*bs2;
620       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
621       Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
622 
623       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
624       Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
625       ierr = PetscLogFlops(bslog*2.0);CHKERRQ(ierr);
626 
627       /* update -U(i,k) */
628       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
629 
630       /* add multiple of row i to k-th row ... */
631       jmin = ili + 1; jmax = bi[i+1];
632       if (jmin < jmax){
633         for (j=jmin; j<jmax; j++) {
634           /* rtmp += -U(i,k)^T * U_bar(i,j) */
635           rtmp_ptr = rtmp + bj[j]*bs2;
636           u = ba + j*bs2;
637           Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
638         }
639         ierr = PetscLogFlops(bslog*(jmax-jmin));CHKERRQ(ierr);
640 
641         /* ... add i to row list for next nonzero entry */
642         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
643         j     = bj[jmin];
644         jl[i] = jl[j]; jl[j] = i; /* update jl */
645       }
646       i = nexti;
647     }
648 
649     /* save nonzero entries in k-th row of U ... */
650 
651     /* invert diagonal block */
652     diag = ba+k*bs2;
653     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
654     ierr = Kernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr);
655 
656     jmin = bi[k]; jmax = bi[k+1];
657     if (jmin < jmax) {
658       for (j=jmin; j<jmax; j++){
659          vj = bj[j];           /* block col. index of U */
660          u   = ba + j*bs2;
661          rtmp_ptr = rtmp + vj*bs2;
662          for (k1=0; k1<bs2; k1++){
663            *u++        = *rtmp_ptr;
664            *rtmp_ptr++ = 0.0;
665          }
666       }
667 
668       /* ... add k to row list for first nonzero entry in k-th row */
669       il[k] = jmin;
670       i     = bj[jmin];
671       jl[k] = jl[i]; jl[i] = k;
672     }
673   }
674 
675   ierr = PetscFree(rtmp);CHKERRQ(ierr);
676   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
677   ierr = PetscFree3(dk,uik,work);CHKERRQ(ierr);
678   ierr = PetscFree(pivots);CHKERRQ(ierr);
679 
680   C->ops->solve          = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
681   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
682   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
683   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
684   C->assembled = PETSC_TRUE;
685   C->preallocated = PETSC_TRUE;
686   ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
687   PetscFunctionReturn(0);
688 }
689 
690 /*
691     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
692     Version for blocks 2 by 2.
693 */
694 #undef __FUNCT__
695 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2"
696 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat C,Mat A,const MatFactorInfo *info)
697 {
698   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
699   IS             perm = b->row;
700   PetscErrorCode ierr;
701   const PetscInt *ai,*aj,*perm_ptr;
702   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
703   PetscInt       *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
704   MatScalar      *ba = b->a,*aa,*ap;
705   MatScalar      *u,*diag,*rtmp,*rtmp_ptr,dk[4],uik[4];
706   PetscReal      shift = info->shiftinblocks;
707 
708   PetscFunctionBegin;
709   /* initialization */
710   /* il and jl record the first nonzero element in each row of the accessing
711      window U(0:k, k:mbs-1).
712      jl:    list of rows to be added to uneliminated rows
713             i>= k: jl(i) is the first row to be added to row i
714             i<  k: jl(i) is the row following row i in some list of rows
715             jl(i) = mbs indicates the end of a list
716      il(i): points to the first nonzero element in columns k,...,mbs-1 of
717             row i of U */
718   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
719   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
720   ierr = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
721   for (i=0; i<mbs; i++) {
722     jl[i] = mbs; il[0] = 0;
723   }
724   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
725 
726   /* check permutation */
727   if (!a->permute){
728     ai = a->i; aj = a->j; aa = a->a;
729   } else {
730     ai   = a->inew; aj = a->jnew;
731     ierr = PetscMalloc(4*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
732     ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
733     ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr);
734     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
735 
736     for (i=0; i<mbs; i++){
737       jmin = ai[i]; jmax = ai[i+1];
738       for (j=jmin; j<jmax; j++){
739         while (a2anew[j] != j){
740           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
741           for (k1=0; k1<4; k1++){
742             dk[k1]       = aa[k*4+k1];
743             aa[k*4+k1] = aa[j*4+k1];
744             aa[j*4+k1] = dk[k1];
745           }
746         }
747         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
748         if (i > aj[j]){
749           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
750           ap = aa + j*4;     /* ptr to the beginning of the block */
751           dk[1] = ap[1];     /* swap ap[1] and ap[2] */
752           ap[1] = ap[2];
753           ap[2] = dk[1];
754         }
755       }
756     }
757     ierr = PetscFree(a2anew);CHKERRQ(ierr);
758   }
759 
760   /* for each row k */
761   for (k = 0; k<mbs; k++){
762 
763     /*initialize k-th row with elements nonzero in row perm(k) of A */
764     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
765     ap = aa + jmin*4;
766     for (j = jmin; j < jmax; j++){
767       vj = perm_ptr[aj[j]];         /* block col. index */
768       rtmp_ptr = rtmp + vj*4;
769       for (i=0; i<4; 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*4,4*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): - ba[ili]*ba[i] */
783       diag = ba + i*4;
784       u    = ba + ili*4;
785       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
786       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
787       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
788       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
789 
790       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
791       dk[0] += uik[0]*u[0] + uik[1]*u[1];
792       dk[1] += uik[2]*u[0] + uik[3]*u[1];
793       dk[2] += uik[0]*u[2] + uik[1]*u[3];
794       dk[3] += uik[2]*u[2] + uik[3]*u[3];
795 
796       ierr = PetscLogFlops(16.0*2.0);CHKERRQ(ierr);
797 
798       /* update -U(i,k): ba[ili] = uik */
799       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
800 
801       /* add multiple of row i to k-th row ... */
802       jmin = ili + 1; jmax = bi[i+1];
803       if (jmin < jmax){
804         for (j=jmin; j<jmax; j++) {
805           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
806           rtmp_ptr = rtmp + bj[j]*4;
807           u = ba + j*4;
808           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
809           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
810           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
811           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
812         }
813         ierr = PetscLogFlops(16.0*(jmax-jmin));CHKERRQ(ierr);
814 
815         /* ... add i to row list for next nonzero entry */
816         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
817         j     = bj[jmin];
818         jl[i] = jl[j]; jl[j] = i; /* update jl */
819       }
820       i = nexti;
821     }
822 
823     /* save nonzero entries in k-th row of U ... */
824 
825     /* invert diagonal block */
826     diag = ba+k*4;
827     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
828     ierr = Kernel_A_gets_inverse_A_2(diag,shift);CHKERRQ(ierr);
829 
830     jmin = bi[k]; jmax = bi[k+1];
831     if (jmin < jmax) {
832       for (j=jmin; j<jmax; j++){
833          vj = bj[j];           /* block col. index of U */
834          u   = ba + j*4;
835          rtmp_ptr = rtmp + vj*4;
836          for (k1=0; k1<4; k1++){
837            *u++        = *rtmp_ptr;
838            *rtmp_ptr++ = 0.0;
839          }
840       }
841 
842       /* ... add k to row list for first nonzero entry in k-th row */
843       il[k] = jmin;
844       i     = bj[jmin];
845       jl[k] = jl[i]; jl[i] = k;
846     }
847   }
848 
849   ierr = PetscFree(rtmp);CHKERRQ(ierr);
850   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
851   if (a->permute) {
852     ierr = PetscFree(aa);CHKERRQ(ierr);
853   }
854   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
855   C->ops->solve          = MatSolve_SeqSBAIJ_2_inplace;
856   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_inplace;
857   C->assembled = PETSC_TRUE;
858   C->preallocated = PETSC_TRUE;
859   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
860   PetscFunctionReturn(0);
861 }
862 
863 /*
864       Version for when blocks are 2 by 2 Using natural ordering
865 */
866 #undef __FUNCT__
867 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering"
868 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
869 {
870   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
871   PetscErrorCode ierr;
872   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
873   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
874   MatScalar      *ba = b->a,*aa,*ap,dk[8],uik[8];
875   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
876   PetscReal      shift = info->shiftinblocks;
877 
878   PetscFunctionBegin;
879   /* initialization */
880   /* il and jl record the first nonzero element in each row of the accessing
881      window U(0:k, k:mbs-1).
882      jl:    list of rows to be added to uneliminated rows
883             i>= k: jl(i) is the first row to be added to row i
884             i<  k: jl(i) is the row following row i in some list of rows
885             jl(i) = mbs indicates the end of a list
886      il(i): points to the first nonzero element in columns k,...,mbs-1 of
887             row i of U */
888   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
889   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
890   ierr = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
891   for (i=0; i<mbs; i++) {
892     jl[i] = mbs; il[0] = 0;
893   }
894   ai = a->i; aj = a->j; aa = a->a;
895 
896   /* for each row k */
897   for (k = 0; k<mbs; k++){
898 
899     /*initialize k-th row with elements nonzero in row k of A */
900     jmin = ai[k]; jmax = ai[k+1];
901     ap = aa + jmin*4;
902     for (j = jmin; j < jmax; j++){
903       vj = aj[j];         /* block col. index */
904       rtmp_ptr = rtmp + vj*4;
905       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
906     }
907 
908     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
909     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
910     i = jl[k]; /* first row to be added to k_th row  */
911 
912     while (i < k){
913       nexti = jl[i]; /* next row to be added to k_th row */
914 
915       /* compute multiplier */
916       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
917 
918       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
919       diag = ba + i*4;
920       u    = ba + ili*4;
921       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
922       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
923       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
924       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
925 
926       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
927       dk[0] += uik[0]*u[0] + uik[1]*u[1];
928       dk[1] += uik[2]*u[0] + uik[3]*u[1];
929       dk[2] += uik[0]*u[2] + uik[1]*u[3];
930       dk[3] += uik[2]*u[2] + uik[3]*u[3];
931 
932       ierr = PetscLogFlops(16.0*2.0);CHKERRQ(ierr);
933 
934       /* update -U(i,k): ba[ili] = uik */
935       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
936 
937       /* add multiple of row i to k-th row ... */
938       jmin = ili + 1; jmax = bi[i+1];
939       if (jmin < jmax){
940         for (j=jmin; j<jmax; j++) {
941           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
942           rtmp_ptr = rtmp + bj[j]*4;
943           u = ba + j*4;
944           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
945           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
946           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
947           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
948         }
949         ierr = PetscLogFlops(16.0*(jmax-jmin));CHKERRQ(ierr);
950 
951         /* ... add i to row list for next nonzero entry */
952         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
953         j     = bj[jmin];
954         jl[i] = jl[j]; jl[j] = i; /* update jl */
955       }
956       i = nexti;
957     }
958 
959     /* save nonzero entries in k-th row of U ... */
960 
961     /* invert diagonal block */
962     diag = ba+k*4;
963     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
964     ierr = Kernel_A_gets_inverse_A_2(diag,shift);CHKERRQ(ierr);
965 
966     jmin = bi[k]; jmax = bi[k+1];
967     if (jmin < jmax) {
968       for (j=jmin; j<jmax; j++){
969          vj = bj[j];           /* block col. index of U */
970          u   = ba + j*4;
971          rtmp_ptr = rtmp + vj*4;
972          for (k1=0; k1<4; k1++){
973            *u++        = *rtmp_ptr;
974            *rtmp_ptr++ = 0.0;
975          }
976       }
977 
978       /* ... add k to row list for first nonzero entry in k-th row */
979       il[k] = jmin;
980       i     = bj[jmin];
981       jl[k] = jl[i]; jl[i] = k;
982     }
983   }
984 
985   ierr = PetscFree(rtmp);CHKERRQ(ierr);
986   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
987 
988   C->ops->solve          = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
989   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
990   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
991   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
992   C->assembled = PETSC_TRUE;
993   C->preallocated = PETSC_TRUE;
994   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
995   PetscFunctionReturn(0);
996 }
997 
998 /*
999     Numeric U^T*D*U factorization for SBAIJ format.
1000     Version for blocks are 1 by 1.
1001 */
1002 #undef __FUNCT__
1003 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1"
1004 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat C,Mat A,const MatFactorInfo *info)
1005 {
1006   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data;
1007   IS             ip=b->row;
1008   PetscErrorCode ierr;
1009   const PetscInt *ai,*aj,*rip;
1010   PetscInt       *a2anew,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j,*bcol;
1011   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1012   MatScalar      *rtmp,*ba=b->a,*bval,*aa,dk,uikdi;
1013   PetscReal      zeropivot,rs,shiftnz;
1014   PetscReal      shiftpd;
1015   ChShift_Ctx    sctx;
1016   PetscInt       newshift;
1017 
1018   PetscFunctionBegin;
1019   /* initialization */
1020   shiftnz   = info->shiftnz;
1021   shiftpd   = info->shiftpd;
1022   zeropivot = info->zeropivot;
1023 
1024   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
1025   if (!a->permute){
1026     ai = a->i; aj = a->j; aa = a->a;
1027   } else {
1028     ai = a->inew; aj = a->jnew;
1029     nz = ai[mbs];
1030     ierr = PetscMalloc(nz*sizeof(MatScalar),&aa);CHKERRQ(ierr);
1031     a2anew = a->a2anew;
1032     bval   = a->a;
1033     for (j=0; j<nz; j++){
1034       aa[a2anew[j]] = *(bval++);
1035     }
1036   }
1037 
1038   /* initialization */
1039   /* il and jl record the first nonzero element in each row of the accessing
1040      window U(0:k, k:mbs-1).
1041      jl:    list of rows to be added to uneliminated rows
1042             i>= k: jl(i) is the first row to be added to row i
1043             i<  k: jl(i) is the row following row i in some list of rows
1044             jl(i) = mbs indicates the end of a list
1045      il(i): points to the first nonzero element in columns k,...,mbs-1 of
1046             row i of U */
1047   ierr = PetscMalloc3(mbs,MatScalar,&rtmp,mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
1048 
1049   sctx.shift_amount = 0;
1050   sctx.nshift       = 0;
1051   do {
1052     sctx.chshift = PETSC_FALSE;
1053     for (i=0; i<mbs; i++) {
1054       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1055     }
1056 
1057     for (k = 0; k<mbs; k++){
1058       /*initialize k-th row by the perm[k]-th row of A */
1059       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
1060       bval = ba + bi[k];
1061       for (j = jmin; j < jmax; j++){
1062         col = rip[aj[j]];
1063         rtmp[col] = aa[j];
1064         *bval++  = 0.0; /* for in-place factorization */
1065       }
1066 
1067       /* shift the diagonal of the matrix */
1068       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
1069 
1070       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1071       dk = rtmp[k];
1072       i = jl[k]; /* first row to be added to k_th row  */
1073 
1074       while (i < k){
1075         nexti = jl[i]; /* next row to be added to k_th row */
1076 
1077         /* compute multiplier, update diag(k) and U(i,k) */
1078         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1079         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
1080         dk += uikdi*ba[ili];
1081         ba[ili] = uikdi; /* -U(i,k) */
1082 
1083         /* add multiple of row i to k-th row */
1084         jmin = ili + 1; jmax = bi[i+1];
1085         if (jmin < jmax){
1086           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1087           ierr = PetscLogFlops(2.0*(jmax-jmin));CHKERRQ(ierr);
1088 
1089           /* update il and jl for row i */
1090           il[i] = jmin;
1091           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1092         }
1093         i = nexti;
1094       }
1095 
1096       /* shift the diagonals when zero pivot is detected */
1097       /* compute rs=sum of abs(off-diagonal) */
1098       rs   = 0.0;
1099       jmin = bi[k]+1;
1100       nz   = bi[k+1] - jmin;
1101       if (nz){
1102         bcol = bj + jmin;
1103         while (nz--){
1104           rs += PetscAbsScalar(rtmp[*bcol]);
1105           bcol++;
1106         }
1107       }
1108 
1109       sctx.rs = rs;
1110       sctx.pv = dk;
1111       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
1112       if (newshift == 1) break;    /* sctx.shift_amount is updated */
1113 
1114       /* copy data into U(k,:) */
1115       ba[bi[k]] = 1.0/dk; /* U(k,k) */
1116       jmin = bi[k]+1; jmax = bi[k+1];
1117       if (jmin < jmax) {
1118         for (j=jmin; j<jmax; j++){
1119           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
1120         }
1121         /* add the k-th row into il and jl */
1122         il[k] = jmin;
1123         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
1124       }
1125     }
1126   } while (sctx.chshift);
1127   ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr);
1128   if (a->permute){ierr = PetscFree(aa);CHKERRQ(ierr);}
1129 
1130   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
1131   C->ops->solve          = MatSolve_SeqSBAIJ_1_inplace;
1132   C->ops->solves         = MatSolves_SeqSBAIJ_1_inplace;
1133   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace;
1134   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_inplace;
1135   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_inplace;
1136   C->assembled    = PETSC_TRUE;
1137   C->preallocated = PETSC_TRUE;
1138   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
1139   if (sctx.nshift){
1140     if (shiftnz) {
1141       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1142     } else if (shiftpd) {
1143       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1144     }
1145   }
1146   PetscFunctionReturn(0);
1147 }
1148 
1149 /*
1150   Version for when blocks are 1 by 1 Using natural ordering
1151 */
1152 #undef __FUNCT__
1153 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering"
1154 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
1155 {
1156   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data;
1157   PetscErrorCode ierr;
1158   PetscInt       i,j,mbs = a->mbs;
1159   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
1160   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
1161   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
1162   PetscReal      zeropivot,rs,shiftnz;
1163   PetscReal      shiftpd;
1164   ChShift_Ctx    sctx;
1165   PetscInt       newshift;
1166 
1167   PetscFunctionBegin;
1168   /* initialization */
1169   shiftnz   = info->shiftnz;
1170   shiftpd   = info->shiftpd;
1171   zeropivot = info->zeropivot;
1172 
1173   /* il and jl record the first nonzero element in each row of the accessing
1174      window U(0:k, k:mbs-1).
1175      jl:    list of rows to be added to uneliminated rows
1176             i>= k: jl(i) is the first row to be added to row i
1177             i<  k: jl(i) is the row following row i in some list of rows
1178             jl(i) = mbs indicates the end of a list
1179      il(i): points to the first nonzero element in U(i,k:mbs-1)
1180   */
1181   ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
1182   ierr = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
1183 
1184   sctx.shift_amount = 0;
1185   sctx.nshift       = 0;
1186   do {
1187     sctx.chshift = PETSC_FALSE;
1188     for (i=0; i<mbs; i++) {
1189       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1190     }
1191 
1192     for (k = 0; k<mbs; k++){
1193       /*initialize k-th row with elements nonzero in row perm(k) of A */
1194       nz   = ai[k+1] - ai[k];
1195       acol = aj + ai[k];
1196       aval = aa + ai[k];
1197       bval = ba + bi[k];
1198       while (nz -- ){
1199         rtmp[*acol++] = *aval++;
1200         *bval++       = 0.0; /* for in-place factorization */
1201       }
1202 
1203       /* shift the diagonal of the matrix */
1204       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
1205 
1206       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1207       dk = rtmp[k];
1208       i  = jl[k]; /* first row to be added to k_th row  */
1209 
1210       while (i < k){
1211         nexti = jl[i]; /* next row to be added to k_th row */
1212         /* compute multiplier, update D(k) and U(i,k) */
1213         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1214         uikdi = - ba[ili]*ba[bi[i]];
1215         dk   += uikdi*ba[ili];
1216         ba[ili] = uikdi; /* -U(i,k) */
1217 
1218         /* add multiple of row i to k-th row ... */
1219         jmin = ili + 1;
1220         nz   = bi[i+1] - jmin;
1221         if (nz > 0){
1222           bcol = bj + jmin;
1223           bval = ba + jmin;
1224           ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr);
1225           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
1226 
1227           /* update il and jl for i-th row */
1228           il[i] = jmin;
1229           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1230         }
1231         i = nexti;
1232       }
1233 
1234       /* shift the diagonals when zero pivot is detected */
1235       /* compute rs=sum of abs(off-diagonal) */
1236       rs   = 0.0;
1237       jmin = bi[k]+1;
1238       nz   = bi[k+1] - jmin;
1239       if (nz){
1240         bcol = bj + jmin;
1241         while (nz--){
1242           rs += PetscAbsScalar(rtmp[*bcol]);
1243           bcol++;
1244         }
1245       }
1246 
1247       sctx.rs = rs;
1248       sctx.pv = dk;
1249       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
1250       if (newshift == 1) break;    /* sctx.shift_amount is updated */
1251 
1252       /* copy data into U(k,:) */
1253       ba[bi[k]] = 1.0/dk;
1254       jmin      = bi[k]+1;
1255       nz        = bi[k+1] - jmin;
1256       if (nz){
1257         bcol = bj + jmin;
1258         bval = ba + jmin;
1259         while (nz--){
1260           *bval++       = rtmp[*bcol];
1261           rtmp[*bcol++] = 0.0;
1262         }
1263         /* add k-th row into il and jl */
1264         il[k] = jmin;
1265         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
1266       }
1267     } /* end of for (k = 0; k<mbs; k++) */
1268   } while (sctx.chshift);
1269   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1270   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
1271 
1272   C->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1273   C->ops->solves         = MatSolves_SeqSBAIJ_1_inplace;
1274   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1275   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1276   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1277 
1278   C->assembled    = PETSC_TRUE;
1279   C->preallocated = PETSC_TRUE;
1280   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
1281   if (sctx.nshift){
1282     if (shiftnz) {
1283       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1284     } else if (shiftpd) {
1285       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1286     }
1287   }
1288   PetscFunctionReturn(0);
1289 }
1290 
1291 #undef __FUNCT__
1292 #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ"
1293 PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,const MatFactorInfo *info)
1294 {
1295   PetscErrorCode ierr;
1296   Mat            C;
1297 
1298   PetscFunctionBegin;
1299   ierr = MatGetFactor(A,"petsc",MAT_FACTOR_CHOLESKY,&C);CHKERRQ(ierr);
1300   ierr = MatCholeskyFactorSymbolic(C,A,perm,info);CHKERRQ(ierr);
1301   ierr = MatCholeskyFactorNumeric(C,A,info);CHKERRQ(ierr);
1302   A->ops->solve            = C->ops->solve;
1303   A->ops->solvetranspose   = C->ops->solvetranspose;
1304   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
1305   PetscFunctionReturn(0);
1306 }
1307 
1308 
1309