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