xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact.c (revision 5fb5fb7741ad1a6b2a484b31956b9f64eccaa44a)
1 /*$Id: sbaijfact.c,v 1.61 2001/08/06 21:15:47 bsmith Exp $*/
2 /* Using Modified Sparse Row (MSR) storage.
3 See page 85, "Iterative Methods ..." by Saad. */
4 
5 /*
6     Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP.
7 */
8 #include "src/mat/impls/baij/seq/baij.h"
9 #include "src/mat/impls/sbaij/seq/sbaij.h"
10 #include "src/vec/vecimpl.h"
11 #include "src/inline/ilu.h"
12 #include "include/petscis.h"
13 
14 /* Use Modified Sparse Row storage for u and ju, see Sasd pp.85 */
15 #undef __FUNCT__
16 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ"
17 int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat A,IS perm,PetscReal f,Mat *B)
18 {
19   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b;
20   int          *rip,ierr,i,mbs = a->mbs,*ai,*aj;
21   int          *jutmp,bs = a->bs,bs2=a->bs2;
22   int          m,realloc = 0,prow;
23   int          *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd;
24   PetscTruth   perm_identity;
25 
26   PetscFunctionBegin;
27 
28   /* check whether perm is the identity mapping */
29   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
30   if (!perm_identity) a->permute = PETSC_TRUE;
31 
32   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
33 
34   if (perm_identity){ /* without permutation */
35     ai = a->i; aj = a->j;
36   } else {            /* non-trivial permutation */
37     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
38     ai = a->inew; aj = a->jnew;
39   }
40 
41   /* initialization */
42   ierr  = PetscMalloc((mbs+1)*sizeof(int),&iu);CHKERRQ(ierr);
43   umax  = (int)(f*ai[mbs] + 1); umax += mbs + 1;
44   ierr  = PetscMalloc(umax*sizeof(int),&ju);CHKERRQ(ierr);
45   iu[0] = mbs+1;
46   juidx = mbs + 1; /* index for ju */
47   ierr  = PetscMalloc(2*mbs*sizeof(int),&jl);CHKERRQ(ierr); /* linked list for pivot row */
48   q     = jl + mbs;   /* linked list for col index */
49   for (i=0; i<mbs; i++){
50     jl[i] = mbs;
51     q[i] = 0;
52   }
53 
54   /* for each row k */
55   for (k=0; k<mbs; k++){
56     nzk  = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */
57     q[k] = mbs;
58     /* initialize nonzero structure of k-th row to row rip[k] of A */
59     jmin = ai[rip[k]];
60     jmax = ai[rip[k]+1];
61     for (j=jmin; j<jmax; j++){
62       vj = rip[aj[j]]; /* col. value */
63       if(vj > k){
64         qm = k;
65         do {
66           m  = qm; qm = q[m];
67         } while(qm < vj);
68         if (qm == vj) {
69           SETERRQ(1," error: duplicate entry in A\n");
70         }
71         nzk++;
72         q[m]  = vj;
73         q[vj] = qm;
74       } /* if(vj > k) */
75     } /* for (j=jmin; j<jmax; j++) */
76 
77     /* modify nonzero structure of k-th row by computing fill-in
78        for each row i to be merged in */
79     prow = k;
80     prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */
81 
82     while (prow < k){
83       /* merge row prow into k-th row */
84       jmin = iu[prow] + 1; jmax = iu[prow+1];
85       qm = k;
86       for (j=jmin; j<jmax; j++){
87         vj = ju[j];
88         do {
89           m = qm; qm = q[m];
90         } while (qm < vj);
91         if (qm != vj){
92          nzk++; q[m] = vj; q[vj] = qm; qm = vj;
93         }
94       }
95       prow = jl[prow]; /* next pivot row */
96     }
97 
98     /* add k to row list for first nonzero element in k-th row */
99     if (nzk > 0){
100       i = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */
101       jl[k] = jl[i]; jl[i] = k;
102     }
103     iu[k+1] = iu[k] + nzk;
104 
105     /* allocate more space to ju if needed */
106     if (iu[k+1] > umax) {
107       /* estimate how much additional space we will need */
108       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
109       /* just double the memory each time */
110       maxadd = umax;
111       if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2;
112       umax += maxadd;
113 
114       /* allocate a longer ju */
115       ierr = PetscMalloc(umax*sizeof(int),&jutmp);CHKERRQ(ierr);
116       ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(int));CHKERRQ(ierr);
117       ierr = PetscFree(ju);CHKERRQ(ierr);
118       ju   = jutmp;
119       realloc++; /* count how many times we realloc */
120     }
121 
122     /* save nonzero structure of k-th row in ju */
123     i=k;
124     while (nzk --) {
125       i           = q[i];
126       ju[juidx++] = i;
127     }
128   }
129 
130   if (ai[mbs] != 0) {
131     PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
132     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af);
133     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Run with -pc_cholesky_fill %g or use \n",af);
134     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:PCCholeskySetFill(pc,%g);\n",af);
135     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:for best performance.\n");
136   } else {
137      PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Empty matrix.\n");
138   }
139 
140   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
141   /* ierr = PetscFree(q);CHKERRQ(ierr); */
142   ierr = PetscFree(jl);CHKERRQ(ierr);
143 
144   /* put together the new matrix */
145   ierr = MatCreateSeqSBAIJ(A->comm,bs,bs*mbs,bs*mbs,0,PETSC_NULL,B);CHKERRQ(ierr);
146   /* PetscLogObjectParent(*B,iperm); */
147   b = (Mat_SeqSBAIJ*)(*B)->data;
148   ierr = PetscFree(b->imax);CHKERRQ(ierr);
149   b->singlemalloc = PETSC_FALSE;
150   /* the next line frees the default space generated by the Create() */
151   ierr = PetscFree(b->a);CHKERRQ(ierr);
152   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
153   ierr = PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
154   b->j    = ju;
155   b->i    = iu;
156   b->diag = 0;
157   b->ilen = 0;
158   b->imax = 0;
159   b->row  = perm;
160   b->pivotinblocks = PETSC_FALSE; /* need to get from MatCholeskyInfo */
161   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
162   b->icol = perm;
163   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
164   ierr    = PetscMalloc((bs*mbs+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
165   /* In b structure:  Free imax, ilen, old a, old j.
166      Allocate idnew, solve_work, new a, new j */
167   PetscLogObjectMemory(*B,(iu[mbs]-mbs)*(sizeof(int)+sizeof(MatScalar)));
168   b->s_maxnz = b->s_nz = iu[mbs];
169 
170   (*B)->factor                 = FACTOR_CHOLESKY;
171   (*B)->info.factor_mallocs    = realloc;
172   (*B)->info.fill_ratio_given  = f;
173   if (ai[mbs] != 0) {
174     (*B)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
175   } else {
176     (*B)->info.fill_ratio_needed = 0.0;
177   }
178 
179   if (perm_identity){
180     switch (bs) {
181       case 1:
182         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering;
183         (*B)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
184         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=1\n");
185         break;
186       case 2:
187         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
188         (*B)->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
189         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
190         break;
191       case 3:
192         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
193         (*B)->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
194         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:sing special in-place natural ordering factor and solve BS=3\n");
195         break;
196       case 4:
197         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
198         (*B)->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
199         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
200         break;
201       case 5:
202         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
203         (*B)->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
204         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
205         break;
206       case 6:
207         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
208         (*B)->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
209         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
210         break;
211       case 7:
212         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
213         (*B)->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
214         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
215       break;
216       default:
217         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering;
218         (*B)->ops->solve           = MatSolve_SeqSBAIJ_N_NaturalOrdering;
219         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS>7\n");
220       break;
221     }
222   }
223 
224   PetscFunctionReturn(0);
225 }
226 
227 
228 #undef __FUNCT__
229 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N"
230 int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat A,Mat *B)
231 {
232   Mat                C = *B;
233   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
234   IS                 perm = b->row;
235   int                *perm_ptr,ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
236   int                *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
237   int                bs=a->bs,bs2 = a->bs2;
238   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
239   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
240   MatScalar          *work;
241   int                *pivots;
242 
243   PetscFunctionBegin;
244 
245   /* initialization */
246   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
247   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
248   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
249   jl   = il + mbs;
250   for (i=0; i<mbs; i++) {
251     jl[i] = mbs; il[0] = 0;
252   }
253   ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr);
254   uik  = dk + bs2;
255   work = uik + bs2;
256   ierr = PetscMalloc(bs*sizeof(int),&pivots);CHKERRQ(ierr);
257 
258   ierr  = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
259 
260   /* check permutation */
261   if (!a->permute){
262     ai = a->i; aj = a->j; aa = a->a;
263   } else {
264     ai   = a->inew; aj = a->jnew;
265     ierr = PetscMalloc(bs2*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
266     ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
267     ierr = PetscMalloc(ai[mbs]*sizeof(int),&a2anew);CHKERRQ(ierr);
268     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr);
269 
270     for (i=0; i<mbs; i++){
271       jmin = ai[i]; jmax = ai[i+1];
272       for (j=jmin; j<jmax; j++){
273         while (a2anew[j] != j){
274           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
275           for (k1=0; k1<bs2; k1++){
276             dk[k1]       = aa[k*bs2+k1];
277             aa[k*bs2+k1] = aa[j*bs2+k1];
278             aa[j*bs2+k1] = dk[k1];
279           }
280         }
281         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
282         if (i > aj[j]){
283           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
284           ap = aa + j*bs2;                     /* ptr to the beginning of j-th block of aa */
285           for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
286           for (k=0; k<bs; k++){               /* j-th block of aa <- dk^T */
287             for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1];
288           }
289         }
290       }
291     }
292     ierr = PetscFree(a2anew);CHKERRQ(ierr);
293   }
294 
295   /* for each row k */
296   for (k = 0; k<mbs; k++){
297 
298     /*initialize k-th row with elements nonzero in row perm(k) of A */
299     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
300     if (jmin < jmax) {
301       ap = aa + jmin*bs2;
302       for (j = jmin; j < jmax; j++){
303         vj = perm_ptr[aj[j]];         /* block col. index */
304         rtmp_ptr = rtmp + vj*bs2;
305         for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
306       }
307     }
308 
309     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
310     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
311     i = jl[k]; /* first row to be added to k_th row  */
312 
313     while (i < mbs){
314       nexti = jl[i]; /* next row to be added to k_th row */
315 
316       /* compute multiplier */
317       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
318 
319       /* uik = -inv(Di)*U_bar(i,k) */
320       diag = ba + i*bs2;
321       u    = ba + ili*bs2;
322       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
323       Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
324 
325       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
326       Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
327 
328       /* update -U(i,k) */
329       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
330 
331       /* add multiple of row i to k-th row ... */
332       jmin = ili + 1; jmax = bi[i+1];
333       if (jmin < jmax){
334         for (j=jmin; j<jmax; j++) {
335           /* rtmp += -U(i,k)^T * U_bar(i,j) */
336           rtmp_ptr = rtmp + bj[j]*bs2;
337           u = ba + j*bs2;
338           Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
339         }
340 
341         /* ... add i to row list for next nonzero entry */
342         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
343         j     = bj[jmin];
344         jl[i] = jl[j]; jl[j] = i; /* update jl */
345       }
346       i = nexti;
347     }
348 
349     /* save nonzero entries in k-th row of U ... */
350 
351     /* invert diagonal block */
352     diag = ba+k*bs2;
353     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
354     Kernel_A_gets_inverse_A(bs,diag,pivots,work);
355 
356     jmin = bi[k]; jmax = bi[k+1];
357     if (jmin < jmax) {
358       for (j=jmin; j<jmax; j++){
359          vj = bj[j];           /* block col. index of U */
360          u   = ba + j*bs2;
361          rtmp_ptr = rtmp + vj*bs2;
362          for (k1=0; k1<bs2; k1++){
363            *u++        = *rtmp_ptr;
364            *rtmp_ptr++ = 0.0;
365          }
366       }
367 
368       /* ... add k to row list for first nonzero entry in k-th row */
369       il[k] = jmin;
370       i     = bj[jmin];
371       jl[k] = jl[i]; jl[i] = k;
372     }
373   }
374 
375   ierr = PetscFree(rtmp);CHKERRQ(ierr);
376   ierr = PetscFree(il);CHKERRQ(ierr);
377   ierr = PetscFree(dk);CHKERRQ(ierr);
378   ierr = PetscFree(pivots);CHKERRQ(ierr);
379   if (a->permute){
380     ierr = PetscFree(aa);CHKERRQ(ierr);
381   }
382 
383   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
384   C->factor    = FACTOR_CHOLESKY;
385   C->assembled = PETSC_TRUE;
386   C->preallocated = PETSC_TRUE;
387   PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */
388   PetscFunctionReturn(0);
389 }
390 
391 #undef __FUNCT__
392 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering"
393 int MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat A,Mat *B)
394 {
395   Mat                C = *B;
396   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
397   int                ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
398   int                *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
399   int                bs=a->bs,bs2 = a->bs2;
400   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
401   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
402   MatScalar          *work;
403   int                *pivots;
404 
405   PetscFunctionBegin;
406 
407   /* initialization */
408 
409   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
410   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
411   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
412   jl   = il + mbs;
413   for (i=0; i<mbs; i++) {
414     jl[i] = mbs; il[0] = 0;
415   }
416   ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr);
417   uik  = dk + bs2;
418   work = uik + bs2;
419   ierr = PetscMalloc(bs*sizeof(int),&pivots);CHKERRQ(ierr);
420 
421   ai = a->i; aj = a->j; aa = a->a;
422 
423   /* for each row k */
424   for (k = 0; k<mbs; k++){
425 
426     /*initialize k-th row with elements nonzero in row k of A */
427     jmin = ai[k]; jmax = ai[k+1];
428     if (jmin < jmax) {
429       ap = aa + jmin*bs2;
430       for (j = jmin; j < jmax; j++){
431         vj = aj[j];         /* block col. index */
432         rtmp_ptr = rtmp + vj*bs2;
433         for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
434       }
435     }
436 
437     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
438     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
439     i = jl[k]; /* first row to be added to k_th row  */
440 
441     while (i < mbs){
442       nexti = jl[i]; /* next row to be added to k_th row */
443 
444       /* compute multiplier */
445       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
446 
447       /* uik = -inv(Di)*U_bar(i,k) */
448       diag = ba + i*bs2;
449       u    = ba + ili*bs2;
450       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
451       Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
452 
453       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
454       Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
455 
456       /* update -U(i,k) */
457       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
458 
459       /* add multiple of row i to k-th row ... */
460       jmin = ili + 1; jmax = bi[i+1];
461       if (jmin < jmax){
462         for (j=jmin; j<jmax; j++) {
463           /* rtmp += -U(i,k)^T * U_bar(i,j) */
464           rtmp_ptr = rtmp + bj[j]*bs2;
465           u = ba + j*bs2;
466           Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
467         }
468 
469         /* ... add i to row list for next nonzero entry */
470         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
471         j     = bj[jmin];
472         jl[i] = jl[j]; jl[j] = i; /* update jl */
473       }
474       i = nexti;
475     }
476 
477     /* save nonzero entries in k-th row of U ... */
478 
479     /* invert diagonal block */
480     diag = ba+k*bs2;
481     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
482     Kernel_A_gets_inverse_A(bs,diag,pivots,work);
483 
484     jmin = bi[k]; jmax = bi[k+1];
485     if (jmin < jmax) {
486       for (j=jmin; j<jmax; j++){
487          vj = bj[j];           /* block col. index of U */
488          u   = ba + j*bs2;
489          rtmp_ptr = rtmp + vj*bs2;
490          for (k1=0; k1<bs2; k1++){
491            *u++        = *rtmp_ptr;
492            *rtmp_ptr++ = 0.0;
493          }
494       }
495 
496       /* ... add k to row list for first nonzero entry in k-th row */
497       il[k] = jmin;
498       i     = bj[jmin];
499       jl[k] = jl[i]; jl[i] = k;
500     }
501   }
502 
503   ierr = PetscFree(rtmp);CHKERRQ(ierr);
504   ierr = PetscFree(il);CHKERRQ(ierr);
505   ierr = PetscFree(dk);CHKERRQ(ierr);
506   ierr = PetscFree(pivots);CHKERRQ(ierr);
507 
508   C->factor    = FACTOR_CHOLESKY;
509   C->assembled = PETSC_TRUE;
510   C->preallocated = PETSC_TRUE;
511   PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */
512   PetscFunctionReturn(0);
513 }
514 
515 /*
516     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
517     Version for blocks 2 by 2.
518 */
519 #undef __FUNCT__
520 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2"
521 int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat A,Mat *B)
522 {
523   Mat                C = *B;
524   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
525   IS                 perm = b->row;
526   int                *perm_ptr,ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
527   int                *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
528   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
529   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
530 
531   PetscFunctionBegin;
532 
533   /* initialization */
534   /* il and jl record the first nonzero element in each row of the accessing
535      window U(0:k, k:mbs-1).
536      jl:    list of rows to be added to uneliminated rows
537             i>= k: jl(i) is the first row to be added to row i
538             i<  k: jl(i) is the row following row i in some list of rows
539             jl(i) = mbs indicates the end of a list
540      il(i): points to the first nonzero element in columns k,...,mbs-1 of
541             row i of U */
542   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
543   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
544   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
545   jl   = il + mbs;
546   for (i=0; i<mbs; i++) {
547     jl[i] = mbs; il[0] = 0;
548   }
549   ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr);
550   uik  = dk + 4;
551   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
552 
553   /* check permutation */
554   if (!a->permute){
555     ai = a->i; aj = a->j; aa = a->a;
556   } else {
557     ai   = a->inew; aj = a->jnew;
558     ierr = PetscMalloc(4*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
559     ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
560     ierr = PetscMalloc(ai[mbs]*sizeof(int),&a2anew);CHKERRQ(ierr);
561     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr);
562 
563     for (i=0; i<mbs; i++){
564       jmin = ai[i]; jmax = ai[i+1];
565       for (j=jmin; j<jmax; j++){
566         while (a2anew[j] != j){
567           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
568           for (k1=0; k1<4; k1++){
569             dk[k1]       = aa[k*4+k1];
570             aa[k*4+k1] = aa[j*4+k1];
571             aa[j*4+k1] = dk[k1];
572           }
573         }
574         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
575         if (i > aj[j]){
576           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
577           ap = aa + j*4;     /* ptr to the beginning of the block */
578           dk[1] = ap[1];     /* swap ap[1] and ap[2] */
579           ap[1] = ap[2];
580           ap[2] = dk[1];
581         }
582       }
583     }
584     ierr = PetscFree(a2anew);CHKERRQ(ierr);
585   }
586 
587   /* for each row k */
588   for (k = 0; k<mbs; k++){
589 
590     /*initialize k-th row with elements nonzero in row perm(k) of A */
591     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
592     if (jmin < jmax) {
593       ap = aa + jmin*4;
594       for (j = jmin; j < jmax; j++){
595         vj = perm_ptr[aj[j]];         /* block col. index */
596         rtmp_ptr = rtmp + vj*4;
597         for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
598       }
599     }
600 
601     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
602     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
603     i = jl[k]; /* first row to be added to k_th row  */
604 
605     while (i < mbs){
606       nexti = jl[i]; /* next row to be added to k_th row */
607 
608       /* compute multiplier */
609       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
610 
611       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
612       diag = ba + i*4;
613       u    = ba + ili*4;
614       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
615       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
616       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
617       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
618 
619       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
620       dk[0] += uik[0]*u[0] + uik[1]*u[1];
621       dk[1] += uik[2]*u[0] + uik[3]*u[1];
622       dk[2] += uik[0]*u[2] + uik[1]*u[3];
623       dk[3] += uik[2]*u[2] + uik[3]*u[3];
624 
625       /* update -U(i,k): ba[ili] = uik */
626       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
627 
628       /* add multiple of row i to k-th row ... */
629       jmin = ili + 1; jmax = bi[i+1];
630       if (jmin < jmax){
631         for (j=jmin; j<jmax; j++) {
632           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
633           rtmp_ptr = rtmp + bj[j]*4;
634           u = ba + j*4;
635           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
636           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
637           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
638           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
639         }
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*4;
653     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
654     ierr = Kernel_A_gets_inverse_A_2(diag);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*4;
661          rtmp_ptr = rtmp + vj*4;
662          for (k1=0; k1<4; 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 = PetscFree(il);CHKERRQ(ierr);
677   ierr = PetscFree(dk);CHKERRQ(ierr);
678   if (a->permute) {
679     ierr = PetscFree(aa);CHKERRQ(ierr);
680   }
681   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
682   C->factor    = FACTOR_CHOLESKY;
683   C->assembled = PETSC_TRUE;
684   C->preallocated = PETSC_TRUE;
685   PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */
686   PetscFunctionReturn(0);
687 }
688 
689 /*
690       Version for when blocks are 2 by 2 Using natural ordering
691 */
692 #undef __FUNCT__
693 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering"
694 int MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat A,Mat *B)
695 {
696   Mat                C = *B;
697   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
698   int                ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
699   int                *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
700   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
701   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
702 
703   PetscFunctionBegin;
704 
705   /* initialization */
706   /* il and jl record the first nonzero element in each row of the accessing
707      window U(0:k, k:mbs-1).
708      jl:    list of rows to be added to uneliminated rows
709             i>= k: jl(i) is the first row to be added to row i
710             i<  k: jl(i) is the row following row i in some list of rows
711             jl(i) = mbs indicates the end of a list
712      il(i): points to the first nonzero element in columns k,...,mbs-1 of
713             row i of U */
714   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
715   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
716   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
717   jl   = il + mbs;
718   for (i=0; i<mbs; i++) {
719     jl[i] = mbs; il[0] = 0;
720   }
721   ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr);
722   uik  = dk + 4;
723 
724   ai = a->i; aj = a->j; aa = a->a;
725 
726   /* for each row k */
727   for (k = 0; k<mbs; k++){
728 
729     /*initialize k-th row with elements nonzero in row k of A */
730     jmin = ai[k]; jmax = ai[k+1];
731     if (jmin < jmax) {
732       ap = aa + jmin*4;
733       for (j = jmin; j < jmax; j++){
734         vj = aj[j];         /* block col. index */
735         rtmp_ptr = rtmp + vj*4;
736         for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
737       }
738     }
739 
740     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
741     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
742     i = jl[k]; /* first row to be added to k_th row  */
743 
744     while (i < mbs){
745       nexti = jl[i]; /* next row to be added to k_th row */
746 
747       /* compute multiplier */
748       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
749 
750       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
751       diag = ba + i*4;
752       u    = ba + ili*4;
753       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
754       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
755       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
756       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
757 
758       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
759       dk[0] += uik[0]*u[0] + uik[1]*u[1];
760       dk[1] += uik[2]*u[0] + uik[3]*u[1];
761       dk[2] += uik[0]*u[2] + uik[1]*u[3];
762       dk[3] += uik[2]*u[2] + uik[3]*u[3];
763 
764       /* update -U(i,k): ba[ili] = uik */
765       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
766 
767       /* add multiple of row i to k-th row ... */
768       jmin = ili + 1; jmax = bi[i+1];
769       if (jmin < jmax){
770         for (j=jmin; j<jmax; j++) {
771           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
772           rtmp_ptr = rtmp + bj[j]*4;
773           u = ba + j*4;
774           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
775           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
776           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
777           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
778         }
779 
780         /* ... add i to row list for next nonzero entry */
781         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
782         j     = bj[jmin];
783         jl[i] = jl[j]; jl[j] = i; /* update jl */
784       }
785       i = nexti;
786     }
787 
788     /* save nonzero entries in k-th row of U ... */
789 
790     /* invert diagonal block */
791     diag = ba+k*4;
792     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
793     ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr);
794 
795     jmin = bi[k]; jmax = bi[k+1];
796     if (jmin < jmax) {
797       for (j=jmin; j<jmax; j++){
798          vj = bj[j];           /* block col. index of U */
799          u   = ba + j*4;
800          rtmp_ptr = rtmp + vj*4;
801          for (k1=0; k1<4; k1++){
802            *u++        = *rtmp_ptr;
803            *rtmp_ptr++ = 0.0;
804          }
805       }
806 
807       /* ... add k to row list for first nonzero entry in k-th row */
808       il[k] = jmin;
809       i     = bj[jmin];
810       jl[k] = jl[i]; jl[i] = k;
811     }
812   }
813 
814   ierr = PetscFree(rtmp);CHKERRQ(ierr);
815   ierr = PetscFree(il);CHKERRQ(ierr);
816   ierr = PetscFree(dk);CHKERRQ(ierr);
817 
818   C->factor    = FACTOR_CHOLESKY;
819   C->assembled = PETSC_TRUE;
820   C->preallocated = PETSC_TRUE;
821   PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */
822   PetscFunctionReturn(0);
823 }
824 
825 /*
826     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
827     Version for blocks are 1 by 1.
828 */
829 #undef __FUNCT__
830 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1"
831 int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat A,Mat *B)
832 {
833   Mat                C = *B;
834   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
835   IS                 ip = b->row;
836   int                *rip,ierr,i,j,mbs = a->mbs,*bi = b->i,*bj = b->j;
837   int                *ai,*aj,*r;
838   int                k,jmin,jmax,*jl,*il,vj,nexti,ili;
839   MatScalar          *rtmp;
840   MatScalar          *ba = b->a,*aa,ak;
841   MatScalar          dk,uikdi;
842 
843   PetscFunctionBegin;
844   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
845   if (!a->permute){
846     ai = a->i; aj = a->j; aa = a->a;
847   } else {
848     ai = a->inew; aj = a->jnew;
849     ierr = PetscMalloc(ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
850     ierr = PetscMemcpy(aa,a->a,ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
851     ierr = PetscMalloc(ai[mbs]*sizeof(int),&r);CHKERRQ(ierr);
852     ierr= PetscMemcpy(r,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr);
853 
854     jmin = ai[0]; jmax = ai[mbs];
855     for (j=jmin; j<jmax; j++){
856       while (r[j] != j){
857         k = r[j]; r[j] = r[k]; r[k] = k;
858         ak = aa[k]; aa[k] = aa[j]; aa[j] = ak;
859       }
860     }
861     ierr = PetscFree(r);CHKERRQ(ierr);
862   }
863 
864   /* initialization */
865   /* il and jl record the first nonzero element in each row of the accessing
866      window U(0:k, k:mbs-1).
867      jl:    list of rows to be added to uneliminated rows
868             i>= k: jl(i) is the first row to be added to row i
869             i<  k: jl(i) is the row following row i in some list of rows
870             jl(i) = mbs indicates the end of a list
871      il(i): points to the first nonzero element in columns k,...,mbs-1 of
872             row i of U */
873   ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
874   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
875   jl   = il + mbs;
876   for (i=0; i<mbs; i++) {
877     rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
878   }
879 
880   /* for each row k */
881   for (k = 0; k<mbs; k++){
882 
883     /*initialize k-th row with elements nonzero in row perm(k) of A */
884     jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
885     if (jmin < jmax) {
886       for (j = jmin; j < jmax; j++){
887         vj = rip[aj[j]];
888         rtmp[vj] = aa[j];
889       }
890     }
891 
892     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
893     dk = rtmp[k];
894     i = jl[k]; /* first row to be added to k_th row  */
895 
896     while (i < mbs){
897       nexti = jl[i]; /* next row to be added to k_th row */
898 
899       /* compute multiplier, update D(k) and U(i,k) */
900       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
901       uikdi = - ba[ili]*ba[i];
902       dk += uikdi*ba[ili];
903       ba[ili] = uikdi; /* -U(i,k) */
904 
905       /* add multiple of row i to k-th row ... */
906       jmin = ili + 1; jmax = bi[i+1];
907       if (jmin < jmax){
908         for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
909         /* ... add i to row list for next nonzero entry */
910         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
911         j     = bj[jmin];
912         jl[i] = jl[j]; jl[j] = i; /* update jl */
913       }
914       i = nexti;
915     }
916 
917     /* check for zero pivot and save diagoanl element */
918     if (dk == 0.0){
919       SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot");
920       /*
921     } else if (PetscRealPart(dk) < 0.0){
922       SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Negative pivot: d[%d] = %g\n",k,dk);
923       */
924     }
925 
926     /* save nonzero entries in k-th row of U ... */
927     ba[k] = 1.0/dk;
928     jmin = bi[k]; jmax = bi[k+1];
929     if (jmin < jmax) {
930       for (j=jmin; j<jmax; j++){
931          vj = bj[j]; ba[j] = rtmp[vj]; rtmp[vj] = 0.0;
932       }
933       /* ... add k to row list for first nonzero entry in k-th row */
934       il[k] = jmin;
935       i     = bj[jmin];
936       jl[k] = jl[i]; jl[i] = k;
937     }
938   }
939 
940   ierr = PetscFree(rtmp);CHKERRQ(ierr);
941   ierr = PetscFree(il);CHKERRQ(ierr);
942   if (a->permute){
943     ierr = PetscFree(aa);CHKERRQ(ierr);
944   }
945 
946   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
947   C->factor    = FACTOR_CHOLESKY;
948   C->assembled = PETSC_TRUE;
949   C->preallocated = PETSC_TRUE;
950   PetscLogFlops(b->mbs);
951   PetscFunctionReturn(0);
952 }
953 
954 /*
955   Version for when blocks are 1 by 1 Using natural ordering
956 */
957 #undef __FUNCT__
958 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering"
959 int MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat A,Mat *B)
960 {
961   Mat                C = *B;
962   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
963   int                ierr,i,j,mbs = a->mbs,*bi = b->i,*bj = b->j;
964   int                *ai,*aj;
965   int                k,jmin,jmax,*jl,*il,vj,nexti,ili;
966   MatScalar          *rtmp,*ba = b->a,*aa,dk,uikdi;
967 
968   PetscFunctionBegin;
969 
970   /* initialization */
971   /* il and jl record the first nonzero element in each row of the accessing
972      window U(0:k, k:mbs-1).
973      jl:    list of rows to be added to uneliminated rows
974             i>= k: jl(i) is the first row to be added to row i
975             i<  k: jl(i) is the row following row i in some list of rows
976             jl(i) = mbs indicates the end of a list
977      il(i): points to the first nonzero element in columns k,...,mbs-1 of
978             row i of U */
979   ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
980   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
981   jl   = il + mbs;
982   for (i=0; i<mbs; i++) {
983     rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
984   }
985 
986   ai = a->i; aj = a->j; aa = a->a;
987 
988   /* for each row k */
989   for (k = 0; k<mbs; k++){
990 
991     /*initialize k-th row with elements nonzero in row perm(k) of A */
992     jmin = ai[k]; jmax = ai[k+1];
993     if (jmin < jmax) {
994       for (j = jmin; j < jmax; j++){
995         vj = aj[j];
996         rtmp[vj] = aa[j];
997       }
998     }
999 
1000     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
1001     dk = rtmp[k];
1002     i = jl[k]; /* first row to be added to k_th row  */
1003 
1004     while (i < mbs){
1005       nexti = jl[i]; /* next row to be added to k_th row */
1006 
1007       /* compute multiplier, update D(k) and U(i,k) */
1008       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1009       uikdi = - ba[ili]*ba[i];
1010       dk += uikdi*ba[ili];
1011       ba[ili] = uikdi; /* -U(i,k) */
1012 
1013       /* add multiple of row i to k-th row ... */
1014       jmin = ili + 1; jmax = bi[i+1];
1015       if (jmin < jmax){
1016         for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1017         /* ... add i to row list for next nonzero entry */
1018         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
1019         j     = bj[jmin];
1020         jl[i] = jl[j]; jl[j] = i; /* update jl */
1021       }
1022       i = nexti;
1023     }
1024 
1025     /* check for zero pivot and save diagoanl element */
1026     if (dk == 0.0){
1027       SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot");
1028       /*
1029     } else if (PetscRealPart(dk) < 0){
1030       SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Negative pivot: d[%d] = %g\n",k,dk);
1031       */
1032     }
1033 
1034     /* save nonzero entries in k-th row of U ... */
1035     ba[k] = 1.0/dk;
1036     jmin = bi[k]; jmax = bi[k+1];
1037     if (jmin < jmax) {
1038       for (j=jmin; j<jmax; j++){
1039          vj = bj[j]; ba[j] = rtmp[vj]; rtmp[vj] = 0.0;
1040       }
1041       /* ... add k to row list for first nonzero entry in k-th row */
1042       il[k] = jmin;
1043       i     = bj[jmin];
1044       jl[k] = jl[i]; jl[i] = k;
1045     }
1046   }
1047 
1048   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1049   ierr = PetscFree(il);CHKERRQ(ierr);
1050 
1051   C->factor    = FACTOR_CHOLESKY;
1052   C->assembled = PETSC_TRUE;
1053   C->preallocated = PETSC_TRUE;
1054   PetscLogFlops(b->mbs);
1055   PetscFunctionReturn(0);
1056 }
1057 
1058 #undef __FUNCT__
1059 #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ"
1060 int MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,PetscReal f)
1061 {
1062   int ierr;
1063   Mat C;
1064 
1065   PetscFunctionBegin;
1066   ierr = MatCholeskyFactorSymbolic(A,perm,f,&C);CHKERRQ(ierr);
1067   ierr = MatCholeskyFactorNumeric(A,&C);CHKERRQ(ierr);
1068   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
1069   PetscFunctionReturn(0);
1070 }
1071 
1072 
1073