xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact.c (revision 3a7fca6b988f1faea172e46abdf950477dbfaf76)
1*3a7fca6bSBarry Smith /*$Id: sbaijfact.c,v 1.59 2001/04/10 22:35:39 balay Exp bsmith $*/
249b5e25fSSatish Balay /* Using Modified Sparse Row (MSR) storage.
349b5e25fSSatish Balay See page 85, "Iterative Methods ..." by Saad. */
449b5e25fSSatish Balay 
549b5e25fSSatish Balay /*
6effa298cSHong Zhang     Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP.
749b5e25fSSatish Balay */
849b5e25fSSatish Balay #include "src/mat/impls/baij/seq/baij.h"
9*3a7fca6bSBarry Smith #include "src/mat/impls/sbaij/seq/sbaij.h"
1049b5e25fSSatish Balay #include "src/vec/vecimpl.h"
1149b5e25fSSatish Balay #include "src/inline/ilu.h"
1249a6740bSHong Zhang #include "include/petscis.h"
1349b5e25fSSatish Balay 
144a2ae208SSatish Balay #undef __FUNCT__
154a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ"
1677f73120SHong Zhang int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat A,IS perm,PetscReal f,Mat *B)
1749b5e25fSSatish Balay {
1849b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b;
19cb718733SHong Zhang   int         *rip,ierr,i,mbs = a->mbs,*ai,*aj;
20066653e3SSatish Balay   int         *jutmp,bs = a->bs,bs2=a->bs2;
21066653e3SSatish Balay   int         m,nzi,realloc = 0;
22066653e3SSatish Balay   int         *jl,*q,jumin,jmin,jmax,juptr,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd;
23671cb588SHong Zhang   PetscTruth  perm_identity;
2449b5e25fSSatish Balay 
2549b5e25fSSatish Balay   PetscFunctionBegin;
262d007305SHong Zhang 
27cb718733SHong Zhang   /* check whether perm is the identity mapping */
28671cb588SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
29671cb588SHong Zhang   if (!perm_identity) a->permute = PETSC_TRUE;
302d007305SHong Zhang 
31671cb588SHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
32671cb588SHong Zhang 
33671cb588SHong Zhang   if (perm_identity){ /* without permutation */
342d007305SHong Zhang     ai = a->i; aj = a->j;
352d007305SHong Zhang   } else {            /* non-trivial permutation */
36323b833fSBarry Smith     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
372d007305SHong Zhang     ai = a->inew; aj = a->jnew;
382d007305SHong Zhang   }
392d007305SHong Zhang 
4049b5e25fSSatish Balay   /* initialization */
4149b5e25fSSatish Balay   /* Don't know how many column pointers are needed so estimate.
4249b5e25fSSatish Balay      Use Modified Sparse Row storage for u and ju, see Sasd pp.85 */
43b0a32e0cSBarry Smith   ierr  = PetscMalloc((mbs+1)*sizeof(int),&iu);CHKERRQ(ierr);
4449b5e25fSSatish Balay   umax  = (int)(f*ai[mbs] + 1); umax += mbs + 1;
4582502324SSatish Balay   ierr  = PetscMalloc(umax*sizeof(int),&ju);CHKERRQ(ierr);
4649b5e25fSSatish Balay   iu[0] = mbs+1;
4749b5e25fSSatish Balay   juptr = mbs;
4882502324SSatish Balay   ierr  = PetscMalloc(mbs*sizeof(int),&jl);CHKERRQ(ierr);
4982502324SSatish Balay   ierr  = PetscMalloc(mbs*sizeof(int),&q);CHKERRQ(ierr);
5049b5e25fSSatish Balay   for (i=0; i<mbs; i++){
5149b5e25fSSatish Balay     jl[i] = mbs; q[i] = 0;
5249b5e25fSSatish Balay   }
5349b5e25fSSatish Balay 
5449b5e25fSSatish Balay   /* for each row k */
5549b5e25fSSatish Balay   for (k=0; k<mbs; k++){
5649b5e25fSSatish Balay     nzk = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */
5749b5e25fSSatish Balay     q[k] = mbs;
5849b5e25fSSatish Balay     /* initialize nonzero structure of k-th row to row rip[k] of A */
5949b5e25fSSatish Balay     jmin = ai[rip[k]];
6049b5e25fSSatish Balay     jmax = ai[rip[k]+1];
6149b5e25fSSatish Balay     for (j=jmin; j<jmax; j++){
62cb718733SHong Zhang       vj = rip[aj[j]]; /* col. value */
6349b5e25fSSatish Balay       if(vj > k){
6449b5e25fSSatish Balay         qm = k;
6549b5e25fSSatish Balay         do {
6649b5e25fSSatish Balay           m  = qm; qm = q[m];
6749b5e25fSSatish Balay         } while(qm < vj);
6849b5e25fSSatish Balay         if (qm == vj) {
6949b5e25fSSatish Balay           printf(" error: duplicate entry in A\n"); break;
7049b5e25fSSatish Balay         }
7149b5e25fSSatish Balay         nzk++;
7249b5e25fSSatish Balay         q[m] = vj;
7349b5e25fSSatish Balay         q[vj] = qm;
7449b5e25fSSatish Balay       } /* if(vj > k) */
7549b5e25fSSatish Balay     } /* for (j=jmin; j<jmax; j++) */
7649b5e25fSSatish Balay 
7749b5e25fSSatish Balay     /* modify nonzero structure of k-th row by computing fill-in
7849b5e25fSSatish Balay        for each row i to be merged in */
7949b5e25fSSatish Balay     i = k;
8049b5e25fSSatish Balay     i = jl[i]; /* next pivot row (== mbs for symbolic factorization) */
8149b5e25fSSatish Balay     /* printf(" next pivot row i=%d\n",i); */
8249b5e25fSSatish Balay     while (i < mbs){
8349b5e25fSSatish Balay       /* merge row i into k-th row */
8449b5e25fSSatish Balay       nzi = iu[i+1] - (iu[i]+1);
8549b5e25fSSatish Balay       jmin = iu[i] + 1; jmax = iu[i] + nzi;
8649b5e25fSSatish Balay       qm = k;
8749b5e25fSSatish Balay       for (j=jmin; j<jmax+1; j++){
8849b5e25fSSatish Balay         vj = ju[j];
8949b5e25fSSatish Balay         do {
9049b5e25fSSatish Balay           m = qm; qm = q[m];
9149b5e25fSSatish Balay         } while (qm < vj);
9249b5e25fSSatish Balay         if (qm != vj){
9349b5e25fSSatish Balay          nzk++; q[m] = vj; q[vj] = qm; qm = vj;
9449b5e25fSSatish Balay         }
9549b5e25fSSatish Balay       }
9649b5e25fSSatish Balay       i = jl[i]; /* next pivot row */
9749b5e25fSSatish Balay     }
9849b5e25fSSatish Balay 
9949b5e25fSSatish Balay     /* add k to row list for first nonzero element in k-th row */
10049b5e25fSSatish Balay     if (nzk > 0){
10149b5e25fSSatish Balay       i = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */
10249b5e25fSSatish Balay       jl[k] = jl[i]; jl[i] = k;
10349b5e25fSSatish Balay     }
10449b5e25fSSatish Balay     iu[k+1] = iu[k] + nzk;   /* printf(" iu[%d]=%d, umax=%d\n", k+1, iu[k+1],umax);*/
10549b5e25fSSatish Balay 
10649b5e25fSSatish Balay     /* allocate more space to ju if needed */
10749b5e25fSSatish Balay     if (iu[k+1] > umax) { printf("allocate more space, iu[%d]=%d > umax=%d\n",k+1, iu[k+1],umax);
10849b5e25fSSatish Balay       /* estimate how much additional space we will need */
10949b5e25fSSatish Balay       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
11049b5e25fSSatish Balay       /* just double the memory each time */
11149b5e25fSSatish Balay       maxadd = umax;
11249b5e25fSSatish Balay       if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2;
11349b5e25fSSatish Balay       umax += maxadd;
11449b5e25fSSatish Balay 
1159f9cb213SHong Zhang       /* allocate a longer ju */
11682502324SSatish Balay       ierr = PetscMalloc(umax*sizeof(int),&jutmp);CHKERRQ(ierr);
11749b5e25fSSatish Balay       ierr  = PetscMemcpy(jutmp,ju,iu[k]*sizeof(int));CHKERRQ(ierr);
11849b5e25fSSatish Balay       ierr  = PetscFree(ju);CHKERRQ(ierr);
1199f9cb213SHong Zhang       ju    = jutmp;
12049b5e25fSSatish Balay       realloc++; /* count how many times we realloc */
12149b5e25fSSatish Balay     }
12249b5e25fSSatish Balay 
12349b5e25fSSatish Balay     /* save nonzero structure of k-th row in ju */
12449b5e25fSSatish Balay     i=k;
12549b5e25fSSatish Balay     jumin = juptr + 1; juptr += nzk;
12649b5e25fSSatish Balay     for (j=jumin; j<juptr+1; j++){
12749b5e25fSSatish Balay       i=q[i];
12849b5e25fSSatish Balay       ju[j]=i;
12949b5e25fSSatish Balay     }
13077f73120SHong Zhang   }
13149b5e25fSSatish Balay 
13249b5e25fSSatish Balay   if (ai[mbs] != 0) {
13349b5e25fSSatish Balay     PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
134b0a32e0cSBarry Smith     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af);
135b0a32e0cSBarry Smith     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Run with -pc_lu_fill %g or use \n",af);
136b0a32e0cSBarry Smith     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:PCLUSetFill(pc,%g);\n",af);
137b0a32e0cSBarry Smith     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:for best performance.\n");
13849b5e25fSSatish Balay   } else {
139b0a32e0cSBarry Smith      PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Empty matrix.\n");
14049b5e25fSSatish Balay   }
14149b5e25fSSatish Balay 
14277f73120SHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
14349b5e25fSSatish Balay   ierr = PetscFree(q);CHKERRQ(ierr);
14449b5e25fSSatish Balay   ierr = PetscFree(jl);CHKERRQ(ierr);
14549b5e25fSSatish Balay 
14649b5e25fSSatish Balay   /* put together the new matrix */
14749b5e25fSSatish Balay   ierr = MatCreateSeqSBAIJ(A->comm,bs,bs*mbs,bs*mbs,0,PETSC_NULL,B);CHKERRQ(ierr);
148b0a32e0cSBarry Smith   /* PetscLogObjectParent(*B,iperm); */
14949b5e25fSSatish Balay   b = (Mat_SeqSBAIJ*)(*B)->data;
15049b5e25fSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
15149b5e25fSSatish Balay   b->singlemalloc = PETSC_FALSE;
15249b5e25fSSatish Balay   /* the next line frees the default space generated by the Create() */
15349b5e25fSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
15449b5e25fSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
15582502324SSatish Balay   ierr = PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
15649b5e25fSSatish Balay   b->j    = ju;
15749b5e25fSSatish Balay   b->i    = iu;
15849b5e25fSSatish Balay   b->diag = 0;
15949b5e25fSSatish Balay   b->ilen = 0;
16049b5e25fSSatish Balay   b->imax = 0;
16177f73120SHong Zhang   b->row  = perm;
16277f73120SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
163cb718733SHong Zhang   b->icol = perm;
164cb718733SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
165b0a32e0cSBarry Smith   ierr    = PetscMalloc((bs*mbs+bs)*sizeof(Scalar),&b->solve_work);CHKERRQ(ierr);
16649b5e25fSSatish Balay   /* In b structure:  Free imax, ilen, old a, old j.
16749b5e25fSSatish Balay      Allocate idnew, solve_work, new a, new j */
168b0a32e0cSBarry Smith   PetscLogObjectMemory(*B,(iu[mbs]-mbs)*(sizeof(int)+sizeof(MatScalar)));
16949b5e25fSSatish Balay   b->s_maxnz = b->s_nz = iu[mbs];
17049b5e25fSSatish Balay 
1715c0bcdfcSHong Zhang   (*B)->factor                 = FACTOR_CHOLESKY;
17249b5e25fSSatish Balay   (*B)->info.factor_mallocs    = realloc;
17349b5e25fSSatish Balay   (*B)->info.fill_ratio_given  = f;
17449b5e25fSSatish Balay   if (ai[mbs] != 0) {
17549b5e25fSSatish Balay     (*B)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
17649b5e25fSSatish Balay   } else {
17749b5e25fSSatish Balay     (*B)->info.fill_ratio_needed = 0.0;
17849b5e25fSSatish Balay   }
1790aa0ce06SHong Zhang 
180671cb588SHong Zhang   if (perm_identity){
181671cb588SHong Zhang     switch (bs) {
182671cb588SHong Zhang       case 1:
183671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering;
184671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1854d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=1\n");
186671cb588SHong Zhang         break;
187671cb588SHong Zhang       case 2:
188671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
189671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
1904d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
191671cb588SHong Zhang         break;
192671cb588SHong Zhang       case 3:
193671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
194671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
1954d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:sing special in-place natural ordering factor and solve BS=3\n");
196671cb588SHong Zhang         break;
197671cb588SHong Zhang       case 4:
198671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
199671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
2004d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
201671cb588SHong Zhang         break;
202671cb588SHong Zhang       case 5:
203671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
204671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
2054d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
206671cb588SHong Zhang         break;
207671cb588SHong Zhang       case 6:
208671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
209671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
2104d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
211671cb588SHong Zhang         break;
212671cb588SHong Zhang       case 7:
213671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
214671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
2154d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
216671cb588SHong Zhang       break;
217671cb588SHong Zhang       default:
218671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering;
219671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_N_NaturalOrdering;
2204d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS>7\n");
221671cb588SHong Zhang       break;
222671cb588SHong Zhang     }
223671cb588SHong Zhang   }
224671cb588SHong Zhang 
22549b5e25fSSatish Balay   PetscFunctionReturn(0);
22649b5e25fSSatish Balay }
22749b5e25fSSatish Balay 
22835d8aa7fSBarry Smith 
2294a2ae208SSatish Balay #undef __FUNCT__
2304a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N"
2316f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat A,Mat *B)
23249b5e25fSSatish Balay {
23349b5e25fSSatish Balay   Mat                C = *B;
2344c16a6a6SHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
2354c16a6a6SHong Zhang   IS                 perm = b->row;
2364c16a6a6SHong Zhang   int                *perm_ptr,ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
2374c16a6a6SHong Zhang   int                *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
2384c16a6a6SHong Zhang   int                bs=a->bs,bs2 = a->bs2;
2394c16a6a6SHong Zhang   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
2404c16a6a6SHong Zhang   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
24128de702eSHong Zhang   MatScalar          *work;
2424c16a6a6SHong Zhang   int                *pivots;
2434c16a6a6SHong Zhang 
2444c16a6a6SHong Zhang   PetscFunctionBegin;
2459706f043SHong Zhang 
2464c16a6a6SHong Zhang   /* initialization */
24782502324SSatish Balay   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
2484c16a6a6SHong Zhang   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
24982502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
25028de702eSHong Zhang   jl   = il + mbs;
2514c16a6a6SHong Zhang   for (i=0; i<mbs; i++) {
2524c16a6a6SHong Zhang     jl[i] = mbs; il[0] = 0;
2534c16a6a6SHong Zhang   }
254b0a32e0cSBarry Smith   ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr);
25528de702eSHong Zhang   uik  = dk + bs2;
25628de702eSHong Zhang   work = uik + bs2;
25782502324SSatish Balay   ierr = PetscMalloc(bs*sizeof(int),&pivots);CHKERRQ(ierr);
2584c16a6a6SHong Zhang 
2594c16a6a6SHong Zhang   ierr  = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
2604c16a6a6SHong Zhang 
2614c16a6a6SHong Zhang   /* check permutation */
2624c16a6a6SHong Zhang   if (!a->permute){
2634c16a6a6SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
2644c16a6a6SHong Zhang   } else {
2654c16a6a6SHong Zhang     ai   = a->inew; aj = a->jnew;
26682502324SSatish Balay     ierr = PetscMalloc(bs2*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
2674c16a6a6SHong Zhang     ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
26882502324SSatish Balay     ierr = PetscMalloc(ai[mbs]*sizeof(int),&a2anew);CHKERRQ(ierr);
2694c16a6a6SHong Zhang     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr);
2704c16a6a6SHong Zhang 
2714c16a6a6SHong Zhang     for (i=0; i<mbs; i++){
2724c16a6a6SHong Zhang       jmin = ai[i]; jmax = ai[i+1];
2734c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++){
2744c16a6a6SHong Zhang         while (a2anew[j] != j){
2754c16a6a6SHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
2764c16a6a6SHong Zhang           for (k1=0; k1<bs2; k1++){
2774c16a6a6SHong Zhang             dk[k1]       = aa[k*bs2+k1];
2784c16a6a6SHong Zhang             aa[k*bs2+k1] = aa[j*bs2+k1];
2794c16a6a6SHong Zhang             aa[j*bs2+k1] = dk[k1];
2804c16a6a6SHong Zhang           }
2814c16a6a6SHong Zhang         }
2824c16a6a6SHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
2834c16a6a6SHong Zhang         if (i > aj[j]){
2844c16a6a6SHong Zhang           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
2854c16a6a6SHong Zhang           ap = aa + j*bs2;                     /* ptr to the beginning of j-th block of aa */
2864c16a6a6SHong Zhang           for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
2874c16a6a6SHong Zhang           for (k=0; k<bs; k++){               /* j-th block of aa <- dk^T */
2884c16a6a6SHong Zhang             for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1];
2894c16a6a6SHong Zhang           }
2904c16a6a6SHong Zhang         }
2914c16a6a6SHong Zhang       }
2924c16a6a6SHong Zhang     }
293323b833fSBarry Smith     ierr = PetscFree(a2anew);CHKERRQ(ierr);
2944c16a6a6SHong Zhang   }
2954c16a6a6SHong Zhang 
2964c16a6a6SHong Zhang   /* for each row k */
2974c16a6a6SHong Zhang   for (k = 0; k<mbs; k++){
2984c16a6a6SHong Zhang 
2994c16a6a6SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
3004c16a6a6SHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
3014c16a6a6SHong Zhang     if (jmin < jmax) {
3024c16a6a6SHong Zhang       ap = aa + jmin*bs2;
3034c16a6a6SHong Zhang       for (j = jmin; j < jmax; j++){
3044c16a6a6SHong Zhang         vj = perm_ptr[aj[j]];         /* block col. index */
3054c16a6a6SHong Zhang         rtmp_ptr = rtmp + vj*bs2;
3064c16a6a6SHong Zhang         for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
3074c16a6a6SHong Zhang       }
3084c16a6a6SHong Zhang     }
3094c16a6a6SHong Zhang 
3104c16a6a6SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
3114c16a6a6SHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
3124c16a6a6SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
3134c16a6a6SHong Zhang 
3144c16a6a6SHong Zhang     while (i < mbs){
3154c16a6a6SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
3164c16a6a6SHong Zhang 
3174c16a6a6SHong Zhang       /* compute multiplier */
3184c16a6a6SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
3194c16a6a6SHong Zhang 
3204c16a6a6SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
3214c16a6a6SHong Zhang       diag = ba + i*bs2;
3224c16a6a6SHong Zhang       u    = ba + ili*bs2;
3234c16a6a6SHong Zhang       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
3244c16a6a6SHong Zhang       Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
3254c16a6a6SHong Zhang 
3264c16a6a6SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
3274c16a6a6SHong Zhang       Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
3284c16a6a6SHong Zhang 
3294c16a6a6SHong Zhang       /* update -U(i,k) */
3304c16a6a6SHong Zhang       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
3314c16a6a6SHong Zhang 
3324c16a6a6SHong Zhang       /* add multiple of row i to k-th row ... */
3334c16a6a6SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
3344c16a6a6SHong Zhang       if (jmin < jmax){
3354c16a6a6SHong Zhang         for (j=jmin; j<jmax; j++) {
3364c16a6a6SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
3374c16a6a6SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
3384c16a6a6SHong Zhang           u = ba + j*bs2;
3394c16a6a6SHong Zhang           Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
3404c16a6a6SHong Zhang         }
3414c16a6a6SHong Zhang 
3424c16a6a6SHong Zhang         /* ... add i to row list for next nonzero entry */
3434c16a6a6SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
3444c16a6a6SHong Zhang         j     = bj[jmin];
3454c16a6a6SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
3464c16a6a6SHong Zhang       }
3474c16a6a6SHong Zhang       i = nexti;
3484c16a6a6SHong Zhang     }
3494c16a6a6SHong Zhang 
3504c16a6a6SHong Zhang     /* save nonzero entries in k-th row of U ... */
3514c16a6a6SHong Zhang 
3524c16a6a6SHong Zhang     /* invert diagonal block */
3534c16a6a6SHong Zhang     diag = ba+k*bs2;
3544c16a6a6SHong Zhang     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
3554c16a6a6SHong Zhang     Kernel_A_gets_inverse_A(bs,diag,pivots,work);
3564c16a6a6SHong Zhang 
3574c16a6a6SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
3584c16a6a6SHong Zhang     if (jmin < jmax) {
3594c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++){
3604c16a6a6SHong Zhang          vj = bj[j];           /* block col. index of U */
3614c16a6a6SHong Zhang          u   = ba + j*bs2;
3624c16a6a6SHong Zhang          rtmp_ptr = rtmp + vj*bs2;
3634c16a6a6SHong Zhang          for (k1=0; k1<bs2; k1++){
3644c16a6a6SHong Zhang            *u++        = *rtmp_ptr;
3654c16a6a6SHong Zhang            *rtmp_ptr++ = 0.0;
3664c16a6a6SHong Zhang          }
3674c16a6a6SHong Zhang       }
3684c16a6a6SHong Zhang 
3694c16a6a6SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
3704c16a6a6SHong Zhang       il[k] = jmin;
3714c16a6a6SHong Zhang       i     = bj[jmin];
3724c16a6a6SHong Zhang       jl[k] = jl[i]; jl[i] = k;
3734c16a6a6SHong Zhang     }
3744c16a6a6SHong Zhang   }
3754c16a6a6SHong Zhang 
3764c16a6a6SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
3774c16a6a6SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
3784c16a6a6SHong Zhang   ierr = PetscFree(dk);CHKERRQ(ierr);
3794c16a6a6SHong Zhang   ierr = PetscFree(pivots);CHKERRQ(ierr);
3804c16a6a6SHong Zhang   if (a->permute){
3814c16a6a6SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
3824c16a6a6SHong Zhang   }
3834c16a6a6SHong Zhang 
3844c16a6a6SHong Zhang   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
3854c16a6a6SHong Zhang   C->factor    = FACTOR_CHOLESKY;
3864c16a6a6SHong Zhang   C->assembled = PETSC_TRUE;
3874c16a6a6SHong Zhang   C->preallocated = PETSC_TRUE;
388b0a32e0cSBarry Smith   PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */
3894c16a6a6SHong Zhang   PetscFunctionReturn(0);
3904c16a6a6SHong Zhang }
391d4edadd4SHong Zhang 
3924a2ae208SSatish Balay #undef __FUNCT__
3934a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering"
394671cb588SHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat A,Mat *B)
395671cb588SHong Zhang {
396671cb588SHong Zhang   Mat                C = *B;
397671cb588SHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
398671cb588SHong Zhang   int                ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
3991ad1de41SHong Zhang   int                *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
400671cb588SHong Zhang   int                bs=a->bs,bs2 = a->bs2;
401671cb588SHong Zhang   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
402671cb588SHong Zhang   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
40328de702eSHong Zhang   MatScalar          *work;
404671cb588SHong Zhang   int                *pivots;
405671cb588SHong Zhang 
406671cb588SHong Zhang   PetscFunctionBegin;
407671cb588SHong Zhang 
408671cb588SHong Zhang   /* initialization */
409671cb588SHong Zhang 
41082502324SSatish Balay   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
411671cb588SHong Zhang   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
41282502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
41328de702eSHong Zhang   jl   = il + mbs;
414671cb588SHong Zhang   for (i=0; i<mbs; i++) {
415671cb588SHong Zhang     jl[i] = mbs; il[0] = 0;
416671cb588SHong Zhang   }
417b0a32e0cSBarry Smith   ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr);
41828de702eSHong Zhang   uik  = dk + bs2;
41928de702eSHong Zhang   work = uik + bs2;
42082502324SSatish Balay   ierr = PetscMalloc(bs*sizeof(int),&pivots);CHKERRQ(ierr);
421671cb588SHong Zhang 
422671cb588SHong Zhang   ai = a->i; aj = a->j; aa = a->a;
423671cb588SHong Zhang 
424671cb588SHong Zhang   /* for each row k */
425671cb588SHong Zhang   for (k = 0; k<mbs; k++){
426671cb588SHong Zhang 
427671cb588SHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
428671cb588SHong Zhang     jmin = ai[k]; jmax = ai[k+1];
429671cb588SHong Zhang     if (jmin < jmax) {
430671cb588SHong Zhang       ap = aa + jmin*bs2;
431671cb588SHong Zhang       for (j = jmin; j < jmax; j++){
432671cb588SHong Zhang         vj = aj[j];         /* block col. index */
433671cb588SHong Zhang         rtmp_ptr = rtmp + vj*bs2;
434671cb588SHong Zhang         for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
435671cb588SHong Zhang       }
436671cb588SHong Zhang     }
437671cb588SHong Zhang 
438671cb588SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
439671cb588SHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
440671cb588SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
441671cb588SHong Zhang 
442671cb588SHong Zhang     while (i < mbs){
443671cb588SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
444671cb588SHong Zhang 
445671cb588SHong Zhang       /* compute multiplier */
446671cb588SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
447671cb588SHong Zhang 
448671cb588SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
449671cb588SHong Zhang       diag = ba + i*bs2;
450671cb588SHong Zhang       u    = ba + ili*bs2;
451671cb588SHong Zhang       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
452671cb588SHong Zhang       Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
453671cb588SHong Zhang 
454671cb588SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
455671cb588SHong Zhang       Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
456671cb588SHong Zhang 
457671cb588SHong Zhang       /* update -U(i,k) */
458671cb588SHong Zhang       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
459671cb588SHong Zhang 
460671cb588SHong Zhang       /* add multiple of row i to k-th row ... */
461671cb588SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
462671cb588SHong Zhang       if (jmin < jmax){
463671cb588SHong Zhang         for (j=jmin; j<jmax; j++) {
464671cb588SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
465671cb588SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
466671cb588SHong Zhang           u = ba + j*bs2;
467671cb588SHong Zhang           Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
468671cb588SHong Zhang         }
469671cb588SHong Zhang 
470671cb588SHong Zhang         /* ... add i to row list for next nonzero entry */
471671cb588SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
472671cb588SHong Zhang         j     = bj[jmin];
473671cb588SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
474671cb588SHong Zhang       }
475671cb588SHong Zhang       i = nexti;
476671cb588SHong Zhang     }
477671cb588SHong Zhang 
478671cb588SHong Zhang     /* save nonzero entries in k-th row of U ... */
479671cb588SHong Zhang 
480671cb588SHong Zhang     /* invert diagonal block */
481671cb588SHong Zhang     diag = ba+k*bs2;
482671cb588SHong Zhang     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
483671cb588SHong Zhang     Kernel_A_gets_inverse_A(bs,diag,pivots,work);
484671cb588SHong Zhang 
485671cb588SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
486671cb588SHong Zhang     if (jmin < jmax) {
487671cb588SHong Zhang       for (j=jmin; j<jmax; j++){
488671cb588SHong Zhang          vj = bj[j];           /* block col. index of U */
489671cb588SHong Zhang          u   = ba + j*bs2;
490671cb588SHong Zhang          rtmp_ptr = rtmp + vj*bs2;
491671cb588SHong Zhang          for (k1=0; k1<bs2; k1++){
492671cb588SHong Zhang            *u++        = *rtmp_ptr;
493671cb588SHong Zhang            *rtmp_ptr++ = 0.0;
494671cb588SHong Zhang          }
495671cb588SHong Zhang       }
496671cb588SHong Zhang 
497671cb588SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
498671cb588SHong Zhang       il[k] = jmin;
499671cb588SHong Zhang       i     = bj[jmin];
500671cb588SHong Zhang       jl[k] = jl[i]; jl[i] = k;
501671cb588SHong Zhang     }
502671cb588SHong Zhang   }
503671cb588SHong Zhang 
504671cb588SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
505671cb588SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
506671cb588SHong Zhang   ierr = PetscFree(dk);CHKERRQ(ierr);
507671cb588SHong Zhang   ierr = PetscFree(pivots);CHKERRQ(ierr);
508671cb588SHong Zhang 
509671cb588SHong Zhang   C->factor    = FACTOR_CHOLESKY;
510671cb588SHong Zhang   C->assembled = PETSC_TRUE;
511671cb588SHong Zhang   C->preallocated = PETSC_TRUE;
512b0a32e0cSBarry Smith   PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */
513671cb588SHong Zhang   PetscFunctionReturn(0);
514671cb588SHong Zhang }
515671cb588SHong Zhang 
51649b5e25fSSatish Balay /*
517fcf159c0SHong Zhang     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
518cc0c071aSHong Zhang     Version for blocks 2 by 2.
51949b5e25fSSatish Balay */
5204a2ae208SSatish Balay #undef __FUNCT__
5214a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2"
5226f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat A,Mat *B)
52349b5e25fSSatish Balay {
52449b5e25fSSatish Balay   Mat                C = *B;
52591602c66SHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
526cc0c071aSHong Zhang   IS                 perm = b->row;
527cc0c071aSHong Zhang   int                *perm_ptr,ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
528cc0c071aSHong Zhang   int                *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
529a1723e09SHong Zhang   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
530cc0c071aSHong Zhang   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
53149b5e25fSSatish Balay 
53249b5e25fSSatish Balay   PetscFunctionBegin;
533671cb588SHong Zhang 
53491602c66SHong Zhang   /* initialization */
53591602c66SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
53691602c66SHong Zhang      window U(0:k, k:mbs-1).
53791602c66SHong Zhang      jl:    list of rows to be added to uneliminated rows
53891602c66SHong Zhang             i>= k: jl(i) is the first row to be added to row i
53991602c66SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
54091602c66SHong Zhang             jl(i) = mbs indicates the end of a list
54191602c66SHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
54291602c66SHong Zhang             row i of U */
54382502324SSatish Balay   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
544cc0c071aSHong Zhang   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
54582502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
54628de702eSHong Zhang   jl   = il + mbs;
54791602c66SHong Zhang   for (i=0; i<mbs; i++) {
5483845f261SHong Zhang     jl[i] = mbs; il[0] = 0;
54991602c66SHong Zhang   }
55082502324SSatish Balay   ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr);
55128de702eSHong Zhang   uik  = dk + 4;
552cc0c071aSHong Zhang   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
553cc0c071aSHong Zhang 
554cc0c071aSHong Zhang   /* check permutation */
555cc0c071aSHong Zhang   if (!a->permute){
556cc0c071aSHong Zhang     ai = a->i; aj = a->j; aa = a->a;
557cc0c071aSHong Zhang   } else {
558cc0c071aSHong Zhang     ai   = a->inew; aj = a->jnew;
55982502324SSatish Balay     ierr = PetscMalloc(4*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
560cc0c071aSHong Zhang     ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
56182502324SSatish Balay     ierr = PetscMalloc(ai[mbs]*sizeof(int),&a2anew);CHKERRQ(ierr);
562cc0c071aSHong Zhang     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr);
563cc0c071aSHong Zhang 
564cc0c071aSHong Zhang     for (i=0; i<mbs; i++){
565cc0c071aSHong Zhang       jmin = ai[i]; jmax = ai[i+1];
566cc0c071aSHong Zhang       for (j=jmin; j<jmax; j++){
567cc0c071aSHong Zhang         while (a2anew[j] != j){
568cc0c071aSHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
569cc0c071aSHong Zhang           for (k1=0; k1<4; k1++){
570cc0c071aSHong Zhang             dk[k1]       = aa[k*4+k1];
571cc0c071aSHong Zhang             aa[k*4+k1] = aa[j*4+k1];
572cc0c071aSHong Zhang             aa[j*4+k1] = dk[k1];
573cc0c071aSHong Zhang           }
574cc0c071aSHong Zhang         }
575cc0c071aSHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
576cc0c071aSHong Zhang         if (i > aj[j]){
577a1723e09SHong Zhang           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
578cc0c071aSHong Zhang           ap = aa + j*4;     /* ptr to the beginning of the block */
579cc0c071aSHong Zhang           dk[1] = ap[1];     /* swap ap[1] and ap[2] */
580cc0c071aSHong Zhang           ap[1] = ap[2];
581cc0c071aSHong Zhang           ap[2] = dk[1];
582cc0c071aSHong Zhang         }
583cc0c071aSHong Zhang       }
584cc0c071aSHong Zhang     }
585ac355199SBarry Smith     ierr = PetscFree(a2anew);CHKERRQ(ierr);
586cc0c071aSHong Zhang   }
5873845f261SHong Zhang 
58891602c66SHong Zhang   /* for each row k */
58991602c66SHong Zhang   for (k = 0; k<mbs; k++){
59091602c66SHong Zhang 
59191602c66SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
592cc0c071aSHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
59391602c66SHong Zhang     if (jmin < jmax) {
594cc0c071aSHong Zhang       ap = aa + jmin*4;
59591602c66SHong Zhang       for (j = jmin; j < jmax; j++){
596cc0c071aSHong Zhang         vj = perm_ptr[aj[j]];         /* block col. index */
597cc0c071aSHong Zhang         rtmp_ptr = rtmp + vj*4;
598cc0c071aSHong Zhang         for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
59991602c66SHong Zhang       }
60091602c66SHong Zhang     }
60191602c66SHong Zhang 
60291602c66SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
603cc0c071aSHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
60491602c66SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
60591602c66SHong Zhang 
60691602c66SHong Zhang     while (i < mbs){
60791602c66SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
60891602c66SHong Zhang 
6093845f261SHong Zhang       /* compute multiplier */
61091602c66SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
6113845f261SHong Zhang 
6123845f261SHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
613cc0c071aSHong Zhang       diag = ba + i*4;
614cc0c071aSHong Zhang       u    = ba + ili*4;
615cc0c071aSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
616cc0c071aSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
617cc0c071aSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
618cc0c071aSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
6193845f261SHong Zhang 
6203845f261SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
621cc0c071aSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
622cc0c071aSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
623cc0c071aSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
624cc0c071aSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
6253845f261SHong Zhang 
6263845f261SHong Zhang       /* update -U(i,k): ba[ili] = uik */
627cc0c071aSHong Zhang       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
62891602c66SHong Zhang 
62991602c66SHong Zhang       /* add multiple of row i to k-th row ... */
63091602c66SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
63191602c66SHong Zhang       if (jmin < jmax){
6323845f261SHong Zhang         for (j=jmin; j<jmax; j++) {
6333845f261SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
634cc0c071aSHong Zhang           rtmp_ptr = rtmp + bj[j]*4;
635cc0c071aSHong Zhang           u = ba + j*4;
636cc0c071aSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
637cc0c071aSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
638cc0c071aSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
639cc0c071aSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
6403845f261SHong Zhang         }
6413845f261SHong Zhang 
64291602c66SHong Zhang         /* ... add i to row list for next nonzero entry */
64391602c66SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
64491602c66SHong Zhang         j     = bj[jmin];
64591602c66SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
64691602c66SHong Zhang       }
647a1723e09SHong Zhang       i = nexti;
64891602c66SHong Zhang     }
649cc0c071aSHong Zhang 
65091602c66SHong Zhang     /* save nonzero entries in k-th row of U ... */
6513845f261SHong Zhang 
652cc0c071aSHong Zhang     /* invert diagonal block */
653cc0c071aSHong Zhang     diag = ba+k*4;
654cc0c071aSHong Zhang     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
6553845f261SHong Zhang     ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr);
6563845f261SHong Zhang 
65791602c66SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
65891602c66SHong Zhang     if (jmin < jmax) {
65991602c66SHong Zhang       for (j=jmin; j<jmax; j++){
660cc0c071aSHong Zhang          vj = bj[j];           /* block col. index of U */
661cc0c071aSHong Zhang          u   = ba + j*4;
662cc0c071aSHong Zhang          rtmp_ptr = rtmp + vj*4;
663cc0c071aSHong Zhang          for (k1=0; k1<4; k1++){
664cc0c071aSHong Zhang            *u++        = *rtmp_ptr;
665cc0c071aSHong Zhang            *rtmp_ptr++ = 0.0;
666cc0c071aSHong Zhang          }
667cc0c071aSHong Zhang       }
6683845f261SHong Zhang 
66991602c66SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
67091602c66SHong Zhang       il[k] = jmin;
67191602c66SHong Zhang       i     = bj[jmin];
67291602c66SHong Zhang       jl[k] = jl[i]; jl[i] = k;
67391602c66SHong Zhang     }
67491602c66SHong Zhang   }
6753845f261SHong Zhang 
67649b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
67791602c66SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
6783845f261SHong Zhang   ierr = PetscFree(dk);CHKERRQ(ierr);
67991602c66SHong Zhang   if (a->permute) {
68091602c66SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
68191602c66SHong Zhang   }
682cc0c071aSHong Zhang   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
68391602c66SHong Zhang   C->factor    = FACTOR_CHOLESKY;
68449b5e25fSSatish Balay   C->assembled = PETSC_TRUE;
6855c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
686b0a32e0cSBarry Smith   PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */
68749b5e25fSSatish Balay   PetscFunctionReturn(0);
68849b5e25fSSatish Balay }
68991602c66SHong Zhang 
69049b5e25fSSatish Balay /*
69149b5e25fSSatish Balay       Version for when blocks are 2 by 2 Using natural ordering
69249b5e25fSSatish Balay */
6934a2ae208SSatish Balay #undef __FUNCT__
6944a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering"
6956f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat A,Mat *B)
69649b5e25fSSatish Balay {
69749b5e25fSSatish Balay   Mat                C = *B;
698ab27746eSHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
699ab27746eSHong Zhang   int                ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
70071149678SHong Zhang   int                *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
701ab27746eSHong Zhang   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
702ab27746eSHong Zhang   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
70349b5e25fSSatish Balay 
70449b5e25fSSatish Balay   PetscFunctionBegin;
705671cb588SHong Zhang 
706ab27746eSHong Zhang   /* initialization */
707ab27746eSHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
708ab27746eSHong Zhang      window U(0:k, k:mbs-1).
709ab27746eSHong Zhang      jl:    list of rows to be added to uneliminated rows
710ab27746eSHong Zhang             i>= k: jl(i) is the first row to be added to row i
711ab27746eSHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
712ab27746eSHong Zhang             jl(i) = mbs indicates the end of a list
713ab27746eSHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
714ab27746eSHong Zhang             row i of U */
71582502324SSatish Balay   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
716ab27746eSHong Zhang   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
71782502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
71828de702eSHong Zhang   jl   = il + mbs;
719ab27746eSHong Zhang   for (i=0; i<mbs; i++) {
720ab27746eSHong Zhang     jl[i] = mbs; il[0] = 0;
72149b5e25fSSatish Balay   }
72282502324SSatish Balay   ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr);
72328de702eSHong Zhang   uik  = dk + 4;
724ab27746eSHong Zhang 
725ab27746eSHong Zhang   ai = a->i; aj = a->j; aa = a->a;
726ab27746eSHong Zhang 
727ab27746eSHong Zhang   /* for each row k */
728ab27746eSHong Zhang   for (k = 0; k<mbs; k++){
729ab27746eSHong Zhang 
730ab27746eSHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
731ab27746eSHong Zhang     jmin = ai[k]; jmax = ai[k+1];
732ab27746eSHong Zhang     if (jmin < jmax) {
733ab27746eSHong Zhang       ap = aa + jmin*4;
734ab27746eSHong Zhang       for (j = jmin; j < jmax; j++){
735ab27746eSHong Zhang         vj = aj[j];         /* block col. index */
736ab27746eSHong Zhang         rtmp_ptr = rtmp + vj*4;
737ab27746eSHong Zhang         for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
73849b5e25fSSatish Balay       }
73949b5e25fSSatish Balay     }
740ab27746eSHong Zhang 
741ab27746eSHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
742ab27746eSHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
743ab27746eSHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
744ab27746eSHong Zhang 
745ab27746eSHong Zhang     while (i < mbs){
746ab27746eSHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
747ab27746eSHong Zhang 
748ab27746eSHong Zhang       /* compute multiplier */
749ab27746eSHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
750ab27746eSHong Zhang 
751ab27746eSHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
752ab27746eSHong Zhang       diag = ba + i*4;
753ab27746eSHong Zhang       u    = ba + ili*4;
754ab27746eSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
755ab27746eSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
756ab27746eSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
757ab27746eSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
758ab27746eSHong Zhang 
759ab27746eSHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
760ab27746eSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
761ab27746eSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
762ab27746eSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
763ab27746eSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
764ab27746eSHong Zhang 
765ab27746eSHong Zhang       /* update -U(i,k): ba[ili] = uik */
766ab27746eSHong Zhang       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
767ab27746eSHong Zhang 
768ab27746eSHong Zhang       /* add multiple of row i to k-th row ... */
769ab27746eSHong Zhang       jmin = ili + 1; jmax = bi[i+1];
770ab27746eSHong Zhang       if (jmin < jmax){
771ab27746eSHong Zhang         for (j=jmin; j<jmax; j++) {
772ab27746eSHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
773ab27746eSHong Zhang           rtmp_ptr = rtmp + bj[j]*4;
774ab27746eSHong Zhang           u = ba + j*4;
775ab27746eSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
776ab27746eSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
777ab27746eSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
778ab27746eSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
77949b5e25fSSatish Balay         }
780ab27746eSHong Zhang 
781ab27746eSHong Zhang         /* ... add i to row list for next nonzero entry */
782ab27746eSHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
783ab27746eSHong Zhang         j     = bj[jmin];
784ab27746eSHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
78549b5e25fSSatish Balay       }
786ab27746eSHong Zhang       i = nexti;
78749b5e25fSSatish Balay     }
788ab27746eSHong Zhang 
789ab27746eSHong Zhang     /* save nonzero entries in k-th row of U ... */
790ab27746eSHong Zhang 
79149b5e25fSSatish Balay     /* invert diagonal block */
792ab27746eSHong Zhang     diag = ba+k*4;
793ab27746eSHong Zhang     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
794ab27746eSHong Zhang     ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr);
795ab27746eSHong Zhang 
796ab27746eSHong Zhang     jmin = bi[k]; jmax = bi[k+1];
797ab27746eSHong Zhang     if (jmin < jmax) {
798ab27746eSHong Zhang       for (j=jmin; j<jmax; j++){
799ab27746eSHong Zhang          vj = bj[j];           /* block col. index of U */
800ab27746eSHong Zhang          u   = ba + j*4;
801ab27746eSHong Zhang          rtmp_ptr = rtmp + vj*4;
802ab27746eSHong Zhang          for (k1=0; k1<4; k1++){
803ab27746eSHong Zhang            *u++        = *rtmp_ptr;
804ab27746eSHong Zhang            *rtmp_ptr++ = 0.0;
805ab27746eSHong Zhang          }
806ab27746eSHong Zhang       }
807ab27746eSHong Zhang 
808ab27746eSHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
809ab27746eSHong Zhang       il[k] = jmin;
810ab27746eSHong Zhang       i     = bj[jmin];
811ab27746eSHong Zhang       jl[k] = jl[i]; jl[i] = k;
812ab27746eSHong Zhang     }
81349b5e25fSSatish Balay   }
81449b5e25fSSatish Balay 
81549b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
816ab27746eSHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
817ab27746eSHong Zhang   ierr = PetscFree(dk);CHKERRQ(ierr);
818ab27746eSHong Zhang 
819ab27746eSHong Zhang   C->factor    = FACTOR_CHOLESKY;
82049b5e25fSSatish Balay   C->assembled = PETSC_TRUE;
8215c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
822b0a32e0cSBarry Smith   PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */
82349b5e25fSSatish Balay   PetscFunctionReturn(0);
82449b5e25fSSatish Balay }
82549b5e25fSSatish Balay 
82649b5e25fSSatish Balay /*
8275c0bcdfcSHong Zhang     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
82891602c66SHong Zhang     Version for blocks are 1 by 1.
82949b5e25fSSatish Balay */
8304a2ae208SSatish Balay #undef __FUNCT__
8314a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1"
8326f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat A,Mat *B)
83349b5e25fSSatish Balay {
83449b5e25fSSatish Balay   Mat                C = *B;
83549b5e25fSSatish Balay   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
83649b5e25fSSatish Balay   IS                 ip = b->row;
837351d0355SHong Zhang   int                *rip,ierr,i,j,mbs = a->mbs,*bi = b->i,*bj = b->j;
838cb718733SHong Zhang   int                *ai,*aj,*r;
839ab27746eSHong Zhang   int                k,jmin,jmax,*jl,*il,vj,nexti,ili;
840066653e3SSatish Balay   MatScalar          *rtmp;
8412d007305SHong Zhang   MatScalar          *ba = b->a,*aa,ak;
84249b5e25fSSatish Balay   MatScalar          dk,uikdi;
84349b5e25fSSatish Balay 
84449b5e25fSSatish Balay   PetscFunctionBegin;
84549b5e25fSSatish Balay   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
846cb718733SHong Zhang   if (!a->permute){
8472d007305SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
8482d007305SHong Zhang   } else {
8492d007305SHong Zhang     ai = a->inew; aj = a->jnew;
85082502324SSatish Balay     ierr = PetscMalloc(ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
851cb718733SHong Zhang     ierr = PetscMemcpy(aa,a->a,ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
85282502324SSatish Balay     ierr = PetscMalloc(ai[mbs]*sizeof(int),&r);CHKERRQ(ierr);
8532d007305SHong Zhang     ierr= PetscMemcpy(r,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr);
8542d007305SHong Zhang 
8552d007305SHong Zhang     jmin = ai[0]; jmax = ai[mbs];
8562d007305SHong Zhang     for (j=jmin; j<jmax; j++){
857cb718733SHong Zhang       while (r[j] != j){
8582d007305SHong Zhang         k = r[j]; r[j] = r[k]; r[k] = k;
8592d007305SHong Zhang         ak = aa[k]; aa[k] = aa[j]; aa[j] = ak;
8602d007305SHong Zhang       }
8612d007305SHong Zhang     }
862ac355199SBarry Smith     ierr = PetscFree(r);CHKERRQ(ierr);
8632d007305SHong Zhang   }
8642d007305SHong Zhang 
86591602c66SHong Zhang   /* initialization */
86649b5e25fSSatish Balay   /* il and jl record the first nonzero element in each row of the accessing
86749b5e25fSSatish Balay      window U(0:k, k:mbs-1).
86849b5e25fSSatish Balay      jl:    list of rows to be added to uneliminated rows
86949b5e25fSSatish Balay             i>= k: jl(i) is the first row to be added to row i
87049b5e25fSSatish Balay             i<  k: jl(i) is the row following row i in some list of rows
87149b5e25fSSatish Balay             jl(i) = mbs indicates the end of a list
87249b5e25fSSatish Balay      il(i): points to the first nonzero element in columns k,...,mbs-1 of
87349b5e25fSSatish Balay             row i of U */
87482502324SSatish Balay   ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
87582502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
87628de702eSHong Zhang   jl   = il + mbs;
87749b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
87849b5e25fSSatish Balay     rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
87949b5e25fSSatish Balay   }
88049b5e25fSSatish Balay 
88191602c66SHong Zhang   /* for each row k */
88249b5e25fSSatish Balay   for (k = 0; k<mbs; k++){
88349b5e25fSSatish Balay 
88491602c66SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
88549b5e25fSSatish Balay     jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
88649b5e25fSSatish Balay     if (jmin < jmax) {
88749b5e25fSSatish Balay       for (j = jmin; j < jmax; j++){
888351d0355SHong Zhang         vj = rip[aj[j]];
889ab27746eSHong Zhang         rtmp[vj] = aa[j];
89049b5e25fSSatish Balay       }
89149b5e25fSSatish Balay     }
89249b5e25fSSatish Balay 
89391602c66SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
89449b5e25fSSatish Balay     dk = rtmp[k];
89549b5e25fSSatish Balay     i = jl[k]; /* first row to be added to k_th row  */
89649b5e25fSSatish Balay 
89749b5e25fSSatish Balay     while (i < mbs){
89849b5e25fSSatish Balay       nexti = jl[i]; /* next row to be added to k_th row */
89949b5e25fSSatish Balay 
90091602c66SHong Zhang       /* compute multiplier, update D(k) and U(i,k) */
90149b5e25fSSatish Balay       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
90249b5e25fSSatish Balay       uikdi = - ba[ili]*ba[i];
90349b5e25fSSatish Balay       dk += uikdi*ba[ili];
904658e7b3eSHong Zhang       ba[ili] = uikdi; /* -U(i,k) */
90549b5e25fSSatish Balay 
90691602c66SHong Zhang       /* add multiple of row i to k-th row ... */
90749b5e25fSSatish Balay       jmin = ili + 1; jmax = bi[i+1];
90849b5e25fSSatish Balay       if (jmin < jmax){
90949b5e25fSSatish Balay         for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
91091602c66SHong Zhang         /* ... add i to row list for next nonzero entry */
91149b5e25fSSatish Balay         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
91249b5e25fSSatish Balay         j     = bj[jmin];
91349b5e25fSSatish Balay         jl[i] = jl[j]; jl[j] = i; /* update jl */
91449b5e25fSSatish Balay       }
915ab27746eSHong Zhang       i = nexti;
91649b5e25fSSatish Balay     }
91749b5e25fSSatish Balay 
91891602c66SHong Zhang     /* check for zero pivot and save diagoanl element */
91949b5e25fSSatish Balay     if (dk == 0.0){
92029bbc08cSBarry Smith       SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot");
921671cb588SHong Zhang       /*
9228be8c0c7SHong Zhang     }else if (PetscRealPart(dk) < 0){
9238be8c0c7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Negative pivot: d[%d] = %g\n",k,dk);
924671cb588SHong Zhang       */
92549b5e25fSSatish Balay     }
92649b5e25fSSatish Balay 
92791602c66SHong Zhang     /* save nonzero entries in k-th row of U ... */
928f3dd7b2dSKris Buschelman     ba[k] = 1.0/dk;
92949b5e25fSSatish Balay     jmin = bi[k]; jmax = bi[k+1];
93049b5e25fSSatish Balay     if (jmin < jmax) {
93149b5e25fSSatish Balay       for (j=jmin; j<jmax; j++){
932cc0c071aSHong Zhang          vj = bj[j]; ba[j] = rtmp[vj]; rtmp[vj] = 0.0;
93349b5e25fSSatish Balay       }
93491602c66SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
93549b5e25fSSatish Balay       il[k] = jmin;
93649b5e25fSSatish Balay       i     = bj[jmin];
93749b5e25fSSatish Balay       jl[k] = jl[i]; jl[i] = k;
93849b5e25fSSatish Balay     }
93949b5e25fSSatish Balay   }
94049b5e25fSSatish Balay 
94149b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
94249b5e25fSSatish Balay   ierr = PetscFree(il);CHKERRQ(ierr);
943cb718733SHong Zhang   if (a->permute){
944cb718733SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
945cb718733SHong Zhang   }
94649b5e25fSSatish Balay 
94749b5e25fSSatish Balay   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
9488be8c0c7SHong Zhang   C->factor    = FACTOR_CHOLESKY;
94949b5e25fSSatish Balay   C->assembled = PETSC_TRUE;
9505c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
951b0a32e0cSBarry Smith   PetscLogFlops(b->mbs);
95249b5e25fSSatish Balay   PetscFunctionReturn(0);
95349b5e25fSSatish Balay }
95449b5e25fSSatish Balay 
955671cb588SHong Zhang /*
956671cb588SHong Zhang   Version for when blocks are 1 by 1 Using natural ordering
957671cb588SHong Zhang */
9584a2ae208SSatish Balay #undef __FUNCT__
9594a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering"
960671cb588SHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat A,Mat *B)
961671cb588SHong Zhang {
962671cb588SHong Zhang   Mat                C = *B;
963671cb588SHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
964671cb588SHong Zhang   int                ierr,i,j,mbs = a->mbs,*bi = b->i,*bj = b->j;
9651ad1de41SHong Zhang   int                *ai,*aj;
966671cb588SHong Zhang   int                k,jmin,jmax,*jl,*il,vj,nexti,ili;
9671ad1de41SHong Zhang   MatScalar          *rtmp,*ba = b->a,*aa,dk,uikdi;
968671cb588SHong Zhang 
969671cb588SHong Zhang   PetscFunctionBegin;
970671cb588SHong Zhang 
971671cb588SHong Zhang   /* initialization */
972671cb588SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
973671cb588SHong Zhang      window U(0:k, k:mbs-1).
974671cb588SHong Zhang      jl:    list of rows to be added to uneliminated rows
975671cb588SHong Zhang             i>= k: jl(i) is the first row to be added to row i
976671cb588SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
977671cb588SHong Zhang             jl(i) = mbs indicates the end of a list
978671cb588SHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
979671cb588SHong Zhang             row i of U */
98082502324SSatish Balay   ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
98182502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
98228de702eSHong Zhang   jl   = il + mbs;
983671cb588SHong Zhang   for (i=0; i<mbs; i++) {
984671cb588SHong Zhang     rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
985671cb588SHong Zhang   }
986671cb588SHong Zhang 
987671cb588SHong Zhang   ai = a->i; aj = a->j; aa = a->a;
988671cb588SHong Zhang 
989671cb588SHong Zhang   /* for each row k */
990671cb588SHong Zhang   for (k = 0; k<mbs; k++){
991671cb588SHong Zhang 
992671cb588SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
993671cb588SHong Zhang     jmin = ai[k]; jmax = ai[k+1];
994671cb588SHong Zhang     if (jmin < jmax) {
995671cb588SHong Zhang       for (j = jmin; j < jmax; j++){
996671cb588SHong Zhang         vj = aj[j];
997671cb588SHong Zhang         rtmp[vj] = aa[j];
998671cb588SHong Zhang       }
999671cb588SHong Zhang     }
1000671cb588SHong Zhang 
1001671cb588SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
1002671cb588SHong Zhang     dk = rtmp[k];
1003671cb588SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
1004671cb588SHong Zhang 
1005671cb588SHong Zhang     while (i < mbs){
1006671cb588SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
1007671cb588SHong Zhang 
1008671cb588SHong Zhang       /* compute multiplier, update D(k) and U(i,k) */
1009671cb588SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1010671cb588SHong Zhang       uikdi = - ba[ili]*ba[i];
1011671cb588SHong Zhang       dk += uikdi*ba[ili];
1012671cb588SHong Zhang       ba[ili] = uikdi; /* -U(i,k) */
1013671cb588SHong Zhang 
1014671cb588SHong Zhang       /* add multiple of row i to k-th row ... */
1015671cb588SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
1016671cb588SHong Zhang       if (jmin < jmax){
1017671cb588SHong Zhang         for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1018671cb588SHong Zhang         /* ... add i to row list for next nonzero entry */
1019671cb588SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
1020671cb588SHong Zhang         j     = bj[jmin];
1021671cb588SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
1022671cb588SHong Zhang       }
1023671cb588SHong Zhang       i = nexti;
1024671cb588SHong Zhang     }
1025671cb588SHong Zhang 
1026671cb588SHong Zhang     /* check for zero pivot and save diagoanl element */
1027671cb588SHong Zhang     if (dk == 0.0){
1028671cb588SHong Zhang       SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot");
1029671cb588SHong Zhang       /*
1030671cb588SHong Zhang     }else if (PetscRealPart(dk) < 0){
1031671cb588SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Negative pivot: d[%d] = %g\n",k,dk);
1032671cb588SHong Zhang       */
1033671cb588SHong Zhang     }
1034671cb588SHong Zhang 
1035671cb588SHong Zhang     /* save nonzero entries in k-th row of U ... */
1036671cb588SHong Zhang     ba[k] = 1.0/dk;
1037671cb588SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
1038671cb588SHong Zhang     if (jmin < jmax) {
1039671cb588SHong Zhang       for (j=jmin; j<jmax; j++){
1040671cb588SHong Zhang          vj = bj[j]; ba[j] = rtmp[vj]; rtmp[vj] = 0.0;
1041671cb588SHong Zhang       }
1042671cb588SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
1043671cb588SHong Zhang       il[k] = jmin;
1044671cb588SHong Zhang       i     = bj[jmin];
1045671cb588SHong Zhang       jl[k] = jl[i]; jl[i] = k;
1046671cb588SHong Zhang     }
1047671cb588SHong Zhang   }
1048671cb588SHong Zhang 
1049671cb588SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1050671cb588SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
1051671cb588SHong Zhang 
1052671cb588SHong Zhang   C->factor    = FACTOR_CHOLESKY;
1053671cb588SHong Zhang   C->assembled = PETSC_TRUE;
1054671cb588SHong Zhang   C->preallocated = PETSC_TRUE;
1055b0a32e0cSBarry Smith   PetscLogFlops(b->mbs);
1056671cb588SHong Zhang   PetscFunctionReturn(0);
1057671cb588SHong Zhang }
1058671cb588SHong Zhang 
10594a2ae208SSatish Balay #undef __FUNCT__
10604a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ"
1061b8b23757SHong Zhang int MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,PetscReal f)
106249b5e25fSSatish Balay {
10634445ddedSHong Zhang   int ierr;
106449b5e25fSSatish Balay   Mat C;
106549b5e25fSSatish Balay 
106649b5e25fSSatish Balay   PetscFunctionBegin;
1067b8b23757SHong Zhang   ierr = MatCholeskyFactorSymbolic(A,perm,f,&C);CHKERRQ(ierr);
1068a4ada70bSHong Zhang   ierr = MatCholeskyFactorNumeric(A,&C);CHKERRQ(ierr);
10694445ddedSHong Zhang   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
107049b5e25fSSatish Balay   PetscFunctionReturn(0);
107149b5e25fSSatish Balay }
107249b5e25fSSatish Balay 
107349b5e25fSSatish Balay 
1074