xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact.c (revision 3078e140be7f8729cec247bba71411da195a6363)
173f4d377SMatthew Knepley /*$Id: sbaijfact.c,v 1.61 2001/08/06 21:15:47 bsmith Exp $*/
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"
93a7fca6bSBarry 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) {
69*3078e140SBarry Smith           SETERRQ(1," error: duplicate entry in A\n");
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 */
107*3078e140SBarry Smith     if (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);
135*3078e140SBarry Smith     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Run with -pc_cholesky_fill %g or use \n",af);
136*3078e140SBarry Smith     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:PCCholeskySetFill(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;
162bcd9e38bSBarry Smith   b->pivotinblocks = PETSC_FALSE; /* need to get from MatCholeskyInfo */
16377f73120SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
164cb718733SHong Zhang   b->icol = perm;
165cb718733SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
16687828ca2SBarry Smith   ierr    = PetscMalloc((bs*mbs+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
16749b5e25fSSatish Balay   /* In b structure:  Free imax, ilen, old a, old j.
16849b5e25fSSatish Balay      Allocate idnew, solve_work, new a, new j */
169b0a32e0cSBarry Smith   PetscLogObjectMemory(*B,(iu[mbs]-mbs)*(sizeof(int)+sizeof(MatScalar)));
17049b5e25fSSatish Balay   b->s_maxnz = b->s_nz = iu[mbs];
17149b5e25fSSatish Balay 
1725c0bcdfcSHong Zhang   (*B)->factor                 = FACTOR_CHOLESKY;
17349b5e25fSSatish Balay   (*B)->info.factor_mallocs    = realloc;
17449b5e25fSSatish Balay   (*B)->info.fill_ratio_given  = f;
17549b5e25fSSatish Balay   if (ai[mbs] != 0) {
17649b5e25fSSatish Balay     (*B)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
17749b5e25fSSatish Balay   } else {
17849b5e25fSSatish Balay     (*B)->info.fill_ratio_needed = 0.0;
17949b5e25fSSatish Balay   }
1800aa0ce06SHong Zhang 
181671cb588SHong Zhang   if (perm_identity){
182671cb588SHong Zhang     switch (bs) {
183671cb588SHong Zhang       case 1:
184671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering;
185671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1864d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=1\n");
187671cb588SHong Zhang         break;
188671cb588SHong Zhang       case 2:
189671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
190671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
1914d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
192671cb588SHong Zhang         break;
193671cb588SHong Zhang       case 3:
194671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
195671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
1964d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:sing special in-place natural ordering factor and solve BS=3\n");
197671cb588SHong Zhang         break;
198671cb588SHong Zhang       case 4:
199671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
200671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
2014d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
202671cb588SHong Zhang         break;
203671cb588SHong Zhang       case 5:
204671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
205671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
2064d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
207671cb588SHong Zhang         break;
208671cb588SHong Zhang       case 6:
209671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
210671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
2114d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
212671cb588SHong Zhang         break;
213671cb588SHong Zhang       case 7:
214671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
215671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
2164d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
217671cb588SHong Zhang       break;
218671cb588SHong Zhang       default:
219671cb588SHong Zhang         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering;
220671cb588SHong Zhang         (*B)->ops->solve           = MatSolve_SeqSBAIJ_N_NaturalOrdering;
2214d101231SSatish Balay         PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS>7\n");
222671cb588SHong Zhang       break;
223671cb588SHong Zhang     }
224671cb588SHong Zhang   }
225671cb588SHong Zhang 
22649b5e25fSSatish Balay   PetscFunctionReturn(0);
22749b5e25fSSatish Balay }
22849b5e25fSSatish Balay 
22935d8aa7fSBarry Smith 
2304a2ae208SSatish Balay #undef __FUNCT__
2314a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N"
2326f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat A,Mat *B)
23349b5e25fSSatish Balay {
23449b5e25fSSatish Balay   Mat                C = *B;
2354c16a6a6SHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
2364c16a6a6SHong Zhang   IS                 perm = b->row;
2374c16a6a6SHong Zhang   int                *perm_ptr,ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
2384c16a6a6SHong Zhang   int                *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
2394c16a6a6SHong Zhang   int                bs=a->bs,bs2 = a->bs2;
2404c16a6a6SHong Zhang   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
2414c16a6a6SHong Zhang   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
24228de702eSHong Zhang   MatScalar          *work;
2434c16a6a6SHong Zhang   int                *pivots;
2444c16a6a6SHong Zhang 
2454c16a6a6SHong Zhang   PetscFunctionBegin;
2469706f043SHong Zhang 
2474c16a6a6SHong Zhang   /* initialization */
24882502324SSatish Balay   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
2494c16a6a6SHong Zhang   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
25082502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
25128de702eSHong Zhang   jl   = il + mbs;
2524c16a6a6SHong Zhang   for (i=0; i<mbs; i++) {
2534c16a6a6SHong Zhang     jl[i] = mbs; il[0] = 0;
2544c16a6a6SHong Zhang   }
255b0a32e0cSBarry Smith   ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr);
25628de702eSHong Zhang   uik  = dk + bs2;
25728de702eSHong Zhang   work = uik + bs2;
25882502324SSatish Balay   ierr = PetscMalloc(bs*sizeof(int),&pivots);CHKERRQ(ierr);
2594c16a6a6SHong Zhang 
2604c16a6a6SHong Zhang   ierr  = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
2614c16a6a6SHong Zhang 
2624c16a6a6SHong Zhang   /* check permutation */
2634c16a6a6SHong Zhang   if (!a->permute){
2644c16a6a6SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
2654c16a6a6SHong Zhang   } else {
2664c16a6a6SHong Zhang     ai   = a->inew; aj = a->jnew;
26782502324SSatish Balay     ierr = PetscMalloc(bs2*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
2684c16a6a6SHong Zhang     ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
26982502324SSatish Balay     ierr = PetscMalloc(ai[mbs]*sizeof(int),&a2anew);CHKERRQ(ierr);
2704c16a6a6SHong Zhang     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr);
2714c16a6a6SHong Zhang 
2724c16a6a6SHong Zhang     for (i=0; i<mbs; i++){
2734c16a6a6SHong Zhang       jmin = ai[i]; jmax = ai[i+1];
2744c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++){
2754c16a6a6SHong Zhang         while (a2anew[j] != j){
2764c16a6a6SHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
2774c16a6a6SHong Zhang           for (k1=0; k1<bs2; k1++){
2784c16a6a6SHong Zhang             dk[k1]       = aa[k*bs2+k1];
2794c16a6a6SHong Zhang             aa[k*bs2+k1] = aa[j*bs2+k1];
2804c16a6a6SHong Zhang             aa[j*bs2+k1] = dk[k1];
2814c16a6a6SHong Zhang           }
2824c16a6a6SHong Zhang         }
2834c16a6a6SHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
2844c16a6a6SHong Zhang         if (i > aj[j]){
2854c16a6a6SHong Zhang           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
2864c16a6a6SHong Zhang           ap = aa + j*bs2;                     /* ptr to the beginning of j-th block of aa */
2874c16a6a6SHong Zhang           for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
2884c16a6a6SHong Zhang           for (k=0; k<bs; k++){               /* j-th block of aa <- dk^T */
2894c16a6a6SHong Zhang             for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1];
2904c16a6a6SHong Zhang           }
2914c16a6a6SHong Zhang         }
2924c16a6a6SHong Zhang       }
2934c16a6a6SHong Zhang     }
294323b833fSBarry Smith     ierr = PetscFree(a2anew);CHKERRQ(ierr);
2954c16a6a6SHong Zhang   }
2964c16a6a6SHong Zhang 
2974c16a6a6SHong Zhang   /* for each row k */
2984c16a6a6SHong Zhang   for (k = 0; k<mbs; k++){
2994c16a6a6SHong Zhang 
3004c16a6a6SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
3014c16a6a6SHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
3024c16a6a6SHong Zhang     if (jmin < jmax) {
3034c16a6a6SHong Zhang       ap = aa + jmin*bs2;
3044c16a6a6SHong Zhang       for (j = jmin; j < jmax; j++){
3054c16a6a6SHong Zhang         vj = perm_ptr[aj[j]];         /* block col. index */
3064c16a6a6SHong Zhang         rtmp_ptr = rtmp + vj*bs2;
3074c16a6a6SHong Zhang         for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
3084c16a6a6SHong Zhang       }
3094c16a6a6SHong Zhang     }
3104c16a6a6SHong Zhang 
3114c16a6a6SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
3124c16a6a6SHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
3134c16a6a6SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
3144c16a6a6SHong Zhang 
3154c16a6a6SHong Zhang     while (i < mbs){
3164c16a6a6SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
3174c16a6a6SHong Zhang 
3184c16a6a6SHong Zhang       /* compute multiplier */
3194c16a6a6SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
3204c16a6a6SHong Zhang 
3214c16a6a6SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
3224c16a6a6SHong Zhang       diag = ba + i*bs2;
3234c16a6a6SHong Zhang       u    = ba + ili*bs2;
3244c16a6a6SHong Zhang       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
3254c16a6a6SHong Zhang       Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
3264c16a6a6SHong Zhang 
3274c16a6a6SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
3284c16a6a6SHong Zhang       Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
3294c16a6a6SHong Zhang 
3304c16a6a6SHong Zhang       /* update -U(i,k) */
3314c16a6a6SHong Zhang       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
3324c16a6a6SHong Zhang 
3334c16a6a6SHong Zhang       /* add multiple of row i to k-th row ... */
3344c16a6a6SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
3354c16a6a6SHong Zhang       if (jmin < jmax){
3364c16a6a6SHong Zhang         for (j=jmin; j<jmax; j++) {
3374c16a6a6SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
3384c16a6a6SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
3394c16a6a6SHong Zhang           u = ba + j*bs2;
3404c16a6a6SHong Zhang           Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
3414c16a6a6SHong Zhang         }
3424c16a6a6SHong Zhang 
3434c16a6a6SHong Zhang         /* ... add i to row list for next nonzero entry */
3444c16a6a6SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
3454c16a6a6SHong Zhang         j     = bj[jmin];
3464c16a6a6SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
3474c16a6a6SHong Zhang       }
3484c16a6a6SHong Zhang       i = nexti;
3494c16a6a6SHong Zhang     }
3504c16a6a6SHong Zhang 
3514c16a6a6SHong Zhang     /* save nonzero entries in k-th row of U ... */
3524c16a6a6SHong Zhang 
3534c16a6a6SHong Zhang     /* invert diagonal block */
3544c16a6a6SHong Zhang     diag = ba+k*bs2;
3554c16a6a6SHong Zhang     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
3564c16a6a6SHong Zhang     Kernel_A_gets_inverse_A(bs,diag,pivots,work);
3574c16a6a6SHong Zhang 
3584c16a6a6SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
3594c16a6a6SHong Zhang     if (jmin < jmax) {
3604c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++){
3614c16a6a6SHong Zhang          vj = bj[j];           /* block col. index of U */
3624c16a6a6SHong Zhang          u   = ba + j*bs2;
3634c16a6a6SHong Zhang          rtmp_ptr = rtmp + vj*bs2;
3644c16a6a6SHong Zhang          for (k1=0; k1<bs2; k1++){
3654c16a6a6SHong Zhang            *u++        = *rtmp_ptr;
3664c16a6a6SHong Zhang            *rtmp_ptr++ = 0.0;
3674c16a6a6SHong Zhang          }
3684c16a6a6SHong Zhang       }
3694c16a6a6SHong Zhang 
3704c16a6a6SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
3714c16a6a6SHong Zhang       il[k] = jmin;
3724c16a6a6SHong Zhang       i     = bj[jmin];
3734c16a6a6SHong Zhang       jl[k] = jl[i]; jl[i] = k;
3744c16a6a6SHong Zhang     }
3754c16a6a6SHong Zhang   }
3764c16a6a6SHong Zhang 
3774c16a6a6SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
3784c16a6a6SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
3794c16a6a6SHong Zhang   ierr = PetscFree(dk);CHKERRQ(ierr);
3804c16a6a6SHong Zhang   ierr = PetscFree(pivots);CHKERRQ(ierr);
3814c16a6a6SHong Zhang   if (a->permute){
3824c16a6a6SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
3834c16a6a6SHong Zhang   }
3844c16a6a6SHong Zhang 
3854c16a6a6SHong Zhang   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
3864c16a6a6SHong Zhang   C->factor    = FACTOR_CHOLESKY;
3874c16a6a6SHong Zhang   C->assembled = PETSC_TRUE;
3884c16a6a6SHong Zhang   C->preallocated = PETSC_TRUE;
389b0a32e0cSBarry Smith   PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */
3904c16a6a6SHong Zhang   PetscFunctionReturn(0);
3914c16a6a6SHong Zhang }
392d4edadd4SHong Zhang 
3934a2ae208SSatish Balay #undef __FUNCT__
3944a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering"
395671cb588SHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat A,Mat *B)
396671cb588SHong Zhang {
397671cb588SHong Zhang   Mat                C = *B;
398671cb588SHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
399671cb588SHong Zhang   int                ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
4001ad1de41SHong Zhang   int                *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
401671cb588SHong Zhang   int                bs=a->bs,bs2 = a->bs2;
402671cb588SHong Zhang   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
403671cb588SHong Zhang   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
40428de702eSHong Zhang   MatScalar          *work;
405671cb588SHong Zhang   int                *pivots;
406671cb588SHong Zhang 
407671cb588SHong Zhang   PetscFunctionBegin;
408671cb588SHong Zhang 
409671cb588SHong Zhang   /* initialization */
410671cb588SHong Zhang 
41182502324SSatish Balay   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
412671cb588SHong Zhang   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
41382502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
41428de702eSHong Zhang   jl   = il + mbs;
415671cb588SHong Zhang   for (i=0; i<mbs; i++) {
416671cb588SHong Zhang     jl[i] = mbs; il[0] = 0;
417671cb588SHong Zhang   }
418b0a32e0cSBarry Smith   ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr);
41928de702eSHong Zhang   uik  = dk + bs2;
42028de702eSHong Zhang   work = uik + bs2;
42182502324SSatish Balay   ierr = PetscMalloc(bs*sizeof(int),&pivots);CHKERRQ(ierr);
422671cb588SHong Zhang 
423671cb588SHong Zhang   ai = a->i; aj = a->j; aa = a->a;
424671cb588SHong Zhang 
425671cb588SHong Zhang   /* for each row k */
426671cb588SHong Zhang   for (k = 0; k<mbs; k++){
427671cb588SHong Zhang 
428671cb588SHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
429671cb588SHong Zhang     jmin = ai[k]; jmax = ai[k+1];
430671cb588SHong Zhang     if (jmin < jmax) {
431671cb588SHong Zhang       ap = aa + jmin*bs2;
432671cb588SHong Zhang       for (j = jmin; j < jmax; j++){
433671cb588SHong Zhang         vj = aj[j];         /* block col. index */
434671cb588SHong Zhang         rtmp_ptr = rtmp + vj*bs2;
435671cb588SHong Zhang         for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
436671cb588SHong Zhang       }
437671cb588SHong Zhang     }
438671cb588SHong Zhang 
439671cb588SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
440671cb588SHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
441671cb588SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
442671cb588SHong Zhang 
443671cb588SHong Zhang     while (i < mbs){
444671cb588SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
445671cb588SHong Zhang 
446671cb588SHong Zhang       /* compute multiplier */
447671cb588SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
448671cb588SHong Zhang 
449671cb588SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
450671cb588SHong Zhang       diag = ba + i*bs2;
451671cb588SHong Zhang       u    = ba + ili*bs2;
452671cb588SHong Zhang       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
453671cb588SHong Zhang       Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
454671cb588SHong Zhang 
455671cb588SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
456671cb588SHong Zhang       Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
457671cb588SHong Zhang 
458671cb588SHong Zhang       /* update -U(i,k) */
459671cb588SHong Zhang       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
460671cb588SHong Zhang 
461671cb588SHong Zhang       /* add multiple of row i to k-th row ... */
462671cb588SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
463671cb588SHong Zhang       if (jmin < jmax){
464671cb588SHong Zhang         for (j=jmin; j<jmax; j++) {
465671cb588SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
466671cb588SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
467671cb588SHong Zhang           u = ba + j*bs2;
468671cb588SHong Zhang           Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
469671cb588SHong Zhang         }
470671cb588SHong Zhang 
471671cb588SHong Zhang         /* ... add i to row list for next nonzero entry */
472671cb588SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
473671cb588SHong Zhang         j     = bj[jmin];
474671cb588SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
475671cb588SHong Zhang       }
476671cb588SHong Zhang       i = nexti;
477671cb588SHong Zhang     }
478671cb588SHong Zhang 
479671cb588SHong Zhang     /* save nonzero entries in k-th row of U ... */
480671cb588SHong Zhang 
481671cb588SHong Zhang     /* invert diagonal block */
482671cb588SHong Zhang     diag = ba+k*bs2;
483671cb588SHong Zhang     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
484671cb588SHong Zhang     Kernel_A_gets_inverse_A(bs,diag,pivots,work);
485671cb588SHong Zhang 
486671cb588SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
487671cb588SHong Zhang     if (jmin < jmax) {
488671cb588SHong Zhang       for (j=jmin; j<jmax; j++){
489671cb588SHong Zhang          vj = bj[j];           /* block col. index of U */
490671cb588SHong Zhang          u   = ba + j*bs2;
491671cb588SHong Zhang          rtmp_ptr = rtmp + vj*bs2;
492671cb588SHong Zhang          for (k1=0; k1<bs2; k1++){
493671cb588SHong Zhang            *u++        = *rtmp_ptr;
494671cb588SHong Zhang            *rtmp_ptr++ = 0.0;
495671cb588SHong Zhang          }
496671cb588SHong Zhang       }
497671cb588SHong Zhang 
498671cb588SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
499671cb588SHong Zhang       il[k] = jmin;
500671cb588SHong Zhang       i     = bj[jmin];
501671cb588SHong Zhang       jl[k] = jl[i]; jl[i] = k;
502671cb588SHong Zhang     }
503671cb588SHong Zhang   }
504671cb588SHong Zhang 
505671cb588SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
506671cb588SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
507671cb588SHong Zhang   ierr = PetscFree(dk);CHKERRQ(ierr);
508671cb588SHong Zhang   ierr = PetscFree(pivots);CHKERRQ(ierr);
509671cb588SHong Zhang 
510671cb588SHong Zhang   C->factor    = FACTOR_CHOLESKY;
511671cb588SHong Zhang   C->assembled = PETSC_TRUE;
512671cb588SHong Zhang   C->preallocated = PETSC_TRUE;
513b0a32e0cSBarry Smith   PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */
514671cb588SHong Zhang   PetscFunctionReturn(0);
515671cb588SHong Zhang }
516671cb588SHong Zhang 
51749b5e25fSSatish Balay /*
518fcf159c0SHong Zhang     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
519cc0c071aSHong Zhang     Version for blocks 2 by 2.
52049b5e25fSSatish Balay */
5214a2ae208SSatish Balay #undef __FUNCT__
5224a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2"
5236f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat A,Mat *B)
52449b5e25fSSatish Balay {
52549b5e25fSSatish Balay   Mat                C = *B;
52691602c66SHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
527cc0c071aSHong Zhang   IS                 perm = b->row;
528cc0c071aSHong Zhang   int                *perm_ptr,ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
529cc0c071aSHong Zhang   int                *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
530a1723e09SHong Zhang   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
531cc0c071aSHong Zhang   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
53249b5e25fSSatish Balay 
53349b5e25fSSatish Balay   PetscFunctionBegin;
534671cb588SHong Zhang 
53591602c66SHong Zhang   /* initialization */
53691602c66SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
53791602c66SHong Zhang      window U(0:k, k:mbs-1).
53891602c66SHong Zhang      jl:    list of rows to be added to uneliminated rows
53991602c66SHong Zhang             i>= k: jl(i) is the first row to be added to row i
54091602c66SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
54191602c66SHong Zhang             jl(i) = mbs indicates the end of a list
54291602c66SHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
54391602c66SHong Zhang             row i of U */
54482502324SSatish Balay   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
545cc0c071aSHong Zhang   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
54682502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
54728de702eSHong Zhang   jl   = il + mbs;
54891602c66SHong Zhang   for (i=0; i<mbs; i++) {
5493845f261SHong Zhang     jl[i] = mbs; il[0] = 0;
55091602c66SHong Zhang   }
55182502324SSatish Balay   ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr);
55228de702eSHong Zhang   uik  = dk + 4;
553cc0c071aSHong Zhang   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
554cc0c071aSHong Zhang 
555cc0c071aSHong Zhang   /* check permutation */
556cc0c071aSHong Zhang   if (!a->permute){
557cc0c071aSHong Zhang     ai = a->i; aj = a->j; aa = a->a;
558cc0c071aSHong Zhang   } else {
559cc0c071aSHong Zhang     ai   = a->inew; aj = a->jnew;
56082502324SSatish Balay     ierr = PetscMalloc(4*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
561cc0c071aSHong Zhang     ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
56282502324SSatish Balay     ierr = PetscMalloc(ai[mbs]*sizeof(int),&a2anew);CHKERRQ(ierr);
563cc0c071aSHong Zhang     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr);
564cc0c071aSHong Zhang 
565cc0c071aSHong Zhang     for (i=0; i<mbs; i++){
566cc0c071aSHong Zhang       jmin = ai[i]; jmax = ai[i+1];
567cc0c071aSHong Zhang       for (j=jmin; j<jmax; j++){
568cc0c071aSHong Zhang         while (a2anew[j] != j){
569cc0c071aSHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
570cc0c071aSHong Zhang           for (k1=0; k1<4; k1++){
571cc0c071aSHong Zhang             dk[k1]       = aa[k*4+k1];
572cc0c071aSHong Zhang             aa[k*4+k1] = aa[j*4+k1];
573cc0c071aSHong Zhang             aa[j*4+k1] = dk[k1];
574cc0c071aSHong Zhang           }
575cc0c071aSHong Zhang         }
576cc0c071aSHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
577cc0c071aSHong Zhang         if (i > aj[j]){
578a1723e09SHong Zhang           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
579cc0c071aSHong Zhang           ap = aa + j*4;     /* ptr to the beginning of the block */
580cc0c071aSHong Zhang           dk[1] = ap[1];     /* swap ap[1] and ap[2] */
581cc0c071aSHong Zhang           ap[1] = ap[2];
582cc0c071aSHong Zhang           ap[2] = dk[1];
583cc0c071aSHong Zhang         }
584cc0c071aSHong Zhang       }
585cc0c071aSHong Zhang     }
586ac355199SBarry Smith     ierr = PetscFree(a2anew);CHKERRQ(ierr);
587cc0c071aSHong Zhang   }
5883845f261SHong Zhang 
58991602c66SHong Zhang   /* for each row k */
59091602c66SHong Zhang   for (k = 0; k<mbs; k++){
59191602c66SHong Zhang 
59291602c66SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
593cc0c071aSHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
59491602c66SHong Zhang     if (jmin < jmax) {
595cc0c071aSHong Zhang       ap = aa + jmin*4;
59691602c66SHong Zhang       for (j = jmin; j < jmax; j++){
597cc0c071aSHong Zhang         vj = perm_ptr[aj[j]];         /* block col. index */
598cc0c071aSHong Zhang         rtmp_ptr = rtmp + vj*4;
599cc0c071aSHong Zhang         for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
60091602c66SHong Zhang       }
60191602c66SHong Zhang     }
60291602c66SHong Zhang 
60391602c66SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
604cc0c071aSHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
60591602c66SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
60691602c66SHong Zhang 
60791602c66SHong Zhang     while (i < mbs){
60891602c66SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
60991602c66SHong Zhang 
6103845f261SHong Zhang       /* compute multiplier */
61191602c66SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
6123845f261SHong Zhang 
6133845f261SHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
614cc0c071aSHong Zhang       diag = ba + i*4;
615cc0c071aSHong Zhang       u    = ba + ili*4;
616cc0c071aSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
617cc0c071aSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
618cc0c071aSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
619cc0c071aSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
6203845f261SHong Zhang 
6213845f261SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
622cc0c071aSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
623cc0c071aSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
624cc0c071aSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
625cc0c071aSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
6263845f261SHong Zhang 
6273845f261SHong Zhang       /* update -U(i,k): ba[ili] = uik */
628cc0c071aSHong Zhang       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
62991602c66SHong Zhang 
63091602c66SHong Zhang       /* add multiple of row i to k-th row ... */
63191602c66SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
63291602c66SHong Zhang       if (jmin < jmax){
6333845f261SHong Zhang         for (j=jmin; j<jmax; j++) {
6343845f261SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
635cc0c071aSHong Zhang           rtmp_ptr = rtmp + bj[j]*4;
636cc0c071aSHong Zhang           u = ba + j*4;
637cc0c071aSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
638cc0c071aSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
639cc0c071aSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
640cc0c071aSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
6413845f261SHong Zhang         }
6423845f261SHong Zhang 
64391602c66SHong Zhang         /* ... add i to row list for next nonzero entry */
64491602c66SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
64591602c66SHong Zhang         j     = bj[jmin];
64691602c66SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
64791602c66SHong Zhang       }
648a1723e09SHong Zhang       i = nexti;
64991602c66SHong Zhang     }
650cc0c071aSHong Zhang 
65191602c66SHong Zhang     /* save nonzero entries in k-th row of U ... */
6523845f261SHong Zhang 
653cc0c071aSHong Zhang     /* invert diagonal block */
654cc0c071aSHong Zhang     diag = ba+k*4;
655cc0c071aSHong Zhang     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
6563845f261SHong Zhang     ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr);
6573845f261SHong Zhang 
65891602c66SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
65991602c66SHong Zhang     if (jmin < jmax) {
66091602c66SHong Zhang       for (j=jmin; j<jmax; j++){
661cc0c071aSHong Zhang          vj = bj[j];           /* block col. index of U */
662cc0c071aSHong Zhang          u   = ba + j*4;
663cc0c071aSHong Zhang          rtmp_ptr = rtmp + vj*4;
664cc0c071aSHong Zhang          for (k1=0; k1<4; k1++){
665cc0c071aSHong Zhang            *u++        = *rtmp_ptr;
666cc0c071aSHong Zhang            *rtmp_ptr++ = 0.0;
667cc0c071aSHong Zhang          }
668cc0c071aSHong Zhang       }
6693845f261SHong Zhang 
67091602c66SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
67191602c66SHong Zhang       il[k] = jmin;
67291602c66SHong Zhang       i     = bj[jmin];
67391602c66SHong Zhang       jl[k] = jl[i]; jl[i] = k;
67491602c66SHong Zhang     }
67591602c66SHong Zhang   }
6763845f261SHong Zhang 
67749b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
67891602c66SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
6793845f261SHong Zhang   ierr = PetscFree(dk);CHKERRQ(ierr);
68091602c66SHong Zhang   if (a->permute) {
68191602c66SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
68291602c66SHong Zhang   }
683cc0c071aSHong Zhang   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
68491602c66SHong Zhang   C->factor    = FACTOR_CHOLESKY;
68549b5e25fSSatish Balay   C->assembled = PETSC_TRUE;
6865c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
687b0a32e0cSBarry Smith   PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */
68849b5e25fSSatish Balay   PetscFunctionReturn(0);
68949b5e25fSSatish Balay }
69091602c66SHong Zhang 
69149b5e25fSSatish Balay /*
69249b5e25fSSatish Balay       Version for when blocks are 2 by 2 Using natural ordering
69349b5e25fSSatish Balay */
6944a2ae208SSatish Balay #undef __FUNCT__
6954a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering"
6966f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat A,Mat *B)
69749b5e25fSSatish Balay {
69849b5e25fSSatish Balay   Mat                C = *B;
699ab27746eSHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
700ab27746eSHong Zhang   int                ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
70171149678SHong Zhang   int                *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
702ab27746eSHong Zhang   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
703ab27746eSHong Zhang   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
70449b5e25fSSatish Balay 
70549b5e25fSSatish Balay   PetscFunctionBegin;
706671cb588SHong Zhang 
707ab27746eSHong Zhang   /* initialization */
708ab27746eSHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
709ab27746eSHong Zhang      window U(0:k, k:mbs-1).
710ab27746eSHong Zhang      jl:    list of rows to be added to uneliminated rows
711ab27746eSHong Zhang             i>= k: jl(i) is the first row to be added to row i
712ab27746eSHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
713ab27746eSHong Zhang             jl(i) = mbs indicates the end of a list
714ab27746eSHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
715ab27746eSHong Zhang             row i of U */
71682502324SSatish Balay   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
717ab27746eSHong Zhang   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
71882502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
71928de702eSHong Zhang   jl   = il + mbs;
720ab27746eSHong Zhang   for (i=0; i<mbs; i++) {
721ab27746eSHong Zhang     jl[i] = mbs; il[0] = 0;
72249b5e25fSSatish Balay   }
72382502324SSatish Balay   ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr);
72428de702eSHong Zhang   uik  = dk + 4;
725ab27746eSHong Zhang 
726ab27746eSHong Zhang   ai = a->i; aj = a->j; aa = a->a;
727ab27746eSHong Zhang 
728ab27746eSHong Zhang   /* for each row k */
729ab27746eSHong Zhang   for (k = 0; k<mbs; k++){
730ab27746eSHong Zhang 
731ab27746eSHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
732ab27746eSHong Zhang     jmin = ai[k]; jmax = ai[k+1];
733ab27746eSHong Zhang     if (jmin < jmax) {
734ab27746eSHong Zhang       ap = aa + jmin*4;
735ab27746eSHong Zhang       for (j = jmin; j < jmax; j++){
736ab27746eSHong Zhang         vj = aj[j];         /* block col. index */
737ab27746eSHong Zhang         rtmp_ptr = rtmp + vj*4;
738ab27746eSHong Zhang         for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
73949b5e25fSSatish Balay       }
74049b5e25fSSatish Balay     }
741ab27746eSHong Zhang 
742ab27746eSHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
743ab27746eSHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
744ab27746eSHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
745ab27746eSHong Zhang 
746ab27746eSHong Zhang     while (i < mbs){
747ab27746eSHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
748ab27746eSHong Zhang 
749ab27746eSHong Zhang       /* compute multiplier */
750ab27746eSHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
751ab27746eSHong Zhang 
752ab27746eSHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
753ab27746eSHong Zhang       diag = ba + i*4;
754ab27746eSHong Zhang       u    = ba + ili*4;
755ab27746eSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
756ab27746eSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
757ab27746eSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
758ab27746eSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
759ab27746eSHong Zhang 
760ab27746eSHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
761ab27746eSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
762ab27746eSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
763ab27746eSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
764ab27746eSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
765ab27746eSHong Zhang 
766ab27746eSHong Zhang       /* update -U(i,k): ba[ili] = uik */
767ab27746eSHong Zhang       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
768ab27746eSHong Zhang 
769ab27746eSHong Zhang       /* add multiple of row i to k-th row ... */
770ab27746eSHong Zhang       jmin = ili + 1; jmax = bi[i+1];
771ab27746eSHong Zhang       if (jmin < jmax){
772ab27746eSHong Zhang         for (j=jmin; j<jmax; j++) {
773ab27746eSHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
774ab27746eSHong Zhang           rtmp_ptr = rtmp + bj[j]*4;
775ab27746eSHong Zhang           u = ba + j*4;
776ab27746eSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
777ab27746eSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
778ab27746eSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
779ab27746eSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
78049b5e25fSSatish Balay         }
781ab27746eSHong Zhang 
782ab27746eSHong Zhang         /* ... add i to row list for next nonzero entry */
783ab27746eSHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
784ab27746eSHong Zhang         j     = bj[jmin];
785ab27746eSHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
78649b5e25fSSatish Balay       }
787ab27746eSHong Zhang       i = nexti;
78849b5e25fSSatish Balay     }
789ab27746eSHong Zhang 
790ab27746eSHong Zhang     /* save nonzero entries in k-th row of U ... */
791ab27746eSHong Zhang 
79249b5e25fSSatish Balay     /* invert diagonal block */
793ab27746eSHong Zhang     diag = ba+k*4;
794ab27746eSHong Zhang     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
795ab27746eSHong Zhang     ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr);
796ab27746eSHong Zhang 
797ab27746eSHong Zhang     jmin = bi[k]; jmax = bi[k+1];
798ab27746eSHong Zhang     if (jmin < jmax) {
799ab27746eSHong Zhang       for (j=jmin; j<jmax; j++){
800ab27746eSHong Zhang          vj = bj[j];           /* block col. index of U */
801ab27746eSHong Zhang          u   = ba + j*4;
802ab27746eSHong Zhang          rtmp_ptr = rtmp + vj*4;
803ab27746eSHong Zhang          for (k1=0; k1<4; k1++){
804ab27746eSHong Zhang            *u++        = *rtmp_ptr;
805ab27746eSHong Zhang            *rtmp_ptr++ = 0.0;
806ab27746eSHong Zhang          }
807ab27746eSHong Zhang       }
808ab27746eSHong Zhang 
809ab27746eSHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
810ab27746eSHong Zhang       il[k] = jmin;
811ab27746eSHong Zhang       i     = bj[jmin];
812ab27746eSHong Zhang       jl[k] = jl[i]; jl[i] = k;
813ab27746eSHong Zhang     }
81449b5e25fSSatish Balay   }
81549b5e25fSSatish Balay 
81649b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
817ab27746eSHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
818ab27746eSHong Zhang   ierr = PetscFree(dk);CHKERRQ(ierr);
819ab27746eSHong Zhang 
820ab27746eSHong Zhang   C->factor    = FACTOR_CHOLESKY;
82149b5e25fSSatish Balay   C->assembled = PETSC_TRUE;
8225c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
823b0a32e0cSBarry Smith   PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */
82449b5e25fSSatish Balay   PetscFunctionReturn(0);
82549b5e25fSSatish Balay }
82649b5e25fSSatish Balay 
82749b5e25fSSatish Balay /*
8285c0bcdfcSHong Zhang     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
82991602c66SHong Zhang     Version for blocks are 1 by 1.
83049b5e25fSSatish Balay */
8314a2ae208SSatish Balay #undef __FUNCT__
8324a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1"
8336f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat A,Mat *B)
83449b5e25fSSatish Balay {
83549b5e25fSSatish Balay   Mat                C = *B;
83649b5e25fSSatish Balay   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
83749b5e25fSSatish Balay   IS                 ip = b->row;
838351d0355SHong Zhang   int                *rip,ierr,i,j,mbs = a->mbs,*bi = b->i,*bj = b->j;
839cb718733SHong Zhang   int                *ai,*aj,*r;
840ab27746eSHong Zhang   int                k,jmin,jmax,*jl,*il,vj,nexti,ili;
841066653e3SSatish Balay   MatScalar          *rtmp;
8422d007305SHong Zhang   MatScalar          *ba = b->a,*aa,ak;
84349b5e25fSSatish Balay   MatScalar          dk,uikdi;
84449b5e25fSSatish Balay 
84549b5e25fSSatish Balay   PetscFunctionBegin;
84649b5e25fSSatish Balay   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
847cb718733SHong Zhang   if (!a->permute){
8482d007305SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
8492d007305SHong Zhang   } else {
8502d007305SHong Zhang     ai = a->inew; aj = a->jnew;
85182502324SSatish Balay     ierr = PetscMalloc(ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
852cb718733SHong Zhang     ierr = PetscMemcpy(aa,a->a,ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
85382502324SSatish Balay     ierr = PetscMalloc(ai[mbs]*sizeof(int),&r);CHKERRQ(ierr);
8542d007305SHong Zhang     ierr= PetscMemcpy(r,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr);
8552d007305SHong Zhang 
8562d007305SHong Zhang     jmin = ai[0]; jmax = ai[mbs];
8572d007305SHong Zhang     for (j=jmin; j<jmax; j++){
858cb718733SHong Zhang       while (r[j] != j){
8592d007305SHong Zhang         k = r[j]; r[j] = r[k]; r[k] = k;
8602d007305SHong Zhang         ak = aa[k]; aa[k] = aa[j]; aa[j] = ak;
8612d007305SHong Zhang       }
8622d007305SHong Zhang     }
863ac355199SBarry Smith     ierr = PetscFree(r);CHKERRQ(ierr);
8642d007305SHong Zhang   }
8652d007305SHong Zhang 
86691602c66SHong Zhang   /* initialization */
86749b5e25fSSatish Balay   /* il and jl record the first nonzero element in each row of the accessing
86849b5e25fSSatish Balay      window U(0:k, k:mbs-1).
86949b5e25fSSatish Balay      jl:    list of rows to be added to uneliminated rows
87049b5e25fSSatish Balay             i>= k: jl(i) is the first row to be added to row i
87149b5e25fSSatish Balay             i<  k: jl(i) is the row following row i in some list of rows
87249b5e25fSSatish Balay             jl(i) = mbs indicates the end of a list
87349b5e25fSSatish Balay      il(i): points to the first nonzero element in columns k,...,mbs-1 of
87449b5e25fSSatish Balay             row i of U */
87582502324SSatish Balay   ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
87682502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
87728de702eSHong Zhang   jl   = il + mbs;
87849b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
87949b5e25fSSatish Balay     rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
88049b5e25fSSatish Balay   }
88149b5e25fSSatish Balay 
88291602c66SHong Zhang   /* for each row k */
88349b5e25fSSatish Balay   for (k = 0; k<mbs; k++){
88449b5e25fSSatish Balay 
88591602c66SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
88649b5e25fSSatish Balay     jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
88749b5e25fSSatish Balay     if (jmin < jmax) {
88849b5e25fSSatish Balay       for (j = jmin; j < jmax; j++){
889351d0355SHong Zhang         vj = rip[aj[j]];
890ab27746eSHong Zhang         rtmp[vj] = aa[j];
89149b5e25fSSatish Balay       }
89249b5e25fSSatish Balay     }
89349b5e25fSSatish Balay 
89491602c66SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
89549b5e25fSSatish Balay     dk = rtmp[k];
89649b5e25fSSatish Balay     i = jl[k]; /* first row to be added to k_th row  */
89749b5e25fSSatish Balay 
89849b5e25fSSatish Balay     while (i < mbs){
89949b5e25fSSatish Balay       nexti = jl[i]; /* next row to be added to k_th row */
90049b5e25fSSatish Balay 
90191602c66SHong Zhang       /* compute multiplier, update D(k) and U(i,k) */
90249b5e25fSSatish Balay       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
90349b5e25fSSatish Balay       uikdi = - ba[ili]*ba[i];
90449b5e25fSSatish Balay       dk += uikdi*ba[ili];
905658e7b3eSHong Zhang       ba[ili] = uikdi; /* -U(i,k) */
90649b5e25fSSatish Balay 
90791602c66SHong Zhang       /* add multiple of row i to k-th row ... */
90849b5e25fSSatish Balay       jmin = ili + 1; jmax = bi[i+1];
90949b5e25fSSatish Balay       if (jmin < jmax){
91049b5e25fSSatish Balay         for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
91191602c66SHong Zhang         /* ... add i to row list for next nonzero entry */
91249b5e25fSSatish Balay         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
91349b5e25fSSatish Balay         j     = bj[jmin];
91449b5e25fSSatish Balay         jl[i] = jl[j]; jl[j] = i; /* update jl */
91549b5e25fSSatish Balay       }
916ab27746eSHong Zhang       i = nexti;
91749b5e25fSSatish Balay     }
91849b5e25fSSatish Balay 
91991602c66SHong Zhang     /* check for zero pivot and save diagoanl element */
92049b5e25fSSatish Balay     if (dk == 0.0){
92129bbc08cSBarry Smith       SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot");
922671cb588SHong Zhang       /*
9238be8c0c7SHong Zhang     }else if (PetscRealPart(dk) < 0){
9248be8c0c7SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Negative pivot: d[%d] = %g\n",k,dk);
925671cb588SHong Zhang       */
92649b5e25fSSatish Balay     }
92749b5e25fSSatish Balay 
92891602c66SHong Zhang     /* save nonzero entries in k-th row of U ... */
929f3dd7b2dSKris Buschelman     ba[k] = 1.0/dk;
93049b5e25fSSatish Balay     jmin = bi[k]; jmax = bi[k+1];
93149b5e25fSSatish Balay     if (jmin < jmax) {
93249b5e25fSSatish Balay       for (j=jmin; j<jmax; j++){
933cc0c071aSHong Zhang          vj = bj[j]; ba[j] = rtmp[vj]; rtmp[vj] = 0.0;
93449b5e25fSSatish Balay       }
93591602c66SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
93649b5e25fSSatish Balay       il[k] = jmin;
93749b5e25fSSatish Balay       i     = bj[jmin];
93849b5e25fSSatish Balay       jl[k] = jl[i]; jl[i] = k;
93949b5e25fSSatish Balay     }
94049b5e25fSSatish Balay   }
94149b5e25fSSatish Balay 
94249b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
94349b5e25fSSatish Balay   ierr = PetscFree(il);CHKERRQ(ierr);
944cb718733SHong Zhang   if (a->permute){
945cb718733SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
946cb718733SHong Zhang   }
94749b5e25fSSatish Balay 
94849b5e25fSSatish Balay   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
9498be8c0c7SHong Zhang   C->factor    = FACTOR_CHOLESKY;
95049b5e25fSSatish Balay   C->assembled = PETSC_TRUE;
9515c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
952b0a32e0cSBarry Smith   PetscLogFlops(b->mbs);
95349b5e25fSSatish Balay   PetscFunctionReturn(0);
95449b5e25fSSatish Balay }
95549b5e25fSSatish Balay 
956671cb588SHong Zhang /*
957671cb588SHong Zhang   Version for when blocks are 1 by 1 Using natural ordering
958671cb588SHong Zhang */
9594a2ae208SSatish Balay #undef __FUNCT__
9604a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering"
961671cb588SHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat A,Mat *B)
962671cb588SHong Zhang {
963671cb588SHong Zhang   Mat                C = *B;
964671cb588SHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
965671cb588SHong Zhang   int                ierr,i,j,mbs = a->mbs,*bi = b->i,*bj = b->j;
9661ad1de41SHong Zhang   int                *ai,*aj;
967671cb588SHong Zhang   int                k,jmin,jmax,*jl,*il,vj,nexti,ili;
9681ad1de41SHong Zhang   MatScalar          *rtmp,*ba = b->a,*aa,dk,uikdi;
969671cb588SHong Zhang 
970671cb588SHong Zhang   PetscFunctionBegin;
971671cb588SHong Zhang 
972671cb588SHong Zhang   /* initialization */
973671cb588SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
974671cb588SHong Zhang      window U(0:k, k:mbs-1).
975671cb588SHong Zhang      jl:    list of rows to be added to uneliminated rows
976671cb588SHong Zhang             i>= k: jl(i) is the first row to be added to row i
977671cb588SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
978671cb588SHong Zhang             jl(i) = mbs indicates the end of a list
979671cb588SHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
980671cb588SHong Zhang             row i of U */
98182502324SSatish Balay   ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
98282502324SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
98328de702eSHong Zhang   jl   = il + mbs;
984671cb588SHong Zhang   for (i=0; i<mbs; i++) {
985671cb588SHong Zhang     rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
986671cb588SHong Zhang   }
987671cb588SHong Zhang 
988671cb588SHong Zhang   ai = a->i; aj = a->j; aa = a->a;
989671cb588SHong Zhang 
990671cb588SHong Zhang   /* for each row k */
991671cb588SHong Zhang   for (k = 0; k<mbs; k++){
992671cb588SHong Zhang 
993671cb588SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
994671cb588SHong Zhang     jmin = ai[k]; jmax = ai[k+1];
995671cb588SHong Zhang     if (jmin < jmax) {
996671cb588SHong Zhang       for (j = jmin; j < jmax; j++){
997671cb588SHong Zhang         vj = aj[j];
998671cb588SHong Zhang         rtmp[vj] = aa[j];
999671cb588SHong Zhang       }
1000671cb588SHong Zhang     }
1001671cb588SHong Zhang 
1002671cb588SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
1003671cb588SHong Zhang     dk = rtmp[k];
1004671cb588SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
1005671cb588SHong Zhang 
1006671cb588SHong Zhang     while (i < mbs){
1007671cb588SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
1008671cb588SHong Zhang 
1009671cb588SHong Zhang       /* compute multiplier, update D(k) and U(i,k) */
1010671cb588SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1011671cb588SHong Zhang       uikdi = - ba[ili]*ba[i];
1012671cb588SHong Zhang       dk += uikdi*ba[ili];
1013671cb588SHong Zhang       ba[ili] = uikdi; /* -U(i,k) */
1014671cb588SHong Zhang 
1015671cb588SHong Zhang       /* add multiple of row i to k-th row ... */
1016671cb588SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
1017671cb588SHong Zhang       if (jmin < jmax){
1018671cb588SHong Zhang         for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1019671cb588SHong Zhang         /* ... add i to row list for next nonzero entry */
1020671cb588SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
1021671cb588SHong Zhang         j     = bj[jmin];
1022671cb588SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
1023671cb588SHong Zhang       }
1024671cb588SHong Zhang       i = nexti;
1025671cb588SHong Zhang     }
1026671cb588SHong Zhang 
1027671cb588SHong Zhang     /* check for zero pivot and save diagoanl element */
1028671cb588SHong Zhang     if (dk == 0.0){
1029671cb588SHong Zhang       SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot");
1030671cb588SHong Zhang       /*
1031671cb588SHong Zhang     }else if (PetscRealPart(dk) < 0){
1032671cb588SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Negative pivot: d[%d] = %g\n",k,dk);
1033671cb588SHong Zhang       */
1034671cb588SHong Zhang     }
1035671cb588SHong Zhang 
1036671cb588SHong Zhang     /* save nonzero entries in k-th row of U ... */
1037671cb588SHong Zhang     ba[k] = 1.0/dk;
1038671cb588SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
1039671cb588SHong Zhang     if (jmin < jmax) {
1040671cb588SHong Zhang       for (j=jmin; j<jmax; j++){
1041671cb588SHong Zhang          vj = bj[j]; ba[j] = rtmp[vj]; rtmp[vj] = 0.0;
1042671cb588SHong Zhang       }
1043671cb588SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
1044671cb588SHong Zhang       il[k] = jmin;
1045671cb588SHong Zhang       i     = bj[jmin];
1046671cb588SHong Zhang       jl[k] = jl[i]; jl[i] = k;
1047671cb588SHong Zhang     }
1048671cb588SHong Zhang   }
1049671cb588SHong Zhang 
1050671cb588SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1051671cb588SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
1052671cb588SHong Zhang 
1053671cb588SHong Zhang   C->factor    = FACTOR_CHOLESKY;
1054671cb588SHong Zhang   C->assembled = PETSC_TRUE;
1055671cb588SHong Zhang   C->preallocated = PETSC_TRUE;
1056b0a32e0cSBarry Smith   PetscLogFlops(b->mbs);
1057671cb588SHong Zhang   PetscFunctionReturn(0);
1058671cb588SHong Zhang }
1059671cb588SHong Zhang 
10604a2ae208SSatish Balay #undef __FUNCT__
10614a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ"
1062b8b23757SHong Zhang int MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,PetscReal f)
106349b5e25fSSatish Balay {
10644445ddedSHong Zhang   int ierr;
106549b5e25fSSatish Balay   Mat C;
106649b5e25fSSatish Balay 
106749b5e25fSSatish Balay   PetscFunctionBegin;
1068b8b23757SHong Zhang   ierr = MatCholeskyFactorSymbolic(A,perm,f,&C);CHKERRQ(ierr);
1069a4ada70bSHong Zhang   ierr = MatCholeskyFactorNumeric(A,&C);CHKERRQ(ierr);
10704445ddedSHong Zhang   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
107149b5e25fSSatish Balay   PetscFunctionReturn(0);
107249b5e25fSSatish Balay }
107349b5e25fSSatish Balay 
107449b5e25fSSatish Balay 
1075