xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact.c (revision eeeff2ec8a6dcd65162de47cb107fc47d432e393)
173f4d377SMatthew Knepley /*$Id: sbaijfact.c,v 1.61 2001/08/06 21:15:47 bsmith Exp $*/
249b5e25fSSatish Balay 
349b5e25fSSatish Balay #include "src/mat/impls/baij/seq/baij.h"
43a7fca6bSBarry Smith #include "src/mat/impls/sbaij/seq/sbaij.h"
549b5e25fSSatish Balay #include "src/vec/vecimpl.h"
649b5e25fSSatish Balay #include "src/inline/ilu.h"
749a6740bSHong Zhang #include "include/petscis.h"
849b5e25fSSatish Balay 
98dc9108eSHong Zhang #if !defined(PETSC_USE_COMPLEX)
105f9f512dSHong Zhang /*
115f9f512dSHong Zhang   input:
12c037c3f7SHong Zhang    F -- numeric factor
135f9f512dSHong Zhang   output:
14c037c3f7SHong Zhang    nneg, nzero, npos: matrix inertia
155f9f512dSHong Zhang */
165f9f512dSHong Zhang 
175f9f512dSHong Zhang #undef __FUNCT__
185f9f512dSHong Zhang #define __FUNCT__ "MatGetInertia_SeqSBAIJ"
19c037c3f7SHong Zhang int MatGetInertia_SeqSBAIJ(Mat F,int *nneig,int *nzero,int *npos)
205f9f512dSHong Zhang {
21638f5ce0SDinesh Kaushik   Mat_SeqSBAIJ *fact_ptr = (Mat_SeqSBAIJ*)F->data;
223e0d88b5SBarry Smith   PetscScalar  *dd = fact_ptr->a;
23*eeeff2ecSHong Zhang   int          mbs=fact_ptr->mbs,bs=fact_ptr->bs,i,nneig_tmp,npos_tmp,
24*eeeff2ecSHong Zhang                *fi = fact_ptr->i;
255f9f512dSHong Zhang 
265f9f512dSHong Zhang   PetscFunctionBegin;
27*eeeff2ecSHong Zhang   if (bs != 1) SETERRQ1(PETSC_ERR_SUP,"No support for bs: %d >1 yet",bs);
28*eeeff2ecSHong Zhang   nneig_tmp = 0; npos_tmp = 0;
29*eeeff2ecSHong Zhang   for (i=0; i<mbs; i++){
30*eeeff2ecSHong Zhang     if (PetscRealPart(dd[*fi]) > 0.0){
31*eeeff2ecSHong Zhang       npos_tmp++;
32*eeeff2ecSHong Zhang     } else if (PetscRealPart(dd[*fi]) < 0.0){
33*eeeff2ecSHong Zhang       nneig_tmp++;
345f9f512dSHong Zhang     }
35*eeeff2ecSHong Zhang     fi++;
363e0d88b5SBarry Smith   }
37*eeeff2ecSHong Zhang   if (nneig) *nneig = nneig_tmp;
38*eeeff2ecSHong Zhang   if (npos)  *npos  = npos_tmp;
39*eeeff2ecSHong Zhang   if (nzero) *nzero = mbs - nneig_tmp - npos_tmp;
40*eeeff2ecSHong Zhang 
415f9f512dSHong Zhang   PetscFunctionReturn(0);
425f9f512dSHong Zhang }
438dc9108eSHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */
445f9f512dSHong Zhang 
455f9f512dSHong Zhang /* Using Modified Sparse Row (MSR) storage.
465f9f512dSHong Zhang See page 85, "Iterative Methods ..." by Saad. */
475f9f512dSHong Zhang /*
485f9f512dSHong Zhang     Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP.
495f9f512dSHong Zhang */
502d9a3abdSHong Zhang /* Use Modified Sparse Row storage for u and ju, see Saad pp.85 */
514a2ae208SSatish Balay #undef __FUNCT__
524a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ"
5315e8a5b3SHong Zhang int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *B)
5449b5e25fSSatish Balay {
5549b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b;
56cb718733SHong Zhang   int          *rip,ierr,i,mbs = a->mbs,*ai,*aj;
57066653e3SSatish Balay   int          *jutmp,bs = a->bs,bs2=a->bs2;
5887fe1012SHong Zhang   int          m,realloc = 0,prow;
59d60e5362SHong Zhang   int          *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd;
60064503c5SHong Zhang   int          *il,ili,nextprow;
6115e8a5b3SHong Zhang   PetscReal    f = info->fill;
62671cb588SHong Zhang   PetscTruth   perm_identity;
6349b5e25fSSatish Balay 
6449b5e25fSSatish Balay   PetscFunctionBegin;
65cb718733SHong Zhang   /* check whether perm is the identity mapping */
66671cb588SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
67064503c5SHong Zhang 
68064503c5SHong Zhang   /* -- inplace factorization, i.e., use sbaij for *B -- */
69064503c5SHong Zhang   if (perm_identity && bs==1 ){
70064503c5SHong Zhang     if (!perm_identity) a->permute = PETSC_TRUE;
71064503c5SHong Zhang 
72064503c5SHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
73064503c5SHong Zhang 
74064503c5SHong Zhang   if (perm_identity){ /* without permutation */
75064503c5SHong Zhang     ai = a->i; aj = a->j;
76064503c5SHong Zhang   } else {            /* non-trivial permutation */
77064503c5SHong Zhang     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
78064503c5SHong Zhang     ai = a->inew; aj = a->jnew;
79064503c5SHong Zhang   }
80064503c5SHong Zhang 
81064503c5SHong Zhang   /* initialization */
82064503c5SHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(int),&iu);CHKERRQ(ierr);
83064503c5SHong Zhang   umax  = (int)(f*ai[mbs] + 1);
84064503c5SHong Zhang   ierr  = PetscMalloc(umax*sizeof(int),&ju);CHKERRQ(ierr);
85064503c5SHong Zhang   iu[0] = 0;
86064503c5SHong Zhang   juidx = 0; /* index for ju */
87064503c5SHong Zhang   ierr  = PetscMalloc((3*mbs+1)*sizeof(int),&jl);CHKERRQ(ierr); /* linked list for getting pivot row */
88064503c5SHong Zhang   q     = jl + mbs;   /* linked list for col index of active row */
89064503c5SHong Zhang   il    = q  + mbs;
90064503c5SHong Zhang   for (i=0; i<mbs; i++){
91064503c5SHong Zhang     jl[i] = mbs;
92064503c5SHong Zhang     q[i]  = 0;
93064503c5SHong Zhang     il[i] = 0;
94064503c5SHong Zhang   }
95064503c5SHong Zhang 
96064503c5SHong Zhang   /* for each row k */
97064503c5SHong Zhang   for (k=0; k<mbs; k++){
98064503c5SHong Zhang     nzk  = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */
99064503c5SHong Zhang     q[k] = mbs;
100064503c5SHong Zhang     /* initialize nonzero structure of k-th row to row rip[k] of A */
101064503c5SHong Zhang     jmin = ai[rip[k]] +1; /* exclude diag[k] */
102064503c5SHong Zhang     jmax = ai[rip[k]+1];
103064503c5SHong Zhang     for (j=jmin; j<jmax; j++){
104064503c5SHong Zhang       vj = rip[aj[j]]; /* col. value */
105064503c5SHong Zhang       if(vj > k){
106064503c5SHong Zhang         qm = k;
107064503c5SHong Zhang         do {
108064503c5SHong Zhang           m  = qm; qm = q[m];
109064503c5SHong Zhang         } while(qm < vj);
110064503c5SHong Zhang         if (qm == vj) {
111064503c5SHong Zhang           SETERRQ(1," error: duplicate entry in A\n");
112064503c5SHong Zhang         }
113064503c5SHong Zhang         nzk++;
114064503c5SHong Zhang         q[m]  = vj;
115064503c5SHong Zhang         q[vj] = qm;
116064503c5SHong Zhang       } /* if(vj > k) */
117064503c5SHong Zhang     } /* for (j=jmin; j<jmax; j++) */
118064503c5SHong Zhang 
119064503c5SHong Zhang     /* modify nonzero structure of k-th row by computing fill-in
120064503c5SHong Zhang        for each row i to be merged in */
121064503c5SHong Zhang     prow = k;
122064503c5SHong Zhang     prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */
123064503c5SHong Zhang 
124064503c5SHong Zhang     while (prow < k){
125064503c5SHong Zhang       nextprow = jl[prow];
126064503c5SHong Zhang 
127064503c5SHong Zhang       /* merge row prow into k-th row */
128064503c5SHong Zhang       ili = il[prow];
129064503c5SHong Zhang       jmin = ili + 1;  /* points to 2nd nzero entry in U(prow,k:mbs-1) */
130064503c5SHong Zhang       jmax = iu[prow+1];
131064503c5SHong Zhang       qm = k;
132064503c5SHong Zhang       for (j=jmin; j<jmax; j++){
133064503c5SHong Zhang         vj = ju[j];
134064503c5SHong Zhang         do {
135064503c5SHong Zhang           m = qm; qm = q[m];
136064503c5SHong Zhang         } while (qm < vj);
137064503c5SHong Zhang         if (qm != vj){  /* a fill */
138064503c5SHong Zhang           nzk++; q[m] = vj; q[vj] = qm; qm = vj;
139064503c5SHong Zhang         }
140064503c5SHong Zhang       } /* end of for (j=jmin; j<jmax; j++) */
141064503c5SHong Zhang       if (jmin < jmax){
142064503c5SHong Zhang         il[prow] = jmin;
143064503c5SHong Zhang         j = ju[jmin];
144064503c5SHong Zhang         jl[prow] = jl[j]; jl[j] = prow;  /* update jl */
145064503c5SHong Zhang       }
146064503c5SHong Zhang       prow = nextprow;
147064503c5SHong Zhang     }
148064503c5SHong Zhang 
149064503c5SHong Zhang     /* update il and jl */
150064503c5SHong Zhang     if (nzk > 0){
151064503c5SHong Zhang       i = q[k]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
152064503c5SHong Zhang       jl[k] = jl[i]; jl[i] = k;
153064503c5SHong Zhang       il[k] = iu[k] + 1;
154064503c5SHong Zhang     }
155064503c5SHong Zhang     iu[k+1] = iu[k] + nzk + 1;  /* include diag[k] */
156064503c5SHong Zhang 
157064503c5SHong Zhang     /* allocate more space to ju if needed */
158064503c5SHong Zhang     if (iu[k+1] > umax) {
159064503c5SHong Zhang       /* estimate how much additional space we will need */
160064503c5SHong Zhang       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
161064503c5SHong Zhang       /* just double the memory each time */
162064503c5SHong Zhang       maxadd = umax;
163064503c5SHong Zhang       if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2;
164064503c5SHong Zhang       umax += maxadd;
165064503c5SHong Zhang 
166064503c5SHong Zhang       /* allocate a longer ju */
167064503c5SHong Zhang       ierr = PetscMalloc(umax*sizeof(int),&jutmp);CHKERRQ(ierr);
168064503c5SHong Zhang       ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(int));CHKERRQ(ierr);
169064503c5SHong Zhang       ierr = PetscFree(ju);CHKERRQ(ierr);
170064503c5SHong Zhang       ju   = jutmp;
171064503c5SHong Zhang       realloc++; /* count how many times we realloc */
172064503c5SHong Zhang     }
173064503c5SHong Zhang 
174064503c5SHong Zhang     /* save nonzero structure of k-th row in ju */
175064503c5SHong Zhang     ju[juidx++] = k; /* diag[k] */
176064503c5SHong Zhang     i = k;
177064503c5SHong Zhang     while (nzk --) {
178064503c5SHong Zhang       i           = q[i];
179064503c5SHong Zhang       ju[juidx++] = i;
180064503c5SHong Zhang     }
181064503c5SHong Zhang   }
182064503c5SHong Zhang 
183064503c5SHong Zhang   if (ai[mbs] != 0) {
184064503c5SHong Zhang     PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
185064503c5SHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af);
186064503c5SHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Run with -pc_cholesky_fill %g or use \n",af);
187064503c5SHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:PCCholeskySetFill(pc,%g);\n",af);
188064503c5SHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:for best performance.\n");
189064503c5SHong Zhang   } else {
190064503c5SHong Zhang      PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Empty matrix.\n");
191064503c5SHong Zhang   }
192064503c5SHong Zhang 
193064503c5SHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
194064503c5SHong Zhang   /* ierr = PetscFree(q);CHKERRQ(ierr); */
195064503c5SHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
196064503c5SHong Zhang 
197064503c5SHong Zhang   /* put together the new matrix */
198064503c5SHong Zhang   ierr = MatCreateSeqSBAIJ(A->comm,bs,bs*mbs,bs*mbs,0,PETSC_NULL,B);CHKERRQ(ierr);
199064503c5SHong Zhang   /* PetscLogObjectParent(*B,iperm); */
200064503c5SHong Zhang   b = (Mat_SeqSBAIJ*)(*B)->data;
201064503c5SHong Zhang   ierr = PetscFree(b->imax);CHKERRQ(ierr);
202064503c5SHong Zhang   b->singlemalloc = PETSC_FALSE;
203064503c5SHong Zhang   /* the next line frees the default space generated by the Create() */
204064503c5SHong Zhang   ierr = PetscFree(b->a);CHKERRQ(ierr);
205064503c5SHong Zhang   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
206064503c5SHong Zhang   ierr = PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
207064503c5SHong Zhang   b->j    = ju;
208064503c5SHong Zhang   b->i    = iu;
209064503c5SHong Zhang   b->diag = 0;
210064503c5SHong Zhang   b->ilen = 0;
211064503c5SHong Zhang   b->imax = 0;
212064503c5SHong Zhang   b->row  = perm;
21315e8a5b3SHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
214064503c5SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
215064503c5SHong Zhang   b->icol = perm;
216064503c5SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
217064503c5SHong Zhang   ierr    = PetscMalloc((bs*mbs+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
218064503c5SHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.
219064503c5SHong Zhang      Allocate idnew, solve_work, new a, new j */
220064503c5SHong Zhang   PetscLogObjectMemory(*B,(iu[mbs]-mbs)*(sizeof(int)+sizeof(MatScalar)));
221064503c5SHong Zhang   b->s_maxnz = b->s_nz = iu[mbs];
222064503c5SHong Zhang 
223064503c5SHong Zhang   (*B)->factor                 = FACTOR_CHOLESKY;
224064503c5SHong Zhang   (*B)->info.factor_mallocs    = realloc;
225064503c5SHong Zhang   (*B)->info.fill_ratio_given  = f;
226064503c5SHong Zhang   if (ai[mbs] != 0) {
227064503c5SHong Zhang     (*B)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
228064503c5SHong Zhang   } else {
229064503c5SHong Zhang     (*B)->info.fill_ratio_needed = 0.0;
230064503c5SHong Zhang   }
231064503c5SHong Zhang 
232064503c5SHong Zhang 
233b45a75daSHong Zhang   (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering;
234b45a75daSHong Zhang   (*B)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
235064503c5SHong Zhang   PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=1\n");
236064503c5SHong Zhang 
237064503c5SHong Zhang   PetscFunctionReturn(0);
238064503c5SHong Zhang   }
239064503c5SHong Zhang   /* -----------  end of new code --------------------*/
240064503c5SHong Zhang 
241064503c5SHong Zhang 
242671cb588SHong Zhang   if (!perm_identity) a->permute = PETSC_TRUE;
2432d007305SHong Zhang 
244671cb588SHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
245671cb588SHong Zhang 
246671cb588SHong Zhang   if (perm_identity){ /* without permutation */
2472d007305SHong Zhang     ai = a->i; aj = a->j;
2482d007305SHong Zhang   } else {            /* non-trivial permutation */
249323b833fSBarry Smith     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
2502d007305SHong Zhang     ai = a->inew; aj = a->jnew;
2512d007305SHong Zhang   }
2522d007305SHong Zhang 
25349b5e25fSSatish Balay   /* initialization */
254b0a32e0cSBarry Smith   ierr  = PetscMalloc((mbs+1)*sizeof(int),&iu);CHKERRQ(ierr);
25549b5e25fSSatish Balay   umax  = (int)(f*ai[mbs] + 1); umax += mbs + 1;
25682502324SSatish Balay   ierr  = PetscMalloc(umax*sizeof(int),&ju);CHKERRQ(ierr);
25749b5e25fSSatish Balay   iu[0] = mbs+1;
25887fe1012SHong Zhang   juidx = mbs + 1; /* index for ju */
25987fe1012SHong Zhang   ierr  = PetscMalloc(2*mbs*sizeof(int),&jl);CHKERRQ(ierr); /* linked list for pivot row */
26087fe1012SHong Zhang   q     = jl + mbs;   /* linked list for col index */
26149b5e25fSSatish Balay   for (i=0; i<mbs; i++){
26287fe1012SHong Zhang     jl[i] = mbs;
26387fe1012SHong Zhang     q[i] = 0;
26449b5e25fSSatish Balay   }
26549b5e25fSSatish Balay 
26649b5e25fSSatish Balay   /* for each row k */
26749b5e25fSSatish Balay   for (k=0; k<mbs; k++){
268064503c5SHong Zhang     for (i=0; i<mbs; i++) q[i] = 0;  /* to be removed! */
26949b5e25fSSatish Balay     nzk  = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */
27049b5e25fSSatish Balay     q[k] = mbs;
27149b5e25fSSatish Balay     /* initialize nonzero structure of k-th row to row rip[k] of A */
272064503c5SHong Zhang     jmin = ai[rip[k]] +1; /* exclude diag[k] */
27349b5e25fSSatish Balay     jmax = ai[rip[k]+1];
27449b5e25fSSatish Balay     for (j=jmin; j<jmax; j++){
275cb718733SHong Zhang       vj = rip[aj[j]]; /* col. value */
27649b5e25fSSatish Balay       if(vj > k){
27749b5e25fSSatish Balay         qm = k;
27849b5e25fSSatish Balay         do {
27949b5e25fSSatish Balay           m  = qm; qm = q[m];
28049b5e25fSSatish Balay         } while(qm < vj);
28149b5e25fSSatish Balay         if (qm == vj) {
2823078e140SBarry Smith           SETERRQ(1," error: duplicate entry in A\n");
28349b5e25fSSatish Balay         }
28449b5e25fSSatish Balay         nzk++;
28549b5e25fSSatish Balay         q[m]  = vj;
28649b5e25fSSatish Balay         q[vj] = qm;
28749b5e25fSSatish Balay       } /* if(vj > k) */
28849b5e25fSSatish Balay     } /* for (j=jmin; j<jmax; j++) */
28949b5e25fSSatish Balay 
29049b5e25fSSatish Balay     /* modify nonzero structure of k-th row by computing fill-in
29149b5e25fSSatish Balay        for each row i to be merged in */
29287fe1012SHong Zhang     prow = k;
29387fe1012SHong Zhang     prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */
29487fe1012SHong Zhang 
29587fe1012SHong Zhang     while (prow < k){
29687fe1012SHong Zhang       /* merge row prow into k-th row */
29787fe1012SHong Zhang       jmin = iu[prow] + 1; jmax = iu[prow+1];
29849b5e25fSSatish Balay       qm = k;
29987fe1012SHong Zhang       for (j=jmin; j<jmax; j++){
30049b5e25fSSatish Balay         vj = ju[j];
30149b5e25fSSatish Balay         do {
30249b5e25fSSatish Balay           m = qm; qm = q[m];
30349b5e25fSSatish Balay         } while (qm < vj);
30449b5e25fSSatish Balay         if (qm != vj){
30549b5e25fSSatish Balay          nzk++; q[m] = vj; q[vj] = qm; qm = vj;
30649b5e25fSSatish Balay         }
30749b5e25fSSatish Balay       }
30887fe1012SHong Zhang       prow = jl[prow]; /* next pivot row */
30949b5e25fSSatish Balay     }
31049b5e25fSSatish Balay 
31149b5e25fSSatish Balay     /* add k to row list for first nonzero element in k-th row */
31249b5e25fSSatish Balay     if (nzk > 0){
31349b5e25fSSatish Balay       i = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */
31449b5e25fSSatish Balay       jl[k] = jl[i]; jl[i] = k;
31549b5e25fSSatish Balay     }
31687fe1012SHong Zhang     iu[k+1] = iu[k] + nzk;
31749b5e25fSSatish Balay 
31849b5e25fSSatish Balay     /* allocate more space to ju if needed */
3193078e140SBarry Smith     if (iu[k+1] > umax) {
32049b5e25fSSatish Balay       /* estimate how much additional space we will need */
32149b5e25fSSatish Balay       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
32249b5e25fSSatish Balay       /* just double the memory each time */
32349b5e25fSSatish Balay       maxadd = umax;
32449b5e25fSSatish Balay       if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2;
32549b5e25fSSatish Balay       umax += maxadd;
32649b5e25fSSatish Balay 
3279f9cb213SHong Zhang       /* allocate a longer ju */
32882502324SSatish Balay       ierr = PetscMalloc(umax*sizeof(int),&jutmp);CHKERRQ(ierr);
32949b5e25fSSatish Balay       ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(int));CHKERRQ(ierr);
33049b5e25fSSatish Balay       ierr = PetscFree(ju);CHKERRQ(ierr);
3319f9cb213SHong Zhang       ju   = jutmp;
33249b5e25fSSatish Balay       realloc++; /* count how many times we realloc */
33349b5e25fSSatish Balay     }
33449b5e25fSSatish Balay 
33549b5e25fSSatish Balay     /* save nonzero structure of k-th row in ju */
33649b5e25fSSatish Balay     i=k;
33787fe1012SHong Zhang     while (nzk --) {
33849b5e25fSSatish Balay       i           = q[i];
33987fe1012SHong Zhang       ju[juidx++] = i;
34049b5e25fSSatish Balay     }
34177f73120SHong Zhang   }
34249b5e25fSSatish Balay 
34349b5e25fSSatish Balay   if (ai[mbs] != 0) {
34449b5e25fSSatish Balay     PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
345b0a32e0cSBarry Smith     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af);
3463078e140SBarry Smith     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Run with -pc_cholesky_fill %g or use \n",af);
3473078e140SBarry Smith     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:PCCholeskySetFill(pc,%g);\n",af);
348b0a32e0cSBarry Smith     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:for best performance.\n");
34949b5e25fSSatish Balay   } else {
350b0a32e0cSBarry Smith      PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Empty matrix.\n");
35149b5e25fSSatish Balay   }
35249b5e25fSSatish Balay 
35377f73120SHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
35487fe1012SHong Zhang   /* ierr = PetscFree(q);CHKERRQ(ierr); */
35549b5e25fSSatish Balay   ierr = PetscFree(jl);CHKERRQ(ierr);
35649b5e25fSSatish Balay 
35749b5e25fSSatish Balay   /* put together the new matrix */
35849b5e25fSSatish Balay   ierr = MatCreateSeqSBAIJ(A->comm,bs,bs*mbs,bs*mbs,0,PETSC_NULL,B);CHKERRQ(ierr);
359b0a32e0cSBarry Smith   /* PetscLogObjectParent(*B,iperm); */
36049b5e25fSSatish Balay   b = (Mat_SeqSBAIJ*)(*B)->data;
36149b5e25fSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
36249b5e25fSSatish Balay   b->singlemalloc = PETSC_FALSE;
36349b5e25fSSatish Balay   /* the next line frees the default space generated by the Create() */
36449b5e25fSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
36549b5e25fSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
36682502324SSatish Balay   ierr = PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
36749b5e25fSSatish Balay   b->j    = ju;
36849b5e25fSSatish Balay   b->i    = iu;
36949b5e25fSSatish Balay   b->diag = 0;
37049b5e25fSSatish Balay   b->ilen = 0;
37149b5e25fSSatish Balay   b->imax = 0;
37277f73120SHong Zhang   b->row  = perm;
37315e8a5b3SHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
37477f73120SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
375cb718733SHong Zhang   b->icol = perm;
376cb718733SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
37787828ca2SBarry Smith   ierr    = PetscMalloc((bs*mbs+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
37849b5e25fSSatish Balay   /* In b structure:  Free imax, ilen, old a, old j.
37949b5e25fSSatish Balay      Allocate idnew, solve_work, new a, new j */
380b0a32e0cSBarry Smith   PetscLogObjectMemory(*B,(iu[mbs]-mbs)*(sizeof(int)+sizeof(MatScalar)));
38149b5e25fSSatish Balay   b->s_maxnz = b->s_nz = iu[mbs];
38249b5e25fSSatish Balay 
3835c0bcdfcSHong Zhang   (*B)->factor                 = FACTOR_CHOLESKY;
38449b5e25fSSatish Balay   (*B)->info.factor_mallocs    = realloc;
38549b5e25fSSatish Balay   (*B)->info.fill_ratio_given  = f;
38649b5e25fSSatish Balay   if (ai[mbs] != 0) {
38749b5e25fSSatish Balay     (*B)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
38849b5e25fSSatish Balay   } else {
38949b5e25fSSatish Balay     (*B)->info.fill_ratio_needed = 0.0;
39049b5e25fSSatish Balay   }
3910aa0ce06SHong Zhang 
392671cb588SHong Zhang   if (perm_identity){
393671cb588SHong Zhang     switch (bs) {
394671cb588SHong Zhang       case 1:
395671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering;
396671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
3974d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=1\n");
398671cb588SHong Zhang         break;
399671cb588SHong Zhang       case 2:
400671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
401671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
4024d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
403671cb588SHong Zhang         break;
404671cb588SHong Zhang       case 3:
405671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
406671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
4074d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:sing special in-place natural ordering factor and solve BS=3\n");
408671cb588SHong Zhang         break;
409671cb588SHong Zhang       case 4:
410671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
411671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
4124d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
413671cb588SHong Zhang         break;
414671cb588SHong Zhang       case 5:
415671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
416671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
4174d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
418671cb588SHong Zhang         break;
419671cb588SHong Zhang       case 6:
420671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
421671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
4224d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
423671cb588SHong Zhang         break;
424671cb588SHong Zhang       case 7:
425671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
426671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
4274d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
428671cb588SHong Zhang       break;
429671cb588SHong Zhang       default:
430671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering;
431671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_N_NaturalOrdering;
4324d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS>7\n");
433671cb588SHong Zhang       break;
434671cb588SHong Zhang     }
435671cb588SHong Zhang   }
436671cb588SHong Zhang 
43749b5e25fSSatish Balay   PetscFunctionReturn(0);
43849b5e25fSSatish Balay }
43949b5e25fSSatish Balay 
44035d8aa7fSBarry Smith 
4414a2ae208SSatish Balay #undef __FUNCT__
4424a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N"
4436f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat A,Mat *B)
44449b5e25fSSatish Balay {
44549b5e25fSSatish Balay   Mat                C = *B;
4464c16a6a6SHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
4474c16a6a6SHong Zhang   IS                 perm = b->row;
4484c16a6a6SHong Zhang   int                *perm_ptr,ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
4494c16a6a6SHong Zhang   int                *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
4504c16a6a6SHong Zhang   int                bs=a->bs,bs2 = a->bs2;
4514c16a6a6SHong Zhang   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
4524c16a6a6SHong Zhang   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
45328de702eSHong Zhang   MatScalar          *work;
4544c16a6a6SHong Zhang   int                *pivots;
4554c16a6a6SHong Zhang 
4564c16a6a6SHong Zhang   PetscFunctionBegin;
4579706f043SHong Zhang 
4584c16a6a6SHong Zhang   /* initialization */
45982502324SSatish Balay   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
4604c16a6a6SHong Zhang   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
46182502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
46228de702eSHong Zhang   jl   = il + mbs;
4634c16a6a6SHong Zhang   for (i=0; i<mbs; i++) {
4644c16a6a6SHong Zhang     jl[i] = mbs; il[0] = 0;
4654c16a6a6SHong Zhang   }
466b0a32e0cSBarry Smith   ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr);
46728de702eSHong Zhang   uik  = dk + bs2;
46828de702eSHong Zhang   work = uik + bs2;
46982502324SSatish Balay   ierr = PetscMalloc(bs*sizeof(int),&pivots);CHKERRQ(ierr);
4704c16a6a6SHong Zhang 
4714c16a6a6SHong Zhang   ierr  = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
4724c16a6a6SHong Zhang 
4734c16a6a6SHong Zhang   /* check permutation */
4744c16a6a6SHong Zhang   if (!a->permute){
4754c16a6a6SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
4764c16a6a6SHong Zhang   } else {
4774c16a6a6SHong Zhang     ai   = a->inew; aj = a->jnew;
47882502324SSatish Balay     ierr = PetscMalloc(bs2*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
4794c16a6a6SHong Zhang     ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
48082502324SSatish Balay     ierr = PetscMalloc(ai[mbs]*sizeof(int),&a2anew);CHKERRQ(ierr);
4814c16a6a6SHong Zhang     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr);
4824c16a6a6SHong Zhang 
4834c16a6a6SHong Zhang     for (i=0; i<mbs; i++){
4844c16a6a6SHong Zhang       jmin = ai[i]; jmax = ai[i+1];
4854c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++){
4864c16a6a6SHong Zhang         while (a2anew[j] != j){
4874c16a6a6SHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
4884c16a6a6SHong Zhang           for (k1=0; k1<bs2; k1++){
4894c16a6a6SHong Zhang             dk[k1]       = aa[k*bs2+k1];
4904c16a6a6SHong Zhang             aa[k*bs2+k1] = aa[j*bs2+k1];
4914c16a6a6SHong Zhang             aa[j*bs2+k1] = dk[k1];
4924c16a6a6SHong Zhang           }
4934c16a6a6SHong Zhang         }
4944c16a6a6SHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
4954c16a6a6SHong Zhang         if (i > aj[j]){
4964c16a6a6SHong Zhang           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
4974c16a6a6SHong Zhang           ap = aa + j*bs2;                     /* ptr to the beginning of j-th block of aa */
4984c16a6a6SHong Zhang           for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
4994c16a6a6SHong Zhang           for (k=0; k<bs; k++){               /* j-th block of aa <- dk^T */
5004c16a6a6SHong Zhang             for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1];
5014c16a6a6SHong Zhang           }
5024c16a6a6SHong Zhang         }
5034c16a6a6SHong Zhang       }
5044c16a6a6SHong Zhang     }
505323b833fSBarry Smith     ierr = PetscFree(a2anew);CHKERRQ(ierr);
5064c16a6a6SHong Zhang   }
5074c16a6a6SHong Zhang 
5084c16a6a6SHong Zhang   /* for each row k */
5094c16a6a6SHong Zhang   for (k = 0; k<mbs; k++){
5104c16a6a6SHong Zhang 
5114c16a6a6SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
5124c16a6a6SHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
513057f5ba7SHong Zhang 
5144c16a6a6SHong Zhang     ap = aa + jmin*bs2;
5154c16a6a6SHong Zhang     for (j = jmin; j < jmax; j++){
5164c16a6a6SHong Zhang       vj = perm_ptr[aj[j]];         /* block col. index */
5174c16a6a6SHong Zhang       rtmp_ptr = rtmp + vj*bs2;
5184c16a6a6SHong Zhang       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
5194c16a6a6SHong Zhang     }
5204c16a6a6SHong Zhang 
5214c16a6a6SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
5224c16a6a6SHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
5234c16a6a6SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
5244c16a6a6SHong Zhang 
525057f5ba7SHong Zhang     while (i < k){
5264c16a6a6SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
5274c16a6a6SHong Zhang 
5284c16a6a6SHong Zhang       /* compute multiplier */
5294c16a6a6SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
5304c16a6a6SHong Zhang 
5314c16a6a6SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
5324c16a6a6SHong Zhang       diag = ba + i*bs2;
5334c16a6a6SHong Zhang       u    = ba + ili*bs2;
5344c16a6a6SHong Zhang       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
5354c16a6a6SHong Zhang       Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
5364c16a6a6SHong Zhang 
5374c16a6a6SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
5384c16a6a6SHong Zhang       Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
5394c16a6a6SHong Zhang 
5404c16a6a6SHong Zhang       /* update -U(i,k) */
5414c16a6a6SHong Zhang       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
5424c16a6a6SHong Zhang 
5434c16a6a6SHong Zhang       /* add multiple of row i to k-th row ... */
5444c16a6a6SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
5454c16a6a6SHong Zhang       if (jmin < jmax){
5464c16a6a6SHong Zhang         for (j=jmin; j<jmax; j++) {
5474c16a6a6SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
5484c16a6a6SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
5494c16a6a6SHong Zhang           u = ba + j*bs2;
5504c16a6a6SHong Zhang           Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
5514c16a6a6SHong Zhang         }
5524c16a6a6SHong Zhang 
5534c16a6a6SHong Zhang         /* ... add i to row list for next nonzero entry */
5544c16a6a6SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
5554c16a6a6SHong Zhang         j     = bj[jmin];
5564c16a6a6SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
5574c16a6a6SHong Zhang       }
5584c16a6a6SHong Zhang       i = nexti;
5594c16a6a6SHong Zhang     }
5604c16a6a6SHong Zhang 
5614c16a6a6SHong Zhang     /* save nonzero entries in k-th row of U ... */
5624c16a6a6SHong Zhang 
5634c16a6a6SHong Zhang     /* invert diagonal block */
5644c16a6a6SHong Zhang     diag = ba+k*bs2;
5654c16a6a6SHong Zhang     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
5664c16a6a6SHong Zhang     Kernel_A_gets_inverse_A(bs,diag,pivots,work);
5674c16a6a6SHong Zhang 
5684c16a6a6SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
5694c16a6a6SHong Zhang     if (jmin < jmax) {
5704c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++){
5714c16a6a6SHong Zhang          vj = bj[j];           /* block col. index of U */
5724c16a6a6SHong Zhang          u   = ba + j*bs2;
5734c16a6a6SHong Zhang          rtmp_ptr = rtmp + vj*bs2;
5744c16a6a6SHong Zhang          for (k1=0; k1<bs2; k1++){
5754c16a6a6SHong Zhang            *u++        = *rtmp_ptr;
5764c16a6a6SHong Zhang            *rtmp_ptr++ = 0.0;
5774c16a6a6SHong Zhang          }
5784c16a6a6SHong Zhang       }
5794c16a6a6SHong Zhang 
5804c16a6a6SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
5814c16a6a6SHong Zhang       il[k] = jmin;
5824c16a6a6SHong Zhang       i     = bj[jmin];
5834c16a6a6SHong Zhang       jl[k] = jl[i]; jl[i] = k;
5844c16a6a6SHong Zhang     }
5854c16a6a6SHong Zhang   }
5864c16a6a6SHong Zhang 
5874c16a6a6SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
5884c16a6a6SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
5894c16a6a6SHong Zhang   ierr = PetscFree(dk);CHKERRQ(ierr);
5904c16a6a6SHong Zhang   ierr = PetscFree(pivots);CHKERRQ(ierr);
5914c16a6a6SHong Zhang   if (a->permute){
5924c16a6a6SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
5934c16a6a6SHong Zhang   }
5944c16a6a6SHong Zhang 
5954c16a6a6SHong Zhang   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
5964c16a6a6SHong Zhang   C->factor       = FACTOR_CHOLESKY;
5974c16a6a6SHong Zhang   C->assembled    = PETSC_TRUE;
5984c16a6a6SHong Zhang   C->preallocated = PETSC_TRUE;
599b0a32e0cSBarry Smith   PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */
6004c16a6a6SHong Zhang   PetscFunctionReturn(0);
6014c16a6a6SHong Zhang }
602d4edadd4SHong Zhang 
6034a2ae208SSatish Balay #undef __FUNCT__
6044a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering"
605671cb588SHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat A,Mat *B)
606671cb588SHong Zhang {
607671cb588SHong Zhang   Mat                C = *B;
608671cb588SHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
609671cb588SHong Zhang   int                ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
6101ad1de41SHong Zhang   int                *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
611671cb588SHong Zhang   int                bs=a->bs,bs2 = a->bs2;
612671cb588SHong Zhang   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
613671cb588SHong Zhang   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
61428de702eSHong Zhang   MatScalar          *work;
615671cb588SHong Zhang   int                *pivots;
616671cb588SHong Zhang 
617671cb588SHong Zhang   PetscFunctionBegin;
618671cb588SHong Zhang 
619671cb588SHong Zhang   /* initialization */
620671cb588SHong Zhang 
62182502324SSatish Balay   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
622671cb588SHong Zhang   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
62382502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
62428de702eSHong Zhang   jl   = il + mbs;
625671cb588SHong Zhang   for (i=0; i<mbs; i++) {
626671cb588SHong Zhang     jl[i] = mbs; il[0] = 0;
627671cb588SHong Zhang   }
628b0a32e0cSBarry Smith   ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr);
62928de702eSHong Zhang   uik  = dk + bs2;
63028de702eSHong Zhang   work = uik + bs2;
63182502324SSatish Balay   ierr = PetscMalloc(bs*sizeof(int),&pivots);CHKERRQ(ierr);
632671cb588SHong Zhang 
633671cb588SHong Zhang   ai = a->i; aj = a->j; aa = a->a;
634671cb588SHong Zhang 
635671cb588SHong Zhang   /* for each row k */
636671cb588SHong Zhang   for (k = 0; k<mbs; k++){
637671cb588SHong Zhang 
638671cb588SHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
639671cb588SHong Zhang     jmin = ai[k]; jmax = ai[k+1];
640671cb588SHong Zhang     ap = aa + jmin*bs2;
641671cb588SHong Zhang     for (j = jmin; j < jmax; j++){
642671cb588SHong Zhang       vj = aj[j];         /* block col. index */
643671cb588SHong Zhang       rtmp_ptr = rtmp + vj*bs2;
644671cb588SHong Zhang       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
645671cb588SHong Zhang     }
646671cb588SHong Zhang 
647671cb588SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
648671cb588SHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
649671cb588SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
650671cb588SHong Zhang 
651057f5ba7SHong Zhang     while (i < k){
652671cb588SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
653671cb588SHong Zhang 
654671cb588SHong Zhang       /* compute multiplier */
655671cb588SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
656671cb588SHong Zhang 
657671cb588SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
658671cb588SHong Zhang       diag = ba + i*bs2;
659671cb588SHong Zhang       u    = ba + ili*bs2;
660671cb588SHong Zhang       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
661671cb588SHong Zhang       Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
662671cb588SHong Zhang 
663671cb588SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
664671cb588SHong Zhang       Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
665671cb588SHong Zhang 
666671cb588SHong Zhang       /* update -U(i,k) */
667671cb588SHong Zhang       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
668671cb588SHong Zhang 
669671cb588SHong Zhang       /* add multiple of row i to k-th row ... */
670671cb588SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
671671cb588SHong Zhang       if (jmin < jmax){
672671cb588SHong Zhang         for (j=jmin; j<jmax; j++) {
673671cb588SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
674671cb588SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
675671cb588SHong Zhang           u = ba + j*bs2;
676671cb588SHong Zhang           Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
677671cb588SHong Zhang         }
678671cb588SHong Zhang 
679671cb588SHong Zhang         /* ... add i to row list for next nonzero entry */
680671cb588SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
681671cb588SHong Zhang         j     = bj[jmin];
682671cb588SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
683671cb588SHong Zhang       }
684671cb588SHong Zhang       i = nexti;
685671cb588SHong Zhang     }
686671cb588SHong Zhang 
687671cb588SHong Zhang     /* save nonzero entries in k-th row of U ... */
688671cb588SHong Zhang 
689671cb588SHong Zhang     /* invert diagonal block */
690671cb588SHong Zhang     diag = ba+k*bs2;
691671cb588SHong Zhang     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
692671cb588SHong Zhang     Kernel_A_gets_inverse_A(bs,diag,pivots,work);
693671cb588SHong Zhang 
694671cb588SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
695671cb588SHong Zhang     if (jmin < jmax) {
696671cb588SHong Zhang       for (j=jmin; j<jmax; j++){
697671cb588SHong Zhang          vj = bj[j];           /* block col. index of U */
698671cb588SHong Zhang          u   = ba + j*bs2;
699671cb588SHong Zhang          rtmp_ptr = rtmp + vj*bs2;
700671cb588SHong Zhang          for (k1=0; k1<bs2; k1++){
701671cb588SHong Zhang            *u++        = *rtmp_ptr;
702671cb588SHong Zhang            *rtmp_ptr++ = 0.0;
703671cb588SHong Zhang          }
704671cb588SHong Zhang       }
705671cb588SHong Zhang 
706671cb588SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
707671cb588SHong Zhang       il[k] = jmin;
708671cb588SHong Zhang       i     = bj[jmin];
709671cb588SHong Zhang       jl[k] = jl[i]; jl[i] = k;
710671cb588SHong Zhang     }
711671cb588SHong Zhang   }
712671cb588SHong Zhang 
713671cb588SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
714671cb588SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
715671cb588SHong Zhang   ierr = PetscFree(dk);CHKERRQ(ierr);
716671cb588SHong Zhang   ierr = PetscFree(pivots);CHKERRQ(ierr);
717671cb588SHong Zhang 
718671cb588SHong Zhang   C->factor    = FACTOR_CHOLESKY;
719671cb588SHong Zhang   C->assembled = PETSC_TRUE;
720671cb588SHong Zhang   C->preallocated = PETSC_TRUE;
721b0a32e0cSBarry Smith   PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */
722671cb588SHong Zhang   PetscFunctionReturn(0);
723671cb588SHong Zhang }
724671cb588SHong Zhang 
72549b5e25fSSatish Balay /*
726fcf159c0SHong Zhang     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
727cc0c071aSHong Zhang     Version for blocks 2 by 2.
72849b5e25fSSatish Balay */
7294a2ae208SSatish Balay #undef __FUNCT__
7304a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2"
7316f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat A,Mat *B)
73249b5e25fSSatish Balay {
73349b5e25fSSatish Balay   Mat                C = *B;
73491602c66SHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
735cc0c071aSHong Zhang   IS                 perm = b->row;
736cc0c071aSHong Zhang   int                *perm_ptr,ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
737cc0c071aSHong Zhang   int                *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
738a1723e09SHong Zhang   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
739cc0c071aSHong Zhang   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
74049b5e25fSSatish Balay 
74149b5e25fSSatish Balay   PetscFunctionBegin;
742671cb588SHong Zhang 
74391602c66SHong Zhang   /* initialization */
74491602c66SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
74591602c66SHong Zhang      window U(0:k, k:mbs-1).
74691602c66SHong Zhang      jl:    list of rows to be added to uneliminated rows
74791602c66SHong Zhang             i>= k: jl(i) is the first row to be added to row i
74891602c66SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
74991602c66SHong Zhang             jl(i) = mbs indicates the end of a list
75091602c66SHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
75191602c66SHong Zhang             row i of U */
75282502324SSatish Balay   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
753cc0c071aSHong Zhang   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
75482502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
75528de702eSHong Zhang   jl   = il + mbs;
75691602c66SHong Zhang   for (i=0; i<mbs; i++) {
7573845f261SHong Zhang     jl[i] = mbs; il[0] = 0;
75891602c66SHong Zhang   }
75982502324SSatish Balay   ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr);
76028de702eSHong Zhang   uik  = dk + 4;
761cc0c071aSHong Zhang   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
762cc0c071aSHong Zhang 
763cc0c071aSHong Zhang   /* check permutation */
764cc0c071aSHong Zhang   if (!a->permute){
765cc0c071aSHong Zhang     ai = a->i; aj = a->j; aa = a->a;
766cc0c071aSHong Zhang   } else {
767cc0c071aSHong Zhang     ai   = a->inew; aj = a->jnew;
76882502324SSatish Balay     ierr = PetscMalloc(4*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
769cc0c071aSHong Zhang     ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
77082502324SSatish Balay     ierr = PetscMalloc(ai[mbs]*sizeof(int),&a2anew);CHKERRQ(ierr);
771cc0c071aSHong Zhang     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr);
772cc0c071aSHong Zhang 
773cc0c071aSHong Zhang     for (i=0; i<mbs; i++){
774cc0c071aSHong Zhang       jmin = ai[i]; jmax = ai[i+1];
775cc0c071aSHong Zhang       for (j=jmin; j<jmax; j++){
776cc0c071aSHong Zhang         while (a2anew[j] != j){
777cc0c071aSHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
778cc0c071aSHong Zhang           for (k1=0; k1<4; k1++){
779cc0c071aSHong Zhang             dk[k1]       = aa[k*4+k1];
780cc0c071aSHong Zhang             aa[k*4+k1] = aa[j*4+k1];
781cc0c071aSHong Zhang             aa[j*4+k1] = dk[k1];
782cc0c071aSHong Zhang           }
783cc0c071aSHong Zhang         }
784cc0c071aSHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
785cc0c071aSHong Zhang         if (i > aj[j]){
786a1723e09SHong Zhang           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
787cc0c071aSHong Zhang           ap = aa + j*4;     /* ptr to the beginning of the block */
788cc0c071aSHong Zhang           dk[1] = ap[1];     /* swap ap[1] and ap[2] */
789cc0c071aSHong Zhang           ap[1] = ap[2];
790cc0c071aSHong Zhang           ap[2] = dk[1];
791cc0c071aSHong Zhang         }
792cc0c071aSHong Zhang       }
793cc0c071aSHong Zhang     }
794ac355199SBarry Smith     ierr = PetscFree(a2anew);CHKERRQ(ierr);
795cc0c071aSHong Zhang   }
7963845f261SHong Zhang 
79791602c66SHong Zhang   /* for each row k */
79891602c66SHong Zhang   for (k = 0; k<mbs; k++){
79991602c66SHong Zhang 
80091602c66SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
801cc0c071aSHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
802cc0c071aSHong Zhang     ap = aa + jmin*4;
80391602c66SHong Zhang     for (j = jmin; j < jmax; j++){
804cc0c071aSHong Zhang       vj = perm_ptr[aj[j]];         /* block col. index */
805cc0c071aSHong Zhang       rtmp_ptr = rtmp + vj*4;
806cc0c071aSHong Zhang       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
80791602c66SHong Zhang     }
80891602c66SHong Zhang 
80991602c66SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
810cc0c071aSHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
81191602c66SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
81291602c66SHong Zhang 
813057f5ba7SHong Zhang     while (i < k){
81491602c66SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
81591602c66SHong Zhang 
8163845f261SHong Zhang       /* compute multiplier */
81791602c66SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
8183845f261SHong Zhang 
8193845f261SHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
820cc0c071aSHong Zhang       diag = ba + i*4;
821cc0c071aSHong Zhang       u    = ba + ili*4;
822cc0c071aSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
823cc0c071aSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
824cc0c071aSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
825cc0c071aSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
8263845f261SHong Zhang 
8273845f261SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
828cc0c071aSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
829cc0c071aSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
830cc0c071aSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
831cc0c071aSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
8323845f261SHong Zhang 
8333845f261SHong Zhang       /* update -U(i,k): ba[ili] = uik */
834cc0c071aSHong Zhang       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
83591602c66SHong Zhang 
83691602c66SHong Zhang       /* add multiple of row i to k-th row ... */
83791602c66SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
83891602c66SHong Zhang       if (jmin < jmax){
8393845f261SHong Zhang         for (j=jmin; j<jmax; j++) {
8403845f261SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
841cc0c071aSHong Zhang           rtmp_ptr = rtmp + bj[j]*4;
842cc0c071aSHong Zhang           u = ba + j*4;
843cc0c071aSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
844cc0c071aSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
845cc0c071aSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
846cc0c071aSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
8473845f261SHong Zhang         }
8483845f261SHong Zhang 
84991602c66SHong Zhang         /* ... add i to row list for next nonzero entry */
85091602c66SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
85191602c66SHong Zhang         j     = bj[jmin];
85291602c66SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
85391602c66SHong Zhang       }
854a1723e09SHong Zhang       i = nexti;
85591602c66SHong Zhang     }
856cc0c071aSHong Zhang 
85791602c66SHong Zhang     /* save nonzero entries in k-th row of U ... */
8583845f261SHong Zhang 
859cc0c071aSHong Zhang     /* invert diagonal block */
860cc0c071aSHong Zhang     diag = ba+k*4;
861cc0c071aSHong Zhang     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
8623845f261SHong Zhang     ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr);
8633845f261SHong Zhang 
86491602c66SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
86591602c66SHong Zhang     if (jmin < jmax) {
86691602c66SHong Zhang       for (j=jmin; j<jmax; j++){
867cc0c071aSHong Zhang          vj = bj[j];           /* block col. index of U */
868cc0c071aSHong Zhang          u   = ba + j*4;
869cc0c071aSHong Zhang          rtmp_ptr = rtmp + vj*4;
870cc0c071aSHong Zhang          for (k1=0; k1<4; k1++){
871cc0c071aSHong Zhang            *u++        = *rtmp_ptr;
872cc0c071aSHong Zhang            *rtmp_ptr++ = 0.0;
873cc0c071aSHong Zhang          }
874cc0c071aSHong Zhang       }
8753845f261SHong Zhang 
87691602c66SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
87791602c66SHong Zhang       il[k] = jmin;
87891602c66SHong Zhang       i     = bj[jmin];
87991602c66SHong Zhang       jl[k] = jl[i]; jl[i] = k;
88091602c66SHong Zhang     }
88191602c66SHong Zhang   }
8823845f261SHong Zhang 
88349b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
88491602c66SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
8853845f261SHong Zhang   ierr = PetscFree(dk);CHKERRQ(ierr);
88691602c66SHong Zhang   if (a->permute) {
88791602c66SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
88891602c66SHong Zhang   }
889cc0c071aSHong Zhang   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
89091602c66SHong Zhang   C->factor    = FACTOR_CHOLESKY;
89149b5e25fSSatish Balay   C->assembled = PETSC_TRUE;
8925c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
893b0a32e0cSBarry Smith   PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */
89449b5e25fSSatish Balay   PetscFunctionReturn(0);
89549b5e25fSSatish Balay }
89691602c66SHong Zhang 
89749b5e25fSSatish Balay /*
89849b5e25fSSatish Balay       Version for when blocks are 2 by 2 Using natural ordering
89949b5e25fSSatish Balay */
9004a2ae208SSatish Balay #undef __FUNCT__
9014a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering"
9026f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat A,Mat *B)
90349b5e25fSSatish Balay {
90449b5e25fSSatish Balay   Mat                C = *B;
905ab27746eSHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
906ab27746eSHong Zhang   int                ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
90771149678SHong Zhang   int                *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
908ab27746eSHong Zhang   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
909ab27746eSHong Zhang   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
91049b5e25fSSatish Balay 
91149b5e25fSSatish Balay   PetscFunctionBegin;
912671cb588SHong Zhang 
913ab27746eSHong Zhang   /* initialization */
914ab27746eSHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
915ab27746eSHong Zhang      window U(0:k, k:mbs-1).
916ab27746eSHong Zhang      jl:    list of rows to be added to uneliminated rows
917ab27746eSHong Zhang             i>= k: jl(i) is the first row to be added to row i
918ab27746eSHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
919ab27746eSHong Zhang             jl(i) = mbs indicates the end of a list
920ab27746eSHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
921ab27746eSHong Zhang             row i of U */
92282502324SSatish Balay   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
923ab27746eSHong Zhang   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
92482502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
92528de702eSHong Zhang   jl   = il + mbs;
926ab27746eSHong Zhang   for (i=0; i<mbs; i++) {
927ab27746eSHong Zhang     jl[i] = mbs; il[0] = 0;
92849b5e25fSSatish Balay   }
92982502324SSatish Balay   ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr);
93028de702eSHong Zhang   uik  = dk + 4;
931ab27746eSHong Zhang 
932ab27746eSHong Zhang   ai = a->i; aj = a->j; aa = a->a;
933ab27746eSHong Zhang 
934ab27746eSHong Zhang   /* for each row k */
935ab27746eSHong Zhang   for (k = 0; k<mbs; k++){
936ab27746eSHong Zhang 
937ab27746eSHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
938ab27746eSHong Zhang     jmin = ai[k]; jmax = ai[k+1];
939ab27746eSHong Zhang     ap = aa + jmin*4;
940ab27746eSHong Zhang     for (j = jmin; j < jmax; j++){
941ab27746eSHong Zhang       vj = aj[j];         /* block col. index */
942ab27746eSHong Zhang       rtmp_ptr = rtmp + vj*4;
943ab27746eSHong Zhang       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
94449b5e25fSSatish Balay     }
945ab27746eSHong Zhang 
946ab27746eSHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
947ab27746eSHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
948ab27746eSHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
949ab27746eSHong Zhang 
950057f5ba7SHong Zhang     while (i < k){
951ab27746eSHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
952ab27746eSHong Zhang 
953ab27746eSHong Zhang       /* compute multiplier */
954ab27746eSHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
955ab27746eSHong Zhang 
956ab27746eSHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
957ab27746eSHong Zhang       diag = ba + i*4;
958ab27746eSHong Zhang       u    = ba + ili*4;
959ab27746eSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
960ab27746eSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
961ab27746eSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
962ab27746eSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
963ab27746eSHong Zhang 
964ab27746eSHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
965ab27746eSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
966ab27746eSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
967ab27746eSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
968ab27746eSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
969ab27746eSHong Zhang 
970ab27746eSHong Zhang       /* update -U(i,k): ba[ili] = uik */
971ab27746eSHong Zhang       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
972ab27746eSHong Zhang 
973ab27746eSHong Zhang       /* add multiple of row i to k-th row ... */
974ab27746eSHong Zhang       jmin = ili + 1; jmax = bi[i+1];
975ab27746eSHong Zhang       if (jmin < jmax){
976ab27746eSHong Zhang         for (j=jmin; j<jmax; j++) {
977ab27746eSHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
978ab27746eSHong Zhang           rtmp_ptr = rtmp + bj[j]*4;
979ab27746eSHong Zhang           u = ba + j*4;
980ab27746eSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
981ab27746eSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
982ab27746eSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
983ab27746eSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
98449b5e25fSSatish Balay         }
985ab27746eSHong Zhang 
986ab27746eSHong Zhang         /* ... add i to row list for next nonzero entry */
987ab27746eSHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
988ab27746eSHong Zhang         j     = bj[jmin];
989ab27746eSHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
99049b5e25fSSatish Balay       }
991ab27746eSHong Zhang       i = nexti;
99249b5e25fSSatish Balay     }
993ab27746eSHong Zhang 
994ab27746eSHong Zhang     /* save nonzero entries in k-th row of U ... */
995ab27746eSHong Zhang 
99649b5e25fSSatish Balay     /* invert diagonal block */
997ab27746eSHong Zhang     diag = ba+k*4;
998ab27746eSHong Zhang     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
999ab27746eSHong Zhang     ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr);
1000ab27746eSHong Zhang 
1001ab27746eSHong Zhang     jmin = bi[k]; jmax = bi[k+1];
1002ab27746eSHong Zhang     if (jmin < jmax) {
1003ab27746eSHong Zhang       for (j=jmin; j<jmax; j++){
1004ab27746eSHong Zhang          vj = bj[j];           /* block col. index of U */
1005ab27746eSHong Zhang          u   = ba + j*4;
1006ab27746eSHong Zhang          rtmp_ptr = rtmp + vj*4;
1007ab27746eSHong Zhang          for (k1=0; k1<4; k1++){
1008ab27746eSHong Zhang            *u++        = *rtmp_ptr;
1009ab27746eSHong Zhang            *rtmp_ptr++ = 0.0;
1010ab27746eSHong Zhang          }
1011ab27746eSHong Zhang       }
1012ab27746eSHong Zhang 
1013ab27746eSHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
1014ab27746eSHong Zhang       il[k] = jmin;
1015ab27746eSHong Zhang       i     = bj[jmin];
1016ab27746eSHong Zhang       jl[k] = jl[i]; jl[i] = k;
1017ab27746eSHong Zhang     }
101849b5e25fSSatish Balay   }
101949b5e25fSSatish Balay 
102049b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1021ab27746eSHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
1022ab27746eSHong Zhang   ierr = PetscFree(dk);CHKERRQ(ierr);
1023ab27746eSHong Zhang 
1024ab27746eSHong Zhang   C->factor    = FACTOR_CHOLESKY;
102549b5e25fSSatish Balay   C->assembled = PETSC_TRUE;
10265c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
1027b0a32e0cSBarry Smith   PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */
102849b5e25fSSatish Balay   PetscFunctionReturn(0);
102949b5e25fSSatish Balay }
103049b5e25fSSatish Balay 
103149b5e25fSSatish Balay /*
10325c0bcdfcSHong Zhang     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
103391602c66SHong Zhang     Version for blocks are 1 by 1.
103449b5e25fSSatish Balay */
10354a2ae208SSatish Balay #undef __FUNCT__
10364a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1"
10376f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat A,Mat *B)
103849b5e25fSSatish Balay {
103949b5e25fSSatish Balay   Mat                C = *B;
104049b5e25fSSatish Balay   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
104149b5e25fSSatish Balay   IS                 ip = b->row;
1042351d0355SHong Zhang   int                *rip,ierr,i,j,mbs = a->mbs,*bi = b->i,*bj = b->j;
1043cb718733SHong Zhang   int                *ai,*aj,*r;
1044ab27746eSHong Zhang   int                k,jmin,jmax,*jl,*il,vj,nexti,ili;
1045066653e3SSatish Balay   MatScalar          *rtmp;
10462d007305SHong Zhang   MatScalar          *ba = b->a,*aa,ak;
104749b5e25fSSatish Balay   MatScalar          dk,uikdi;
104849b5e25fSSatish Balay 
104949b5e25fSSatish Balay   PetscFunctionBegin;
105049b5e25fSSatish Balay   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
1051cb718733SHong Zhang   if (!a->permute){
10522d007305SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
10532d007305SHong Zhang   } else {
10542d007305SHong Zhang     ai = a->inew; aj = a->jnew;
105582502324SSatish Balay     ierr = PetscMalloc(ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
1056cb718733SHong Zhang     ierr = PetscMemcpy(aa,a->a,ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
105782502324SSatish Balay     ierr = PetscMalloc(ai[mbs]*sizeof(int),&r);CHKERRQ(ierr);
10582d007305SHong Zhang     ierr= PetscMemcpy(r,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr);
10592d007305SHong Zhang 
10602d007305SHong Zhang     jmin = ai[0]; jmax = ai[mbs];
10612d007305SHong Zhang     for (j=jmin; j<jmax; j++){
1062cb718733SHong Zhang       while (r[j] != j){
10632d007305SHong Zhang         k = r[j]; r[j] = r[k]; r[k] = k;
10642d007305SHong Zhang         ak = aa[k]; aa[k] = aa[j]; aa[j] = ak;
10652d007305SHong Zhang       }
10662d007305SHong Zhang     }
1067ac355199SBarry Smith     ierr = PetscFree(r);CHKERRQ(ierr);
10682d007305SHong Zhang   }
10692d007305SHong Zhang 
107091602c66SHong Zhang   /* initialization */
107149b5e25fSSatish Balay   /* il and jl record the first nonzero element in each row of the accessing
107249b5e25fSSatish Balay      window U(0:k, k:mbs-1).
107349b5e25fSSatish Balay      jl:    list of rows to be added to uneliminated rows
107449b5e25fSSatish Balay             i>= k: jl(i) is the first row to be added to row i
107549b5e25fSSatish Balay             i<  k: jl(i) is the row following row i in some list of rows
107649b5e25fSSatish Balay             jl(i) = mbs indicates the end of a list
107749b5e25fSSatish Balay      il(i): points to the first nonzero element in columns k,...,mbs-1 of
107849b5e25fSSatish Balay             row i of U */
107982502324SSatish Balay   ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
108082502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
108128de702eSHong Zhang   jl   = il + mbs;
108249b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
108349b5e25fSSatish Balay     rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
108449b5e25fSSatish Balay   }
108549b5e25fSSatish Balay 
108691602c66SHong Zhang   /* for each row k */
108749b5e25fSSatish Balay   for (k = 0; k<mbs; k++){
108849b5e25fSSatish Balay 
108991602c66SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
109049b5e25fSSatish Balay     jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
1091057f5ba7SHong Zhang 
109249b5e25fSSatish Balay     for (j = jmin; j < jmax; j++){
1093351d0355SHong Zhang       vj = rip[aj[j]];
1094ab27746eSHong Zhang       rtmp[vj] = aa[j];
109549b5e25fSSatish Balay     }
109649b5e25fSSatish Balay 
109791602c66SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
109849b5e25fSSatish Balay     dk = rtmp[k];
109949b5e25fSSatish Balay     i = jl[k]; /* first row to be added to k_th row  */
110049b5e25fSSatish Balay 
1101057f5ba7SHong Zhang     while (i < k){
110249b5e25fSSatish Balay       nexti = jl[i]; /* next row to be added to k_th row */
110349b5e25fSSatish Balay 
110491602c66SHong Zhang       /* compute multiplier, update D(k) and U(i,k) */
110549b5e25fSSatish Balay       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
110649b5e25fSSatish Balay       uikdi = - ba[ili]*ba[i];
110749b5e25fSSatish Balay       dk += uikdi*ba[ili];
1108658e7b3eSHong Zhang       ba[ili] = uikdi; /* -U(i,k) */
110949b5e25fSSatish Balay 
111091602c66SHong Zhang       /* add multiple of row i to k-th row ... */
111149b5e25fSSatish Balay       jmin = ili + 1; jmax = bi[i+1];
111249b5e25fSSatish Balay       if (jmin < jmax){
111349b5e25fSSatish Balay         for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
111491602c66SHong Zhang         /* ... add i to row list for next nonzero entry */
111549b5e25fSSatish Balay         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
111649b5e25fSSatish Balay         j     = bj[jmin];
111749b5e25fSSatish Balay         jl[i] = jl[j]; jl[j] = i; /* update jl */
111849b5e25fSSatish Balay       }
1119ab27746eSHong Zhang       i = nexti;
112049b5e25fSSatish Balay     }
112149b5e25fSSatish Balay 
112291602c66SHong Zhang     /* check for zero pivot and save diagoanl element */
112349b5e25fSSatish Balay     if (dk == 0.0){
112429bbc08cSBarry Smith       SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot");
1125671cb588SHong Zhang       /*
11261223c01bSHong Zhang     } else if (PetscRealPart(dk) < 0.0){
11271223c01bSHong Zhang       SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Negative pivot: d[%d] = %g\n",k,dk);
1128671cb588SHong Zhang       */
112949b5e25fSSatish Balay     }
113049b5e25fSSatish Balay 
113191602c66SHong Zhang     /* save nonzero entries in k-th row of U ... */
1132f3dd7b2dSKris Buschelman     ba[k] = 1.0/dk;
113349b5e25fSSatish Balay     jmin = bi[k]; jmax = bi[k+1];
113449b5e25fSSatish Balay     if (jmin < jmax) {
113549b5e25fSSatish Balay       for (j=jmin; j<jmax; j++){
1136cc0c071aSHong Zhang          vj = bj[j]; ba[j] = rtmp[vj]; rtmp[vj] = 0.0;
113749b5e25fSSatish Balay       }
113891602c66SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
113949b5e25fSSatish Balay       il[k] = jmin;
114049b5e25fSSatish Balay       i     = bj[jmin];
114149b5e25fSSatish Balay       jl[k] = jl[i]; jl[i] = k;
114249b5e25fSSatish Balay     }
114349b5e25fSSatish Balay   }
114449b5e25fSSatish Balay 
114549b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
114649b5e25fSSatish Balay   ierr = PetscFree(il);CHKERRQ(ierr);
1147cb718733SHong Zhang   if (a->permute){
1148cb718733SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
1149cb718733SHong Zhang   }
115049b5e25fSSatish Balay 
115149b5e25fSSatish Balay   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
11528be8c0c7SHong Zhang   C->factor    = FACTOR_CHOLESKY;
115349b5e25fSSatish Balay   C->assembled = PETSC_TRUE;
11545c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
1155b0a32e0cSBarry Smith   PetscLogFlops(b->mbs);
115649b5e25fSSatish Balay   PetscFunctionReturn(0);
115749b5e25fSSatish Balay }
115849b5e25fSSatish Balay 
1159671cb588SHong Zhang /*
1160671cb588SHong Zhang   Version for when blocks are 1 by 1 Using natural ordering
1161671cb588SHong Zhang */
11624a2ae208SSatish Balay #undef __FUNCT__
11634a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering"
1164671cb588SHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat A,Mat *B)
1165671cb588SHong Zhang {
1166671cb588SHong Zhang   Mat                C = *B;
1167671cb588SHong Zhang   Mat_SeqSBAIJ       *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data;
1168653a6975SHong Zhang   int                ierr,i,j,mbs = a->mbs;
1169653a6975SHong Zhang   int                *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
1170b00f7748SHong Zhang   int                k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz,ndamp = 0;
1171653a6975SHong Zhang   MatScalar          *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
117221c26570Svictorle   PetscReal          damping=b->factor_damping, zeropivot=b->factor_zeropivot,shift_amount;
117321c26570Svictorle   PetscTruth         damp,chshift;
117421c26570Svictorle   int                nshift=0;
1175653a6975SHong Zhang 
1176653a6975SHong Zhang   PetscFunctionBegin;
1177653a6975SHong Zhang   /* initialization */
1178653a6975SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
1179653a6975SHong Zhang      window U(0:k, k:mbs-1).
1180653a6975SHong Zhang      jl:    list of rows to be added to uneliminated rows
1181653a6975SHong Zhang             i>= k: jl(i) is the first row to be added to row i
1182653a6975SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
1183653a6975SHong Zhang             jl(i) = mbs indicates the end of a list
1184653a6975SHong Zhang      il(i): points to the first nonzero element in U(i,k:mbs-1)
1185653a6975SHong Zhang   */
1186653a6975SHong Zhang   ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
1187653a6975SHong Zhang   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
1188653a6975SHong Zhang   jl   = il + mbs;
1189b00f7748SHong Zhang 
119021c26570Svictorle   shift_amount = 0;
1191b00f7748SHong Zhang   do {
1192b00f7748SHong Zhang     damp = PETSC_FALSE;
119321c26570Svictorle     chshift = PETSC_FALSE;
1194653a6975SHong Zhang     for (i=0; i<mbs; i++) {
1195653a6975SHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1196653a6975SHong Zhang     }
1197653a6975SHong Zhang 
1198b00f7748SHong Zhang     for (k = 0; k<mbs; k++){ /* row k */
1199653a6975SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
1200653a6975SHong Zhang       nz   = ai[k+1] - ai[k];
1201653a6975SHong Zhang       acol = aj + ai[k];
1202653a6975SHong Zhang       aval = aa + ai[k];
1203653a6975SHong Zhang       bval = ba + bi[k];
1204653a6975SHong Zhang       while (nz -- ){
1205653a6975SHong Zhang         rtmp[*acol++] = *aval++;
1206653a6975SHong Zhang         *bval++       = 0.0; /* for in-place factorization */
1207653a6975SHong Zhang       }
1208b00f7748SHong Zhang       /* damp the diagonal of the matrix */
120921c26570Svictorle       if (ndamp||nshift) rtmp[k] += damping+shift_amount;
1210653a6975SHong Zhang 
1211653a6975SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k) != 0 */
1212653a6975SHong Zhang       dk = rtmp[k];
1213653a6975SHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
1214653a6975SHong Zhang 
1215653a6975SHong Zhang       while (i < k){
1216653a6975SHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
1217653a6975SHong Zhang 
1218653a6975SHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
1219653a6975SHong Zhang         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1220653a6975SHong Zhang         uikdi = - ba[ili]*ba[bi[i]];
1221653a6975SHong Zhang         dk   += uikdi*ba[ili];
1222653a6975SHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
1223653a6975SHong Zhang 
1224653a6975SHong Zhang         /* add multiple of row i to k-th row ... */
1225653a6975SHong Zhang         jmin = ili + 1;
1226653a6975SHong Zhang         nz   = bi[i+1] - jmin;
1227653a6975SHong Zhang         if (nz > 0){
1228653a6975SHong Zhang           bcol = bj + jmin;
1229653a6975SHong Zhang           bval = ba + jmin;
1230653a6975SHong Zhang           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
1231653a6975SHong Zhang           /* ... add i to row list for next nonzero entry */
1232653a6975SHong Zhang           il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
1233653a6975SHong Zhang           j     = bj[jmin];
1234653a6975SHong Zhang           jl[i] = jl[j]; jl[j] = i; /* update jl */
1235653a6975SHong Zhang         }
1236653a6975SHong Zhang         i = nexti;
1237653a6975SHong Zhang       }
1238653a6975SHong Zhang 
123921c26570Svictorle       /* check for zero pivot and save diagonal element */
124021c26570Svictorle       if (PetscRealPart(dk) < zeropivot && b->factor_shift){
124121c26570Svictorle 	PetscReal rs = -PetscRealPart(dk);
124221c26570Svictorle 	jmin      = bi[k]+1;
124321c26570Svictorle 	nz        = bi[k+1] - jmin;
124421c26570Svictorle 	if (nz){
124521c26570Svictorle 	  bcol = bj + jmin;
124621c26570Svictorle 	  bval = ba + jmin;
124721c26570Svictorle 	  while (nz--){
124821c26570Svictorle 	    rs += PetscAbsScalar(rtmp[*bcol++]);
124921c26570Svictorle 	  }
125021c26570Svictorle 	}
125121c26570Svictorle 	shift_amount = rs;
125221c26570Svictorle 	chshift  = PETSC_TRUE;
125321c26570Svictorle 	nshift++;
125421c26570Svictorle 	break;
125521c26570Svictorle       }
1256ffa0b812SHong Zhang       if (PetscRealPart(dk) < zeropivot){
1257ffa0b812SHong Zhang         if (damping == (PetscReal) PETSC_DECIDE) damping = -PetscRealPart(dk)/(k+1);
1258b00f7748SHong Zhang         if (damping > 0.0) {
1259b00f7748SHong Zhang           if (ndamp) damping *= 2.0;
1260b00f7748SHong Zhang           damp = PETSC_TRUE;
1261b00f7748SHong Zhang           ndamp++;
1262b00f7748SHong Zhang           break;
1263481c2c92SHong Zhang         } else if (PetscAbsScalar(dk) < zeropivot){
1264ffa0b812SHong Zhang           SETERRQ3(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %d value %g tolerance %g",k,PetscRealPart(dk),zeropivot);
1265481c2c92SHong Zhang         } else {
1266ffa0b812SHong Zhang           PetscLogInfo((PetscObject)A,"Negative pivot %g in row %d of Cholesky factorization\n",PetscRealPart(dk),k);
1267b00f7748SHong Zhang         }
1268064503c5SHong Zhang       }
1269653a6975SHong Zhang 
1270653a6975SHong Zhang       /* save nonzero entries in k-th row of U ... */
127182799104SHong Zhang       /* printf("%d, dk: %g, 1/dk: %g\n",k,dk,1/dk); */
1272653a6975SHong Zhang       ba[bi[k]] = 1.0/dk;
1273653a6975SHong Zhang       jmin      = bi[k]+1;
1274653a6975SHong Zhang       nz        = bi[k+1] - jmin;
1275653a6975SHong Zhang       if (nz){
1276653a6975SHong Zhang         bcol = bj + jmin;
1277653a6975SHong Zhang         bval = ba + jmin;
1278653a6975SHong Zhang         while (nz--){
1279653a6975SHong Zhang           *bval++       = rtmp[*bcol];
1280653a6975SHong Zhang           rtmp[*bcol++] = 0.0;
1281653a6975SHong Zhang         }
1282653a6975SHong Zhang         /* ... add k to row list for first nonzero entry in k-th row */
1283653a6975SHong Zhang         il[k] = jmin;
1284653a6975SHong Zhang         i     = bj[jmin];
1285653a6975SHong Zhang         jl[k] = jl[i]; jl[i] = k;
1286653a6975SHong Zhang       }
1287b00f7748SHong Zhang     } /* end of for (k = 0; k<mbs; k++) */
128821c26570Svictorle   } while (damp||chshift);
1289653a6975SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1290653a6975SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
1291653a6975SHong Zhang 
1292653a6975SHong Zhang   C->factor       = FACTOR_CHOLESKY;
1293653a6975SHong Zhang   C->assembled    = PETSC_TRUE;
1294653a6975SHong Zhang   C->preallocated = PETSC_TRUE;
1295653a6975SHong Zhang   PetscLogFlops(b->mbs);
1296b00f7748SHong Zhang   if (ndamp) {
1297b45a75daSHong Zhang     PetscLogInfo(0,"MatCholeskyFactorNumerical_SeqSBAIJ_1_NaturalOrdering: number of damping tries %d damping value %g\n",ndamp,damping);
1298b00f7748SHong Zhang   }
1299653a6975SHong Zhang   PetscFunctionReturn(0);
1300653a6975SHong Zhang }
1301653a6975SHong Zhang 
1302653a6975SHong Zhang #undef __FUNCT__
13034a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ"
130415e8a5b3SHong Zhang int MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,MatFactorInfo *info)
130549b5e25fSSatish Balay {
13064445ddedSHong Zhang   int ierr;
130749b5e25fSSatish Balay   Mat C;
130849b5e25fSSatish Balay 
130949b5e25fSSatish Balay   PetscFunctionBegin;
131015e8a5b3SHong Zhang   ierr = MatCholeskyFactorSymbolic(A,perm,info,&C);CHKERRQ(ierr);
1311a4ada70bSHong Zhang   ierr = MatCholeskyFactorNumeric(A,&C);CHKERRQ(ierr);
13124445ddedSHong Zhang   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
131349b5e25fSSatish Balay   PetscFunctionReturn(0);
131449b5e25fSSatish Balay }
131549b5e25fSSatish Balay 
131649b5e25fSSatish Balay 
1317