173f4d377SMatthew Knepley /*$Id: sbaijfact.c,v 1.61 2001/08/06 21:15:47 bsmith Exp $*/ 249b5e25fSSatish Balay 349b5e25fSSatish Balay #include "src/mat/impls/baij/seq/baij.h" 43a7fca6bSBarry Smith #include "src/mat/impls/sbaij/seq/sbaij.h" 549b5e25fSSatish Balay #include "src/vec/vecimpl.h" 649b5e25fSSatish Balay #include "src/inline/ilu.h" 749a6740bSHong Zhang #include "include/petscis.h" 849b5e25fSSatish Balay 98dc9108eSHong Zhang #if !defined(PETSC_USE_COMPLEX) 105f9f512dSHong Zhang /* 115f9f512dSHong Zhang input: 12c037c3f7SHong Zhang F -- numeric factor 135f9f512dSHong Zhang output: 14c037c3f7SHong Zhang nneg, nzero, npos: matrix inertia 155f9f512dSHong Zhang */ 165f9f512dSHong Zhang 175f9f512dSHong Zhang #undef __FUNCT__ 185f9f512dSHong Zhang #define __FUNCT__ "MatGetInertia_SeqSBAIJ" 19c037c3f7SHong Zhang int MatGetInertia_SeqSBAIJ(Mat F,int *nneig,int *nzero,int *npos) 205f9f512dSHong Zhang { 21638f5ce0SDinesh Kaushik Mat_SeqSBAIJ *fact_ptr = (Mat_SeqSBAIJ*)F->data; 223e0d88b5SBarry Smith PetscScalar *dd = fact_ptr->a; 23c037c3f7SHong Zhang int m = F->m,i; 245f9f512dSHong Zhang 255f9f512dSHong Zhang PetscFunctionBegin; 263e0d88b5SBarry Smith if (nneig){ 275f9f512dSHong Zhang *nneig = 0; 285f9f512dSHong Zhang for (i=0; i<m; i++){ 29b36a9721SBarry Smith if (PetscRealPart(dd[i]) < 0.0) (*nneig)++; 305f9f512dSHong Zhang } 313e0d88b5SBarry Smith } 323e0d88b5SBarry Smith if (nzero){ 333e0d88b5SBarry Smith *nzero = 0; 343e0d88b5SBarry Smith for (i=0; i<m; i++){ 353e0d88b5SBarry Smith if (PetscRealPart(dd[i]) == 0.0) (*nzero)++; 363e0d88b5SBarry Smith } 373e0d88b5SBarry Smith } 383e0d88b5SBarry Smith if (npos){ 393e0d88b5SBarry Smith *npos = 0; 403e0d88b5SBarry Smith for (i=0; i<m; i++){ 413e0d88b5SBarry Smith if (PetscRealPart(dd[i]) > 0.0) (*npos)++; 423e0d88b5SBarry Smith } 433e0d88b5SBarry Smith } 445f9f512dSHong Zhang PetscFunctionReturn(0); 455f9f512dSHong Zhang } 468dc9108eSHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */ 475f9f512dSHong Zhang 485f9f512dSHong Zhang /* Using Modified Sparse Row (MSR) storage. 495f9f512dSHong Zhang See page 85, "Iterative Methods ..." by Saad. */ 505f9f512dSHong Zhang /* 515f9f512dSHong Zhang Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP. 525f9f512dSHong Zhang */ 532d9a3abdSHong Zhang /* Use Modified Sparse Row storage for u and ju, see Saad pp.85 */ 544a2ae208SSatish Balay #undef __FUNCT__ 554a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ" 5615e8a5b3SHong Zhang int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *B) 5749b5e25fSSatish Balay { 5849b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b; 59cb718733SHong Zhang int *rip,ierr,i,mbs = a->mbs,*ai,*aj; 60066653e3SSatish Balay int *jutmp,bs = a->bs,bs2=a->bs2; 6187fe1012SHong Zhang int m,realloc = 0,prow; 62d60e5362SHong Zhang int *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd; 63064503c5SHong Zhang int *il,ili,nextprow; 6415e8a5b3SHong Zhang PetscReal f = info->fill; 65671cb588SHong Zhang PetscTruth perm_identity; 6649b5e25fSSatish Balay 6749b5e25fSSatish Balay PetscFunctionBegin; 68cb718733SHong Zhang /* check whether perm is the identity mapping */ 69671cb588SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 70064503c5SHong Zhang 71064503c5SHong Zhang /* -- inplace factorization, i.e., use sbaij for *B -- */ 72064503c5SHong Zhang if (perm_identity && bs==1 ){ 73064503c5SHong Zhang if (!perm_identity) a->permute = PETSC_TRUE; 74064503c5SHong Zhang 75064503c5SHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 76064503c5SHong Zhang 77064503c5SHong Zhang if (perm_identity){ /* without permutation */ 78064503c5SHong Zhang ai = a->i; aj = a->j; 79064503c5SHong Zhang } else { /* non-trivial permutation */ 80064503c5SHong Zhang ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr); 81064503c5SHong Zhang ai = a->inew; aj = a->jnew; 82064503c5SHong Zhang } 83064503c5SHong Zhang 84064503c5SHong Zhang /* initialization */ 85064503c5SHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(int),&iu);CHKERRQ(ierr); 86064503c5SHong Zhang umax = (int)(f*ai[mbs] + 1); 87064503c5SHong Zhang ierr = PetscMalloc(umax*sizeof(int),&ju);CHKERRQ(ierr); 88064503c5SHong Zhang iu[0] = 0; 89064503c5SHong Zhang juidx = 0; /* index for ju */ 90064503c5SHong Zhang ierr = PetscMalloc((3*mbs+1)*sizeof(int),&jl);CHKERRQ(ierr); /* linked list for getting pivot row */ 91064503c5SHong Zhang q = jl + mbs; /* linked list for col index of active row */ 92064503c5SHong Zhang il = q + mbs; 93064503c5SHong Zhang for (i=0; i<mbs; i++){ 94064503c5SHong Zhang jl[i] = mbs; 95064503c5SHong Zhang q[i] = 0; 96064503c5SHong Zhang il[i] = 0; 97064503c5SHong Zhang } 98064503c5SHong Zhang 99064503c5SHong Zhang /* for each row k */ 100064503c5SHong Zhang for (k=0; k<mbs; k++){ 101064503c5SHong Zhang nzk = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */ 102064503c5SHong Zhang q[k] = mbs; 103064503c5SHong Zhang /* initialize nonzero structure of k-th row to row rip[k] of A */ 104064503c5SHong Zhang jmin = ai[rip[k]] +1; /* exclude diag[k] */ 105064503c5SHong Zhang jmax = ai[rip[k]+1]; 106064503c5SHong Zhang for (j=jmin; j<jmax; j++){ 107064503c5SHong Zhang vj = rip[aj[j]]; /* col. value */ 108064503c5SHong Zhang if(vj > k){ 109064503c5SHong Zhang qm = k; 110064503c5SHong Zhang do { 111064503c5SHong Zhang m = qm; qm = q[m]; 112064503c5SHong Zhang } while(qm < vj); 113064503c5SHong Zhang if (qm == vj) { 114064503c5SHong Zhang SETERRQ(1," error: duplicate entry in A\n"); 115064503c5SHong Zhang } 116064503c5SHong Zhang nzk++; 117064503c5SHong Zhang q[m] = vj; 118064503c5SHong Zhang q[vj] = qm; 119064503c5SHong Zhang } /* if(vj > k) */ 120064503c5SHong Zhang } /* for (j=jmin; j<jmax; j++) */ 121064503c5SHong Zhang 122064503c5SHong Zhang /* modify nonzero structure of k-th row by computing fill-in 123064503c5SHong Zhang for each row i to be merged in */ 124064503c5SHong Zhang prow = k; 125064503c5SHong Zhang prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */ 126064503c5SHong Zhang 127064503c5SHong Zhang while (prow < k){ 128064503c5SHong Zhang nextprow = jl[prow]; 129064503c5SHong Zhang 130064503c5SHong Zhang /* merge row prow into k-th row */ 131064503c5SHong Zhang ili = il[prow]; 132064503c5SHong Zhang jmin = ili + 1; /* points to 2nd nzero entry in U(prow,k:mbs-1) */ 133064503c5SHong Zhang jmax = iu[prow+1]; 134064503c5SHong Zhang qm = k; 135064503c5SHong Zhang for (j=jmin; j<jmax; j++){ 136064503c5SHong Zhang vj = ju[j]; 137064503c5SHong Zhang do { 138064503c5SHong Zhang m = qm; qm = q[m]; 139064503c5SHong Zhang } while (qm < vj); 140064503c5SHong Zhang if (qm != vj){ /* a fill */ 141064503c5SHong Zhang nzk++; q[m] = vj; q[vj] = qm; qm = vj; 142064503c5SHong Zhang } 143064503c5SHong Zhang } /* end of for (j=jmin; j<jmax; j++) */ 144064503c5SHong Zhang if (jmin < jmax){ 145064503c5SHong Zhang il[prow] = jmin; 146064503c5SHong Zhang j = ju[jmin]; 147064503c5SHong Zhang jl[prow] = jl[j]; jl[j] = prow; /* update jl */ 148064503c5SHong Zhang } 149064503c5SHong Zhang prow = nextprow; 150064503c5SHong Zhang } 151064503c5SHong Zhang 152064503c5SHong Zhang /* update il and jl */ 153064503c5SHong Zhang if (nzk > 0){ 154064503c5SHong Zhang i = q[k]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */ 155064503c5SHong Zhang jl[k] = jl[i]; jl[i] = k; 156064503c5SHong Zhang il[k] = iu[k] + 1; 157064503c5SHong Zhang } 158064503c5SHong Zhang iu[k+1] = iu[k] + nzk + 1; /* include diag[k] */ 159064503c5SHong Zhang 160064503c5SHong Zhang /* allocate more space to ju if needed */ 161064503c5SHong Zhang if (iu[k+1] > umax) { 162064503c5SHong Zhang /* estimate how much additional space we will need */ 163064503c5SHong Zhang /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 164064503c5SHong Zhang /* just double the memory each time */ 165064503c5SHong Zhang maxadd = umax; 166064503c5SHong Zhang if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2; 167064503c5SHong Zhang umax += maxadd; 168064503c5SHong Zhang 169064503c5SHong Zhang /* allocate a longer ju */ 170064503c5SHong Zhang ierr = PetscMalloc(umax*sizeof(int),&jutmp);CHKERRQ(ierr); 171064503c5SHong Zhang ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(int));CHKERRQ(ierr); 172064503c5SHong Zhang ierr = PetscFree(ju);CHKERRQ(ierr); 173064503c5SHong Zhang ju = jutmp; 174064503c5SHong Zhang realloc++; /* count how many times we realloc */ 175064503c5SHong Zhang } 176064503c5SHong Zhang 177064503c5SHong Zhang /* save nonzero structure of k-th row in ju */ 178064503c5SHong Zhang ju[juidx++] = k; /* diag[k] */ 179064503c5SHong Zhang i = k; 180064503c5SHong Zhang while (nzk --) { 181064503c5SHong Zhang i = q[i]; 182064503c5SHong Zhang ju[juidx++] = i; 183064503c5SHong Zhang } 184064503c5SHong Zhang } 185064503c5SHong Zhang 186064503c5SHong Zhang if (ai[mbs] != 0) { 187064503c5SHong Zhang PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]); 188064503c5SHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af); 189064503c5SHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Run with -pc_cholesky_fill %g or use \n",af); 190064503c5SHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:PCCholeskySetFill(pc,%g);\n",af); 191064503c5SHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:for best performance.\n"); 192064503c5SHong Zhang } else { 193064503c5SHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Empty matrix.\n"); 194064503c5SHong Zhang } 195064503c5SHong Zhang 196064503c5SHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 197064503c5SHong Zhang /* ierr = PetscFree(q);CHKERRQ(ierr); */ 198064503c5SHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 199064503c5SHong Zhang 200064503c5SHong Zhang /* put together the new matrix */ 201064503c5SHong Zhang ierr = MatCreateSeqSBAIJ(A->comm,bs,bs*mbs,bs*mbs,0,PETSC_NULL,B);CHKERRQ(ierr); 202064503c5SHong Zhang /* PetscLogObjectParent(*B,iperm); */ 203064503c5SHong Zhang b = (Mat_SeqSBAIJ*)(*B)->data; 204064503c5SHong Zhang ierr = PetscFree(b->imax);CHKERRQ(ierr); 205064503c5SHong Zhang b->singlemalloc = PETSC_FALSE; 206064503c5SHong Zhang /* the next line frees the default space generated by the Create() */ 207064503c5SHong Zhang ierr = PetscFree(b->a);CHKERRQ(ierr); 208064503c5SHong Zhang ierr = PetscFree(b->ilen);CHKERRQ(ierr); 209064503c5SHong Zhang ierr = PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr); 210064503c5SHong Zhang b->j = ju; 211064503c5SHong Zhang b->i = iu; 212064503c5SHong Zhang b->diag = 0; 213064503c5SHong Zhang b->ilen = 0; 214064503c5SHong Zhang b->imax = 0; 215064503c5SHong Zhang b->row = perm; 21615e8a5b3SHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 217064503c5SHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 218064503c5SHong Zhang b->icol = perm; 219064503c5SHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 220064503c5SHong Zhang ierr = PetscMalloc((bs*mbs+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 221064503c5SHong Zhang /* In b structure: Free imax, ilen, old a, old j. 222064503c5SHong Zhang Allocate idnew, solve_work, new a, new j */ 223064503c5SHong Zhang PetscLogObjectMemory(*B,(iu[mbs]-mbs)*(sizeof(int)+sizeof(MatScalar))); 224064503c5SHong Zhang b->s_maxnz = b->s_nz = iu[mbs]; 225064503c5SHong Zhang 226064503c5SHong Zhang (*B)->factor = FACTOR_CHOLESKY; 227064503c5SHong Zhang (*B)->info.factor_mallocs = realloc; 228064503c5SHong Zhang (*B)->info.fill_ratio_given = f; 229064503c5SHong Zhang if (ai[mbs] != 0) { 230064503c5SHong Zhang (*B)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]); 231064503c5SHong Zhang } else { 232064503c5SHong Zhang (*B)->info.fill_ratio_needed = 0.0; 233064503c5SHong Zhang } 234064503c5SHong Zhang 235064503c5SHong Zhang 236b45a75daSHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering; 237b45a75daSHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 238064503c5SHong Zhang PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=1\n"); 239064503c5SHong Zhang 240064503c5SHong Zhang PetscFunctionReturn(0); 241064503c5SHong Zhang } 242064503c5SHong Zhang /* ----------- end of new code --------------------*/ 243064503c5SHong Zhang 244064503c5SHong Zhang 245671cb588SHong Zhang if (!perm_identity) a->permute = PETSC_TRUE; 2462d007305SHong Zhang 247671cb588SHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 248671cb588SHong Zhang 249671cb588SHong Zhang if (perm_identity){ /* without permutation */ 2502d007305SHong Zhang ai = a->i; aj = a->j; 2512d007305SHong Zhang } else { /* non-trivial permutation */ 252323b833fSBarry Smith ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr); 2532d007305SHong Zhang ai = a->inew; aj = a->jnew; 2542d007305SHong Zhang } 2552d007305SHong Zhang 25649b5e25fSSatish Balay /* initialization */ 257b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&iu);CHKERRQ(ierr); 25849b5e25fSSatish Balay umax = (int)(f*ai[mbs] + 1); umax += mbs + 1; 25982502324SSatish Balay ierr = PetscMalloc(umax*sizeof(int),&ju);CHKERRQ(ierr); 26049b5e25fSSatish Balay iu[0] = mbs+1; 26187fe1012SHong Zhang juidx = mbs + 1; /* index for ju */ 26287fe1012SHong Zhang ierr = PetscMalloc(2*mbs*sizeof(int),&jl);CHKERRQ(ierr); /* linked list for pivot row */ 26387fe1012SHong Zhang q = jl + mbs; /* linked list for col index */ 26449b5e25fSSatish Balay for (i=0; i<mbs; i++){ 26587fe1012SHong Zhang jl[i] = mbs; 26687fe1012SHong Zhang q[i] = 0; 26749b5e25fSSatish Balay } 26849b5e25fSSatish Balay 26949b5e25fSSatish Balay /* for each row k */ 27049b5e25fSSatish Balay for (k=0; k<mbs; k++){ 271064503c5SHong Zhang for (i=0; i<mbs; i++) q[i] = 0; /* to be removed! */ 27249b5e25fSSatish Balay nzk = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */ 27349b5e25fSSatish Balay q[k] = mbs; 27449b5e25fSSatish Balay /* initialize nonzero structure of k-th row to row rip[k] of A */ 275064503c5SHong Zhang jmin = ai[rip[k]] +1; /* exclude diag[k] */ 27649b5e25fSSatish Balay jmax = ai[rip[k]+1]; 27749b5e25fSSatish Balay for (j=jmin; j<jmax; j++){ 278cb718733SHong Zhang vj = rip[aj[j]]; /* col. value */ 27949b5e25fSSatish Balay if(vj > k){ 28049b5e25fSSatish Balay qm = k; 28149b5e25fSSatish Balay do { 28249b5e25fSSatish Balay m = qm; qm = q[m]; 28349b5e25fSSatish Balay } while(qm < vj); 28449b5e25fSSatish Balay if (qm == vj) { 2853078e140SBarry Smith SETERRQ(1," error: duplicate entry in A\n"); 28649b5e25fSSatish Balay } 28749b5e25fSSatish Balay nzk++; 28849b5e25fSSatish Balay q[m] = vj; 28949b5e25fSSatish Balay q[vj] = qm; 29049b5e25fSSatish Balay } /* if(vj > k) */ 29149b5e25fSSatish Balay } /* for (j=jmin; j<jmax; j++) */ 29249b5e25fSSatish Balay 29349b5e25fSSatish Balay /* modify nonzero structure of k-th row by computing fill-in 29449b5e25fSSatish Balay for each row i to be merged in */ 29587fe1012SHong Zhang prow = k; 29687fe1012SHong Zhang prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */ 29787fe1012SHong Zhang 29887fe1012SHong Zhang while (prow < k){ 29987fe1012SHong Zhang /* merge row prow into k-th row */ 30087fe1012SHong Zhang jmin = iu[prow] + 1; jmax = iu[prow+1]; 30149b5e25fSSatish Balay qm = k; 30287fe1012SHong Zhang for (j=jmin; j<jmax; j++){ 30349b5e25fSSatish Balay vj = ju[j]; 30449b5e25fSSatish Balay do { 30549b5e25fSSatish Balay m = qm; qm = q[m]; 30649b5e25fSSatish Balay } while (qm < vj); 30749b5e25fSSatish Balay if (qm != vj){ 30849b5e25fSSatish Balay nzk++; q[m] = vj; q[vj] = qm; qm = vj; 30949b5e25fSSatish Balay } 31049b5e25fSSatish Balay } 31187fe1012SHong Zhang prow = jl[prow]; /* next pivot row */ 31249b5e25fSSatish Balay } 31349b5e25fSSatish Balay 31449b5e25fSSatish Balay /* add k to row list for first nonzero element in k-th row */ 31549b5e25fSSatish Balay if (nzk > 0){ 31649b5e25fSSatish Balay i = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */ 31749b5e25fSSatish Balay jl[k] = jl[i]; jl[i] = k; 31849b5e25fSSatish Balay } 31987fe1012SHong Zhang iu[k+1] = iu[k] + nzk; 32049b5e25fSSatish Balay 32149b5e25fSSatish Balay /* allocate more space to ju if needed */ 3223078e140SBarry Smith if (iu[k+1] > umax) { 32349b5e25fSSatish Balay /* estimate how much additional space we will need */ 32449b5e25fSSatish Balay /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 32549b5e25fSSatish Balay /* just double the memory each time */ 32649b5e25fSSatish Balay maxadd = umax; 32749b5e25fSSatish Balay if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2; 32849b5e25fSSatish Balay umax += maxadd; 32949b5e25fSSatish Balay 3309f9cb213SHong Zhang /* allocate a longer ju */ 33182502324SSatish Balay ierr = PetscMalloc(umax*sizeof(int),&jutmp);CHKERRQ(ierr); 33249b5e25fSSatish Balay ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(int));CHKERRQ(ierr); 33349b5e25fSSatish Balay ierr = PetscFree(ju);CHKERRQ(ierr); 3349f9cb213SHong Zhang ju = jutmp; 33549b5e25fSSatish Balay realloc++; /* count how many times we realloc */ 33649b5e25fSSatish Balay } 33749b5e25fSSatish Balay 33849b5e25fSSatish Balay /* save nonzero structure of k-th row in ju */ 33949b5e25fSSatish Balay i=k; 34087fe1012SHong Zhang while (nzk --) { 34149b5e25fSSatish Balay i = q[i]; 34287fe1012SHong Zhang ju[juidx++] = i; 34349b5e25fSSatish Balay } 34477f73120SHong Zhang } 34549b5e25fSSatish Balay 34649b5e25fSSatish Balay if (ai[mbs] != 0) { 34749b5e25fSSatish Balay PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]); 348b0a32e0cSBarry Smith PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af); 3493078e140SBarry Smith PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Run with -pc_cholesky_fill %g or use \n",af); 3503078e140SBarry Smith PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:PCCholeskySetFill(pc,%g);\n",af); 351b0a32e0cSBarry Smith PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:for best performance.\n"); 35249b5e25fSSatish Balay } else { 353b0a32e0cSBarry Smith PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Empty matrix.\n"); 35449b5e25fSSatish Balay } 35549b5e25fSSatish Balay 35677f73120SHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 35787fe1012SHong Zhang /* ierr = PetscFree(q);CHKERRQ(ierr); */ 35849b5e25fSSatish Balay ierr = PetscFree(jl);CHKERRQ(ierr); 35949b5e25fSSatish Balay 36049b5e25fSSatish Balay /* put together the new matrix */ 36149b5e25fSSatish Balay ierr = MatCreateSeqSBAIJ(A->comm,bs,bs*mbs,bs*mbs,0,PETSC_NULL,B);CHKERRQ(ierr); 362b0a32e0cSBarry Smith /* PetscLogObjectParent(*B,iperm); */ 36349b5e25fSSatish Balay b = (Mat_SeqSBAIJ*)(*B)->data; 36449b5e25fSSatish Balay ierr = PetscFree(b->imax);CHKERRQ(ierr); 36549b5e25fSSatish Balay b->singlemalloc = PETSC_FALSE; 36649b5e25fSSatish Balay /* the next line frees the default space generated by the Create() */ 36749b5e25fSSatish Balay ierr = PetscFree(b->a);CHKERRQ(ierr); 36849b5e25fSSatish Balay ierr = PetscFree(b->ilen);CHKERRQ(ierr); 36982502324SSatish Balay ierr = PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr); 37049b5e25fSSatish Balay b->j = ju; 37149b5e25fSSatish Balay b->i = iu; 37249b5e25fSSatish Balay b->diag = 0; 37349b5e25fSSatish Balay b->ilen = 0; 37449b5e25fSSatish Balay b->imax = 0; 37577f73120SHong Zhang b->row = perm; 37615e8a5b3SHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 37777f73120SHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 378cb718733SHong Zhang b->icol = perm; 379cb718733SHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 38087828ca2SBarry Smith ierr = PetscMalloc((bs*mbs+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 38149b5e25fSSatish Balay /* In b structure: Free imax, ilen, old a, old j. 38249b5e25fSSatish Balay Allocate idnew, solve_work, new a, new j */ 383b0a32e0cSBarry Smith PetscLogObjectMemory(*B,(iu[mbs]-mbs)*(sizeof(int)+sizeof(MatScalar))); 38449b5e25fSSatish Balay b->s_maxnz = b->s_nz = iu[mbs]; 38549b5e25fSSatish Balay 3865c0bcdfcSHong Zhang (*B)->factor = FACTOR_CHOLESKY; 38749b5e25fSSatish Balay (*B)->info.factor_mallocs = realloc; 38849b5e25fSSatish Balay (*B)->info.fill_ratio_given = f; 38949b5e25fSSatish Balay if (ai[mbs] != 0) { 39049b5e25fSSatish Balay (*B)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]); 39149b5e25fSSatish Balay } else { 39249b5e25fSSatish Balay (*B)->info.fill_ratio_needed = 0.0; 39349b5e25fSSatish Balay } 3940aa0ce06SHong Zhang 395671cb588SHong Zhang if (perm_identity){ 396671cb588SHong Zhang switch (bs) { 397671cb588SHong Zhang case 1: 398671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering; 399671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 4004d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=1\n"); 401671cb588SHong Zhang break; 402671cb588SHong Zhang case 2: 403671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 404671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 4054d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 406671cb588SHong Zhang break; 407671cb588SHong Zhang case 3: 408671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 409671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 4104d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:sing special in-place natural ordering factor and solve BS=3\n"); 411671cb588SHong Zhang break; 412671cb588SHong Zhang case 4: 413671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 414671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 4154d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 416671cb588SHong Zhang break; 417671cb588SHong Zhang case 5: 418671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 419671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 4204d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 421671cb588SHong Zhang break; 422671cb588SHong Zhang case 6: 423671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 424671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 4254d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 426671cb588SHong Zhang break; 427671cb588SHong Zhang case 7: 428671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 429671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 4304d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 431671cb588SHong Zhang break; 432671cb588SHong Zhang default: 433671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering; 434671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_N_NaturalOrdering; 4354d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS>7\n"); 436671cb588SHong Zhang break; 437671cb588SHong Zhang } 438671cb588SHong Zhang } 439671cb588SHong Zhang 44049b5e25fSSatish Balay PetscFunctionReturn(0); 44149b5e25fSSatish Balay } 44249b5e25fSSatish Balay 44335d8aa7fSBarry Smith 4444a2ae208SSatish Balay #undef __FUNCT__ 4454a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N" 4466f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat A,Mat *B) 44749b5e25fSSatish Balay { 44849b5e25fSSatish Balay Mat C = *B; 4494c16a6a6SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 4504c16a6a6SHong Zhang IS perm = b->row; 4514c16a6a6SHong Zhang int *perm_ptr,ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 4524c16a6a6SHong Zhang int *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 4534c16a6a6SHong Zhang int bs=a->bs,bs2 = a->bs2; 4544c16a6a6SHong Zhang MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 4554c16a6a6SHong Zhang MatScalar *u,*diag,*rtmp,*rtmp_ptr; 45628de702eSHong Zhang MatScalar *work; 4574c16a6a6SHong Zhang int *pivots; 4584c16a6a6SHong Zhang 4594c16a6a6SHong Zhang PetscFunctionBegin; 4609706f043SHong Zhang 4614c16a6a6SHong Zhang /* initialization */ 46282502324SSatish Balay ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 4634c16a6a6SHong Zhang ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr); 46482502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr); 46528de702eSHong Zhang jl = il + mbs; 4664c16a6a6SHong Zhang for (i=0; i<mbs; i++) { 4674c16a6a6SHong Zhang jl[i] = mbs; il[0] = 0; 4684c16a6a6SHong Zhang } 469b0a32e0cSBarry Smith ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr); 47028de702eSHong Zhang uik = dk + bs2; 47128de702eSHong Zhang work = uik + bs2; 47282502324SSatish Balay ierr = PetscMalloc(bs*sizeof(int),&pivots);CHKERRQ(ierr); 4734c16a6a6SHong Zhang 4744c16a6a6SHong Zhang ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr); 4754c16a6a6SHong Zhang 4764c16a6a6SHong Zhang /* check permutation */ 4774c16a6a6SHong Zhang if (!a->permute){ 4784c16a6a6SHong Zhang ai = a->i; aj = a->j; aa = a->a; 4794c16a6a6SHong Zhang } else { 4804c16a6a6SHong Zhang ai = a->inew; aj = a->jnew; 48182502324SSatish Balay ierr = PetscMalloc(bs2*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); 4824c16a6a6SHong Zhang ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); 48382502324SSatish Balay ierr = PetscMalloc(ai[mbs]*sizeof(int),&a2anew);CHKERRQ(ierr); 4844c16a6a6SHong Zhang ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr); 4854c16a6a6SHong Zhang 4864c16a6a6SHong Zhang for (i=0; i<mbs; i++){ 4874c16a6a6SHong Zhang jmin = ai[i]; jmax = ai[i+1]; 4884c16a6a6SHong Zhang for (j=jmin; j<jmax; j++){ 4894c16a6a6SHong Zhang while (a2anew[j] != j){ 4904c16a6a6SHong Zhang k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k; 4914c16a6a6SHong Zhang for (k1=0; k1<bs2; k1++){ 4924c16a6a6SHong Zhang dk[k1] = aa[k*bs2+k1]; 4934c16a6a6SHong Zhang aa[k*bs2+k1] = aa[j*bs2+k1]; 4944c16a6a6SHong Zhang aa[j*bs2+k1] = dk[k1]; 4954c16a6a6SHong Zhang } 4964c16a6a6SHong Zhang } 4974c16a6a6SHong Zhang /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */ 4984c16a6a6SHong Zhang if (i > aj[j]){ 4994c16a6a6SHong Zhang /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */ 5004c16a6a6SHong Zhang ap = aa + j*bs2; /* ptr to the beginning of j-th block of aa */ 5014c16a6a6SHong Zhang for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */ 5024c16a6a6SHong Zhang for (k=0; k<bs; k++){ /* j-th block of aa <- dk^T */ 5034c16a6a6SHong Zhang for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1]; 5044c16a6a6SHong Zhang } 5054c16a6a6SHong Zhang } 5064c16a6a6SHong Zhang } 5074c16a6a6SHong Zhang } 508323b833fSBarry Smith ierr = PetscFree(a2anew);CHKERRQ(ierr); 5094c16a6a6SHong Zhang } 5104c16a6a6SHong Zhang 5114c16a6a6SHong Zhang /* for each row k */ 5124c16a6a6SHong Zhang for (k = 0; k<mbs; k++){ 5134c16a6a6SHong Zhang 5144c16a6a6SHong Zhang /*initialize k-th row with elements nonzero in row perm(k) of A */ 5154c16a6a6SHong Zhang jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1]; 516057f5ba7SHong Zhang 5174c16a6a6SHong Zhang ap = aa + jmin*bs2; 5184c16a6a6SHong Zhang for (j = jmin; j < jmax; j++){ 5194c16a6a6SHong Zhang vj = perm_ptr[aj[j]]; /* block col. index */ 5204c16a6a6SHong Zhang rtmp_ptr = rtmp + vj*bs2; 5214c16a6a6SHong Zhang for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++; 5224c16a6a6SHong Zhang } 5234c16a6a6SHong Zhang 5244c16a6a6SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 5254c16a6a6SHong Zhang ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 5264c16a6a6SHong Zhang i = jl[k]; /* first row to be added to k_th row */ 5274c16a6a6SHong Zhang 528057f5ba7SHong Zhang while (i < k){ 5294c16a6a6SHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 5304c16a6a6SHong Zhang 5314c16a6a6SHong Zhang /* compute multiplier */ 5324c16a6a6SHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 5334c16a6a6SHong Zhang 5344c16a6a6SHong Zhang /* uik = -inv(Di)*U_bar(i,k) */ 5354c16a6a6SHong Zhang diag = ba + i*bs2; 5364c16a6a6SHong Zhang u = ba + ili*bs2; 5374c16a6a6SHong Zhang ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 5384c16a6a6SHong Zhang Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u); 5394c16a6a6SHong Zhang 5404c16a6a6SHong Zhang /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 5414c16a6a6SHong Zhang Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u); 5424c16a6a6SHong Zhang 5434c16a6a6SHong Zhang /* update -U(i,k) */ 5444c16a6a6SHong Zhang ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 5454c16a6a6SHong Zhang 5464c16a6a6SHong Zhang /* add multiple of row i to k-th row ... */ 5474c16a6a6SHong Zhang jmin = ili + 1; jmax = bi[i+1]; 5484c16a6a6SHong Zhang if (jmin < jmax){ 5494c16a6a6SHong Zhang for (j=jmin; j<jmax; j++) { 5504c16a6a6SHong Zhang /* rtmp += -U(i,k)^T * U_bar(i,j) */ 5514c16a6a6SHong Zhang rtmp_ptr = rtmp + bj[j]*bs2; 5524c16a6a6SHong Zhang u = ba + j*bs2; 5534c16a6a6SHong Zhang Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u); 5544c16a6a6SHong Zhang } 5554c16a6a6SHong Zhang 5564c16a6a6SHong Zhang /* ... add i to row list for next nonzero entry */ 5574c16a6a6SHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 5584c16a6a6SHong Zhang j = bj[jmin]; 5594c16a6a6SHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 5604c16a6a6SHong Zhang } 5614c16a6a6SHong Zhang i = nexti; 5624c16a6a6SHong Zhang } 5634c16a6a6SHong Zhang 5644c16a6a6SHong Zhang /* save nonzero entries in k-th row of U ... */ 5654c16a6a6SHong Zhang 5664c16a6a6SHong Zhang /* invert diagonal block */ 5674c16a6a6SHong Zhang diag = ba+k*bs2; 5684c16a6a6SHong Zhang ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr); 5694c16a6a6SHong Zhang Kernel_A_gets_inverse_A(bs,diag,pivots,work); 5704c16a6a6SHong Zhang 5714c16a6a6SHong Zhang jmin = bi[k]; jmax = bi[k+1]; 5724c16a6a6SHong Zhang if (jmin < jmax) { 5734c16a6a6SHong Zhang for (j=jmin; j<jmax; j++){ 5744c16a6a6SHong Zhang vj = bj[j]; /* block col. index of U */ 5754c16a6a6SHong Zhang u = ba + j*bs2; 5764c16a6a6SHong Zhang rtmp_ptr = rtmp + vj*bs2; 5774c16a6a6SHong Zhang for (k1=0; k1<bs2; k1++){ 5784c16a6a6SHong Zhang *u++ = *rtmp_ptr; 5794c16a6a6SHong Zhang *rtmp_ptr++ = 0.0; 5804c16a6a6SHong Zhang } 5814c16a6a6SHong Zhang } 5824c16a6a6SHong Zhang 5834c16a6a6SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 5844c16a6a6SHong Zhang il[k] = jmin; 5854c16a6a6SHong Zhang i = bj[jmin]; 5864c16a6a6SHong Zhang jl[k] = jl[i]; jl[i] = k; 5874c16a6a6SHong Zhang } 5884c16a6a6SHong Zhang } 5894c16a6a6SHong Zhang 5904c16a6a6SHong Zhang ierr = PetscFree(rtmp);CHKERRQ(ierr); 5914c16a6a6SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 5924c16a6a6SHong Zhang ierr = PetscFree(dk);CHKERRQ(ierr); 5934c16a6a6SHong Zhang ierr = PetscFree(pivots);CHKERRQ(ierr); 5944c16a6a6SHong Zhang if (a->permute){ 5954c16a6a6SHong Zhang ierr = PetscFree(aa);CHKERRQ(ierr); 5964c16a6a6SHong Zhang } 5974c16a6a6SHong Zhang 5984c16a6a6SHong Zhang ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr); 5994c16a6a6SHong Zhang C->factor = FACTOR_CHOLESKY; 6004c16a6a6SHong Zhang C->assembled = PETSC_TRUE; 6014c16a6a6SHong Zhang C->preallocated = PETSC_TRUE; 602b0a32e0cSBarry Smith PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */ 6034c16a6a6SHong Zhang PetscFunctionReturn(0); 6044c16a6a6SHong Zhang } 605d4edadd4SHong Zhang 6064a2ae208SSatish Balay #undef __FUNCT__ 6074a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering" 608671cb588SHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat A,Mat *B) 609671cb588SHong Zhang { 610671cb588SHong Zhang Mat C = *B; 611671cb588SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 612671cb588SHong Zhang int ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 6131ad1de41SHong Zhang int *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 614671cb588SHong Zhang int bs=a->bs,bs2 = a->bs2; 615671cb588SHong Zhang MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 616671cb588SHong Zhang MatScalar *u,*diag,*rtmp,*rtmp_ptr; 61728de702eSHong Zhang MatScalar *work; 618671cb588SHong Zhang int *pivots; 619671cb588SHong Zhang 620671cb588SHong Zhang PetscFunctionBegin; 621671cb588SHong Zhang 622671cb588SHong Zhang /* initialization */ 623671cb588SHong Zhang 62482502324SSatish Balay ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 625671cb588SHong Zhang ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr); 62682502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr); 62728de702eSHong Zhang jl = il + mbs; 628671cb588SHong Zhang for (i=0; i<mbs; i++) { 629671cb588SHong Zhang jl[i] = mbs; il[0] = 0; 630671cb588SHong Zhang } 631b0a32e0cSBarry Smith ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr); 63228de702eSHong Zhang uik = dk + bs2; 63328de702eSHong Zhang work = uik + bs2; 63482502324SSatish Balay ierr = PetscMalloc(bs*sizeof(int),&pivots);CHKERRQ(ierr); 635671cb588SHong Zhang 636671cb588SHong Zhang ai = a->i; aj = a->j; aa = a->a; 637671cb588SHong Zhang 638671cb588SHong Zhang /* for each row k */ 639671cb588SHong Zhang for (k = 0; k<mbs; k++){ 640671cb588SHong Zhang 641671cb588SHong Zhang /*initialize k-th row with elements nonzero in row k of A */ 642671cb588SHong Zhang jmin = ai[k]; jmax = ai[k+1]; 643671cb588SHong Zhang ap = aa + jmin*bs2; 644671cb588SHong Zhang for (j = jmin; j < jmax; j++){ 645671cb588SHong Zhang vj = aj[j]; /* block col. index */ 646671cb588SHong Zhang rtmp_ptr = rtmp + vj*bs2; 647671cb588SHong Zhang for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++; 648671cb588SHong Zhang } 649671cb588SHong Zhang 650671cb588SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 651671cb588SHong Zhang ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 652671cb588SHong Zhang i = jl[k]; /* first row to be added to k_th row */ 653671cb588SHong Zhang 654057f5ba7SHong Zhang while (i < k){ 655671cb588SHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 656671cb588SHong Zhang 657671cb588SHong Zhang /* compute multiplier */ 658671cb588SHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 659671cb588SHong Zhang 660671cb588SHong Zhang /* uik = -inv(Di)*U_bar(i,k) */ 661671cb588SHong Zhang diag = ba + i*bs2; 662671cb588SHong Zhang u = ba + ili*bs2; 663671cb588SHong Zhang ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 664671cb588SHong Zhang Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u); 665671cb588SHong Zhang 666671cb588SHong Zhang /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 667671cb588SHong Zhang Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u); 668671cb588SHong Zhang 669671cb588SHong Zhang /* update -U(i,k) */ 670671cb588SHong Zhang ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 671671cb588SHong Zhang 672671cb588SHong Zhang /* add multiple of row i to k-th row ... */ 673671cb588SHong Zhang jmin = ili + 1; jmax = bi[i+1]; 674671cb588SHong Zhang if (jmin < jmax){ 675671cb588SHong Zhang for (j=jmin; j<jmax; j++) { 676671cb588SHong Zhang /* rtmp += -U(i,k)^T * U_bar(i,j) */ 677671cb588SHong Zhang rtmp_ptr = rtmp + bj[j]*bs2; 678671cb588SHong Zhang u = ba + j*bs2; 679671cb588SHong Zhang Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u); 680671cb588SHong Zhang } 681671cb588SHong Zhang 682671cb588SHong Zhang /* ... add i to row list for next nonzero entry */ 683671cb588SHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 684671cb588SHong Zhang j = bj[jmin]; 685671cb588SHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 686671cb588SHong Zhang } 687671cb588SHong Zhang i = nexti; 688671cb588SHong Zhang } 689671cb588SHong Zhang 690671cb588SHong Zhang /* save nonzero entries in k-th row of U ... */ 691671cb588SHong Zhang 692671cb588SHong Zhang /* invert diagonal block */ 693671cb588SHong Zhang diag = ba+k*bs2; 694671cb588SHong Zhang ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr); 695671cb588SHong Zhang Kernel_A_gets_inverse_A(bs,diag,pivots,work); 696671cb588SHong Zhang 697671cb588SHong Zhang jmin = bi[k]; jmax = bi[k+1]; 698671cb588SHong Zhang if (jmin < jmax) { 699671cb588SHong Zhang for (j=jmin; j<jmax; j++){ 700671cb588SHong Zhang vj = bj[j]; /* block col. index of U */ 701671cb588SHong Zhang u = ba + j*bs2; 702671cb588SHong Zhang rtmp_ptr = rtmp + vj*bs2; 703671cb588SHong Zhang for (k1=0; k1<bs2; k1++){ 704671cb588SHong Zhang *u++ = *rtmp_ptr; 705671cb588SHong Zhang *rtmp_ptr++ = 0.0; 706671cb588SHong Zhang } 707671cb588SHong Zhang } 708671cb588SHong Zhang 709671cb588SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 710671cb588SHong Zhang il[k] = jmin; 711671cb588SHong Zhang i = bj[jmin]; 712671cb588SHong Zhang jl[k] = jl[i]; jl[i] = k; 713671cb588SHong Zhang } 714671cb588SHong Zhang } 715671cb588SHong Zhang 716671cb588SHong Zhang ierr = PetscFree(rtmp);CHKERRQ(ierr); 717671cb588SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 718671cb588SHong Zhang ierr = PetscFree(dk);CHKERRQ(ierr); 719671cb588SHong Zhang ierr = PetscFree(pivots);CHKERRQ(ierr); 720671cb588SHong Zhang 721671cb588SHong Zhang C->factor = FACTOR_CHOLESKY; 722671cb588SHong Zhang C->assembled = PETSC_TRUE; 723671cb588SHong Zhang C->preallocated = PETSC_TRUE; 724b0a32e0cSBarry Smith PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */ 725671cb588SHong Zhang PetscFunctionReturn(0); 726671cb588SHong Zhang } 727671cb588SHong Zhang 72849b5e25fSSatish Balay /* 729fcf159c0SHong Zhang Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP. 730cc0c071aSHong Zhang Version for blocks 2 by 2. 73149b5e25fSSatish Balay */ 7324a2ae208SSatish Balay #undef __FUNCT__ 7334a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2" 7346f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat A,Mat *B) 73549b5e25fSSatish Balay { 73649b5e25fSSatish Balay Mat C = *B; 73791602c66SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 738cc0c071aSHong Zhang IS perm = b->row; 739cc0c071aSHong Zhang int *perm_ptr,ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 740cc0c071aSHong Zhang int *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 741a1723e09SHong Zhang MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 742cc0c071aSHong Zhang MatScalar *u,*diag,*rtmp,*rtmp_ptr; 74349b5e25fSSatish Balay 74449b5e25fSSatish Balay PetscFunctionBegin; 745671cb588SHong Zhang 74691602c66SHong Zhang /* initialization */ 74791602c66SHong Zhang /* il and jl record the first nonzero element in each row of the accessing 74891602c66SHong Zhang window U(0:k, k:mbs-1). 74991602c66SHong Zhang jl: list of rows to be added to uneliminated rows 75091602c66SHong Zhang i>= k: jl(i) is the first row to be added to row i 75191602c66SHong Zhang i< k: jl(i) is the row following row i in some list of rows 75291602c66SHong Zhang jl(i) = mbs indicates the end of a list 75391602c66SHong Zhang il(i): points to the first nonzero element in columns k,...,mbs-1 of 75491602c66SHong Zhang row i of U */ 75582502324SSatish Balay ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 756cc0c071aSHong Zhang ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr); 75782502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr); 75828de702eSHong Zhang jl = il + mbs; 75991602c66SHong Zhang for (i=0; i<mbs; i++) { 7603845f261SHong Zhang jl[i] = mbs; il[0] = 0; 76191602c66SHong Zhang } 76282502324SSatish Balay ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr); 76328de702eSHong Zhang uik = dk + 4; 764cc0c071aSHong Zhang ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr); 765cc0c071aSHong Zhang 766cc0c071aSHong Zhang /* check permutation */ 767cc0c071aSHong Zhang if (!a->permute){ 768cc0c071aSHong Zhang ai = a->i; aj = a->j; aa = a->a; 769cc0c071aSHong Zhang } else { 770cc0c071aSHong Zhang ai = a->inew; aj = a->jnew; 77182502324SSatish Balay ierr = PetscMalloc(4*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); 772cc0c071aSHong Zhang ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); 77382502324SSatish Balay ierr = PetscMalloc(ai[mbs]*sizeof(int),&a2anew);CHKERRQ(ierr); 774cc0c071aSHong Zhang ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr); 775cc0c071aSHong Zhang 776cc0c071aSHong Zhang for (i=0; i<mbs; i++){ 777cc0c071aSHong Zhang jmin = ai[i]; jmax = ai[i+1]; 778cc0c071aSHong Zhang for (j=jmin; j<jmax; j++){ 779cc0c071aSHong Zhang while (a2anew[j] != j){ 780cc0c071aSHong Zhang k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k; 781cc0c071aSHong Zhang for (k1=0; k1<4; k1++){ 782cc0c071aSHong Zhang dk[k1] = aa[k*4+k1]; 783cc0c071aSHong Zhang aa[k*4+k1] = aa[j*4+k1]; 784cc0c071aSHong Zhang aa[j*4+k1] = dk[k1]; 785cc0c071aSHong Zhang } 786cc0c071aSHong Zhang } 787cc0c071aSHong Zhang /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */ 788cc0c071aSHong Zhang if (i > aj[j]){ 789a1723e09SHong Zhang /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */ 790cc0c071aSHong Zhang ap = aa + j*4; /* ptr to the beginning of the block */ 791cc0c071aSHong Zhang dk[1] = ap[1]; /* swap ap[1] and ap[2] */ 792cc0c071aSHong Zhang ap[1] = ap[2]; 793cc0c071aSHong Zhang ap[2] = dk[1]; 794cc0c071aSHong Zhang } 795cc0c071aSHong Zhang } 796cc0c071aSHong Zhang } 797ac355199SBarry Smith ierr = PetscFree(a2anew);CHKERRQ(ierr); 798cc0c071aSHong Zhang } 7993845f261SHong Zhang 80091602c66SHong Zhang /* for each row k */ 80191602c66SHong Zhang for (k = 0; k<mbs; k++){ 80291602c66SHong Zhang 80391602c66SHong Zhang /*initialize k-th row with elements nonzero in row perm(k) of A */ 804cc0c071aSHong Zhang jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1]; 805cc0c071aSHong Zhang ap = aa + jmin*4; 80691602c66SHong Zhang for (j = jmin; j < jmax; j++){ 807cc0c071aSHong Zhang vj = perm_ptr[aj[j]]; /* block col. index */ 808cc0c071aSHong Zhang rtmp_ptr = rtmp + vj*4; 809cc0c071aSHong Zhang for (i=0; i<4; i++) *rtmp_ptr++ = *ap++; 81091602c66SHong Zhang } 81191602c66SHong Zhang 81291602c66SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 813cc0c071aSHong Zhang ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr); 81491602c66SHong Zhang i = jl[k]; /* first row to be added to k_th row */ 81591602c66SHong Zhang 816057f5ba7SHong Zhang while (i < k){ 81791602c66SHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 81891602c66SHong Zhang 8193845f261SHong Zhang /* compute multiplier */ 82091602c66SHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 8213845f261SHong Zhang 8223845f261SHong Zhang /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */ 823cc0c071aSHong Zhang diag = ba + i*4; 824cc0c071aSHong Zhang u = ba + ili*4; 825cc0c071aSHong Zhang uik[0] = -(diag[0]*u[0] + diag[2]*u[1]); 826cc0c071aSHong Zhang uik[1] = -(diag[1]*u[0] + diag[3]*u[1]); 827cc0c071aSHong Zhang uik[2] = -(diag[0]*u[2] + diag[2]*u[3]); 828cc0c071aSHong Zhang uik[3] = -(diag[1]*u[2] + diag[3]*u[3]); 8293845f261SHong Zhang 8303845f261SHong Zhang /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */ 831cc0c071aSHong Zhang dk[0] += uik[0]*u[0] + uik[1]*u[1]; 832cc0c071aSHong Zhang dk[1] += uik[2]*u[0] + uik[3]*u[1]; 833cc0c071aSHong Zhang dk[2] += uik[0]*u[2] + uik[1]*u[3]; 834cc0c071aSHong Zhang dk[3] += uik[2]*u[2] + uik[3]*u[3]; 8353845f261SHong Zhang 8363845f261SHong Zhang /* update -U(i,k): ba[ili] = uik */ 837cc0c071aSHong Zhang ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr); 83891602c66SHong Zhang 83991602c66SHong Zhang /* add multiple of row i to k-th row ... */ 84091602c66SHong Zhang jmin = ili + 1; jmax = bi[i+1]; 84191602c66SHong Zhang if (jmin < jmax){ 8423845f261SHong Zhang for (j=jmin; j<jmax; j++) { 8433845f261SHong Zhang /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */ 844cc0c071aSHong Zhang rtmp_ptr = rtmp + bj[j]*4; 845cc0c071aSHong Zhang u = ba + j*4; 846cc0c071aSHong Zhang rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1]; 847cc0c071aSHong Zhang rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1]; 848cc0c071aSHong Zhang rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3]; 849cc0c071aSHong Zhang rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3]; 8503845f261SHong Zhang } 8513845f261SHong Zhang 85291602c66SHong Zhang /* ... add i to row list for next nonzero entry */ 85391602c66SHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 85491602c66SHong Zhang j = bj[jmin]; 85591602c66SHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 85691602c66SHong Zhang } 857a1723e09SHong Zhang i = nexti; 85891602c66SHong Zhang } 859cc0c071aSHong Zhang 86091602c66SHong Zhang /* save nonzero entries in k-th row of U ... */ 8613845f261SHong Zhang 862cc0c071aSHong Zhang /* invert diagonal block */ 863cc0c071aSHong Zhang diag = ba+k*4; 864cc0c071aSHong Zhang ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr); 8653845f261SHong Zhang ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr); 8663845f261SHong Zhang 86791602c66SHong Zhang jmin = bi[k]; jmax = bi[k+1]; 86891602c66SHong Zhang if (jmin < jmax) { 86991602c66SHong Zhang for (j=jmin; j<jmax; j++){ 870cc0c071aSHong Zhang vj = bj[j]; /* block col. index of U */ 871cc0c071aSHong Zhang u = ba + j*4; 872cc0c071aSHong Zhang rtmp_ptr = rtmp + vj*4; 873cc0c071aSHong Zhang for (k1=0; k1<4; k1++){ 874cc0c071aSHong Zhang *u++ = *rtmp_ptr; 875cc0c071aSHong Zhang *rtmp_ptr++ = 0.0; 876cc0c071aSHong Zhang } 877cc0c071aSHong Zhang } 8783845f261SHong Zhang 87991602c66SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 88091602c66SHong Zhang il[k] = jmin; 88191602c66SHong Zhang i = bj[jmin]; 88291602c66SHong Zhang jl[k] = jl[i]; jl[i] = k; 88391602c66SHong Zhang } 88491602c66SHong Zhang } 8853845f261SHong Zhang 88649b5e25fSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 88791602c66SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 8883845f261SHong Zhang ierr = PetscFree(dk);CHKERRQ(ierr); 88991602c66SHong Zhang if (a->permute) { 89091602c66SHong Zhang ierr = PetscFree(aa);CHKERRQ(ierr); 89191602c66SHong Zhang } 892cc0c071aSHong Zhang ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr); 89391602c66SHong Zhang C->factor = FACTOR_CHOLESKY; 89449b5e25fSSatish Balay C->assembled = PETSC_TRUE; 8955c0bcdfcSHong Zhang C->preallocated = PETSC_TRUE; 896b0a32e0cSBarry Smith PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */ 89749b5e25fSSatish Balay PetscFunctionReturn(0); 89849b5e25fSSatish Balay } 89991602c66SHong Zhang 90049b5e25fSSatish Balay /* 90149b5e25fSSatish Balay Version for when blocks are 2 by 2 Using natural ordering 90249b5e25fSSatish Balay */ 9034a2ae208SSatish Balay #undef __FUNCT__ 9044a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering" 9056f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat A,Mat *B) 90649b5e25fSSatish Balay { 90749b5e25fSSatish Balay Mat C = *B; 908ab27746eSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 909ab27746eSHong Zhang int ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 91071149678SHong Zhang int *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 911ab27746eSHong Zhang MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 912ab27746eSHong Zhang MatScalar *u,*diag,*rtmp,*rtmp_ptr; 91349b5e25fSSatish Balay 91449b5e25fSSatish Balay PetscFunctionBegin; 915671cb588SHong Zhang 916ab27746eSHong Zhang /* initialization */ 917ab27746eSHong Zhang /* il and jl record the first nonzero element in each row of the accessing 918ab27746eSHong Zhang window U(0:k, k:mbs-1). 919ab27746eSHong Zhang jl: list of rows to be added to uneliminated rows 920ab27746eSHong Zhang i>= k: jl(i) is the first row to be added to row i 921ab27746eSHong Zhang i< k: jl(i) is the row following row i in some list of rows 922ab27746eSHong Zhang jl(i) = mbs indicates the end of a list 923ab27746eSHong Zhang il(i): points to the first nonzero element in columns k,...,mbs-1 of 924ab27746eSHong Zhang row i of U */ 92582502324SSatish Balay ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 926ab27746eSHong Zhang ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr); 92782502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr); 92828de702eSHong Zhang jl = il + mbs; 929ab27746eSHong Zhang for (i=0; i<mbs; i++) { 930ab27746eSHong Zhang jl[i] = mbs; il[0] = 0; 93149b5e25fSSatish Balay } 93282502324SSatish Balay ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr); 93328de702eSHong Zhang uik = dk + 4; 934ab27746eSHong Zhang 935ab27746eSHong Zhang ai = a->i; aj = a->j; aa = a->a; 936ab27746eSHong Zhang 937ab27746eSHong Zhang /* for each row k */ 938ab27746eSHong Zhang for (k = 0; k<mbs; k++){ 939ab27746eSHong Zhang 940ab27746eSHong Zhang /*initialize k-th row with elements nonzero in row k of A */ 941ab27746eSHong Zhang jmin = ai[k]; jmax = ai[k+1]; 942ab27746eSHong Zhang ap = aa + jmin*4; 943ab27746eSHong Zhang for (j = jmin; j < jmax; j++){ 944ab27746eSHong Zhang vj = aj[j]; /* block col. index */ 945ab27746eSHong Zhang rtmp_ptr = rtmp + vj*4; 946ab27746eSHong Zhang for (i=0; i<4; i++) *rtmp_ptr++ = *ap++; 94749b5e25fSSatish Balay } 948ab27746eSHong Zhang 949ab27746eSHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 950ab27746eSHong Zhang ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr); 951ab27746eSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 952ab27746eSHong Zhang 953057f5ba7SHong Zhang while (i < k){ 954ab27746eSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 955ab27746eSHong Zhang 956ab27746eSHong Zhang /* compute multiplier */ 957ab27746eSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 958ab27746eSHong Zhang 959ab27746eSHong Zhang /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */ 960ab27746eSHong Zhang diag = ba + i*4; 961ab27746eSHong Zhang u = ba + ili*4; 962ab27746eSHong Zhang uik[0] = -(diag[0]*u[0] + diag[2]*u[1]); 963ab27746eSHong Zhang uik[1] = -(diag[1]*u[0] + diag[3]*u[1]); 964ab27746eSHong Zhang uik[2] = -(diag[0]*u[2] + diag[2]*u[3]); 965ab27746eSHong Zhang uik[3] = -(diag[1]*u[2] + diag[3]*u[3]); 966ab27746eSHong Zhang 967ab27746eSHong Zhang /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */ 968ab27746eSHong Zhang dk[0] += uik[0]*u[0] + uik[1]*u[1]; 969ab27746eSHong Zhang dk[1] += uik[2]*u[0] + uik[3]*u[1]; 970ab27746eSHong Zhang dk[2] += uik[0]*u[2] + uik[1]*u[3]; 971ab27746eSHong Zhang dk[3] += uik[2]*u[2] + uik[3]*u[3]; 972ab27746eSHong Zhang 973ab27746eSHong Zhang /* update -U(i,k): ba[ili] = uik */ 974ab27746eSHong Zhang ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr); 975ab27746eSHong Zhang 976ab27746eSHong Zhang /* add multiple of row i to k-th row ... */ 977ab27746eSHong Zhang jmin = ili + 1; jmax = bi[i+1]; 978ab27746eSHong Zhang if (jmin < jmax){ 979ab27746eSHong Zhang for (j=jmin; j<jmax; j++) { 980ab27746eSHong Zhang /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */ 981ab27746eSHong Zhang rtmp_ptr = rtmp + bj[j]*4; 982ab27746eSHong Zhang u = ba + j*4; 983ab27746eSHong Zhang rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1]; 984ab27746eSHong Zhang rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1]; 985ab27746eSHong Zhang rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3]; 986ab27746eSHong Zhang rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3]; 98749b5e25fSSatish Balay } 988ab27746eSHong Zhang 989ab27746eSHong Zhang /* ... add i to row list for next nonzero entry */ 990ab27746eSHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 991ab27746eSHong Zhang j = bj[jmin]; 992ab27746eSHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 99349b5e25fSSatish Balay } 994ab27746eSHong Zhang i = nexti; 99549b5e25fSSatish Balay } 996ab27746eSHong Zhang 997ab27746eSHong Zhang /* save nonzero entries in k-th row of U ... */ 998ab27746eSHong Zhang 99949b5e25fSSatish Balay /* invert diagonal block */ 1000ab27746eSHong Zhang diag = ba+k*4; 1001ab27746eSHong Zhang ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr); 1002ab27746eSHong Zhang ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr); 1003ab27746eSHong Zhang 1004ab27746eSHong Zhang jmin = bi[k]; jmax = bi[k+1]; 1005ab27746eSHong Zhang if (jmin < jmax) { 1006ab27746eSHong Zhang for (j=jmin; j<jmax; j++){ 1007ab27746eSHong Zhang vj = bj[j]; /* block col. index of U */ 1008ab27746eSHong Zhang u = ba + j*4; 1009ab27746eSHong Zhang rtmp_ptr = rtmp + vj*4; 1010ab27746eSHong Zhang for (k1=0; k1<4; k1++){ 1011ab27746eSHong Zhang *u++ = *rtmp_ptr; 1012ab27746eSHong Zhang *rtmp_ptr++ = 0.0; 1013ab27746eSHong Zhang } 1014ab27746eSHong Zhang } 1015ab27746eSHong Zhang 1016ab27746eSHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 1017ab27746eSHong Zhang il[k] = jmin; 1018ab27746eSHong Zhang i = bj[jmin]; 1019ab27746eSHong Zhang jl[k] = jl[i]; jl[i] = k; 1020ab27746eSHong Zhang } 102149b5e25fSSatish Balay } 102249b5e25fSSatish Balay 102349b5e25fSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 1024ab27746eSHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 1025ab27746eSHong Zhang ierr = PetscFree(dk);CHKERRQ(ierr); 1026ab27746eSHong Zhang 1027ab27746eSHong Zhang C->factor = FACTOR_CHOLESKY; 102849b5e25fSSatish Balay C->assembled = PETSC_TRUE; 10295c0bcdfcSHong Zhang C->preallocated = PETSC_TRUE; 1030b0a32e0cSBarry Smith PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */ 103149b5e25fSSatish Balay PetscFunctionReturn(0); 103249b5e25fSSatish Balay } 103349b5e25fSSatish Balay 103449b5e25fSSatish Balay /* 10355c0bcdfcSHong Zhang Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP. 103691602c66SHong Zhang Version for blocks are 1 by 1. 103749b5e25fSSatish Balay */ 10384a2ae208SSatish Balay #undef __FUNCT__ 10394a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1" 10406f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat A,Mat *B) 104149b5e25fSSatish Balay { 104249b5e25fSSatish Balay Mat C = *B; 104349b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 104449b5e25fSSatish Balay IS ip = b->row; 1045351d0355SHong Zhang int *rip,ierr,i,j,mbs = a->mbs,*bi = b->i,*bj = b->j; 1046cb718733SHong Zhang int *ai,*aj,*r; 1047ab27746eSHong Zhang int k,jmin,jmax,*jl,*il,vj,nexti,ili; 1048066653e3SSatish Balay MatScalar *rtmp; 10492d007305SHong Zhang MatScalar *ba = b->a,*aa,ak; 105049b5e25fSSatish Balay MatScalar dk,uikdi; 105149b5e25fSSatish Balay 105249b5e25fSSatish Balay PetscFunctionBegin; 105349b5e25fSSatish Balay ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 1054cb718733SHong Zhang if (!a->permute){ 10552d007305SHong Zhang ai = a->i; aj = a->j; aa = a->a; 10562d007305SHong Zhang } else { 10572d007305SHong Zhang ai = a->inew; aj = a->jnew; 105882502324SSatish Balay ierr = PetscMalloc(ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); 1059cb718733SHong Zhang ierr = PetscMemcpy(aa,a->a,ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); 106082502324SSatish Balay ierr = PetscMalloc(ai[mbs]*sizeof(int),&r);CHKERRQ(ierr); 10612d007305SHong Zhang ierr= PetscMemcpy(r,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr); 10622d007305SHong Zhang 10632d007305SHong Zhang jmin = ai[0]; jmax = ai[mbs]; 10642d007305SHong Zhang for (j=jmin; j<jmax; j++){ 1065cb718733SHong Zhang while (r[j] != j){ 10662d007305SHong Zhang k = r[j]; r[j] = r[k]; r[k] = k; 10672d007305SHong Zhang ak = aa[k]; aa[k] = aa[j]; aa[j] = ak; 10682d007305SHong Zhang } 10692d007305SHong Zhang } 1070ac355199SBarry Smith ierr = PetscFree(r);CHKERRQ(ierr); 10712d007305SHong Zhang } 10722d007305SHong Zhang 107391602c66SHong Zhang /* initialization */ 107449b5e25fSSatish Balay /* il and jl record the first nonzero element in each row of the accessing 107549b5e25fSSatish Balay window U(0:k, k:mbs-1). 107649b5e25fSSatish Balay jl: list of rows to be added to uneliminated rows 107749b5e25fSSatish Balay i>= k: jl(i) is the first row to be added to row i 107849b5e25fSSatish Balay i< k: jl(i) is the row following row i in some list of rows 107949b5e25fSSatish Balay jl(i) = mbs indicates the end of a list 108049b5e25fSSatish Balay il(i): points to the first nonzero element in columns k,...,mbs-1 of 108149b5e25fSSatish Balay row i of U */ 108282502324SSatish Balay ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 108382502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr); 108428de702eSHong Zhang jl = il + mbs; 108549b5e25fSSatish Balay for (i=0; i<mbs; i++) { 108649b5e25fSSatish Balay rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 108749b5e25fSSatish Balay } 108849b5e25fSSatish Balay 108991602c66SHong Zhang /* for each row k */ 109049b5e25fSSatish Balay for (k = 0; k<mbs; k++){ 109149b5e25fSSatish Balay 109291602c66SHong Zhang /*initialize k-th row with elements nonzero in row perm(k) of A */ 109349b5e25fSSatish Balay jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 1094057f5ba7SHong Zhang 109549b5e25fSSatish Balay for (j = jmin; j < jmax; j++){ 1096351d0355SHong Zhang vj = rip[aj[j]]; 1097ab27746eSHong Zhang rtmp[vj] = aa[j]; 109849b5e25fSSatish Balay } 109949b5e25fSSatish Balay 110091602c66SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 110149b5e25fSSatish Balay dk = rtmp[k]; 110249b5e25fSSatish Balay i = jl[k]; /* first row to be added to k_th row */ 110349b5e25fSSatish Balay 1104057f5ba7SHong Zhang while (i < k){ 110549b5e25fSSatish Balay nexti = jl[i]; /* next row to be added to k_th row */ 110649b5e25fSSatish Balay 110791602c66SHong Zhang /* compute multiplier, update D(k) and U(i,k) */ 110849b5e25fSSatish Balay ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 110949b5e25fSSatish Balay uikdi = - ba[ili]*ba[i]; 111049b5e25fSSatish Balay dk += uikdi*ba[ili]; 1111658e7b3eSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 111249b5e25fSSatish Balay 111391602c66SHong Zhang /* add multiple of row i to k-th row ... */ 111449b5e25fSSatish Balay jmin = ili + 1; jmax = bi[i+1]; 111549b5e25fSSatish Balay if (jmin < jmax){ 111649b5e25fSSatish Balay for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 111791602c66SHong Zhang /* ... add i to row list for next nonzero entry */ 111849b5e25fSSatish Balay il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 111949b5e25fSSatish Balay j = bj[jmin]; 112049b5e25fSSatish Balay jl[i] = jl[j]; jl[j] = i; /* update jl */ 112149b5e25fSSatish Balay } 1122ab27746eSHong Zhang i = nexti; 112349b5e25fSSatish Balay } 112449b5e25fSSatish Balay 112591602c66SHong Zhang /* check for zero pivot and save diagoanl element */ 112649b5e25fSSatish Balay if (dk == 0.0){ 112729bbc08cSBarry Smith SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot"); 1128671cb588SHong Zhang /* 11291223c01bSHong Zhang } else if (PetscRealPart(dk) < 0.0){ 11301223c01bSHong Zhang SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Negative pivot: d[%d] = %g\n",k,dk); 1131671cb588SHong Zhang */ 113249b5e25fSSatish Balay } 113349b5e25fSSatish Balay 113491602c66SHong Zhang /* save nonzero entries in k-th row of U ... */ 1135f3dd7b2dSKris Buschelman ba[k] = 1.0/dk; 113649b5e25fSSatish Balay jmin = bi[k]; jmax = bi[k+1]; 113749b5e25fSSatish Balay if (jmin < jmax) { 113849b5e25fSSatish Balay for (j=jmin; j<jmax; j++){ 1139cc0c071aSHong Zhang vj = bj[j]; ba[j] = rtmp[vj]; rtmp[vj] = 0.0; 114049b5e25fSSatish Balay } 114191602c66SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 114249b5e25fSSatish Balay il[k] = jmin; 114349b5e25fSSatish Balay i = bj[jmin]; 114449b5e25fSSatish Balay jl[k] = jl[i]; jl[i] = k; 114549b5e25fSSatish Balay } 114649b5e25fSSatish Balay } 114749b5e25fSSatish Balay 114849b5e25fSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 114949b5e25fSSatish Balay ierr = PetscFree(il);CHKERRQ(ierr); 1150cb718733SHong Zhang if (a->permute){ 1151cb718733SHong Zhang ierr = PetscFree(aa);CHKERRQ(ierr); 1152cb718733SHong Zhang } 115349b5e25fSSatish Balay 115449b5e25fSSatish Balay ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 11558be8c0c7SHong Zhang C->factor = FACTOR_CHOLESKY; 115649b5e25fSSatish Balay C->assembled = PETSC_TRUE; 11575c0bcdfcSHong Zhang C->preallocated = PETSC_TRUE; 1158b0a32e0cSBarry Smith PetscLogFlops(b->mbs); 115949b5e25fSSatish Balay PetscFunctionReturn(0); 116049b5e25fSSatish Balay } 116149b5e25fSSatish Balay 1162671cb588SHong Zhang /* 1163671cb588SHong Zhang Version for when blocks are 1 by 1 Using natural ordering 1164671cb588SHong Zhang */ 11654a2ae208SSatish Balay #undef __FUNCT__ 11664a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering" 1167671cb588SHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat A,Mat *B) 1168671cb588SHong Zhang { 1169671cb588SHong Zhang Mat C = *B; 1170671cb588SHong Zhang Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data; 1171653a6975SHong Zhang int ierr,i,j,mbs = a->mbs; 1172653a6975SHong Zhang int *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 1173b00f7748SHong Zhang int k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz,ndamp = 0; 1174653a6975SHong Zhang MatScalar *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval; 117521c26570Svictorle PetscReal damping=b->factor_damping, zeropivot=b->factor_zeropivot,shift_amount; 117621c26570Svictorle PetscTruth damp,chshift; 117721c26570Svictorle int nshift=0; 1178653a6975SHong Zhang 1179653a6975SHong Zhang PetscFunctionBegin; 1180653a6975SHong Zhang /* initialization */ 1181653a6975SHong Zhang /* il and jl record the first nonzero element in each row of the accessing 1182653a6975SHong Zhang window U(0:k, k:mbs-1). 1183653a6975SHong Zhang jl: list of rows to be added to uneliminated rows 1184653a6975SHong Zhang i>= k: jl(i) is the first row to be added to row i 1185653a6975SHong Zhang i< k: jl(i) is the row following row i in some list of rows 1186653a6975SHong Zhang jl(i) = mbs indicates the end of a list 1187653a6975SHong Zhang il(i): points to the first nonzero element in U(i,k:mbs-1) 1188653a6975SHong Zhang */ 1189653a6975SHong Zhang ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 1190653a6975SHong Zhang ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr); 1191653a6975SHong Zhang jl = il + mbs; 1192b00f7748SHong Zhang 119321c26570Svictorle shift_amount = 0; 1194b00f7748SHong Zhang do { 1195b00f7748SHong Zhang damp = PETSC_FALSE; 119621c26570Svictorle chshift = PETSC_FALSE; 1197653a6975SHong Zhang for (i=0; i<mbs; i++) { 1198653a6975SHong Zhang rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 1199653a6975SHong Zhang } 1200653a6975SHong Zhang 1201b00f7748SHong Zhang for (k = 0; k<mbs; k++){ /* row k */ 1202653a6975SHong Zhang /*initialize k-th row with elements nonzero in row perm(k) of A */ 1203653a6975SHong Zhang nz = ai[k+1] - ai[k]; 1204653a6975SHong Zhang acol = aj + ai[k]; 1205653a6975SHong Zhang aval = aa + ai[k]; 1206653a6975SHong Zhang bval = ba + bi[k]; 1207653a6975SHong Zhang while (nz -- ){ 1208653a6975SHong Zhang rtmp[*acol++] = *aval++; 1209653a6975SHong Zhang *bval++ = 0.0; /* for in-place factorization */ 1210653a6975SHong Zhang } 1211b00f7748SHong Zhang /* damp the diagonal of the matrix */ 121221c26570Svictorle if (ndamp||nshift) rtmp[k] += damping+shift_amount; 1213653a6975SHong Zhang 1214653a6975SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 1215653a6975SHong Zhang dk = rtmp[k]; 1216653a6975SHong Zhang i = jl[k]; /* first row to be added to k_th row */ 1217653a6975SHong Zhang 1218653a6975SHong Zhang while (i < k){ 1219653a6975SHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 1220653a6975SHong Zhang 1221653a6975SHong Zhang /* compute multiplier, update D(k) and U(i,k) */ 1222653a6975SHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 1223653a6975SHong Zhang uikdi = - ba[ili]*ba[bi[i]]; 1224653a6975SHong Zhang dk += uikdi*ba[ili]; 1225653a6975SHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 1226653a6975SHong Zhang 1227653a6975SHong Zhang /* add multiple of row i to k-th row ... */ 1228653a6975SHong Zhang jmin = ili + 1; 1229653a6975SHong Zhang nz = bi[i+1] - jmin; 1230653a6975SHong Zhang if (nz > 0){ 1231653a6975SHong Zhang bcol = bj + jmin; 1232653a6975SHong Zhang bval = ba + jmin; 1233653a6975SHong Zhang while (nz --) rtmp[*bcol++] += uikdi*(*bval++); 1234653a6975SHong Zhang /* ... add i to row list for next nonzero entry */ 1235653a6975SHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 1236653a6975SHong Zhang j = bj[jmin]; 1237653a6975SHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 1238653a6975SHong Zhang } 1239653a6975SHong Zhang i = nexti; 1240653a6975SHong Zhang } 1241653a6975SHong Zhang 124221c26570Svictorle /* check for zero pivot and save diagonal element */ 124321c26570Svictorle if (PetscRealPart(dk) < zeropivot && b->factor_shift){ 124421c26570Svictorle PetscReal rs = -PetscRealPart(dk); 124521c26570Svictorle jmin = bi[k]+1; 124621c26570Svictorle nz = bi[k+1] - jmin; 124721c26570Svictorle if (nz){ 124821c26570Svictorle bcol = bj + jmin; 124921c26570Svictorle bval = ba + jmin; 125021c26570Svictorle while (nz--){ 125121c26570Svictorle rs += PetscAbsScalar(rtmp[*bcol++]); 125221c26570Svictorle } 125321c26570Svictorle } 125421c26570Svictorle shift_amount = rs; 125521c26570Svictorle chshift = PETSC_TRUE; 125621c26570Svictorle nshift++; 125721c26570Svictorle break; 125821c26570Svictorle } 1259ffa0b812SHong Zhang if (PetscRealPart(dk) < zeropivot){ 1260ffa0b812SHong Zhang if (damping == (PetscReal) PETSC_DECIDE) damping = -PetscRealPart(dk)/(k+1); 1261b00f7748SHong Zhang if (damping > 0.0) { 1262b00f7748SHong Zhang if (ndamp) damping *= 2.0; 1263b00f7748SHong Zhang damp = PETSC_TRUE; 1264b00f7748SHong Zhang ndamp++; 1265b00f7748SHong Zhang break; 1266481c2c92SHong Zhang } else if (PetscAbsScalar(dk) < zeropivot){ 1267ffa0b812SHong Zhang SETERRQ3(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %d value %g tolerance %g",k,PetscRealPart(dk),zeropivot); 1268481c2c92SHong Zhang } else { 1269ffa0b812SHong Zhang PetscLogInfo((PetscObject)A,"Negative pivot %g in row %d of Cholesky factorization\n",PetscRealPart(dk),k); 1270b00f7748SHong Zhang } 1271064503c5SHong Zhang } 1272653a6975SHong Zhang 1273653a6975SHong Zhang /* save nonzero entries in k-th row of U ... */ 127482799104SHong Zhang /* printf("%d, dk: %g, 1/dk: %g\n",k,dk,1/dk); */ 1275653a6975SHong Zhang ba[bi[k]] = 1.0/dk; 1276653a6975SHong Zhang jmin = bi[k]+1; 1277653a6975SHong Zhang nz = bi[k+1] - jmin; 1278653a6975SHong Zhang if (nz){ 1279653a6975SHong Zhang bcol = bj + jmin; 1280653a6975SHong Zhang bval = ba + jmin; 1281653a6975SHong Zhang while (nz--){ 1282653a6975SHong Zhang *bval++ = rtmp[*bcol]; 1283653a6975SHong Zhang rtmp[*bcol++] = 0.0; 1284653a6975SHong Zhang } 1285653a6975SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 1286653a6975SHong Zhang il[k] = jmin; 1287653a6975SHong Zhang i = bj[jmin]; 1288653a6975SHong Zhang jl[k] = jl[i]; jl[i] = k; 1289653a6975SHong Zhang } 1290b00f7748SHong Zhang } /* end of for (k = 0; k<mbs; k++) */ 129121c26570Svictorle } while (damp||chshift); 1292653a6975SHong Zhang ierr = PetscFree(rtmp);CHKERRQ(ierr); 1293653a6975SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 1294653a6975SHong Zhang 1295653a6975SHong Zhang C->factor = FACTOR_CHOLESKY; 1296653a6975SHong Zhang C->assembled = PETSC_TRUE; 1297653a6975SHong Zhang C->preallocated = PETSC_TRUE; 1298653a6975SHong Zhang PetscLogFlops(b->mbs); 1299b00f7748SHong Zhang if (ndamp) { 1300b45a75daSHong Zhang PetscLogInfo(0,"MatCholeskyFactorNumerical_SeqSBAIJ_1_NaturalOrdering: number of damping tries %d damping value %g\n",ndamp,damping); 1301b00f7748SHong Zhang } 1302*fb10cecfSBarry Smith if (nshift) { 1303*fb10cecfSBarry Smith PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering diagonal shifted %d shifts\n",nshift); 1304*fb10cecfSBarry Smith } 1305*fb10cecfSBarry Smith 1306653a6975SHong Zhang PetscFunctionReturn(0); 1307653a6975SHong Zhang } 1308653a6975SHong Zhang 1309653a6975SHong Zhang #undef __FUNCT__ 13104a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ" 131115e8a5b3SHong Zhang int MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,MatFactorInfo *info) 131249b5e25fSSatish Balay { 13134445ddedSHong Zhang int ierr; 131449b5e25fSSatish Balay Mat C; 131549b5e25fSSatish Balay 131649b5e25fSSatish Balay PetscFunctionBegin; 131715e8a5b3SHong Zhang ierr = MatCholeskyFactorSymbolic(A,perm,info,&C);CHKERRQ(ierr); 1318a4ada70bSHong Zhang ierr = MatCholeskyFactorNumeric(A,&C);CHKERRQ(ierr); 13194445ddedSHong Zhang ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 132049b5e25fSSatish Balay PetscFunctionReturn(0); 132149b5e25fSSatish Balay } 132249b5e25fSSatish Balay 132349b5e25fSSatish Balay 1324