xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact.c (revision 77431f27ef6a796c710ee53da3b02f01a7247aeb)
149b5e25fSSatish Balay 
249b5e25fSSatish Balay #include "src/mat/impls/baij/seq/baij.h"
33a7fca6bSBarry Smith #include "src/mat/impls/sbaij/seq/sbaij.h"
449b5e25fSSatish Balay #include "src/inline/ilu.h"
549a6740bSHong Zhang #include "include/petscis.h"
649b5e25fSSatish Balay 
78dc9108eSHong Zhang #if !defined(PETSC_USE_COMPLEX)
85f9f512dSHong Zhang /*
95f9f512dSHong Zhang   input:
10c037c3f7SHong Zhang    F -- numeric factor
115f9f512dSHong Zhang   output:
12c037c3f7SHong Zhang    nneg, nzero, npos: matrix inertia
135f9f512dSHong Zhang */
145f9f512dSHong Zhang 
155f9f512dSHong Zhang #undef __FUNCT__
165f9f512dSHong Zhang #define __FUNCT__ "MatGetInertia_SeqSBAIJ"
17dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SeqSBAIJ(Mat F,int *nneig,int *nzero,int *npos)
185f9f512dSHong Zhang {
19638f5ce0SDinesh Kaushik   Mat_SeqSBAIJ *fact_ptr = (Mat_SeqSBAIJ*)F->data;
203e0d88b5SBarry Smith   PetscScalar  *dd = fact_ptr->a;
21eeeff2ecSHong Zhang   int          mbs=fact_ptr->mbs,bs=fact_ptr->bs,i,nneig_tmp,npos_tmp,
22eeeff2ecSHong Zhang                *fi = fact_ptr->i;
235f9f512dSHong Zhang 
245f9f512dSHong Zhang   PetscFunctionBegin;
25*77431f27SBarry Smith   if (bs != 1) SETERRQ1(PETSC_ERR_SUP,"No support for bs: %D >1 yet",bs);
26eeeff2ecSHong Zhang   nneig_tmp = 0; npos_tmp = 0;
27eeeff2ecSHong Zhang   for (i=0; i<mbs; i++){
28eeeff2ecSHong Zhang     if (PetscRealPart(dd[*fi]) > 0.0){
29eeeff2ecSHong Zhang       npos_tmp++;
30eeeff2ecSHong Zhang     } else if (PetscRealPart(dd[*fi]) < 0.0){
31eeeff2ecSHong Zhang       nneig_tmp++;
325f9f512dSHong Zhang     }
33eeeff2ecSHong Zhang     fi++;
343e0d88b5SBarry Smith   }
35eeeff2ecSHong Zhang   if (nneig) *nneig = nneig_tmp;
36eeeff2ecSHong Zhang   if (npos)  *npos  = npos_tmp;
37eeeff2ecSHong Zhang   if (nzero) *nzero = mbs - nneig_tmp - npos_tmp;
38eeeff2ecSHong Zhang 
395f9f512dSHong Zhang   PetscFunctionReturn(0);
405f9f512dSHong Zhang }
418dc9108eSHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */
425f9f512dSHong Zhang 
435f9f512dSHong Zhang /* Using Modified Sparse Row (MSR) storage.
445f9f512dSHong Zhang See page 85, "Iterative Methods ..." by Saad. */
455f9f512dSHong Zhang /*
465f9f512dSHong Zhang     Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP.
475f9f512dSHong Zhang */
482d9a3abdSHong Zhang /* Use Modified Sparse Row storage for u and ju, see Saad pp.85 */
494a2ae208SSatish Balay #undef __FUNCT__
504a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ"
51dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *B)
5249b5e25fSSatish Balay {
5349b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b;
546849ba73SBarry Smith   PetscErrorCode ierr;
556849ba73SBarry Smith   int          *rip,i,mbs = a->mbs,*ai,*aj;
56066653e3SSatish Balay   int          *jutmp,bs = a->bs,bs2=a->bs2;
5787fe1012SHong Zhang   int          m,realloc = 0,prow;
58d60e5362SHong Zhang   int          *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd;
59064503c5SHong Zhang   int          *il,ili,nextprow;
6015e8a5b3SHong Zhang   PetscReal    f = info->fill;
61671cb588SHong Zhang   PetscTruth   perm_identity;
6249b5e25fSSatish Balay 
6349b5e25fSSatish Balay   PetscFunctionBegin;
64cb718733SHong Zhang   /* check whether perm is the identity mapping */
65671cb588SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
66064503c5SHong Zhang 
67064503c5SHong Zhang   /* -- inplace factorization, i.e., use sbaij for *B -- */
68064503c5SHong Zhang   if (perm_identity && bs==1 ){
69064503c5SHong Zhang     if (!perm_identity) a->permute = PETSC_TRUE;
70064503c5SHong Zhang 
71064503c5SHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
72064503c5SHong Zhang 
73064503c5SHong Zhang   if (perm_identity){ /* without permutation */
74064503c5SHong Zhang     ai = a->i; aj = a->j;
75064503c5SHong Zhang   } else {            /* non-trivial permutation */
76064503c5SHong Zhang     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
77064503c5SHong Zhang     ai = a->inew; aj = a->jnew;
78064503c5SHong Zhang   }
79064503c5SHong Zhang 
80064503c5SHong Zhang   /* initialization */
81064503c5SHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(int),&iu);CHKERRQ(ierr);
82064503c5SHong Zhang   umax  = (int)(f*ai[mbs] + 1);
83064503c5SHong Zhang   ierr  = PetscMalloc(umax*sizeof(int),&ju);CHKERRQ(ierr);
84064503c5SHong Zhang   iu[0] = 0;
85064503c5SHong Zhang   juidx = 0; /* index for ju */
86064503c5SHong Zhang   ierr  = PetscMalloc((3*mbs+1)*sizeof(int),&jl);CHKERRQ(ierr); /* linked list for getting pivot row */
87064503c5SHong Zhang   q     = jl + mbs;   /* linked list for col index of active row */
88064503c5SHong Zhang   il    = q  + mbs;
89064503c5SHong Zhang   for (i=0; i<mbs; i++){
90064503c5SHong Zhang     jl[i] = mbs;
91064503c5SHong Zhang     q[i]  = 0;
92064503c5SHong Zhang     il[i] = 0;
93064503c5SHong Zhang   }
94064503c5SHong Zhang 
95064503c5SHong Zhang   /* for each row k */
96064503c5SHong Zhang   for (k=0; k<mbs; k++){
97064503c5SHong Zhang     nzk  = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */
98064503c5SHong Zhang     q[k] = mbs;
99064503c5SHong Zhang     /* initialize nonzero structure of k-th row to row rip[k] of A */
100064503c5SHong Zhang     jmin = ai[rip[k]] +1; /* exclude diag[k] */
101064503c5SHong Zhang     jmax = ai[rip[k]+1];
102064503c5SHong Zhang     for (j=jmin; j<jmax; j++){
103064503c5SHong Zhang       vj = rip[aj[j]]; /* col. value */
104064503c5SHong Zhang       if(vj > k){
105064503c5SHong Zhang         qm = k;
106064503c5SHong Zhang         do {
107064503c5SHong Zhang           m  = qm; qm = q[m];
108064503c5SHong Zhang         } while(qm < vj);
109064503c5SHong Zhang         if (qm == vj) {
110e005ede5SBarry Smith           SETERRQ(PETSC_ERR_PLIB,"Duplicate entry in A\n");
111064503c5SHong Zhang         }
112064503c5SHong Zhang         nzk++;
113064503c5SHong Zhang         q[m]  = vj;
114064503c5SHong Zhang         q[vj] = qm;
115064503c5SHong Zhang       } /* if(vj > k) */
116064503c5SHong Zhang     } /* for (j=jmin; j<jmax; j++) */
117064503c5SHong Zhang 
118064503c5SHong Zhang     /* modify nonzero structure of k-th row by computing fill-in
119064503c5SHong Zhang        for each row i to be merged in */
120064503c5SHong Zhang     prow = k;
121064503c5SHong Zhang     prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */
122064503c5SHong Zhang 
123064503c5SHong Zhang     while (prow < k){
124064503c5SHong Zhang       nextprow = jl[prow];
125064503c5SHong Zhang 
126064503c5SHong Zhang       /* merge row prow into k-th row */
127064503c5SHong Zhang       ili = il[prow];
128064503c5SHong Zhang       jmin = ili + 1;  /* points to 2nd nzero entry in U(prow,k:mbs-1) */
129064503c5SHong Zhang       jmax = iu[prow+1];
130064503c5SHong Zhang       qm = k;
131064503c5SHong Zhang       for (j=jmin; j<jmax; j++){
132064503c5SHong Zhang         vj = ju[j];
133064503c5SHong Zhang         do {
134064503c5SHong Zhang           m = qm; qm = q[m];
135064503c5SHong Zhang         } while (qm < vj);
136064503c5SHong Zhang         if (qm != vj){  /* a fill */
137064503c5SHong Zhang           nzk++; q[m] = vj; q[vj] = qm; qm = vj;
138064503c5SHong Zhang         }
139064503c5SHong Zhang       } /* end of for (j=jmin; j<jmax; j++) */
140064503c5SHong Zhang       if (jmin < jmax){
141064503c5SHong Zhang         il[prow] = jmin;
142064503c5SHong Zhang         j = ju[jmin];
143064503c5SHong Zhang         jl[prow] = jl[j]; jl[j] = prow;  /* update jl */
144064503c5SHong Zhang       }
145064503c5SHong Zhang       prow = nextprow;
146064503c5SHong Zhang     }
147064503c5SHong Zhang 
148064503c5SHong Zhang     /* update il and jl */
149064503c5SHong Zhang     if (nzk > 0){
150064503c5SHong Zhang       i = q[k]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
151064503c5SHong Zhang       jl[k] = jl[i]; jl[i] = k;
152064503c5SHong Zhang       il[k] = iu[k] + 1;
153064503c5SHong Zhang     }
154064503c5SHong Zhang     iu[k+1] = iu[k] + nzk + 1;  /* include diag[k] */
155064503c5SHong Zhang 
156064503c5SHong Zhang     /* allocate more space to ju if needed */
157064503c5SHong Zhang     if (iu[k+1] > umax) {
158064503c5SHong Zhang       /* estimate how much additional space we will need */
159064503c5SHong Zhang       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
160064503c5SHong Zhang       /* just double the memory each time */
161064503c5SHong Zhang       maxadd = umax;
162064503c5SHong Zhang       if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2;
163064503c5SHong Zhang       umax += maxadd;
164064503c5SHong Zhang 
165064503c5SHong Zhang       /* allocate a longer ju */
166064503c5SHong Zhang       ierr = PetscMalloc(umax*sizeof(int),&jutmp);CHKERRQ(ierr);
167064503c5SHong Zhang       ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(int));CHKERRQ(ierr);
168064503c5SHong Zhang       ierr = PetscFree(ju);CHKERRQ(ierr);
169064503c5SHong Zhang       ju   = jutmp;
170064503c5SHong Zhang       realloc++; /* count how many times we realloc */
171064503c5SHong Zhang     }
172064503c5SHong Zhang 
173064503c5SHong Zhang     /* save nonzero structure of k-th row in ju */
174064503c5SHong Zhang     ju[juidx++] = k; /* diag[k] */
175064503c5SHong Zhang     i = k;
176064503c5SHong Zhang     while (nzk --) {
177064503c5SHong Zhang       i           = q[i];
178064503c5SHong Zhang       ju[juidx++] = i;
179064503c5SHong Zhang     }
180064503c5SHong Zhang   }
181064503c5SHong Zhang 
182064503c5SHong Zhang   if (ai[mbs] != 0) {
183064503c5SHong Zhang     PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
184*77431f27SBarry Smith     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Reallocs %D Fill ratio:given %g needed %g\n",realloc,f,af);
185064503c5SHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Run with -pc_cholesky_fill %g or use \n",af);
186064503c5SHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:PCCholeskySetFill(pc,%g);\n",af);
187064503c5SHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:for best performance.\n");
188064503c5SHong Zhang   } else {
189064503c5SHong Zhang      PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Empty matrix.\n");
190064503c5SHong Zhang   }
191064503c5SHong Zhang 
192064503c5SHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
193064503c5SHong Zhang   /* ierr = PetscFree(q);CHKERRQ(ierr); */
194064503c5SHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
195064503c5SHong Zhang 
196064503c5SHong Zhang   /* put together the new matrix */
197be5d1d56SKris Buschelman   ierr = MatCreate(A->comm,bs*mbs,bs*mbs,bs*mbs,bs*mbs,B);CHKERRQ(ierr);
198be5d1d56SKris Buschelman   ierr = MatSetType(*B,A->type_name);CHKERRQ(ierr);
199be5d1d56SKris Buschelman   ierr = MatSeqSBAIJSetPreallocation(*B,bs,0,PETSC_NULL);CHKERRQ(ierr);
200be5d1d56SKris Buschelman 
201064503c5SHong Zhang   /* PetscLogObjectParent(*B,iperm); */
202064503c5SHong Zhang   b = (Mat_SeqSBAIJ*)(*B)->data;
203064503c5SHong Zhang   ierr = PetscFree(b->imax);CHKERRQ(ierr);
204064503c5SHong Zhang   b->singlemalloc = PETSC_FALSE;
205064503c5SHong Zhang   /* the next line frees the default space generated by the Create() */
206064503c5SHong Zhang   ierr = PetscFree(b->a);CHKERRQ(ierr);
207064503c5SHong Zhang   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
208064503c5SHong Zhang   ierr = PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
209064503c5SHong Zhang   b->j    = ju;
210064503c5SHong Zhang   b->i    = iu;
211064503c5SHong Zhang   b->diag = 0;
212064503c5SHong Zhang   b->ilen = 0;
213064503c5SHong Zhang   b->imax = 0;
214064503c5SHong Zhang   b->row  = perm;
21515e8a5b3SHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
216064503c5SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
217064503c5SHong Zhang   b->icol = perm;
218064503c5SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
219064503c5SHong Zhang   ierr    = PetscMalloc((bs*mbs+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
220064503c5SHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.
221064503c5SHong Zhang      Allocate idnew, solve_work, new a, new j */
222064503c5SHong Zhang   PetscLogObjectMemory(*B,(iu[mbs]-mbs)*(sizeof(int)+sizeof(MatScalar)));
2236c6c5352SBarry Smith   b->maxnz = b->nz = iu[mbs];
224064503c5SHong Zhang 
225064503c5SHong Zhang   (*B)->factor                 = FACTOR_CHOLESKY;
226064503c5SHong Zhang   (*B)->info.factor_mallocs    = realloc;
227064503c5SHong Zhang   (*B)->info.fill_ratio_given  = f;
228064503c5SHong Zhang   if (ai[mbs] != 0) {
229064503c5SHong Zhang     (*B)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
230064503c5SHong Zhang   } else {
231064503c5SHong Zhang     (*B)->info.fill_ratio_needed = 0.0;
232064503c5SHong Zhang   }
233064503c5SHong Zhang 
234064503c5SHong Zhang 
235b45a75daSHong Zhang   (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering;
236b45a75daSHong Zhang   (*B)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
237e005ede5SBarry Smith   (*B)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
238064503c5SHong Zhang   PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=1\n");
239064503c5SHong Zhang 
240064503c5SHong Zhang   PetscFunctionReturn(0);
241064503c5SHong Zhang   }
242064503c5SHong Zhang   /* -----------  end of new code --------------------*/
243064503c5SHong Zhang 
244064503c5SHong Zhang 
245671cb588SHong Zhang   if (!perm_identity) a->permute = PETSC_TRUE;
2462d007305SHong Zhang 
247671cb588SHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
248671cb588SHong Zhang 
249671cb588SHong Zhang   if (perm_identity){ /* without permutation */
2502d007305SHong Zhang     ai = a->i; aj = a->j;
2512d007305SHong Zhang   } else {            /* non-trivial permutation */
252323b833fSBarry Smith     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
2532d007305SHong Zhang     ai = a->inew; aj = a->jnew;
2542d007305SHong Zhang   }
2552d007305SHong Zhang 
25649b5e25fSSatish Balay   /* initialization */
257b0a32e0cSBarry Smith   ierr  = PetscMalloc((mbs+1)*sizeof(int),&iu);CHKERRQ(ierr);
25849b5e25fSSatish Balay   umax  = (int)(f*ai[mbs] + 1); umax += mbs + 1;
25982502324SSatish Balay   ierr  = PetscMalloc(umax*sizeof(int),&ju);CHKERRQ(ierr);
26049b5e25fSSatish Balay   iu[0] = mbs+1;
26187fe1012SHong Zhang   juidx = mbs + 1; /* index for ju */
26287fe1012SHong Zhang   ierr  = PetscMalloc(2*mbs*sizeof(int),&jl);CHKERRQ(ierr); /* linked list for pivot row */
26387fe1012SHong Zhang   q     = jl + mbs;   /* linked list for col index */
26449b5e25fSSatish Balay   for (i=0; i<mbs; i++){
26587fe1012SHong Zhang     jl[i] = mbs;
26687fe1012SHong Zhang     q[i] = 0;
26749b5e25fSSatish Balay   }
26849b5e25fSSatish Balay 
26949b5e25fSSatish Balay   /* for each row k */
27049b5e25fSSatish Balay   for (k=0; k<mbs; k++){
271064503c5SHong Zhang     for (i=0; i<mbs; i++) q[i] = 0;  /* to be removed! */
27249b5e25fSSatish Balay     nzk  = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */
27349b5e25fSSatish Balay     q[k] = mbs;
27449b5e25fSSatish Balay     /* initialize nonzero structure of k-th row to row rip[k] of A */
275064503c5SHong Zhang     jmin = ai[rip[k]] +1; /* exclude diag[k] */
27649b5e25fSSatish Balay     jmax = ai[rip[k]+1];
27749b5e25fSSatish Balay     for (j=jmin; j<jmax; j++){
278cb718733SHong Zhang       vj = rip[aj[j]]; /* col. value */
27949b5e25fSSatish Balay       if(vj > k){
28049b5e25fSSatish Balay         qm = k;
28149b5e25fSSatish Balay         do {
28249b5e25fSSatish Balay           m  = qm; qm = q[m];
28349b5e25fSSatish Balay         } while(qm < vj);
28449b5e25fSSatish Balay         if (qm == vj) {
285e005ede5SBarry Smith           SETERRQ(PETSC_ERR_PLIB,"Duplicate entry in A\n");
28649b5e25fSSatish Balay         }
28749b5e25fSSatish Balay         nzk++;
28849b5e25fSSatish Balay         q[m]  = vj;
28949b5e25fSSatish Balay         q[vj] = qm;
29049b5e25fSSatish Balay       } /* if(vj > k) */
29149b5e25fSSatish Balay     } /* for (j=jmin; j<jmax; j++) */
29249b5e25fSSatish Balay 
29349b5e25fSSatish Balay     /* modify nonzero structure of k-th row by computing fill-in
29449b5e25fSSatish Balay        for each row i to be merged in */
29587fe1012SHong Zhang     prow = k;
29687fe1012SHong Zhang     prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */
29787fe1012SHong Zhang 
29887fe1012SHong Zhang     while (prow < k){
29987fe1012SHong Zhang       /* merge row prow into k-th row */
30087fe1012SHong Zhang       jmin = iu[prow] + 1; jmax = iu[prow+1];
30149b5e25fSSatish Balay       qm = k;
30287fe1012SHong Zhang       for (j=jmin; j<jmax; j++){
30349b5e25fSSatish Balay         vj = ju[j];
30449b5e25fSSatish Balay         do {
30549b5e25fSSatish Balay           m = qm; qm = q[m];
30649b5e25fSSatish Balay         } while (qm < vj);
30749b5e25fSSatish Balay         if (qm != vj){
30849b5e25fSSatish Balay          nzk++; q[m] = vj; q[vj] = qm; qm = vj;
30949b5e25fSSatish Balay         }
31049b5e25fSSatish Balay       }
31187fe1012SHong Zhang       prow = jl[prow]; /* next pivot row */
31249b5e25fSSatish Balay     }
31349b5e25fSSatish Balay 
31449b5e25fSSatish Balay     /* add k to row list for first nonzero element in k-th row */
31549b5e25fSSatish Balay     if (nzk > 0){
31649b5e25fSSatish Balay       i = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */
31749b5e25fSSatish Balay       jl[k] = jl[i]; jl[i] = k;
31849b5e25fSSatish Balay     }
31987fe1012SHong Zhang     iu[k+1] = iu[k] + nzk;
32049b5e25fSSatish Balay 
32149b5e25fSSatish Balay     /* allocate more space to ju if needed */
3223078e140SBarry Smith     if (iu[k+1] > umax) {
32349b5e25fSSatish Balay       /* estimate how much additional space we will need */
32449b5e25fSSatish Balay       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
32549b5e25fSSatish Balay       /* just double the memory each time */
32649b5e25fSSatish Balay       maxadd = umax;
32749b5e25fSSatish Balay       if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2;
32849b5e25fSSatish Balay       umax += maxadd;
32949b5e25fSSatish Balay 
3309f9cb213SHong Zhang       /* allocate a longer ju */
33182502324SSatish Balay       ierr = PetscMalloc(umax*sizeof(int),&jutmp);CHKERRQ(ierr);
33249b5e25fSSatish Balay       ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(int));CHKERRQ(ierr);
33349b5e25fSSatish Balay       ierr = PetscFree(ju);CHKERRQ(ierr);
3349f9cb213SHong Zhang       ju   = jutmp;
33549b5e25fSSatish Balay       realloc++; /* count how many times we realloc */
33649b5e25fSSatish Balay     }
33749b5e25fSSatish Balay 
33849b5e25fSSatish Balay     /* save nonzero structure of k-th row in ju */
33949b5e25fSSatish Balay     i=k;
34087fe1012SHong Zhang     while (nzk --) {
34149b5e25fSSatish Balay       i           = q[i];
34287fe1012SHong Zhang       ju[juidx++] = i;
34349b5e25fSSatish Balay     }
34477f73120SHong Zhang   }
34549b5e25fSSatish Balay 
34649b5e25fSSatish Balay   if (ai[mbs] != 0) {
34749b5e25fSSatish Balay     PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
348*77431f27SBarry Smith     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Reallocs %D Fill ratio:given %g needed %g\n",realloc,f,af);
3493078e140SBarry Smith     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Run with -pc_cholesky_fill %g or use \n",af);
3503078e140SBarry Smith     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:PCCholeskySetFill(pc,%g);\n",af);
351b0a32e0cSBarry Smith     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:for best performance.\n");
35249b5e25fSSatish Balay   } else {
353b0a32e0cSBarry Smith      PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Empty matrix.\n");
35449b5e25fSSatish Balay   }
35549b5e25fSSatish Balay 
35677f73120SHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
35787fe1012SHong Zhang   /* ierr = PetscFree(q);CHKERRQ(ierr); */
35849b5e25fSSatish Balay   ierr = PetscFree(jl);CHKERRQ(ierr);
35949b5e25fSSatish Balay 
36049b5e25fSSatish Balay   /* put together the new matrix */
361be5d1d56SKris Buschelman   ierr = MatCreate(A->comm,bs*mbs,bs*mbs,bs*mbs,bs*mbs,B);CHKERRQ(ierr);
362be5d1d56SKris Buschelman   ierr = MatSetType(*B,A->type_name);CHKERRQ(ierr);
363be5d1d56SKris Buschelman   ierr = MatSeqSBAIJSetPreallocation(*B,bs,0,PETSC_NULL);CHKERRQ(ierr);
364be5d1d56SKris Buschelman 
365b0a32e0cSBarry Smith   /* PetscLogObjectParent(*B,iperm); */
36649b5e25fSSatish Balay   b = (Mat_SeqSBAIJ*)(*B)->data;
36749b5e25fSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
36849b5e25fSSatish Balay   b->singlemalloc = PETSC_FALSE;
36949b5e25fSSatish Balay   /* the next line frees the default space generated by the Create() */
37049b5e25fSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
37149b5e25fSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
37282502324SSatish Balay   ierr = PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
37349b5e25fSSatish Balay   b->j    = ju;
37449b5e25fSSatish Balay   b->i    = iu;
37549b5e25fSSatish Balay   b->diag = 0;
37649b5e25fSSatish Balay   b->ilen = 0;
37749b5e25fSSatish Balay   b->imax = 0;
37877f73120SHong Zhang   b->row  = perm;
37915e8a5b3SHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
38077f73120SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
381cb718733SHong Zhang   b->icol = perm;
382cb718733SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
38387828ca2SBarry Smith   ierr    = PetscMalloc((bs*mbs+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
38449b5e25fSSatish Balay   /* In b structure:  Free imax, ilen, old a, old j.
38549b5e25fSSatish Balay      Allocate idnew, solve_work, new a, new j */
386b0a32e0cSBarry Smith   PetscLogObjectMemory(*B,(iu[mbs]-mbs)*(sizeof(int)+sizeof(MatScalar)));
3876c6c5352SBarry Smith   b->maxnz = b->nz = iu[mbs];
38849b5e25fSSatish Balay 
3895c0bcdfcSHong Zhang   (*B)->factor                 = FACTOR_CHOLESKY;
39049b5e25fSSatish Balay   (*B)->info.factor_mallocs    = realloc;
39149b5e25fSSatish Balay   (*B)->info.fill_ratio_given  = f;
39249b5e25fSSatish Balay   if (ai[mbs] != 0) {
39349b5e25fSSatish Balay     (*B)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
39449b5e25fSSatish Balay   } else {
39549b5e25fSSatish Balay     (*B)->info.fill_ratio_needed = 0.0;
39649b5e25fSSatish Balay   }
3970aa0ce06SHong Zhang 
398671cb588SHong Zhang   if (perm_identity){
399671cb588SHong Zhang     switch (bs) {
400671cb588SHong Zhang       case 1:
401671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering;
402671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
4034d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=1\n");
404671cb588SHong Zhang         break;
405671cb588SHong Zhang       case 2:
406671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
407671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
4084d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
409671cb588SHong Zhang         break;
410671cb588SHong Zhang       case 3:
411671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
412671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
4134d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:sing special in-place natural ordering factor and solve BS=3\n");
414671cb588SHong Zhang         break;
415671cb588SHong Zhang       case 4:
416671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
417671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
4184d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
419671cb588SHong Zhang         break;
420671cb588SHong Zhang       case 5:
421671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
422671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
4234d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
424671cb588SHong Zhang         break;
425671cb588SHong Zhang       case 6:
426671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
427671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
4284d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
429671cb588SHong Zhang         break;
430671cb588SHong Zhang       case 7:
431671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
432671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
4334d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
434671cb588SHong Zhang       break;
435671cb588SHong Zhang       default:
436671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering;
437671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_N_NaturalOrdering;
4384d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS>7\n");
439671cb588SHong Zhang       break;
440671cb588SHong Zhang     }
441671cb588SHong Zhang   }
442671cb588SHong Zhang 
44349b5e25fSSatish Balay   PetscFunctionReturn(0);
44449b5e25fSSatish Balay }
44549b5e25fSSatish Balay 
44635d8aa7fSBarry Smith 
4474a2ae208SSatish Balay #undef __FUNCT__
4484a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N"
449dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat A,Mat *B)
45049b5e25fSSatish Balay {
45149b5e25fSSatish Balay   Mat                C = *B;
4524c16a6a6SHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
4534c16a6a6SHong Zhang   IS                 perm = b->row;
4546849ba73SBarry Smith   PetscErrorCode ierr;
4556849ba73SBarry Smith   int                *perm_ptr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
4564c16a6a6SHong Zhang   int                *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
4574c16a6a6SHong Zhang   int                bs=a->bs,bs2 = a->bs2;
4584c16a6a6SHong Zhang   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
4594c16a6a6SHong Zhang   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
46028de702eSHong Zhang   MatScalar          *work;
4614c16a6a6SHong Zhang   int                *pivots;
4624c16a6a6SHong Zhang 
4634c16a6a6SHong Zhang   PetscFunctionBegin;
4649706f043SHong Zhang 
4654c16a6a6SHong Zhang   /* initialization */
46682502324SSatish Balay   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
4674c16a6a6SHong Zhang   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
46882502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
46928de702eSHong Zhang   jl   = il + mbs;
4704c16a6a6SHong Zhang   for (i=0; i<mbs; i++) {
4714c16a6a6SHong Zhang     jl[i] = mbs; il[0] = 0;
4724c16a6a6SHong Zhang   }
473b0a32e0cSBarry Smith   ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr);
47428de702eSHong Zhang   uik  = dk + bs2;
47528de702eSHong Zhang   work = uik + bs2;
47682502324SSatish Balay   ierr = PetscMalloc(bs*sizeof(int),&pivots);CHKERRQ(ierr);
4774c16a6a6SHong Zhang 
4784c16a6a6SHong Zhang   ierr  = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
4794c16a6a6SHong Zhang 
4804c16a6a6SHong Zhang   /* check permutation */
4814c16a6a6SHong Zhang   if (!a->permute){
4824c16a6a6SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
4834c16a6a6SHong Zhang   } else {
4844c16a6a6SHong Zhang     ai   = a->inew; aj = a->jnew;
48582502324SSatish Balay     ierr = PetscMalloc(bs2*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
4864c16a6a6SHong Zhang     ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
48782502324SSatish Balay     ierr = PetscMalloc(ai[mbs]*sizeof(int),&a2anew);CHKERRQ(ierr);
4884c16a6a6SHong Zhang     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr);
4894c16a6a6SHong Zhang 
4904c16a6a6SHong Zhang     for (i=0; i<mbs; i++){
4914c16a6a6SHong Zhang       jmin = ai[i]; jmax = ai[i+1];
4924c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++){
4934c16a6a6SHong Zhang         while (a2anew[j] != j){
4944c16a6a6SHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
4954c16a6a6SHong Zhang           for (k1=0; k1<bs2; k1++){
4964c16a6a6SHong Zhang             dk[k1]       = aa[k*bs2+k1];
4974c16a6a6SHong Zhang             aa[k*bs2+k1] = aa[j*bs2+k1];
4984c16a6a6SHong Zhang             aa[j*bs2+k1] = dk[k1];
4994c16a6a6SHong Zhang           }
5004c16a6a6SHong Zhang         }
5014c16a6a6SHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
5024c16a6a6SHong Zhang         if (i > aj[j]){
5034c16a6a6SHong Zhang           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
5044c16a6a6SHong Zhang           ap = aa + j*bs2;                     /* ptr to the beginning of j-th block of aa */
5054c16a6a6SHong Zhang           for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
5064c16a6a6SHong Zhang           for (k=0; k<bs; k++){               /* j-th block of aa <- dk^T */
5074c16a6a6SHong Zhang             for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1];
5084c16a6a6SHong Zhang           }
5094c16a6a6SHong Zhang         }
5104c16a6a6SHong Zhang       }
5114c16a6a6SHong Zhang     }
512323b833fSBarry Smith     ierr = PetscFree(a2anew);CHKERRQ(ierr);
5134c16a6a6SHong Zhang   }
5144c16a6a6SHong Zhang 
5154c16a6a6SHong Zhang   /* for each row k */
5164c16a6a6SHong Zhang   for (k = 0; k<mbs; k++){
5174c16a6a6SHong Zhang 
5184c16a6a6SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
5194c16a6a6SHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
520057f5ba7SHong Zhang 
5214c16a6a6SHong Zhang     ap = aa + jmin*bs2;
5224c16a6a6SHong Zhang     for (j = jmin; j < jmax; j++){
5234c16a6a6SHong Zhang       vj = perm_ptr[aj[j]];         /* block col. index */
5244c16a6a6SHong Zhang       rtmp_ptr = rtmp + vj*bs2;
5254c16a6a6SHong Zhang       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
5264c16a6a6SHong Zhang     }
5274c16a6a6SHong Zhang 
5284c16a6a6SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
5294c16a6a6SHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
5304c16a6a6SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
5314c16a6a6SHong Zhang 
532057f5ba7SHong Zhang     while (i < k){
5334c16a6a6SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
5344c16a6a6SHong Zhang 
5354c16a6a6SHong Zhang       /* compute multiplier */
5364c16a6a6SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
5374c16a6a6SHong Zhang 
5384c16a6a6SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
5394c16a6a6SHong Zhang       diag = ba + i*bs2;
5404c16a6a6SHong Zhang       u    = ba + ili*bs2;
5414c16a6a6SHong Zhang       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
5424c16a6a6SHong Zhang       Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
5434c16a6a6SHong Zhang 
5444c16a6a6SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
5454c16a6a6SHong Zhang       Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
5464c16a6a6SHong Zhang 
5474c16a6a6SHong Zhang       /* update -U(i,k) */
5484c16a6a6SHong Zhang       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
5494c16a6a6SHong Zhang 
5504c16a6a6SHong Zhang       /* add multiple of row i to k-th row ... */
5514c16a6a6SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
5524c16a6a6SHong Zhang       if (jmin < jmax){
5534c16a6a6SHong Zhang         for (j=jmin; j<jmax; j++) {
5544c16a6a6SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
5554c16a6a6SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
5564c16a6a6SHong Zhang           u = ba + j*bs2;
5574c16a6a6SHong Zhang           Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
5584c16a6a6SHong Zhang         }
5594c16a6a6SHong Zhang 
5604c16a6a6SHong Zhang         /* ... add i to row list for next nonzero entry */
5614c16a6a6SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
5624c16a6a6SHong Zhang         j     = bj[jmin];
5634c16a6a6SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
5644c16a6a6SHong Zhang       }
5654c16a6a6SHong Zhang       i = nexti;
5664c16a6a6SHong Zhang     }
5674c16a6a6SHong Zhang 
5684c16a6a6SHong Zhang     /* save nonzero entries in k-th row of U ... */
5694c16a6a6SHong Zhang 
5704c16a6a6SHong Zhang     /* invert diagonal block */
5714c16a6a6SHong Zhang     diag = ba+k*bs2;
5724c16a6a6SHong Zhang     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
573d230e6fdSBarry Smith     ierr = Kernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr);
5744c16a6a6SHong Zhang 
5754c16a6a6SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
5764c16a6a6SHong Zhang     if (jmin < jmax) {
5774c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++){
5784c16a6a6SHong Zhang          vj = bj[j];           /* block col. index of U */
5794c16a6a6SHong Zhang          u   = ba + j*bs2;
5804c16a6a6SHong Zhang          rtmp_ptr = rtmp + vj*bs2;
5814c16a6a6SHong Zhang          for (k1=0; k1<bs2; k1++){
5824c16a6a6SHong Zhang            *u++        = *rtmp_ptr;
5834c16a6a6SHong Zhang            *rtmp_ptr++ = 0.0;
5844c16a6a6SHong Zhang          }
5854c16a6a6SHong Zhang       }
5864c16a6a6SHong Zhang 
5874c16a6a6SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
5884c16a6a6SHong Zhang       il[k] = jmin;
5894c16a6a6SHong Zhang       i     = bj[jmin];
5904c16a6a6SHong Zhang       jl[k] = jl[i]; jl[i] = k;
5914c16a6a6SHong Zhang     }
5924c16a6a6SHong Zhang   }
5934c16a6a6SHong Zhang 
5944c16a6a6SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
5954c16a6a6SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
5964c16a6a6SHong Zhang   ierr = PetscFree(dk);CHKERRQ(ierr);
5974c16a6a6SHong Zhang   ierr = PetscFree(pivots);CHKERRQ(ierr);
5984c16a6a6SHong Zhang   if (a->permute){
5994c16a6a6SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
6004c16a6a6SHong Zhang   }
6014c16a6a6SHong Zhang 
6024c16a6a6SHong Zhang   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
6034c16a6a6SHong Zhang   C->factor       = FACTOR_CHOLESKY;
6044c16a6a6SHong Zhang   C->assembled    = PETSC_TRUE;
6054c16a6a6SHong Zhang   C->preallocated = PETSC_TRUE;
606b0a32e0cSBarry Smith   PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */
6074c16a6a6SHong Zhang   PetscFunctionReturn(0);
6084c16a6a6SHong Zhang }
609d4edadd4SHong Zhang 
6104a2ae208SSatish Balay #undef __FUNCT__
6114a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering"
612dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat A,Mat *B)
613671cb588SHong Zhang {
614671cb588SHong Zhang   Mat                C = *B;
615671cb588SHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
616dfbe8321SBarry Smith   PetscErrorCode ierr;
617dfbe8321SBarry Smith   int i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
6181ad1de41SHong Zhang   int                *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
619671cb588SHong Zhang   int                bs=a->bs,bs2 = a->bs2;
620671cb588SHong Zhang   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
621671cb588SHong Zhang   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
62228de702eSHong Zhang   MatScalar          *work;
623671cb588SHong Zhang   int                *pivots;
624671cb588SHong Zhang 
625671cb588SHong Zhang   PetscFunctionBegin;
626671cb588SHong Zhang 
627671cb588SHong Zhang   /* initialization */
628671cb588SHong Zhang 
62982502324SSatish Balay   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
630671cb588SHong Zhang   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
63182502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
63228de702eSHong Zhang   jl   = il + mbs;
633671cb588SHong Zhang   for (i=0; i<mbs; i++) {
634671cb588SHong Zhang     jl[i] = mbs; il[0] = 0;
635671cb588SHong Zhang   }
636b0a32e0cSBarry Smith   ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr);
63728de702eSHong Zhang   uik  = dk + bs2;
63828de702eSHong Zhang   work = uik + bs2;
63982502324SSatish Balay   ierr = PetscMalloc(bs*sizeof(int),&pivots);CHKERRQ(ierr);
640671cb588SHong Zhang 
641671cb588SHong Zhang   ai = a->i; aj = a->j; aa = a->a;
642671cb588SHong Zhang 
643671cb588SHong Zhang   /* for each row k */
644671cb588SHong Zhang   for (k = 0; k<mbs; k++){
645671cb588SHong Zhang 
646671cb588SHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
647671cb588SHong Zhang     jmin = ai[k]; jmax = ai[k+1];
648671cb588SHong Zhang     ap = aa + jmin*bs2;
649671cb588SHong Zhang     for (j = jmin; j < jmax; j++){
650671cb588SHong Zhang       vj = aj[j];         /* block col. index */
651671cb588SHong Zhang       rtmp_ptr = rtmp + vj*bs2;
652671cb588SHong Zhang       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
653671cb588SHong Zhang     }
654671cb588SHong Zhang 
655671cb588SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
656671cb588SHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
657671cb588SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
658671cb588SHong Zhang 
659057f5ba7SHong Zhang     while (i < k){
660671cb588SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
661671cb588SHong Zhang 
662671cb588SHong Zhang       /* compute multiplier */
663671cb588SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
664671cb588SHong Zhang 
665671cb588SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
666671cb588SHong Zhang       diag = ba + i*bs2;
667671cb588SHong Zhang       u    = ba + ili*bs2;
668671cb588SHong Zhang       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
669671cb588SHong Zhang       Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
670671cb588SHong Zhang 
671671cb588SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
672671cb588SHong Zhang       Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
673671cb588SHong Zhang 
674671cb588SHong Zhang       /* update -U(i,k) */
675671cb588SHong Zhang       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
676671cb588SHong Zhang 
677671cb588SHong Zhang       /* add multiple of row i to k-th row ... */
678671cb588SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
679671cb588SHong Zhang       if (jmin < jmax){
680671cb588SHong Zhang         for (j=jmin; j<jmax; j++) {
681671cb588SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
682671cb588SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
683671cb588SHong Zhang           u = ba + j*bs2;
684671cb588SHong Zhang           Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
685671cb588SHong Zhang         }
686671cb588SHong Zhang 
687671cb588SHong Zhang         /* ... add i to row list for next nonzero entry */
688671cb588SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
689671cb588SHong Zhang         j     = bj[jmin];
690671cb588SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
691671cb588SHong Zhang       }
692671cb588SHong Zhang       i = nexti;
693671cb588SHong Zhang     }
694671cb588SHong Zhang 
695671cb588SHong Zhang     /* save nonzero entries in k-th row of U ... */
696671cb588SHong Zhang 
697671cb588SHong Zhang     /* invert diagonal block */
698671cb588SHong Zhang     diag = ba+k*bs2;
699671cb588SHong Zhang     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
700d230e6fdSBarry Smith     ierr = Kernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr);
701671cb588SHong Zhang 
702671cb588SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
703671cb588SHong Zhang     if (jmin < jmax) {
704671cb588SHong Zhang       for (j=jmin; j<jmax; j++){
705671cb588SHong Zhang          vj = bj[j];           /* block col. index of U */
706671cb588SHong Zhang          u   = ba + j*bs2;
707671cb588SHong Zhang          rtmp_ptr = rtmp + vj*bs2;
708671cb588SHong Zhang          for (k1=0; k1<bs2; k1++){
709671cb588SHong Zhang            *u++        = *rtmp_ptr;
710671cb588SHong Zhang            *rtmp_ptr++ = 0.0;
711671cb588SHong Zhang          }
712671cb588SHong Zhang       }
713671cb588SHong Zhang 
714671cb588SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
715671cb588SHong Zhang       il[k] = jmin;
716671cb588SHong Zhang       i     = bj[jmin];
717671cb588SHong Zhang       jl[k] = jl[i]; jl[i] = k;
718671cb588SHong Zhang     }
719671cb588SHong Zhang   }
720671cb588SHong Zhang 
721671cb588SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
722671cb588SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
723671cb588SHong Zhang   ierr = PetscFree(dk);CHKERRQ(ierr);
724671cb588SHong Zhang   ierr = PetscFree(pivots);CHKERRQ(ierr);
725671cb588SHong Zhang 
726671cb588SHong Zhang   C->factor    = FACTOR_CHOLESKY;
727671cb588SHong Zhang   C->assembled = PETSC_TRUE;
728671cb588SHong Zhang   C->preallocated = PETSC_TRUE;
729b0a32e0cSBarry Smith   PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */
730671cb588SHong Zhang   PetscFunctionReturn(0);
731671cb588SHong Zhang }
732671cb588SHong Zhang 
73349b5e25fSSatish Balay /*
734fcf159c0SHong Zhang     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
735cc0c071aSHong Zhang     Version for blocks 2 by 2.
73649b5e25fSSatish Balay */
7374a2ae208SSatish Balay #undef __FUNCT__
7384a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2"
739dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat A,Mat *B)
74049b5e25fSSatish Balay {
74149b5e25fSSatish Balay   Mat                C = *B;
74291602c66SHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
743cc0c071aSHong Zhang   IS                 perm = b->row;
7446849ba73SBarry Smith   PetscErrorCode ierr;
7456849ba73SBarry Smith   int                *perm_ptr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
746cc0c071aSHong Zhang   int                *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
747a1723e09SHong Zhang   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
748cc0c071aSHong Zhang   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
74949b5e25fSSatish Balay 
75049b5e25fSSatish Balay   PetscFunctionBegin;
751671cb588SHong Zhang 
75291602c66SHong Zhang   /* initialization */
75391602c66SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
75491602c66SHong Zhang      window U(0:k, k:mbs-1).
75591602c66SHong Zhang      jl:    list of rows to be added to uneliminated rows
75691602c66SHong Zhang             i>= k: jl(i) is the first row to be added to row i
75791602c66SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
75891602c66SHong Zhang             jl(i) = mbs indicates the end of a list
75991602c66SHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
76091602c66SHong Zhang             row i of U */
76182502324SSatish Balay   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
762cc0c071aSHong Zhang   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
76382502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
76428de702eSHong Zhang   jl   = il + mbs;
76591602c66SHong Zhang   for (i=0; i<mbs; i++) {
7663845f261SHong Zhang     jl[i] = mbs; il[0] = 0;
76791602c66SHong Zhang   }
76882502324SSatish Balay   ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr);
76928de702eSHong Zhang   uik  = dk + 4;
770cc0c071aSHong Zhang   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
771cc0c071aSHong Zhang 
772cc0c071aSHong Zhang   /* check permutation */
773cc0c071aSHong Zhang   if (!a->permute){
774cc0c071aSHong Zhang     ai = a->i; aj = a->j; aa = a->a;
775cc0c071aSHong Zhang   } else {
776cc0c071aSHong Zhang     ai   = a->inew; aj = a->jnew;
77782502324SSatish Balay     ierr = PetscMalloc(4*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
778cc0c071aSHong Zhang     ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
77982502324SSatish Balay     ierr = PetscMalloc(ai[mbs]*sizeof(int),&a2anew);CHKERRQ(ierr);
780cc0c071aSHong Zhang     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr);
781cc0c071aSHong Zhang 
782cc0c071aSHong Zhang     for (i=0; i<mbs; i++){
783cc0c071aSHong Zhang       jmin = ai[i]; jmax = ai[i+1];
784cc0c071aSHong Zhang       for (j=jmin; j<jmax; j++){
785cc0c071aSHong Zhang         while (a2anew[j] != j){
786cc0c071aSHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
787cc0c071aSHong Zhang           for (k1=0; k1<4; k1++){
788cc0c071aSHong Zhang             dk[k1]       = aa[k*4+k1];
789cc0c071aSHong Zhang             aa[k*4+k1] = aa[j*4+k1];
790cc0c071aSHong Zhang             aa[j*4+k1] = dk[k1];
791cc0c071aSHong Zhang           }
792cc0c071aSHong Zhang         }
793cc0c071aSHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
794cc0c071aSHong Zhang         if (i > aj[j]){
795a1723e09SHong Zhang           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
796cc0c071aSHong Zhang           ap = aa + j*4;     /* ptr to the beginning of the block */
797cc0c071aSHong Zhang           dk[1] = ap[1];     /* swap ap[1] and ap[2] */
798cc0c071aSHong Zhang           ap[1] = ap[2];
799cc0c071aSHong Zhang           ap[2] = dk[1];
800cc0c071aSHong Zhang         }
801cc0c071aSHong Zhang       }
802cc0c071aSHong Zhang     }
803ac355199SBarry Smith     ierr = PetscFree(a2anew);CHKERRQ(ierr);
804cc0c071aSHong Zhang   }
8053845f261SHong Zhang 
80691602c66SHong Zhang   /* for each row k */
80791602c66SHong Zhang   for (k = 0; k<mbs; k++){
80891602c66SHong Zhang 
80991602c66SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
810cc0c071aSHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
811cc0c071aSHong Zhang     ap = aa + jmin*4;
81291602c66SHong Zhang     for (j = jmin; j < jmax; j++){
813cc0c071aSHong Zhang       vj = perm_ptr[aj[j]];         /* block col. index */
814cc0c071aSHong Zhang       rtmp_ptr = rtmp + vj*4;
815cc0c071aSHong Zhang       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
81691602c66SHong Zhang     }
81791602c66SHong Zhang 
81891602c66SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
819cc0c071aSHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
82091602c66SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
82191602c66SHong Zhang 
822057f5ba7SHong Zhang     while (i < k){
82391602c66SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
82491602c66SHong Zhang 
8253845f261SHong Zhang       /* compute multiplier */
82691602c66SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
8273845f261SHong Zhang 
8283845f261SHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
829cc0c071aSHong Zhang       diag = ba + i*4;
830cc0c071aSHong Zhang       u    = ba + ili*4;
831cc0c071aSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
832cc0c071aSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
833cc0c071aSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
834cc0c071aSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
8353845f261SHong Zhang 
8363845f261SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
837cc0c071aSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
838cc0c071aSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
839cc0c071aSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
840cc0c071aSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
8413845f261SHong Zhang 
8423845f261SHong Zhang       /* update -U(i,k): ba[ili] = uik */
843cc0c071aSHong Zhang       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
84491602c66SHong Zhang 
84591602c66SHong Zhang       /* add multiple of row i to k-th row ... */
84691602c66SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
84791602c66SHong Zhang       if (jmin < jmax){
8483845f261SHong Zhang         for (j=jmin; j<jmax; j++) {
8493845f261SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
850cc0c071aSHong Zhang           rtmp_ptr = rtmp + bj[j]*4;
851cc0c071aSHong Zhang           u = ba + j*4;
852cc0c071aSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
853cc0c071aSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
854cc0c071aSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
855cc0c071aSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
8563845f261SHong Zhang         }
8573845f261SHong Zhang 
85891602c66SHong Zhang         /* ... add i to row list for next nonzero entry */
85991602c66SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
86091602c66SHong Zhang         j     = bj[jmin];
86191602c66SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
86291602c66SHong Zhang       }
863a1723e09SHong Zhang       i = nexti;
86491602c66SHong Zhang     }
865cc0c071aSHong Zhang 
86691602c66SHong Zhang     /* save nonzero entries in k-th row of U ... */
8673845f261SHong Zhang 
868cc0c071aSHong Zhang     /* invert diagonal block */
869cc0c071aSHong Zhang     diag = ba+k*4;
870cc0c071aSHong Zhang     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
8713845f261SHong Zhang     ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr);
8723845f261SHong Zhang 
87391602c66SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
87491602c66SHong Zhang     if (jmin < jmax) {
87591602c66SHong Zhang       for (j=jmin; j<jmax; j++){
876cc0c071aSHong Zhang          vj = bj[j];           /* block col. index of U */
877cc0c071aSHong Zhang          u   = ba + j*4;
878cc0c071aSHong Zhang          rtmp_ptr = rtmp + vj*4;
879cc0c071aSHong Zhang          for (k1=0; k1<4; k1++){
880cc0c071aSHong Zhang            *u++        = *rtmp_ptr;
881cc0c071aSHong Zhang            *rtmp_ptr++ = 0.0;
882cc0c071aSHong Zhang          }
883cc0c071aSHong Zhang       }
8843845f261SHong Zhang 
88591602c66SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
88691602c66SHong Zhang       il[k] = jmin;
88791602c66SHong Zhang       i     = bj[jmin];
88891602c66SHong Zhang       jl[k] = jl[i]; jl[i] = k;
88991602c66SHong Zhang     }
89091602c66SHong Zhang   }
8913845f261SHong Zhang 
89249b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
89391602c66SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
8943845f261SHong Zhang   ierr = PetscFree(dk);CHKERRQ(ierr);
89591602c66SHong Zhang   if (a->permute) {
89691602c66SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
89791602c66SHong Zhang   }
898cc0c071aSHong Zhang   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
89991602c66SHong Zhang   C->factor    = FACTOR_CHOLESKY;
90049b5e25fSSatish Balay   C->assembled = PETSC_TRUE;
9015c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
902b0a32e0cSBarry Smith   PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */
90349b5e25fSSatish Balay   PetscFunctionReturn(0);
90449b5e25fSSatish Balay }
90591602c66SHong Zhang 
90649b5e25fSSatish Balay /*
90749b5e25fSSatish Balay       Version for when blocks are 2 by 2 Using natural ordering
90849b5e25fSSatish Balay */
9094a2ae208SSatish Balay #undef __FUNCT__
9104a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering"
911dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat A,Mat *B)
91249b5e25fSSatish Balay {
91349b5e25fSSatish Balay   Mat                C = *B;
914ab27746eSHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
915dfbe8321SBarry Smith   PetscErrorCode ierr;
916dfbe8321SBarry Smith   int i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
91771149678SHong Zhang   int                *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
918ab27746eSHong Zhang   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
919ab27746eSHong Zhang   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
92049b5e25fSSatish Balay 
92149b5e25fSSatish Balay   PetscFunctionBegin;
922671cb588SHong Zhang 
923ab27746eSHong Zhang   /* initialization */
924ab27746eSHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
925ab27746eSHong Zhang      window U(0:k, k:mbs-1).
926ab27746eSHong Zhang      jl:    list of rows to be added to uneliminated rows
927ab27746eSHong Zhang             i>= k: jl(i) is the first row to be added to row i
928ab27746eSHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
929ab27746eSHong Zhang             jl(i) = mbs indicates the end of a list
930ab27746eSHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
931ab27746eSHong Zhang             row i of U */
93282502324SSatish Balay   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
933ab27746eSHong Zhang   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
93482502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
93528de702eSHong Zhang   jl   = il + mbs;
936ab27746eSHong Zhang   for (i=0; i<mbs; i++) {
937ab27746eSHong Zhang     jl[i] = mbs; il[0] = 0;
93849b5e25fSSatish Balay   }
93982502324SSatish Balay   ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr);
94028de702eSHong Zhang   uik  = dk + 4;
941ab27746eSHong Zhang 
942ab27746eSHong Zhang   ai = a->i; aj = a->j; aa = a->a;
943ab27746eSHong Zhang 
944ab27746eSHong Zhang   /* for each row k */
945ab27746eSHong Zhang   for (k = 0; k<mbs; k++){
946ab27746eSHong Zhang 
947ab27746eSHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
948ab27746eSHong Zhang     jmin = ai[k]; jmax = ai[k+1];
949ab27746eSHong Zhang     ap = aa + jmin*4;
950ab27746eSHong Zhang     for (j = jmin; j < jmax; j++){
951ab27746eSHong Zhang       vj = aj[j];         /* block col. index */
952ab27746eSHong Zhang       rtmp_ptr = rtmp + vj*4;
953ab27746eSHong Zhang       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
95449b5e25fSSatish Balay     }
955ab27746eSHong Zhang 
956ab27746eSHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
957ab27746eSHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
958ab27746eSHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
959ab27746eSHong Zhang 
960057f5ba7SHong Zhang     while (i < k){
961ab27746eSHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
962ab27746eSHong Zhang 
963ab27746eSHong Zhang       /* compute multiplier */
964ab27746eSHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
965ab27746eSHong Zhang 
966ab27746eSHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
967ab27746eSHong Zhang       diag = ba + i*4;
968ab27746eSHong Zhang       u    = ba + ili*4;
969ab27746eSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
970ab27746eSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
971ab27746eSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
972ab27746eSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
973ab27746eSHong Zhang 
974ab27746eSHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
975ab27746eSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
976ab27746eSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
977ab27746eSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
978ab27746eSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
979ab27746eSHong Zhang 
980ab27746eSHong Zhang       /* update -U(i,k): ba[ili] = uik */
981ab27746eSHong Zhang       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
982ab27746eSHong Zhang 
983ab27746eSHong Zhang       /* add multiple of row i to k-th row ... */
984ab27746eSHong Zhang       jmin = ili + 1; jmax = bi[i+1];
985ab27746eSHong Zhang       if (jmin < jmax){
986ab27746eSHong Zhang         for (j=jmin; j<jmax; j++) {
987ab27746eSHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
988ab27746eSHong Zhang           rtmp_ptr = rtmp + bj[j]*4;
989ab27746eSHong Zhang           u = ba + j*4;
990ab27746eSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
991ab27746eSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
992ab27746eSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
993ab27746eSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
99449b5e25fSSatish Balay         }
995ab27746eSHong Zhang 
996ab27746eSHong Zhang         /* ... add i to row list for next nonzero entry */
997ab27746eSHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
998ab27746eSHong Zhang         j     = bj[jmin];
999ab27746eSHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
100049b5e25fSSatish Balay       }
1001ab27746eSHong Zhang       i = nexti;
100249b5e25fSSatish Balay     }
1003ab27746eSHong Zhang 
1004ab27746eSHong Zhang     /* save nonzero entries in k-th row of U ... */
1005ab27746eSHong Zhang 
100649b5e25fSSatish Balay     /* invert diagonal block */
1007ab27746eSHong Zhang     diag = ba+k*4;
1008ab27746eSHong Zhang     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
1009ab27746eSHong Zhang     ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr);
1010ab27746eSHong Zhang 
1011ab27746eSHong Zhang     jmin = bi[k]; jmax = bi[k+1];
1012ab27746eSHong Zhang     if (jmin < jmax) {
1013ab27746eSHong Zhang       for (j=jmin; j<jmax; j++){
1014ab27746eSHong Zhang          vj = bj[j];           /* block col. index of U */
1015ab27746eSHong Zhang          u   = ba + j*4;
1016ab27746eSHong Zhang          rtmp_ptr = rtmp + vj*4;
1017ab27746eSHong Zhang          for (k1=0; k1<4; k1++){
1018ab27746eSHong Zhang            *u++        = *rtmp_ptr;
1019ab27746eSHong Zhang            *rtmp_ptr++ = 0.0;
1020ab27746eSHong Zhang          }
1021ab27746eSHong Zhang       }
1022ab27746eSHong Zhang 
1023ab27746eSHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
1024ab27746eSHong Zhang       il[k] = jmin;
1025ab27746eSHong Zhang       i     = bj[jmin];
1026ab27746eSHong Zhang       jl[k] = jl[i]; jl[i] = k;
1027ab27746eSHong Zhang     }
102849b5e25fSSatish Balay   }
102949b5e25fSSatish Balay 
103049b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1031ab27746eSHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
1032ab27746eSHong Zhang   ierr = PetscFree(dk);CHKERRQ(ierr);
1033ab27746eSHong Zhang 
1034ab27746eSHong Zhang   C->factor    = FACTOR_CHOLESKY;
103549b5e25fSSatish Balay   C->assembled = PETSC_TRUE;
10365c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
1037b0a32e0cSBarry Smith   PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */
103849b5e25fSSatish Balay   PetscFunctionReturn(0);
103949b5e25fSSatish Balay }
104049b5e25fSSatish Balay 
104149b5e25fSSatish Balay /*
10425c0bcdfcSHong Zhang     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
104391602c66SHong Zhang     Version for blocks are 1 by 1.
104449b5e25fSSatish Balay */
10454a2ae208SSatish Balay #undef __FUNCT__
10464a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1"
1047dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat A,Mat *B)
104849b5e25fSSatish Balay {
104949b5e25fSSatish Balay   Mat                C = *B;
105049b5e25fSSatish Balay   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
105149b5e25fSSatish Balay   IS                 ip = b->row;
10526849ba73SBarry Smith   PetscErrorCode ierr;
10536849ba73SBarry Smith   int                *rip,i,j,mbs = a->mbs,*bi = b->i,*bj = b->j;
1054cb718733SHong Zhang   int                *ai,*aj,*r;
1055ab27746eSHong Zhang   int                k,jmin,jmax,*jl,*il,vj,nexti,ili;
1056066653e3SSatish Balay   MatScalar          *rtmp;
10572d007305SHong Zhang   MatScalar          *ba = b->a,*aa,ak;
105849b5e25fSSatish Balay   MatScalar          dk,uikdi;
105949b5e25fSSatish Balay 
106049b5e25fSSatish Balay   PetscFunctionBegin;
106149b5e25fSSatish Balay   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
1062cb718733SHong Zhang   if (!a->permute){
10632d007305SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
10642d007305SHong Zhang   } else {
10652d007305SHong Zhang     ai = a->inew; aj = a->jnew;
106682502324SSatish Balay     ierr = PetscMalloc(ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
1067cb718733SHong Zhang     ierr = PetscMemcpy(aa,a->a,ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
106882502324SSatish Balay     ierr = PetscMalloc(ai[mbs]*sizeof(int),&r);CHKERRQ(ierr);
10692d007305SHong Zhang     ierr= PetscMemcpy(r,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr);
10702d007305SHong Zhang 
10712d007305SHong Zhang     jmin = ai[0]; jmax = ai[mbs];
10722d007305SHong Zhang     for (j=jmin; j<jmax; j++){
1073cb718733SHong Zhang       while (r[j] != j){
10742d007305SHong Zhang         k = r[j]; r[j] = r[k]; r[k] = k;
10752d007305SHong Zhang         ak = aa[k]; aa[k] = aa[j]; aa[j] = ak;
10762d007305SHong Zhang       }
10772d007305SHong Zhang     }
1078ac355199SBarry Smith     ierr = PetscFree(r);CHKERRQ(ierr);
10792d007305SHong Zhang   }
10802d007305SHong Zhang 
108191602c66SHong Zhang   /* initialization */
108249b5e25fSSatish Balay   /* il and jl record the first nonzero element in each row of the accessing
108349b5e25fSSatish Balay      window U(0:k, k:mbs-1).
108449b5e25fSSatish Balay      jl:    list of rows to be added to uneliminated rows
108549b5e25fSSatish Balay             i>= k: jl(i) is the first row to be added to row i
108649b5e25fSSatish Balay             i<  k: jl(i) is the row following row i in some list of rows
108749b5e25fSSatish Balay             jl(i) = mbs indicates the end of a list
108849b5e25fSSatish Balay      il(i): points to the first nonzero element in columns k,...,mbs-1 of
108949b5e25fSSatish Balay             row i of U */
109082502324SSatish Balay   ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
109182502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
109228de702eSHong Zhang   jl   = il + mbs;
109349b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
109449b5e25fSSatish Balay     rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
109549b5e25fSSatish Balay   }
109649b5e25fSSatish Balay 
109791602c66SHong Zhang   /* for each row k */
109849b5e25fSSatish Balay   for (k = 0; k<mbs; k++){
109949b5e25fSSatish Balay 
110091602c66SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
110149b5e25fSSatish Balay     jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
1102057f5ba7SHong Zhang 
110349b5e25fSSatish Balay     for (j = jmin; j < jmax; j++){
1104351d0355SHong Zhang       vj = rip[aj[j]];
1105ab27746eSHong Zhang       rtmp[vj] = aa[j];
110649b5e25fSSatish Balay     }
110749b5e25fSSatish Balay 
110891602c66SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
110949b5e25fSSatish Balay     dk = rtmp[k];
111049b5e25fSSatish Balay     i = jl[k]; /* first row to be added to k_th row  */
111149b5e25fSSatish Balay 
1112057f5ba7SHong Zhang     while (i < k){
111349b5e25fSSatish Balay       nexti = jl[i]; /* next row to be added to k_th row */
111449b5e25fSSatish Balay 
111591602c66SHong Zhang       /* compute multiplier, update D(k) and U(i,k) */
111649b5e25fSSatish Balay       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
111749b5e25fSSatish Balay       uikdi = - ba[ili]*ba[i];
111849b5e25fSSatish Balay       dk += uikdi*ba[ili];
1119658e7b3eSHong Zhang       ba[ili] = uikdi; /* -U(i,k) */
112049b5e25fSSatish Balay 
112191602c66SHong Zhang       /* add multiple of row i to k-th row ... */
112249b5e25fSSatish Balay       jmin = ili + 1; jmax = bi[i+1];
112349b5e25fSSatish Balay       if (jmin < jmax){
112449b5e25fSSatish Balay         for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
112591602c66SHong Zhang         /* ... add i to row list for next nonzero entry */
112649b5e25fSSatish Balay         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
112749b5e25fSSatish Balay         j     = bj[jmin];
112849b5e25fSSatish Balay         jl[i] = jl[j]; jl[j] = i; /* update jl */
112949b5e25fSSatish Balay       }
1130ab27746eSHong Zhang       i = nexti;
113149b5e25fSSatish Balay     }
113249b5e25fSSatish Balay 
113391602c66SHong Zhang     /* check for zero pivot and save diagoanl element */
11345b8514ebSBarry Smith     if (dk == 0.0){
113529bbc08cSBarry Smith       SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot");
1136671cb588SHong Zhang       /*
11371223c01bSHong Zhang     } else if (PetscRealPart(dk) < 0.0){
11381223c01bSHong Zhang       SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Negative pivot: d[%d] = %g\n",k,dk);
1139671cb588SHong Zhang       */
114049b5e25fSSatish Balay     }
114149b5e25fSSatish Balay 
114291602c66SHong Zhang     /* save nonzero entries in k-th row of U ... */
1143f3dd7b2dSKris Buschelman     ba[k] = 1.0/dk;
114449b5e25fSSatish Balay     jmin = bi[k]; jmax = bi[k+1];
114549b5e25fSSatish Balay     if (jmin < jmax) {
114649b5e25fSSatish Balay       for (j=jmin; j<jmax; j++){
1147cc0c071aSHong Zhang          vj = bj[j]; ba[j] = rtmp[vj]; rtmp[vj] = 0.0;
114849b5e25fSSatish Balay       }
114991602c66SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
115049b5e25fSSatish Balay       il[k] = jmin;
115149b5e25fSSatish Balay       i     = bj[jmin];
115249b5e25fSSatish Balay       jl[k] = jl[i]; jl[i] = k;
115349b5e25fSSatish Balay     }
115449b5e25fSSatish Balay   }
115549b5e25fSSatish Balay 
115649b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
115749b5e25fSSatish Balay   ierr = PetscFree(il);CHKERRQ(ierr);
1158cb718733SHong Zhang   if (a->permute){
1159cb718733SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
1160cb718733SHong Zhang   }
116149b5e25fSSatish Balay 
116249b5e25fSSatish Balay   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
11638be8c0c7SHong Zhang   C->factor    = FACTOR_CHOLESKY;
116449b5e25fSSatish Balay   C->assembled = PETSC_TRUE;
11655c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
1166b0a32e0cSBarry Smith   PetscLogFlops(b->mbs);
116749b5e25fSSatish Balay   PetscFunctionReturn(0);
116849b5e25fSSatish Balay }
116949b5e25fSSatish Balay 
1170671cb588SHong Zhang /*
1171671cb588SHong Zhang   Version for when blocks are 1 by 1 Using natural ordering
1172671cb588SHong Zhang */
11734a2ae208SSatish Balay #undef __FUNCT__
11744a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering"
1175dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat A,Mat *B)
1176671cb588SHong Zhang {
1177671cb588SHong Zhang   Mat                C = *B;
1178671cb588SHong Zhang   Mat_SeqSBAIJ       *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data;
1179dfbe8321SBarry Smith   PetscErrorCode ierr;
1180dfbe8321SBarry Smith   int i,j,mbs = a->mbs;
1181653a6975SHong Zhang   int                *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
1182b00f7748SHong Zhang   int                k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz,ndamp = 0;
1183653a6975SHong Zhang   MatScalar          *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
118421c26570Svictorle   PetscReal          damping=b->factor_damping, zeropivot=b->factor_zeropivot,shift_amount;
118521c26570Svictorle   PetscTruth         damp,chshift;
118621c26570Svictorle   int                nshift=0;
1187653a6975SHong Zhang 
1188653a6975SHong Zhang   PetscFunctionBegin;
1189653a6975SHong Zhang   /* initialization */
1190653a6975SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
1191653a6975SHong Zhang      window U(0:k, k:mbs-1).
1192653a6975SHong Zhang      jl:    list of rows to be added to uneliminated rows
1193653a6975SHong Zhang             i>= k: jl(i) is the first row to be added to row i
1194653a6975SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
1195653a6975SHong Zhang             jl(i) = mbs indicates the end of a list
1196653a6975SHong Zhang      il(i): points to the first nonzero element in U(i,k:mbs-1)
1197653a6975SHong Zhang   */
1198653a6975SHong Zhang   ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
1199653a6975SHong Zhang   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
1200653a6975SHong Zhang   jl   = il + mbs;
1201b00f7748SHong Zhang 
120221c26570Svictorle   shift_amount = 0;
1203b00f7748SHong Zhang   do {
1204b00f7748SHong Zhang     damp = PETSC_FALSE;
120521c26570Svictorle     chshift = PETSC_FALSE;
1206653a6975SHong Zhang     for (i=0; i<mbs; i++) {
1207653a6975SHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1208653a6975SHong Zhang     }
1209653a6975SHong Zhang 
1210b00f7748SHong Zhang     for (k = 0; k<mbs; k++){ /* row k */
1211653a6975SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
1212653a6975SHong Zhang       nz   = ai[k+1] - ai[k];
1213653a6975SHong Zhang       acol = aj + ai[k];
1214653a6975SHong Zhang       aval = aa + ai[k];
1215653a6975SHong Zhang       bval = ba + bi[k];
1216653a6975SHong Zhang       while (nz -- ){
1217653a6975SHong Zhang         rtmp[*acol++] = *aval++;
1218653a6975SHong Zhang         *bval++       = 0.0; /* for in-place factorization */
1219653a6975SHong Zhang       }
1220b00f7748SHong Zhang       /* damp the diagonal of the matrix */
122121c26570Svictorle       if (ndamp||nshift) rtmp[k] += damping+shift_amount;
1222653a6975SHong Zhang 
1223653a6975SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k) != 0 */
1224653a6975SHong Zhang       dk = rtmp[k];
1225653a6975SHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
1226653a6975SHong Zhang 
1227653a6975SHong Zhang       while (i < k){
1228653a6975SHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
1229653a6975SHong Zhang 
1230653a6975SHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
1231653a6975SHong Zhang         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1232653a6975SHong Zhang         uikdi = - ba[ili]*ba[bi[i]];
1233653a6975SHong Zhang         dk   += uikdi*ba[ili];
1234653a6975SHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
1235653a6975SHong Zhang 
1236653a6975SHong Zhang         /* add multiple of row i to k-th row ... */
1237653a6975SHong Zhang         jmin = ili + 1;
1238653a6975SHong Zhang         nz   = bi[i+1] - jmin;
1239653a6975SHong Zhang         if (nz > 0){
1240653a6975SHong Zhang           bcol = bj + jmin;
1241653a6975SHong Zhang           bval = ba + jmin;
1242653a6975SHong Zhang           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
1243653a6975SHong Zhang           /* ... add i to row list for next nonzero entry */
1244653a6975SHong Zhang           il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
1245653a6975SHong Zhang           j     = bj[jmin];
1246653a6975SHong Zhang           jl[i] = jl[j]; jl[j] = i; /* update jl */
1247653a6975SHong Zhang         }
1248653a6975SHong Zhang         i = nexti;
1249653a6975SHong Zhang       }
1250653a6975SHong Zhang 
125121c26570Svictorle       if (PetscRealPart(dk) < zeropivot && b->factor_shift){
1252d530a6e7Svictorle 	/* calculate a shift that would make this row diagonally dominant */
125329f69fd4SBarry Smith 	PetscReal rs = PetscAbs(PetscRealPart(dk));
125421c26570Svictorle 	jmin      = bi[k]+1;
125521c26570Svictorle 	nz        = bi[k+1] - jmin;
125621c26570Svictorle 	if (nz){
125721c26570Svictorle 	  bcol = bj + jmin;
125821c26570Svictorle 	  bval = ba + jmin;
125921c26570Svictorle 	  while (nz--){
126021c26570Svictorle 	    rs += PetscAbsScalar(rtmp[*bcol++]);
126121c26570Svictorle 	  }
126221c26570Svictorle 	}
1263d530a6e7Svictorle 	/* if this shift is less than the previous, just up the previous
1264d530a6e7Svictorle 	   one by a bit */
1265d530a6e7Svictorle 	shift_amount = PetscMax(rs,1.1*shift_amount);
126621c26570Svictorle 	chshift  = PETSC_TRUE;
1267d530a6e7Svictorle 	/* Unlike in the ILU case there is no exit condition on nshift:
1268d530a6e7Svictorle 	   we increase the shift until it converges. There is no guarantee that
1269d530a6e7Svictorle 	   this algorithm converges faster or slower, or is better or worse
1270d530a6e7Svictorle 	   than the ILU algorithm. */
127121c26570Svictorle 	nshift++;
127221c26570Svictorle 	break;
127321c26570Svictorle       }
1274ffa0b812SHong Zhang       if (PetscRealPart(dk) < zeropivot){
1275ffa0b812SHong Zhang         if (damping == (PetscReal) PETSC_DECIDE) damping = -PetscRealPart(dk)/(k+1);
1276b00f7748SHong Zhang         if (damping > 0.0) {
1277b00f7748SHong Zhang           if (ndamp) damping *= 2.0;
1278b00f7748SHong Zhang           damp = PETSC_TRUE;
1279b00f7748SHong Zhang           ndamp++;
1280b00f7748SHong Zhang           break;
1281481c2c92SHong Zhang         } else if (PetscAbsScalar(dk) < zeropivot){
1282*77431f27SBarry Smith           SETERRQ3(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g",k,PetscRealPart(dk),zeropivot);
1283481c2c92SHong Zhang         } else {
1284*77431f27SBarry Smith           PetscLogInfo((PetscObject)A,"Negative pivot %g in row %D of Cholesky factorization\n",PetscRealPart(dk),k);
1285b00f7748SHong Zhang         }
1286064503c5SHong Zhang       }
1287653a6975SHong Zhang 
1288653a6975SHong Zhang       /* save nonzero entries in k-th row of U ... */
1289*77431f27SBarry Smith       /* printf("%D, dk: %g, 1/dk: %g\n",k,dk,1/dk); */
1290653a6975SHong Zhang       ba[bi[k]] = 1.0/dk;
1291653a6975SHong Zhang       jmin      = bi[k]+1;
1292653a6975SHong Zhang       nz        = bi[k+1] - jmin;
1293653a6975SHong Zhang       if (nz){
1294653a6975SHong Zhang         bcol = bj + jmin;
1295653a6975SHong Zhang         bval = ba + jmin;
1296653a6975SHong Zhang         while (nz--){
1297653a6975SHong Zhang           *bval++       = rtmp[*bcol];
1298653a6975SHong Zhang           rtmp[*bcol++] = 0.0;
1299653a6975SHong Zhang         }
1300653a6975SHong Zhang         /* ... add k to row list for first nonzero entry in k-th row */
1301653a6975SHong Zhang         il[k] = jmin;
1302653a6975SHong Zhang         i     = bj[jmin];
1303653a6975SHong Zhang         jl[k] = jl[i]; jl[i] = k;
1304653a6975SHong Zhang       }
1305b00f7748SHong Zhang     } /* end of for (k = 0; k<mbs; k++) */
130621c26570Svictorle   } while (damp||chshift);
1307653a6975SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1308653a6975SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
1309653a6975SHong Zhang 
1310653a6975SHong Zhang   C->factor       = FACTOR_CHOLESKY;
1311653a6975SHong Zhang   C->assembled    = PETSC_TRUE;
1312653a6975SHong Zhang   C->preallocated = PETSC_TRUE;
1313653a6975SHong Zhang   PetscLogFlops(b->mbs);
1314b00f7748SHong Zhang   if (ndamp) {
1315*77431f27SBarry Smith     PetscLogInfo(0,"MatCholeskyFactorNumerical_SeqSBAIJ_1_NaturalOrdering: number of damping tries %D damping value %g\n",ndamp,damping);
1316b00f7748SHong Zhang   }
1317fb10cecfSBarry Smith  if (nshift) {
1318*77431f27SBarry Smith     PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering diagonal shifted %D shifts\n",nshift);
1319fb10cecfSBarry Smith   }
1320fb10cecfSBarry Smith 
1321653a6975SHong Zhang   PetscFunctionReturn(0);
1322653a6975SHong Zhang }
1323653a6975SHong Zhang 
1324653a6975SHong Zhang #undef __FUNCT__
13254a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ"
1326dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,MatFactorInfo *info)
132749b5e25fSSatish Balay {
1328dfbe8321SBarry Smith   PetscErrorCode ierr;
132949b5e25fSSatish Balay   Mat C;
133049b5e25fSSatish Balay 
133149b5e25fSSatish Balay   PetscFunctionBegin;
133215e8a5b3SHong Zhang   ierr = MatCholeskyFactorSymbolic(A,perm,info,&C);CHKERRQ(ierr);
1333a4ada70bSHong Zhang   ierr = MatCholeskyFactorNumeric(A,&C);CHKERRQ(ierr);
13344445ddedSHong Zhang   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
133549b5e25fSSatish Balay   PetscFunctionReturn(0);
133649b5e25fSSatish Balay }
133749b5e25fSSatish Balay 
133849b5e25fSSatish Balay 
1339