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/inline/ilu.h" 649a6740bSHong Zhang #include "include/petscis.h" 749b5e25fSSatish Balay 88dc9108eSHong Zhang #if !defined(PETSC_USE_COMPLEX) 95f9f512dSHong Zhang /* 105f9f512dSHong Zhang input: 11c037c3f7SHong Zhang F -- numeric factor 125f9f512dSHong Zhang output: 13c037c3f7SHong Zhang nneg, nzero, npos: matrix inertia 145f9f512dSHong Zhang */ 155f9f512dSHong Zhang 165f9f512dSHong Zhang #undef __FUNCT__ 175f9f512dSHong Zhang #define __FUNCT__ "MatGetInertia_SeqSBAIJ" 18c037c3f7SHong Zhang int MatGetInertia_SeqSBAIJ(Mat F,int *nneig,int *nzero,int *npos) 195f9f512dSHong Zhang { 20638f5ce0SDinesh Kaushik Mat_SeqSBAIJ *fact_ptr = (Mat_SeqSBAIJ*)F->data; 213e0d88b5SBarry Smith PetscScalar *dd = fact_ptr->a; 22eeeff2ecSHong Zhang int mbs=fact_ptr->mbs,bs=fact_ptr->bs,i,nneig_tmp,npos_tmp, 23eeeff2ecSHong Zhang *fi = fact_ptr->i; 245f9f512dSHong Zhang 255f9f512dSHong Zhang PetscFunctionBegin; 26eeeff2ecSHong Zhang if (bs != 1) SETERRQ1(PETSC_ERR_SUP,"No support for bs: %d >1 yet",bs); 27eeeff2ecSHong Zhang nneig_tmp = 0; npos_tmp = 0; 28eeeff2ecSHong Zhang for (i=0; i<mbs; i++){ 29eeeff2ecSHong Zhang if (PetscRealPart(dd[*fi]) > 0.0){ 30eeeff2ecSHong Zhang npos_tmp++; 31eeeff2ecSHong Zhang } else if (PetscRealPart(dd[*fi]) < 0.0){ 32eeeff2ecSHong Zhang nneig_tmp++; 335f9f512dSHong Zhang } 34eeeff2ecSHong Zhang fi++; 353e0d88b5SBarry Smith } 36eeeff2ecSHong Zhang if (nneig) *nneig = nneig_tmp; 37eeeff2ecSHong Zhang if (npos) *npos = npos_tmp; 38eeeff2ecSHong Zhang if (nzero) *nzero = mbs - nneig_tmp - npos_tmp; 39eeeff2ecSHong Zhang 405f9f512dSHong Zhang PetscFunctionReturn(0); 415f9f512dSHong Zhang } 428dc9108eSHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */ 435f9f512dSHong Zhang 445f9f512dSHong Zhang /* Using Modified Sparse Row (MSR) storage. 455f9f512dSHong Zhang See page 85, "Iterative Methods ..." by Saad. */ 465f9f512dSHong Zhang /* 475f9f512dSHong Zhang Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP. 485f9f512dSHong Zhang */ 492d9a3abdSHong Zhang /* Use Modified Sparse Row storage for u and ju, see Saad pp.85 */ 504a2ae208SSatish Balay #undef __FUNCT__ 514a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ" 5215e8a5b3SHong Zhang int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *B) 5349b5e25fSSatish Balay { 5449b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b; 55cb718733SHong Zhang int *rip,ierr,i,mbs = a->mbs,*ai,*aj; 56066653e3SSatish Balay int *jutmp,bs = a->bs,bs2=a->bs2; 5787fe1012SHong Zhang int m,realloc = 0,prow; 58d60e5362SHong Zhang int *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd; 59064503c5SHong Zhang int *il,ili,nextprow; 6015e8a5b3SHong Zhang PetscReal f = info->fill; 61671cb588SHong Zhang PetscTruth perm_identity; 6249b5e25fSSatish Balay 6349b5e25fSSatish Balay PetscFunctionBegin; 64cb718733SHong Zhang /* check whether perm is the identity mapping */ 65671cb588SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 66064503c5SHong Zhang 67064503c5SHong Zhang /* -- inplace factorization, i.e., use sbaij for *B -- */ 68064503c5SHong Zhang if (perm_identity && bs==1 ){ 69064503c5SHong Zhang if (!perm_identity) a->permute = PETSC_TRUE; 70064503c5SHong Zhang 71064503c5SHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 72064503c5SHong Zhang 73064503c5SHong Zhang if (perm_identity){ /* without permutation */ 74064503c5SHong Zhang ai = a->i; aj = a->j; 75064503c5SHong Zhang } else { /* non-trivial permutation */ 76064503c5SHong Zhang ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr); 77064503c5SHong Zhang ai = a->inew; aj = a->jnew; 78064503c5SHong Zhang } 79064503c5SHong Zhang 80064503c5SHong Zhang /* initialization */ 81064503c5SHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(int),&iu);CHKERRQ(ierr); 82064503c5SHong Zhang umax = (int)(f*ai[mbs] + 1); 83064503c5SHong Zhang ierr = PetscMalloc(umax*sizeof(int),&ju);CHKERRQ(ierr); 84064503c5SHong Zhang iu[0] = 0; 85064503c5SHong Zhang juidx = 0; /* index for ju */ 86064503c5SHong Zhang ierr = PetscMalloc((3*mbs+1)*sizeof(int),&jl);CHKERRQ(ierr); /* linked list for getting pivot row */ 87064503c5SHong Zhang q = jl + mbs; /* linked list for col index of active row */ 88064503c5SHong Zhang il = q + mbs; 89064503c5SHong Zhang for (i=0; i<mbs; i++){ 90064503c5SHong Zhang jl[i] = mbs; 91064503c5SHong Zhang q[i] = 0; 92064503c5SHong Zhang il[i] = 0; 93064503c5SHong Zhang } 94064503c5SHong Zhang 95064503c5SHong Zhang /* for each row k */ 96064503c5SHong Zhang for (k=0; k<mbs; k++){ 97064503c5SHong Zhang nzk = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */ 98064503c5SHong Zhang q[k] = mbs; 99064503c5SHong Zhang /* initialize nonzero structure of k-th row to row rip[k] of A */ 100064503c5SHong Zhang jmin = ai[rip[k]] +1; /* exclude diag[k] */ 101064503c5SHong Zhang jmax = ai[rip[k]+1]; 102064503c5SHong Zhang for (j=jmin; j<jmax; j++){ 103064503c5SHong Zhang vj = rip[aj[j]]; /* col. value */ 104064503c5SHong Zhang if(vj > k){ 105064503c5SHong Zhang qm = k; 106064503c5SHong Zhang do { 107064503c5SHong Zhang m = qm; qm = q[m]; 108064503c5SHong Zhang } while(qm < vj); 109064503c5SHong Zhang if (qm == vj) { 110064503c5SHong Zhang SETERRQ(1," error: duplicate entry in A\n"); 111064503c5SHong Zhang } 112064503c5SHong Zhang nzk++; 113064503c5SHong Zhang q[m] = vj; 114064503c5SHong Zhang q[vj] = qm; 115064503c5SHong Zhang } /* if(vj > k) */ 116064503c5SHong Zhang } /* for (j=jmin; j<jmax; j++) */ 117064503c5SHong Zhang 118064503c5SHong Zhang /* modify nonzero structure of k-th row by computing fill-in 119064503c5SHong Zhang for each row i to be merged in */ 120064503c5SHong Zhang prow = k; 121064503c5SHong Zhang prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */ 122064503c5SHong Zhang 123064503c5SHong Zhang while (prow < k){ 124064503c5SHong Zhang nextprow = jl[prow]; 125064503c5SHong Zhang 126064503c5SHong Zhang /* merge row prow into k-th row */ 127064503c5SHong Zhang ili = il[prow]; 128064503c5SHong Zhang jmin = ili + 1; /* points to 2nd nzero entry in U(prow,k:mbs-1) */ 129064503c5SHong Zhang jmax = iu[prow+1]; 130064503c5SHong Zhang qm = k; 131064503c5SHong Zhang for (j=jmin; j<jmax; j++){ 132064503c5SHong Zhang vj = ju[j]; 133064503c5SHong Zhang do { 134064503c5SHong Zhang m = qm; qm = q[m]; 135064503c5SHong Zhang } while (qm < vj); 136064503c5SHong Zhang if (qm != vj){ /* a fill */ 137064503c5SHong Zhang nzk++; q[m] = vj; q[vj] = qm; qm = vj; 138064503c5SHong Zhang } 139064503c5SHong Zhang } /* end of for (j=jmin; j<jmax; j++) */ 140064503c5SHong Zhang if (jmin < jmax){ 141064503c5SHong Zhang il[prow] = jmin; 142064503c5SHong Zhang j = ju[jmin]; 143064503c5SHong Zhang jl[prow] = jl[j]; jl[j] = prow; /* update jl */ 144064503c5SHong Zhang } 145064503c5SHong Zhang prow = nextprow; 146064503c5SHong Zhang } 147064503c5SHong Zhang 148064503c5SHong Zhang /* update il and jl */ 149064503c5SHong Zhang if (nzk > 0){ 150064503c5SHong Zhang i = q[k]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */ 151064503c5SHong Zhang jl[k] = jl[i]; jl[i] = k; 152064503c5SHong Zhang il[k] = iu[k] + 1; 153064503c5SHong Zhang } 154064503c5SHong Zhang iu[k+1] = iu[k] + nzk + 1; /* include diag[k] */ 155064503c5SHong Zhang 156064503c5SHong Zhang /* allocate more space to ju if needed */ 157064503c5SHong Zhang if (iu[k+1] > umax) { 158064503c5SHong Zhang /* estimate how much additional space we will need */ 159064503c5SHong Zhang /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 160064503c5SHong Zhang /* just double the memory each time */ 161064503c5SHong Zhang maxadd = umax; 162064503c5SHong Zhang if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2; 163064503c5SHong Zhang umax += maxadd; 164064503c5SHong Zhang 165064503c5SHong Zhang /* allocate a longer ju */ 166064503c5SHong Zhang ierr = PetscMalloc(umax*sizeof(int),&jutmp);CHKERRQ(ierr); 167064503c5SHong Zhang ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(int));CHKERRQ(ierr); 168064503c5SHong Zhang ierr = PetscFree(ju);CHKERRQ(ierr); 169064503c5SHong Zhang ju = jutmp; 170064503c5SHong Zhang realloc++; /* count how many times we realloc */ 171064503c5SHong Zhang } 172064503c5SHong Zhang 173064503c5SHong Zhang /* save nonzero structure of k-th row in ju */ 174064503c5SHong Zhang ju[juidx++] = k; /* diag[k] */ 175064503c5SHong Zhang i = k; 176064503c5SHong Zhang while (nzk --) { 177064503c5SHong Zhang i = q[i]; 178064503c5SHong Zhang ju[juidx++] = i; 179064503c5SHong Zhang } 180064503c5SHong Zhang } 181064503c5SHong Zhang 182064503c5SHong Zhang if (ai[mbs] != 0) { 183064503c5SHong Zhang PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]); 184064503c5SHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af); 185064503c5SHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Run with -pc_cholesky_fill %g or use \n",af); 186064503c5SHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:PCCholeskySetFill(pc,%g);\n",af); 187064503c5SHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:for best performance.\n"); 188064503c5SHong Zhang } else { 189064503c5SHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Empty matrix.\n"); 190064503c5SHong Zhang } 191064503c5SHong Zhang 192064503c5SHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 193064503c5SHong Zhang /* ierr = PetscFree(q);CHKERRQ(ierr); */ 194064503c5SHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 195064503c5SHong Zhang 196064503c5SHong Zhang /* put together the new matrix */ 197*be5d1d56SKris Buschelman ierr = MatCreate(A->comm,bs*mbs,bs*mbs,bs*mbs,bs*mbs,B);CHKERRQ(ierr); 198*be5d1d56SKris Buschelman ierr = MatSetType(*B,A->type_name);CHKERRQ(ierr); 199*be5d1d56SKris Buschelman ierr = MatSeqSBAIJSetPreallocation(*B,bs,0,PETSC_NULL);CHKERRQ(ierr); 200*be5d1d56SKris Buschelman 201064503c5SHong Zhang /* PetscLogObjectParent(*B,iperm); */ 202064503c5SHong Zhang b = (Mat_SeqSBAIJ*)(*B)->data; 203064503c5SHong Zhang ierr = PetscFree(b->imax);CHKERRQ(ierr); 204064503c5SHong Zhang b->singlemalloc = PETSC_FALSE; 205064503c5SHong Zhang /* the next line frees the default space generated by the Create() */ 206064503c5SHong Zhang ierr = PetscFree(b->a);CHKERRQ(ierr); 207064503c5SHong Zhang ierr = PetscFree(b->ilen);CHKERRQ(ierr); 208064503c5SHong Zhang ierr = PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr); 209064503c5SHong Zhang b->j = ju; 210064503c5SHong Zhang b->i = iu; 211064503c5SHong Zhang b->diag = 0; 212064503c5SHong Zhang b->ilen = 0; 213064503c5SHong Zhang b->imax = 0; 214064503c5SHong Zhang b->row = perm; 21515e8a5b3SHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 216064503c5SHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 217064503c5SHong Zhang b->icol = perm; 218064503c5SHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 219064503c5SHong Zhang ierr = PetscMalloc((bs*mbs+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 220064503c5SHong Zhang /* In b structure: Free imax, ilen, old a, old j. 221064503c5SHong Zhang Allocate idnew, solve_work, new a, new j */ 222064503c5SHong Zhang PetscLogObjectMemory(*B,(iu[mbs]-mbs)*(sizeof(int)+sizeof(MatScalar))); 2236c6c5352SBarry Smith b->maxnz = b->nz = iu[mbs]; 224064503c5SHong Zhang 225064503c5SHong Zhang (*B)->factor = FACTOR_CHOLESKY; 226064503c5SHong Zhang (*B)->info.factor_mallocs = realloc; 227064503c5SHong Zhang (*B)->info.fill_ratio_given = f; 228064503c5SHong Zhang if (ai[mbs] != 0) { 229064503c5SHong Zhang (*B)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]); 230064503c5SHong Zhang } else { 231064503c5SHong Zhang (*B)->info.fill_ratio_needed = 0.0; 232064503c5SHong Zhang } 233064503c5SHong Zhang 234064503c5SHong Zhang 235b45a75daSHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering; 236b45a75daSHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 237064503c5SHong Zhang PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=1\n"); 238064503c5SHong Zhang 239064503c5SHong Zhang PetscFunctionReturn(0); 240064503c5SHong Zhang } 241064503c5SHong Zhang /* ----------- end of new code --------------------*/ 242064503c5SHong Zhang 243064503c5SHong Zhang 244671cb588SHong Zhang if (!perm_identity) a->permute = PETSC_TRUE; 2452d007305SHong Zhang 246671cb588SHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 247671cb588SHong Zhang 248671cb588SHong Zhang if (perm_identity){ /* without permutation */ 2492d007305SHong Zhang ai = a->i; aj = a->j; 2502d007305SHong Zhang } else { /* non-trivial permutation */ 251323b833fSBarry Smith ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr); 2522d007305SHong Zhang ai = a->inew; aj = a->jnew; 2532d007305SHong Zhang } 2542d007305SHong Zhang 25549b5e25fSSatish Balay /* initialization */ 256b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&iu);CHKERRQ(ierr); 25749b5e25fSSatish Balay umax = (int)(f*ai[mbs] + 1); umax += mbs + 1; 25882502324SSatish Balay ierr = PetscMalloc(umax*sizeof(int),&ju);CHKERRQ(ierr); 25949b5e25fSSatish Balay iu[0] = mbs+1; 26087fe1012SHong Zhang juidx = mbs + 1; /* index for ju */ 26187fe1012SHong Zhang ierr = PetscMalloc(2*mbs*sizeof(int),&jl);CHKERRQ(ierr); /* linked list for pivot row */ 26287fe1012SHong Zhang q = jl + mbs; /* linked list for col index */ 26349b5e25fSSatish Balay for (i=0; i<mbs; i++){ 26487fe1012SHong Zhang jl[i] = mbs; 26587fe1012SHong Zhang q[i] = 0; 26649b5e25fSSatish Balay } 26749b5e25fSSatish Balay 26849b5e25fSSatish Balay /* for each row k */ 26949b5e25fSSatish Balay for (k=0; k<mbs; k++){ 270064503c5SHong Zhang for (i=0; i<mbs; i++) q[i] = 0; /* to be removed! */ 27149b5e25fSSatish Balay nzk = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */ 27249b5e25fSSatish Balay q[k] = mbs; 27349b5e25fSSatish Balay /* initialize nonzero structure of k-th row to row rip[k] of A */ 274064503c5SHong Zhang jmin = ai[rip[k]] +1; /* exclude diag[k] */ 27549b5e25fSSatish Balay jmax = ai[rip[k]+1]; 27649b5e25fSSatish Balay for (j=jmin; j<jmax; j++){ 277cb718733SHong Zhang vj = rip[aj[j]]; /* col. value */ 27849b5e25fSSatish Balay if(vj > k){ 27949b5e25fSSatish Balay qm = k; 28049b5e25fSSatish Balay do { 28149b5e25fSSatish Balay m = qm; qm = q[m]; 28249b5e25fSSatish Balay } while(qm < vj); 28349b5e25fSSatish Balay if (qm == vj) { 2843078e140SBarry Smith SETERRQ(1," error: duplicate entry in A\n"); 28549b5e25fSSatish Balay } 28649b5e25fSSatish Balay nzk++; 28749b5e25fSSatish Balay q[m] = vj; 28849b5e25fSSatish Balay q[vj] = qm; 28949b5e25fSSatish Balay } /* if(vj > k) */ 29049b5e25fSSatish Balay } /* for (j=jmin; j<jmax; j++) */ 29149b5e25fSSatish Balay 29249b5e25fSSatish Balay /* modify nonzero structure of k-th row by computing fill-in 29349b5e25fSSatish Balay for each row i to be merged in */ 29487fe1012SHong Zhang prow = k; 29587fe1012SHong Zhang prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */ 29687fe1012SHong Zhang 29787fe1012SHong Zhang while (prow < k){ 29887fe1012SHong Zhang /* merge row prow into k-th row */ 29987fe1012SHong Zhang jmin = iu[prow] + 1; jmax = iu[prow+1]; 30049b5e25fSSatish Balay qm = k; 30187fe1012SHong Zhang for (j=jmin; j<jmax; j++){ 30249b5e25fSSatish Balay vj = ju[j]; 30349b5e25fSSatish Balay do { 30449b5e25fSSatish Balay m = qm; qm = q[m]; 30549b5e25fSSatish Balay } while (qm < vj); 30649b5e25fSSatish Balay if (qm != vj){ 30749b5e25fSSatish Balay nzk++; q[m] = vj; q[vj] = qm; qm = vj; 30849b5e25fSSatish Balay } 30949b5e25fSSatish Balay } 31087fe1012SHong Zhang prow = jl[prow]; /* next pivot row */ 31149b5e25fSSatish Balay } 31249b5e25fSSatish Balay 31349b5e25fSSatish Balay /* add k to row list for first nonzero element in k-th row */ 31449b5e25fSSatish Balay if (nzk > 0){ 31549b5e25fSSatish Balay i = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */ 31649b5e25fSSatish Balay jl[k] = jl[i]; jl[i] = k; 31749b5e25fSSatish Balay } 31887fe1012SHong Zhang iu[k+1] = iu[k] + nzk; 31949b5e25fSSatish Balay 32049b5e25fSSatish Balay /* allocate more space to ju if needed */ 3213078e140SBarry Smith if (iu[k+1] > umax) { 32249b5e25fSSatish Balay /* estimate how much additional space we will need */ 32349b5e25fSSatish Balay /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 32449b5e25fSSatish Balay /* just double the memory each time */ 32549b5e25fSSatish Balay maxadd = umax; 32649b5e25fSSatish Balay if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2; 32749b5e25fSSatish Balay umax += maxadd; 32849b5e25fSSatish Balay 3299f9cb213SHong Zhang /* allocate a longer ju */ 33082502324SSatish Balay ierr = PetscMalloc(umax*sizeof(int),&jutmp);CHKERRQ(ierr); 33149b5e25fSSatish Balay ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(int));CHKERRQ(ierr); 33249b5e25fSSatish Balay ierr = PetscFree(ju);CHKERRQ(ierr); 3339f9cb213SHong Zhang ju = jutmp; 33449b5e25fSSatish Balay realloc++; /* count how many times we realloc */ 33549b5e25fSSatish Balay } 33649b5e25fSSatish Balay 33749b5e25fSSatish Balay /* save nonzero structure of k-th row in ju */ 33849b5e25fSSatish Balay i=k; 33987fe1012SHong Zhang while (nzk --) { 34049b5e25fSSatish Balay i = q[i]; 34187fe1012SHong Zhang ju[juidx++] = i; 34249b5e25fSSatish Balay } 34377f73120SHong Zhang } 34449b5e25fSSatish Balay 34549b5e25fSSatish Balay if (ai[mbs] != 0) { 34649b5e25fSSatish Balay PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]); 347b0a32e0cSBarry Smith PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af); 3483078e140SBarry Smith PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Run with -pc_cholesky_fill %g or use \n",af); 3493078e140SBarry Smith PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:PCCholeskySetFill(pc,%g);\n",af); 350b0a32e0cSBarry Smith PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:for best performance.\n"); 35149b5e25fSSatish Balay } else { 352b0a32e0cSBarry Smith PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Empty matrix.\n"); 35349b5e25fSSatish Balay } 35449b5e25fSSatish Balay 35577f73120SHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 35687fe1012SHong Zhang /* ierr = PetscFree(q);CHKERRQ(ierr); */ 35749b5e25fSSatish Balay ierr = PetscFree(jl);CHKERRQ(ierr); 35849b5e25fSSatish Balay 35949b5e25fSSatish Balay /* put together the new matrix */ 360*be5d1d56SKris Buschelman ierr = MatCreate(A->comm,bs*mbs,bs*mbs,bs*mbs,bs*mbs,B);CHKERRQ(ierr); 361*be5d1d56SKris Buschelman ierr = MatSetType(*B,A->type_name);CHKERRQ(ierr); 362*be5d1d56SKris Buschelman ierr = MatSeqSBAIJSetPreallocation(*B,bs,0,PETSC_NULL);CHKERRQ(ierr); 363*be5d1d56SKris Buschelman 364b0a32e0cSBarry Smith /* PetscLogObjectParent(*B,iperm); */ 36549b5e25fSSatish Balay b = (Mat_SeqSBAIJ*)(*B)->data; 36649b5e25fSSatish Balay ierr = PetscFree(b->imax);CHKERRQ(ierr); 36749b5e25fSSatish Balay b->singlemalloc = PETSC_FALSE; 36849b5e25fSSatish Balay /* the next line frees the default space generated by the Create() */ 36949b5e25fSSatish Balay ierr = PetscFree(b->a);CHKERRQ(ierr); 37049b5e25fSSatish Balay ierr = PetscFree(b->ilen);CHKERRQ(ierr); 37182502324SSatish Balay ierr = PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr); 37249b5e25fSSatish Balay b->j = ju; 37349b5e25fSSatish Balay b->i = iu; 37449b5e25fSSatish Balay b->diag = 0; 37549b5e25fSSatish Balay b->ilen = 0; 37649b5e25fSSatish Balay b->imax = 0; 37777f73120SHong Zhang b->row = perm; 37815e8a5b3SHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 37977f73120SHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 380cb718733SHong Zhang b->icol = perm; 381cb718733SHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 38287828ca2SBarry Smith ierr = PetscMalloc((bs*mbs+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 38349b5e25fSSatish Balay /* In b structure: Free imax, ilen, old a, old j. 38449b5e25fSSatish Balay Allocate idnew, solve_work, new a, new j */ 385b0a32e0cSBarry Smith PetscLogObjectMemory(*B,(iu[mbs]-mbs)*(sizeof(int)+sizeof(MatScalar))); 3866c6c5352SBarry Smith b->maxnz = b->nz = iu[mbs]; 38749b5e25fSSatish Balay 3885c0bcdfcSHong Zhang (*B)->factor = FACTOR_CHOLESKY; 38949b5e25fSSatish Balay (*B)->info.factor_mallocs = realloc; 39049b5e25fSSatish Balay (*B)->info.fill_ratio_given = f; 39149b5e25fSSatish Balay if (ai[mbs] != 0) { 39249b5e25fSSatish Balay (*B)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]); 39349b5e25fSSatish Balay } else { 39449b5e25fSSatish Balay (*B)->info.fill_ratio_needed = 0.0; 39549b5e25fSSatish Balay } 3960aa0ce06SHong Zhang 397671cb588SHong Zhang if (perm_identity){ 398671cb588SHong Zhang switch (bs) { 399671cb588SHong Zhang case 1: 400671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering; 401671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 4024d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=1\n"); 403671cb588SHong Zhang break; 404671cb588SHong Zhang case 2: 405671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 406671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 4074d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 408671cb588SHong Zhang break; 409671cb588SHong Zhang case 3: 410671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 411671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 4124d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:sing special in-place natural ordering factor and solve BS=3\n"); 413671cb588SHong Zhang break; 414671cb588SHong Zhang case 4: 415671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 416671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 4174d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 418671cb588SHong Zhang break; 419671cb588SHong Zhang case 5: 420671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 421671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 4224d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 423671cb588SHong Zhang break; 424671cb588SHong Zhang case 6: 425671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 426671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 4274d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 428671cb588SHong Zhang break; 429671cb588SHong Zhang case 7: 430671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 431671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 4324d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 433671cb588SHong Zhang break; 434671cb588SHong Zhang default: 435671cb588SHong Zhang (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering; 436671cb588SHong Zhang (*B)->ops->solve = MatSolve_SeqSBAIJ_N_NaturalOrdering; 4374d101231SSatish Balay PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS>7\n"); 438671cb588SHong Zhang break; 439671cb588SHong Zhang } 440671cb588SHong Zhang } 441671cb588SHong Zhang 44249b5e25fSSatish Balay PetscFunctionReturn(0); 44349b5e25fSSatish Balay } 44449b5e25fSSatish Balay 44535d8aa7fSBarry Smith 4464a2ae208SSatish Balay #undef __FUNCT__ 4474a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N" 4486f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat A,Mat *B) 44949b5e25fSSatish Balay { 45049b5e25fSSatish Balay Mat C = *B; 4514c16a6a6SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 4524c16a6a6SHong Zhang IS perm = b->row; 4534c16a6a6SHong Zhang int *perm_ptr,ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 4544c16a6a6SHong Zhang int *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 4554c16a6a6SHong Zhang int bs=a->bs,bs2 = a->bs2; 4564c16a6a6SHong Zhang MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 4574c16a6a6SHong Zhang MatScalar *u,*diag,*rtmp,*rtmp_ptr; 45828de702eSHong Zhang MatScalar *work; 4594c16a6a6SHong Zhang int *pivots; 4604c16a6a6SHong Zhang 4614c16a6a6SHong Zhang PetscFunctionBegin; 4629706f043SHong Zhang 4634c16a6a6SHong Zhang /* initialization */ 46482502324SSatish Balay ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 4654c16a6a6SHong Zhang ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr); 46682502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr); 46728de702eSHong Zhang jl = il + mbs; 4684c16a6a6SHong Zhang for (i=0; i<mbs; i++) { 4694c16a6a6SHong Zhang jl[i] = mbs; il[0] = 0; 4704c16a6a6SHong Zhang } 471b0a32e0cSBarry Smith ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr); 47228de702eSHong Zhang uik = dk + bs2; 47328de702eSHong Zhang work = uik + bs2; 47482502324SSatish Balay ierr = PetscMalloc(bs*sizeof(int),&pivots);CHKERRQ(ierr); 4754c16a6a6SHong Zhang 4764c16a6a6SHong Zhang ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr); 4774c16a6a6SHong Zhang 4784c16a6a6SHong Zhang /* check permutation */ 4794c16a6a6SHong Zhang if (!a->permute){ 4804c16a6a6SHong Zhang ai = a->i; aj = a->j; aa = a->a; 4814c16a6a6SHong Zhang } else { 4824c16a6a6SHong Zhang ai = a->inew; aj = a->jnew; 48382502324SSatish Balay ierr = PetscMalloc(bs2*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); 4844c16a6a6SHong Zhang ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); 48582502324SSatish Balay ierr = PetscMalloc(ai[mbs]*sizeof(int),&a2anew);CHKERRQ(ierr); 4864c16a6a6SHong Zhang ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr); 4874c16a6a6SHong Zhang 4884c16a6a6SHong Zhang for (i=0; i<mbs; i++){ 4894c16a6a6SHong Zhang jmin = ai[i]; jmax = ai[i+1]; 4904c16a6a6SHong Zhang for (j=jmin; j<jmax; j++){ 4914c16a6a6SHong Zhang while (a2anew[j] != j){ 4924c16a6a6SHong Zhang k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k; 4934c16a6a6SHong Zhang for (k1=0; k1<bs2; k1++){ 4944c16a6a6SHong Zhang dk[k1] = aa[k*bs2+k1]; 4954c16a6a6SHong Zhang aa[k*bs2+k1] = aa[j*bs2+k1]; 4964c16a6a6SHong Zhang aa[j*bs2+k1] = dk[k1]; 4974c16a6a6SHong Zhang } 4984c16a6a6SHong Zhang } 4994c16a6a6SHong Zhang /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */ 5004c16a6a6SHong Zhang if (i > aj[j]){ 5014c16a6a6SHong Zhang /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */ 5024c16a6a6SHong Zhang ap = aa + j*bs2; /* ptr to the beginning of j-th block of aa */ 5034c16a6a6SHong Zhang for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */ 5044c16a6a6SHong Zhang for (k=0; k<bs; k++){ /* j-th block of aa <- dk^T */ 5054c16a6a6SHong Zhang for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1]; 5064c16a6a6SHong Zhang } 5074c16a6a6SHong Zhang } 5084c16a6a6SHong Zhang } 5094c16a6a6SHong Zhang } 510323b833fSBarry Smith ierr = PetscFree(a2anew);CHKERRQ(ierr); 5114c16a6a6SHong Zhang } 5124c16a6a6SHong Zhang 5134c16a6a6SHong Zhang /* for each row k */ 5144c16a6a6SHong Zhang for (k = 0; k<mbs; k++){ 5154c16a6a6SHong Zhang 5164c16a6a6SHong Zhang /*initialize k-th row with elements nonzero in row perm(k) of A */ 5174c16a6a6SHong Zhang jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1]; 518057f5ba7SHong Zhang 5194c16a6a6SHong Zhang ap = aa + jmin*bs2; 5204c16a6a6SHong Zhang for (j = jmin; j < jmax; j++){ 5214c16a6a6SHong Zhang vj = perm_ptr[aj[j]]; /* block col. index */ 5224c16a6a6SHong Zhang rtmp_ptr = rtmp + vj*bs2; 5234c16a6a6SHong Zhang for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++; 5244c16a6a6SHong Zhang } 5254c16a6a6SHong Zhang 5264c16a6a6SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 5274c16a6a6SHong Zhang ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 5284c16a6a6SHong Zhang i = jl[k]; /* first row to be added to k_th row */ 5294c16a6a6SHong Zhang 530057f5ba7SHong Zhang while (i < k){ 5314c16a6a6SHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 5324c16a6a6SHong Zhang 5334c16a6a6SHong Zhang /* compute multiplier */ 5344c16a6a6SHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 5354c16a6a6SHong Zhang 5364c16a6a6SHong Zhang /* uik = -inv(Di)*U_bar(i,k) */ 5374c16a6a6SHong Zhang diag = ba + i*bs2; 5384c16a6a6SHong Zhang u = ba + ili*bs2; 5394c16a6a6SHong Zhang ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 5404c16a6a6SHong Zhang Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u); 5414c16a6a6SHong Zhang 5424c16a6a6SHong Zhang /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 5434c16a6a6SHong Zhang Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u); 5444c16a6a6SHong Zhang 5454c16a6a6SHong Zhang /* update -U(i,k) */ 5464c16a6a6SHong Zhang ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 5474c16a6a6SHong Zhang 5484c16a6a6SHong Zhang /* add multiple of row i to k-th row ... */ 5494c16a6a6SHong Zhang jmin = ili + 1; jmax = bi[i+1]; 5504c16a6a6SHong Zhang if (jmin < jmax){ 5514c16a6a6SHong Zhang for (j=jmin; j<jmax; j++) { 5524c16a6a6SHong Zhang /* rtmp += -U(i,k)^T * U_bar(i,j) */ 5534c16a6a6SHong Zhang rtmp_ptr = rtmp + bj[j]*bs2; 5544c16a6a6SHong Zhang u = ba + j*bs2; 5554c16a6a6SHong Zhang Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u); 5564c16a6a6SHong Zhang } 5574c16a6a6SHong Zhang 5584c16a6a6SHong Zhang /* ... add i to row list for next nonzero entry */ 5594c16a6a6SHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 5604c16a6a6SHong Zhang j = bj[jmin]; 5614c16a6a6SHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 5624c16a6a6SHong Zhang } 5634c16a6a6SHong Zhang i = nexti; 5644c16a6a6SHong Zhang } 5654c16a6a6SHong Zhang 5664c16a6a6SHong Zhang /* save nonzero entries in k-th row of U ... */ 5674c16a6a6SHong Zhang 5684c16a6a6SHong Zhang /* invert diagonal block */ 5694c16a6a6SHong Zhang diag = ba+k*bs2; 5704c16a6a6SHong Zhang ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr); 5714c16a6a6SHong Zhang Kernel_A_gets_inverse_A(bs,diag,pivots,work); 5724c16a6a6SHong Zhang 5734c16a6a6SHong Zhang jmin = bi[k]; jmax = bi[k+1]; 5744c16a6a6SHong Zhang if (jmin < jmax) { 5754c16a6a6SHong Zhang for (j=jmin; j<jmax; j++){ 5764c16a6a6SHong Zhang vj = bj[j]; /* block col. index of U */ 5774c16a6a6SHong Zhang u = ba + j*bs2; 5784c16a6a6SHong Zhang rtmp_ptr = rtmp + vj*bs2; 5794c16a6a6SHong Zhang for (k1=0; k1<bs2; k1++){ 5804c16a6a6SHong Zhang *u++ = *rtmp_ptr; 5814c16a6a6SHong Zhang *rtmp_ptr++ = 0.0; 5824c16a6a6SHong Zhang } 5834c16a6a6SHong Zhang } 5844c16a6a6SHong Zhang 5854c16a6a6SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 5864c16a6a6SHong Zhang il[k] = jmin; 5874c16a6a6SHong Zhang i = bj[jmin]; 5884c16a6a6SHong Zhang jl[k] = jl[i]; jl[i] = k; 5894c16a6a6SHong Zhang } 5904c16a6a6SHong Zhang } 5914c16a6a6SHong Zhang 5924c16a6a6SHong Zhang ierr = PetscFree(rtmp);CHKERRQ(ierr); 5934c16a6a6SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 5944c16a6a6SHong Zhang ierr = PetscFree(dk);CHKERRQ(ierr); 5954c16a6a6SHong Zhang ierr = PetscFree(pivots);CHKERRQ(ierr); 5964c16a6a6SHong Zhang if (a->permute){ 5974c16a6a6SHong Zhang ierr = PetscFree(aa);CHKERRQ(ierr); 5984c16a6a6SHong Zhang } 5994c16a6a6SHong Zhang 6004c16a6a6SHong Zhang ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr); 6014c16a6a6SHong Zhang C->factor = FACTOR_CHOLESKY; 6024c16a6a6SHong Zhang C->assembled = PETSC_TRUE; 6034c16a6a6SHong Zhang C->preallocated = PETSC_TRUE; 604b0a32e0cSBarry Smith PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */ 6054c16a6a6SHong Zhang PetscFunctionReturn(0); 6064c16a6a6SHong Zhang } 607d4edadd4SHong Zhang 6084a2ae208SSatish Balay #undef __FUNCT__ 6094a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering" 610671cb588SHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat A,Mat *B) 611671cb588SHong Zhang { 612671cb588SHong Zhang Mat C = *B; 613671cb588SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 614671cb588SHong Zhang int ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 6151ad1de41SHong Zhang int *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 616671cb588SHong Zhang int bs=a->bs,bs2 = a->bs2; 617671cb588SHong Zhang MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 618671cb588SHong Zhang MatScalar *u,*diag,*rtmp,*rtmp_ptr; 61928de702eSHong Zhang MatScalar *work; 620671cb588SHong Zhang int *pivots; 621671cb588SHong Zhang 622671cb588SHong Zhang PetscFunctionBegin; 623671cb588SHong Zhang 624671cb588SHong Zhang /* initialization */ 625671cb588SHong Zhang 62682502324SSatish Balay ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 627671cb588SHong Zhang ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr); 62882502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr); 62928de702eSHong Zhang jl = il + mbs; 630671cb588SHong Zhang for (i=0; i<mbs; i++) { 631671cb588SHong Zhang jl[i] = mbs; il[0] = 0; 632671cb588SHong Zhang } 633b0a32e0cSBarry Smith ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr); 63428de702eSHong Zhang uik = dk + bs2; 63528de702eSHong Zhang work = uik + bs2; 63682502324SSatish Balay ierr = PetscMalloc(bs*sizeof(int),&pivots);CHKERRQ(ierr); 637671cb588SHong Zhang 638671cb588SHong Zhang ai = a->i; aj = a->j; aa = a->a; 639671cb588SHong Zhang 640671cb588SHong Zhang /* for each row k */ 641671cb588SHong Zhang for (k = 0; k<mbs; k++){ 642671cb588SHong Zhang 643671cb588SHong Zhang /*initialize k-th row with elements nonzero in row k of A */ 644671cb588SHong Zhang jmin = ai[k]; jmax = ai[k+1]; 645671cb588SHong Zhang ap = aa + jmin*bs2; 646671cb588SHong Zhang for (j = jmin; j < jmax; j++){ 647671cb588SHong Zhang vj = aj[j]; /* block col. index */ 648671cb588SHong Zhang rtmp_ptr = rtmp + vj*bs2; 649671cb588SHong Zhang for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++; 650671cb588SHong Zhang } 651671cb588SHong Zhang 652671cb588SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 653671cb588SHong Zhang ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 654671cb588SHong Zhang i = jl[k]; /* first row to be added to k_th row */ 655671cb588SHong Zhang 656057f5ba7SHong Zhang while (i < k){ 657671cb588SHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 658671cb588SHong Zhang 659671cb588SHong Zhang /* compute multiplier */ 660671cb588SHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 661671cb588SHong Zhang 662671cb588SHong Zhang /* uik = -inv(Di)*U_bar(i,k) */ 663671cb588SHong Zhang diag = ba + i*bs2; 664671cb588SHong Zhang u = ba + ili*bs2; 665671cb588SHong Zhang ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 666671cb588SHong Zhang Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u); 667671cb588SHong Zhang 668671cb588SHong Zhang /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 669671cb588SHong Zhang Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u); 670671cb588SHong Zhang 671671cb588SHong Zhang /* update -U(i,k) */ 672671cb588SHong Zhang ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 673671cb588SHong Zhang 674671cb588SHong Zhang /* add multiple of row i to k-th row ... */ 675671cb588SHong Zhang jmin = ili + 1; jmax = bi[i+1]; 676671cb588SHong Zhang if (jmin < jmax){ 677671cb588SHong Zhang for (j=jmin; j<jmax; j++) { 678671cb588SHong Zhang /* rtmp += -U(i,k)^T * U_bar(i,j) */ 679671cb588SHong Zhang rtmp_ptr = rtmp + bj[j]*bs2; 680671cb588SHong Zhang u = ba + j*bs2; 681671cb588SHong Zhang Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u); 682671cb588SHong Zhang } 683671cb588SHong Zhang 684671cb588SHong Zhang /* ... add i to row list for next nonzero entry */ 685671cb588SHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 686671cb588SHong Zhang j = bj[jmin]; 687671cb588SHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 688671cb588SHong Zhang } 689671cb588SHong Zhang i = nexti; 690671cb588SHong Zhang } 691671cb588SHong Zhang 692671cb588SHong Zhang /* save nonzero entries in k-th row of U ... */ 693671cb588SHong Zhang 694671cb588SHong Zhang /* invert diagonal block */ 695671cb588SHong Zhang diag = ba+k*bs2; 696671cb588SHong Zhang ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr); 697671cb588SHong Zhang Kernel_A_gets_inverse_A(bs,diag,pivots,work); 698671cb588SHong Zhang 699671cb588SHong Zhang jmin = bi[k]; jmax = bi[k+1]; 700671cb588SHong Zhang if (jmin < jmax) { 701671cb588SHong Zhang for (j=jmin; j<jmax; j++){ 702671cb588SHong Zhang vj = bj[j]; /* block col. index of U */ 703671cb588SHong Zhang u = ba + j*bs2; 704671cb588SHong Zhang rtmp_ptr = rtmp + vj*bs2; 705671cb588SHong Zhang for (k1=0; k1<bs2; k1++){ 706671cb588SHong Zhang *u++ = *rtmp_ptr; 707671cb588SHong Zhang *rtmp_ptr++ = 0.0; 708671cb588SHong Zhang } 709671cb588SHong Zhang } 710671cb588SHong Zhang 711671cb588SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 712671cb588SHong Zhang il[k] = jmin; 713671cb588SHong Zhang i = bj[jmin]; 714671cb588SHong Zhang jl[k] = jl[i]; jl[i] = k; 715671cb588SHong Zhang } 716671cb588SHong Zhang } 717671cb588SHong Zhang 718671cb588SHong Zhang ierr = PetscFree(rtmp);CHKERRQ(ierr); 719671cb588SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 720671cb588SHong Zhang ierr = PetscFree(dk);CHKERRQ(ierr); 721671cb588SHong Zhang ierr = PetscFree(pivots);CHKERRQ(ierr); 722671cb588SHong Zhang 723671cb588SHong Zhang C->factor = FACTOR_CHOLESKY; 724671cb588SHong Zhang C->assembled = PETSC_TRUE; 725671cb588SHong Zhang C->preallocated = PETSC_TRUE; 726b0a32e0cSBarry Smith PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */ 727671cb588SHong Zhang PetscFunctionReturn(0); 728671cb588SHong Zhang } 729671cb588SHong Zhang 73049b5e25fSSatish Balay /* 731fcf159c0SHong Zhang Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP. 732cc0c071aSHong Zhang Version for blocks 2 by 2. 73349b5e25fSSatish Balay */ 7344a2ae208SSatish Balay #undef __FUNCT__ 7354a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2" 7366f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat A,Mat *B) 73749b5e25fSSatish Balay { 73849b5e25fSSatish Balay Mat C = *B; 73991602c66SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 740cc0c071aSHong Zhang IS perm = b->row; 741cc0c071aSHong Zhang int *perm_ptr,ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 742cc0c071aSHong Zhang int *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 743a1723e09SHong Zhang MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 744cc0c071aSHong Zhang MatScalar *u,*diag,*rtmp,*rtmp_ptr; 74549b5e25fSSatish Balay 74649b5e25fSSatish Balay PetscFunctionBegin; 747671cb588SHong Zhang 74891602c66SHong Zhang /* initialization */ 74991602c66SHong Zhang /* il and jl record the first nonzero element in each row of the accessing 75091602c66SHong Zhang window U(0:k, k:mbs-1). 75191602c66SHong Zhang jl: list of rows to be added to uneliminated rows 75291602c66SHong Zhang i>= k: jl(i) is the first row to be added to row i 75391602c66SHong Zhang i< k: jl(i) is the row following row i in some list of rows 75491602c66SHong Zhang jl(i) = mbs indicates the end of a list 75591602c66SHong Zhang il(i): points to the first nonzero element in columns k,...,mbs-1 of 75691602c66SHong Zhang row i of U */ 75782502324SSatish Balay ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 758cc0c071aSHong Zhang ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr); 75982502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr); 76028de702eSHong Zhang jl = il + mbs; 76191602c66SHong Zhang for (i=0; i<mbs; i++) { 7623845f261SHong Zhang jl[i] = mbs; il[0] = 0; 76391602c66SHong Zhang } 76482502324SSatish Balay ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr); 76528de702eSHong Zhang uik = dk + 4; 766cc0c071aSHong Zhang ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr); 767cc0c071aSHong Zhang 768cc0c071aSHong Zhang /* check permutation */ 769cc0c071aSHong Zhang if (!a->permute){ 770cc0c071aSHong Zhang ai = a->i; aj = a->j; aa = a->a; 771cc0c071aSHong Zhang } else { 772cc0c071aSHong Zhang ai = a->inew; aj = a->jnew; 77382502324SSatish Balay ierr = PetscMalloc(4*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); 774cc0c071aSHong Zhang ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); 77582502324SSatish Balay ierr = PetscMalloc(ai[mbs]*sizeof(int),&a2anew);CHKERRQ(ierr); 776cc0c071aSHong Zhang ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr); 777cc0c071aSHong Zhang 778cc0c071aSHong Zhang for (i=0; i<mbs; i++){ 779cc0c071aSHong Zhang jmin = ai[i]; jmax = ai[i+1]; 780cc0c071aSHong Zhang for (j=jmin; j<jmax; j++){ 781cc0c071aSHong Zhang while (a2anew[j] != j){ 782cc0c071aSHong Zhang k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k; 783cc0c071aSHong Zhang for (k1=0; k1<4; k1++){ 784cc0c071aSHong Zhang dk[k1] = aa[k*4+k1]; 785cc0c071aSHong Zhang aa[k*4+k1] = aa[j*4+k1]; 786cc0c071aSHong Zhang aa[j*4+k1] = dk[k1]; 787cc0c071aSHong Zhang } 788cc0c071aSHong Zhang } 789cc0c071aSHong Zhang /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */ 790cc0c071aSHong Zhang if (i > aj[j]){ 791a1723e09SHong Zhang /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */ 792cc0c071aSHong Zhang ap = aa + j*4; /* ptr to the beginning of the block */ 793cc0c071aSHong Zhang dk[1] = ap[1]; /* swap ap[1] and ap[2] */ 794cc0c071aSHong Zhang ap[1] = ap[2]; 795cc0c071aSHong Zhang ap[2] = dk[1]; 796cc0c071aSHong Zhang } 797cc0c071aSHong Zhang } 798cc0c071aSHong Zhang } 799ac355199SBarry Smith ierr = PetscFree(a2anew);CHKERRQ(ierr); 800cc0c071aSHong Zhang } 8013845f261SHong Zhang 80291602c66SHong Zhang /* for each row k */ 80391602c66SHong Zhang for (k = 0; k<mbs; k++){ 80491602c66SHong Zhang 80591602c66SHong Zhang /*initialize k-th row with elements nonzero in row perm(k) of A */ 806cc0c071aSHong Zhang jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1]; 807cc0c071aSHong Zhang ap = aa + jmin*4; 80891602c66SHong Zhang for (j = jmin; j < jmax; j++){ 809cc0c071aSHong Zhang vj = perm_ptr[aj[j]]; /* block col. index */ 810cc0c071aSHong Zhang rtmp_ptr = rtmp + vj*4; 811cc0c071aSHong Zhang for (i=0; i<4; i++) *rtmp_ptr++ = *ap++; 81291602c66SHong Zhang } 81391602c66SHong Zhang 81491602c66SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 815cc0c071aSHong Zhang ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr); 81691602c66SHong Zhang i = jl[k]; /* first row to be added to k_th row */ 81791602c66SHong Zhang 818057f5ba7SHong Zhang while (i < k){ 81991602c66SHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 82091602c66SHong Zhang 8213845f261SHong Zhang /* compute multiplier */ 82291602c66SHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 8233845f261SHong Zhang 8243845f261SHong Zhang /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */ 825cc0c071aSHong Zhang diag = ba + i*4; 826cc0c071aSHong Zhang u = ba + ili*4; 827cc0c071aSHong Zhang uik[0] = -(diag[0]*u[0] + diag[2]*u[1]); 828cc0c071aSHong Zhang uik[1] = -(diag[1]*u[0] + diag[3]*u[1]); 829cc0c071aSHong Zhang uik[2] = -(diag[0]*u[2] + diag[2]*u[3]); 830cc0c071aSHong Zhang uik[3] = -(diag[1]*u[2] + diag[3]*u[3]); 8313845f261SHong Zhang 8323845f261SHong Zhang /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */ 833cc0c071aSHong Zhang dk[0] += uik[0]*u[0] + uik[1]*u[1]; 834cc0c071aSHong Zhang dk[1] += uik[2]*u[0] + uik[3]*u[1]; 835cc0c071aSHong Zhang dk[2] += uik[0]*u[2] + uik[1]*u[3]; 836cc0c071aSHong Zhang dk[3] += uik[2]*u[2] + uik[3]*u[3]; 8373845f261SHong Zhang 8383845f261SHong Zhang /* update -U(i,k): ba[ili] = uik */ 839cc0c071aSHong Zhang ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr); 84091602c66SHong Zhang 84191602c66SHong Zhang /* add multiple of row i to k-th row ... */ 84291602c66SHong Zhang jmin = ili + 1; jmax = bi[i+1]; 84391602c66SHong Zhang if (jmin < jmax){ 8443845f261SHong Zhang for (j=jmin; j<jmax; j++) { 8453845f261SHong Zhang /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */ 846cc0c071aSHong Zhang rtmp_ptr = rtmp + bj[j]*4; 847cc0c071aSHong Zhang u = ba + j*4; 848cc0c071aSHong Zhang rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1]; 849cc0c071aSHong Zhang rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1]; 850cc0c071aSHong Zhang rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3]; 851cc0c071aSHong Zhang rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3]; 8523845f261SHong Zhang } 8533845f261SHong Zhang 85491602c66SHong Zhang /* ... add i to row list for next nonzero entry */ 85591602c66SHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 85691602c66SHong Zhang j = bj[jmin]; 85791602c66SHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 85891602c66SHong Zhang } 859a1723e09SHong Zhang i = nexti; 86091602c66SHong Zhang } 861cc0c071aSHong Zhang 86291602c66SHong Zhang /* save nonzero entries in k-th row of U ... */ 8633845f261SHong Zhang 864cc0c071aSHong Zhang /* invert diagonal block */ 865cc0c071aSHong Zhang diag = ba+k*4; 866cc0c071aSHong Zhang ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr); 8673845f261SHong Zhang ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr); 8683845f261SHong Zhang 86991602c66SHong Zhang jmin = bi[k]; jmax = bi[k+1]; 87091602c66SHong Zhang if (jmin < jmax) { 87191602c66SHong Zhang for (j=jmin; j<jmax; j++){ 872cc0c071aSHong Zhang vj = bj[j]; /* block col. index of U */ 873cc0c071aSHong Zhang u = ba + j*4; 874cc0c071aSHong Zhang rtmp_ptr = rtmp + vj*4; 875cc0c071aSHong Zhang for (k1=0; k1<4; k1++){ 876cc0c071aSHong Zhang *u++ = *rtmp_ptr; 877cc0c071aSHong Zhang *rtmp_ptr++ = 0.0; 878cc0c071aSHong Zhang } 879cc0c071aSHong Zhang } 8803845f261SHong Zhang 88191602c66SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 88291602c66SHong Zhang il[k] = jmin; 88391602c66SHong Zhang i = bj[jmin]; 88491602c66SHong Zhang jl[k] = jl[i]; jl[i] = k; 88591602c66SHong Zhang } 88691602c66SHong Zhang } 8873845f261SHong Zhang 88849b5e25fSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 88991602c66SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 8903845f261SHong Zhang ierr = PetscFree(dk);CHKERRQ(ierr); 89191602c66SHong Zhang if (a->permute) { 89291602c66SHong Zhang ierr = PetscFree(aa);CHKERRQ(ierr); 89391602c66SHong Zhang } 894cc0c071aSHong Zhang ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr); 89591602c66SHong Zhang C->factor = FACTOR_CHOLESKY; 89649b5e25fSSatish Balay C->assembled = PETSC_TRUE; 8975c0bcdfcSHong Zhang C->preallocated = PETSC_TRUE; 898b0a32e0cSBarry Smith PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */ 89949b5e25fSSatish Balay PetscFunctionReturn(0); 90049b5e25fSSatish Balay } 90191602c66SHong Zhang 90249b5e25fSSatish Balay /* 90349b5e25fSSatish Balay Version for when blocks are 2 by 2 Using natural ordering 90449b5e25fSSatish Balay */ 9054a2ae208SSatish Balay #undef __FUNCT__ 9064a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering" 9076f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat A,Mat *B) 90849b5e25fSSatish Balay { 90949b5e25fSSatish Balay Mat C = *B; 910ab27746eSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 911ab27746eSHong Zhang int ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 91271149678SHong Zhang int *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 913ab27746eSHong Zhang MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 914ab27746eSHong Zhang MatScalar *u,*diag,*rtmp,*rtmp_ptr; 91549b5e25fSSatish Balay 91649b5e25fSSatish Balay PetscFunctionBegin; 917671cb588SHong Zhang 918ab27746eSHong Zhang /* initialization */ 919ab27746eSHong Zhang /* il and jl record the first nonzero element in each row of the accessing 920ab27746eSHong Zhang window U(0:k, k:mbs-1). 921ab27746eSHong Zhang jl: list of rows to be added to uneliminated rows 922ab27746eSHong Zhang i>= k: jl(i) is the first row to be added to row i 923ab27746eSHong Zhang i< k: jl(i) is the row following row i in some list of rows 924ab27746eSHong Zhang jl(i) = mbs indicates the end of a list 925ab27746eSHong Zhang il(i): points to the first nonzero element in columns k,...,mbs-1 of 926ab27746eSHong Zhang row i of U */ 92782502324SSatish Balay ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 928ab27746eSHong Zhang ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr); 92982502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr); 93028de702eSHong Zhang jl = il + mbs; 931ab27746eSHong Zhang for (i=0; i<mbs; i++) { 932ab27746eSHong Zhang jl[i] = mbs; il[0] = 0; 93349b5e25fSSatish Balay } 93482502324SSatish Balay ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr); 93528de702eSHong Zhang uik = dk + 4; 936ab27746eSHong Zhang 937ab27746eSHong Zhang ai = a->i; aj = a->j; aa = a->a; 938ab27746eSHong Zhang 939ab27746eSHong Zhang /* for each row k */ 940ab27746eSHong Zhang for (k = 0; k<mbs; k++){ 941ab27746eSHong Zhang 942ab27746eSHong Zhang /*initialize k-th row with elements nonzero in row k of A */ 943ab27746eSHong Zhang jmin = ai[k]; jmax = ai[k+1]; 944ab27746eSHong Zhang ap = aa + jmin*4; 945ab27746eSHong Zhang for (j = jmin; j < jmax; j++){ 946ab27746eSHong Zhang vj = aj[j]; /* block col. index */ 947ab27746eSHong Zhang rtmp_ptr = rtmp + vj*4; 948ab27746eSHong Zhang for (i=0; i<4; i++) *rtmp_ptr++ = *ap++; 94949b5e25fSSatish Balay } 950ab27746eSHong Zhang 951ab27746eSHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 952ab27746eSHong Zhang ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr); 953ab27746eSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 954ab27746eSHong Zhang 955057f5ba7SHong Zhang while (i < k){ 956ab27746eSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 957ab27746eSHong Zhang 958ab27746eSHong Zhang /* compute multiplier */ 959ab27746eSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 960ab27746eSHong Zhang 961ab27746eSHong Zhang /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */ 962ab27746eSHong Zhang diag = ba + i*4; 963ab27746eSHong Zhang u = ba + ili*4; 964ab27746eSHong Zhang uik[0] = -(diag[0]*u[0] + diag[2]*u[1]); 965ab27746eSHong Zhang uik[1] = -(diag[1]*u[0] + diag[3]*u[1]); 966ab27746eSHong Zhang uik[2] = -(diag[0]*u[2] + diag[2]*u[3]); 967ab27746eSHong Zhang uik[3] = -(diag[1]*u[2] + diag[3]*u[3]); 968ab27746eSHong Zhang 969ab27746eSHong Zhang /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */ 970ab27746eSHong Zhang dk[0] += uik[0]*u[0] + uik[1]*u[1]; 971ab27746eSHong Zhang dk[1] += uik[2]*u[0] + uik[3]*u[1]; 972ab27746eSHong Zhang dk[2] += uik[0]*u[2] + uik[1]*u[3]; 973ab27746eSHong Zhang dk[3] += uik[2]*u[2] + uik[3]*u[3]; 974ab27746eSHong Zhang 975ab27746eSHong Zhang /* update -U(i,k): ba[ili] = uik */ 976ab27746eSHong Zhang ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr); 977ab27746eSHong Zhang 978ab27746eSHong Zhang /* add multiple of row i to k-th row ... */ 979ab27746eSHong Zhang jmin = ili + 1; jmax = bi[i+1]; 980ab27746eSHong Zhang if (jmin < jmax){ 981ab27746eSHong Zhang for (j=jmin; j<jmax; j++) { 982ab27746eSHong Zhang /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */ 983ab27746eSHong Zhang rtmp_ptr = rtmp + bj[j]*4; 984ab27746eSHong Zhang u = ba + j*4; 985ab27746eSHong Zhang rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1]; 986ab27746eSHong Zhang rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1]; 987ab27746eSHong Zhang rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3]; 988ab27746eSHong Zhang rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3]; 98949b5e25fSSatish Balay } 990ab27746eSHong Zhang 991ab27746eSHong Zhang /* ... add i to row list for next nonzero entry */ 992ab27746eSHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 993ab27746eSHong Zhang j = bj[jmin]; 994ab27746eSHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 99549b5e25fSSatish Balay } 996ab27746eSHong Zhang i = nexti; 99749b5e25fSSatish Balay } 998ab27746eSHong Zhang 999ab27746eSHong Zhang /* save nonzero entries in k-th row of U ... */ 1000ab27746eSHong Zhang 100149b5e25fSSatish Balay /* invert diagonal block */ 1002ab27746eSHong Zhang diag = ba+k*4; 1003ab27746eSHong Zhang ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr); 1004ab27746eSHong Zhang ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr); 1005ab27746eSHong Zhang 1006ab27746eSHong Zhang jmin = bi[k]; jmax = bi[k+1]; 1007ab27746eSHong Zhang if (jmin < jmax) { 1008ab27746eSHong Zhang for (j=jmin; j<jmax; j++){ 1009ab27746eSHong Zhang vj = bj[j]; /* block col. index of U */ 1010ab27746eSHong Zhang u = ba + j*4; 1011ab27746eSHong Zhang rtmp_ptr = rtmp + vj*4; 1012ab27746eSHong Zhang for (k1=0; k1<4; k1++){ 1013ab27746eSHong Zhang *u++ = *rtmp_ptr; 1014ab27746eSHong Zhang *rtmp_ptr++ = 0.0; 1015ab27746eSHong Zhang } 1016ab27746eSHong Zhang } 1017ab27746eSHong Zhang 1018ab27746eSHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 1019ab27746eSHong Zhang il[k] = jmin; 1020ab27746eSHong Zhang i = bj[jmin]; 1021ab27746eSHong Zhang jl[k] = jl[i]; jl[i] = k; 1022ab27746eSHong Zhang } 102349b5e25fSSatish Balay } 102449b5e25fSSatish Balay 102549b5e25fSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 1026ab27746eSHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 1027ab27746eSHong Zhang ierr = PetscFree(dk);CHKERRQ(ierr); 1028ab27746eSHong Zhang 1029ab27746eSHong Zhang C->factor = FACTOR_CHOLESKY; 103049b5e25fSSatish Balay C->assembled = PETSC_TRUE; 10315c0bcdfcSHong Zhang C->preallocated = PETSC_TRUE; 1032b0a32e0cSBarry Smith PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */ 103349b5e25fSSatish Balay PetscFunctionReturn(0); 103449b5e25fSSatish Balay } 103549b5e25fSSatish Balay 103649b5e25fSSatish Balay /* 10375c0bcdfcSHong Zhang Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP. 103891602c66SHong Zhang Version for blocks are 1 by 1. 103949b5e25fSSatish Balay */ 10404a2ae208SSatish Balay #undef __FUNCT__ 10414a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1" 10426f9a564bSHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat A,Mat *B) 104349b5e25fSSatish Balay { 104449b5e25fSSatish Balay Mat C = *B; 104549b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 104649b5e25fSSatish Balay IS ip = b->row; 1047351d0355SHong Zhang int *rip,ierr,i,j,mbs = a->mbs,*bi = b->i,*bj = b->j; 1048cb718733SHong Zhang int *ai,*aj,*r; 1049ab27746eSHong Zhang int k,jmin,jmax,*jl,*il,vj,nexti,ili; 1050066653e3SSatish Balay MatScalar *rtmp; 10512d007305SHong Zhang MatScalar *ba = b->a,*aa,ak; 105249b5e25fSSatish Balay MatScalar dk,uikdi; 105349b5e25fSSatish Balay 105449b5e25fSSatish Balay PetscFunctionBegin; 105549b5e25fSSatish Balay ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 1056cb718733SHong Zhang if (!a->permute){ 10572d007305SHong Zhang ai = a->i; aj = a->j; aa = a->a; 10582d007305SHong Zhang } else { 10592d007305SHong Zhang ai = a->inew; aj = a->jnew; 106082502324SSatish Balay ierr = PetscMalloc(ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); 1061cb718733SHong Zhang ierr = PetscMemcpy(aa,a->a,ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); 106282502324SSatish Balay ierr = PetscMalloc(ai[mbs]*sizeof(int),&r);CHKERRQ(ierr); 10632d007305SHong Zhang ierr= PetscMemcpy(r,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr); 10642d007305SHong Zhang 10652d007305SHong Zhang jmin = ai[0]; jmax = ai[mbs]; 10662d007305SHong Zhang for (j=jmin; j<jmax; j++){ 1067cb718733SHong Zhang while (r[j] != j){ 10682d007305SHong Zhang k = r[j]; r[j] = r[k]; r[k] = k; 10692d007305SHong Zhang ak = aa[k]; aa[k] = aa[j]; aa[j] = ak; 10702d007305SHong Zhang } 10712d007305SHong Zhang } 1072ac355199SBarry Smith ierr = PetscFree(r);CHKERRQ(ierr); 10732d007305SHong Zhang } 10742d007305SHong Zhang 107591602c66SHong Zhang /* initialization */ 107649b5e25fSSatish Balay /* il and jl record the first nonzero element in each row of the accessing 107749b5e25fSSatish Balay window U(0:k, k:mbs-1). 107849b5e25fSSatish Balay jl: list of rows to be added to uneliminated rows 107949b5e25fSSatish Balay i>= k: jl(i) is the first row to be added to row i 108049b5e25fSSatish Balay i< k: jl(i) is the row following row i in some list of rows 108149b5e25fSSatish Balay jl(i) = mbs indicates the end of a list 108249b5e25fSSatish Balay il(i): points to the first nonzero element in columns k,...,mbs-1 of 108349b5e25fSSatish Balay row i of U */ 108482502324SSatish Balay ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 108582502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr); 108628de702eSHong Zhang jl = il + mbs; 108749b5e25fSSatish Balay for (i=0; i<mbs; i++) { 108849b5e25fSSatish Balay rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 108949b5e25fSSatish Balay } 109049b5e25fSSatish Balay 109191602c66SHong Zhang /* for each row k */ 109249b5e25fSSatish Balay for (k = 0; k<mbs; k++){ 109349b5e25fSSatish Balay 109491602c66SHong Zhang /*initialize k-th row with elements nonzero in row perm(k) of A */ 109549b5e25fSSatish Balay jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 1096057f5ba7SHong Zhang 109749b5e25fSSatish Balay for (j = jmin; j < jmax; j++){ 1098351d0355SHong Zhang vj = rip[aj[j]]; 1099ab27746eSHong Zhang rtmp[vj] = aa[j]; 110049b5e25fSSatish Balay } 110149b5e25fSSatish Balay 110291602c66SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 110349b5e25fSSatish Balay dk = rtmp[k]; 110449b5e25fSSatish Balay i = jl[k]; /* first row to be added to k_th row */ 110549b5e25fSSatish Balay 1106057f5ba7SHong Zhang while (i < k){ 110749b5e25fSSatish Balay nexti = jl[i]; /* next row to be added to k_th row */ 110849b5e25fSSatish Balay 110991602c66SHong Zhang /* compute multiplier, update D(k) and U(i,k) */ 111049b5e25fSSatish Balay ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 111149b5e25fSSatish Balay uikdi = - ba[ili]*ba[i]; 111249b5e25fSSatish Balay dk += uikdi*ba[ili]; 1113658e7b3eSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 111449b5e25fSSatish Balay 111591602c66SHong Zhang /* add multiple of row i to k-th row ... */ 111649b5e25fSSatish Balay jmin = ili + 1; jmax = bi[i+1]; 111749b5e25fSSatish Balay if (jmin < jmax){ 111849b5e25fSSatish Balay for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 111991602c66SHong Zhang /* ... add i to row list for next nonzero entry */ 112049b5e25fSSatish Balay il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 112149b5e25fSSatish Balay j = bj[jmin]; 112249b5e25fSSatish Balay jl[i] = jl[j]; jl[j] = i; /* update jl */ 112349b5e25fSSatish Balay } 1124ab27746eSHong Zhang i = nexti; 112549b5e25fSSatish Balay } 112649b5e25fSSatish Balay 112791602c66SHong Zhang /* check for zero pivot and save diagoanl element */ 112849b5e25fSSatish Balay if (dk == 0.0){ 112929bbc08cSBarry Smith SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot"); 1130671cb588SHong Zhang /* 11311223c01bSHong Zhang } else if (PetscRealPart(dk) < 0.0){ 11321223c01bSHong Zhang SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Negative pivot: d[%d] = %g\n",k,dk); 1133671cb588SHong Zhang */ 113449b5e25fSSatish Balay } 113549b5e25fSSatish Balay 113691602c66SHong Zhang /* save nonzero entries in k-th row of U ... */ 1137f3dd7b2dSKris Buschelman ba[k] = 1.0/dk; 113849b5e25fSSatish Balay jmin = bi[k]; jmax = bi[k+1]; 113949b5e25fSSatish Balay if (jmin < jmax) { 114049b5e25fSSatish Balay for (j=jmin; j<jmax; j++){ 1141cc0c071aSHong Zhang vj = bj[j]; ba[j] = rtmp[vj]; rtmp[vj] = 0.0; 114249b5e25fSSatish Balay } 114391602c66SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 114449b5e25fSSatish Balay il[k] = jmin; 114549b5e25fSSatish Balay i = bj[jmin]; 114649b5e25fSSatish Balay jl[k] = jl[i]; jl[i] = k; 114749b5e25fSSatish Balay } 114849b5e25fSSatish Balay } 114949b5e25fSSatish Balay 115049b5e25fSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 115149b5e25fSSatish Balay ierr = PetscFree(il);CHKERRQ(ierr); 1152cb718733SHong Zhang if (a->permute){ 1153cb718733SHong Zhang ierr = PetscFree(aa);CHKERRQ(ierr); 1154cb718733SHong Zhang } 115549b5e25fSSatish Balay 115649b5e25fSSatish Balay ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 11578be8c0c7SHong Zhang C->factor = FACTOR_CHOLESKY; 115849b5e25fSSatish Balay C->assembled = PETSC_TRUE; 11595c0bcdfcSHong Zhang C->preallocated = PETSC_TRUE; 1160b0a32e0cSBarry Smith PetscLogFlops(b->mbs); 116149b5e25fSSatish Balay PetscFunctionReturn(0); 116249b5e25fSSatish Balay } 116349b5e25fSSatish Balay 1164671cb588SHong Zhang /* 1165671cb588SHong Zhang Version for when blocks are 1 by 1 Using natural ordering 1166671cb588SHong Zhang */ 11674a2ae208SSatish Balay #undef __FUNCT__ 11684a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering" 1169671cb588SHong Zhang int MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat A,Mat *B) 1170671cb588SHong Zhang { 1171671cb588SHong Zhang Mat C = *B; 1172671cb588SHong Zhang Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data; 1173653a6975SHong Zhang int ierr,i,j,mbs = a->mbs; 1174653a6975SHong Zhang int *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 1175b00f7748SHong Zhang int k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz,ndamp = 0; 1176653a6975SHong Zhang MatScalar *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval; 117721c26570Svictorle PetscReal damping=b->factor_damping, zeropivot=b->factor_zeropivot,shift_amount; 117821c26570Svictorle PetscTruth damp,chshift; 117921c26570Svictorle int nshift=0; 1180653a6975SHong Zhang 1181653a6975SHong Zhang PetscFunctionBegin; 1182653a6975SHong Zhang /* initialization */ 1183653a6975SHong Zhang /* il and jl record the first nonzero element in each row of the accessing 1184653a6975SHong Zhang window U(0:k, k:mbs-1). 1185653a6975SHong Zhang jl: list of rows to be added to uneliminated rows 1186653a6975SHong Zhang i>= k: jl(i) is the first row to be added to row i 1187653a6975SHong Zhang i< k: jl(i) is the row following row i in some list of rows 1188653a6975SHong Zhang jl(i) = mbs indicates the end of a list 1189653a6975SHong Zhang il(i): points to the first nonzero element in U(i,k:mbs-1) 1190653a6975SHong Zhang */ 1191653a6975SHong Zhang ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 1192653a6975SHong Zhang ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr); 1193653a6975SHong Zhang jl = il + mbs; 1194b00f7748SHong Zhang 119521c26570Svictorle shift_amount = 0; 1196b00f7748SHong Zhang do { 1197b00f7748SHong Zhang damp = PETSC_FALSE; 119821c26570Svictorle chshift = PETSC_FALSE; 1199653a6975SHong Zhang for (i=0; i<mbs; i++) { 1200653a6975SHong Zhang rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 1201653a6975SHong Zhang } 1202653a6975SHong Zhang 1203b00f7748SHong Zhang for (k = 0; k<mbs; k++){ /* row k */ 1204653a6975SHong Zhang /*initialize k-th row with elements nonzero in row perm(k) of A */ 1205653a6975SHong Zhang nz = ai[k+1] - ai[k]; 1206653a6975SHong Zhang acol = aj + ai[k]; 1207653a6975SHong Zhang aval = aa + ai[k]; 1208653a6975SHong Zhang bval = ba + bi[k]; 1209653a6975SHong Zhang while (nz -- ){ 1210653a6975SHong Zhang rtmp[*acol++] = *aval++; 1211653a6975SHong Zhang *bval++ = 0.0; /* for in-place factorization */ 1212653a6975SHong Zhang } 1213b00f7748SHong Zhang /* damp the diagonal of the matrix */ 121421c26570Svictorle if (ndamp||nshift) rtmp[k] += damping+shift_amount; 1215653a6975SHong Zhang 1216653a6975SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 1217653a6975SHong Zhang dk = rtmp[k]; 1218653a6975SHong Zhang i = jl[k]; /* first row to be added to k_th row */ 1219653a6975SHong Zhang 1220653a6975SHong Zhang while (i < k){ 1221653a6975SHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 1222653a6975SHong Zhang 1223653a6975SHong Zhang /* compute multiplier, update D(k) and U(i,k) */ 1224653a6975SHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 1225653a6975SHong Zhang uikdi = - ba[ili]*ba[bi[i]]; 1226653a6975SHong Zhang dk += uikdi*ba[ili]; 1227653a6975SHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 1228653a6975SHong Zhang 1229653a6975SHong Zhang /* add multiple of row i to k-th row ... */ 1230653a6975SHong Zhang jmin = ili + 1; 1231653a6975SHong Zhang nz = bi[i+1] - jmin; 1232653a6975SHong Zhang if (nz > 0){ 1233653a6975SHong Zhang bcol = bj + jmin; 1234653a6975SHong Zhang bval = ba + jmin; 1235653a6975SHong Zhang while (nz --) rtmp[*bcol++] += uikdi*(*bval++); 1236653a6975SHong Zhang /* ... add i to row list for next nonzero entry */ 1237653a6975SHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 1238653a6975SHong Zhang j = bj[jmin]; 1239653a6975SHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 1240653a6975SHong Zhang } 1241653a6975SHong Zhang i = nexti; 1242653a6975SHong Zhang } 1243653a6975SHong Zhang 124421c26570Svictorle if (PetscRealPart(dk) < zeropivot && b->factor_shift){ 1245d530a6e7Svictorle /* calculate a shift that would make this row diagonally dominant */ 124629f69fd4SBarry Smith PetscReal rs = PetscAbs(PetscRealPart(dk)); 124721c26570Svictorle jmin = bi[k]+1; 124821c26570Svictorle nz = bi[k+1] - jmin; 124921c26570Svictorle if (nz){ 125021c26570Svictorle bcol = bj + jmin; 125121c26570Svictorle bval = ba + jmin; 125221c26570Svictorle while (nz--){ 125321c26570Svictorle rs += PetscAbsScalar(rtmp[*bcol++]); 125421c26570Svictorle } 125521c26570Svictorle } 1256d530a6e7Svictorle /* if this shift is less than the previous, just up the previous 1257d530a6e7Svictorle one by a bit */ 1258d530a6e7Svictorle shift_amount = PetscMax(rs,1.1*shift_amount); 125921c26570Svictorle chshift = PETSC_TRUE; 1260d530a6e7Svictorle /* Unlike in the ILU case there is no exit condition on nshift: 1261d530a6e7Svictorle we increase the shift until it converges. There is no guarantee that 1262d530a6e7Svictorle this algorithm converges faster or slower, or is better or worse 1263d530a6e7Svictorle than the ILU algorithm. */ 126421c26570Svictorle nshift++; 126521c26570Svictorle break; 126621c26570Svictorle } 1267ffa0b812SHong Zhang if (PetscRealPart(dk) < zeropivot){ 1268ffa0b812SHong Zhang if (damping == (PetscReal) PETSC_DECIDE) damping = -PetscRealPart(dk)/(k+1); 1269b00f7748SHong Zhang if (damping > 0.0) { 1270b00f7748SHong Zhang if (ndamp) damping *= 2.0; 1271b00f7748SHong Zhang damp = PETSC_TRUE; 1272b00f7748SHong Zhang ndamp++; 1273b00f7748SHong Zhang break; 1274481c2c92SHong Zhang } else if (PetscAbsScalar(dk) < zeropivot){ 1275ffa0b812SHong Zhang SETERRQ3(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %d value %g tolerance %g",k,PetscRealPart(dk),zeropivot); 1276481c2c92SHong Zhang } else { 1277ffa0b812SHong Zhang PetscLogInfo((PetscObject)A,"Negative pivot %g in row %d of Cholesky factorization\n",PetscRealPart(dk),k); 1278b00f7748SHong Zhang } 1279064503c5SHong Zhang } 1280653a6975SHong Zhang 1281653a6975SHong Zhang /* save nonzero entries in k-th row of U ... */ 128282799104SHong Zhang /* printf("%d, dk: %g, 1/dk: %g\n",k,dk,1/dk); */ 1283653a6975SHong Zhang ba[bi[k]] = 1.0/dk; 1284653a6975SHong Zhang jmin = bi[k]+1; 1285653a6975SHong Zhang nz = bi[k+1] - jmin; 1286653a6975SHong Zhang if (nz){ 1287653a6975SHong Zhang bcol = bj + jmin; 1288653a6975SHong Zhang bval = ba + jmin; 1289653a6975SHong Zhang while (nz--){ 1290653a6975SHong Zhang *bval++ = rtmp[*bcol]; 1291653a6975SHong Zhang rtmp[*bcol++] = 0.0; 1292653a6975SHong Zhang } 1293653a6975SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 1294653a6975SHong Zhang il[k] = jmin; 1295653a6975SHong Zhang i = bj[jmin]; 1296653a6975SHong Zhang jl[k] = jl[i]; jl[i] = k; 1297653a6975SHong Zhang } 1298b00f7748SHong Zhang } /* end of for (k = 0; k<mbs; k++) */ 129921c26570Svictorle } while (damp||chshift); 1300653a6975SHong Zhang ierr = PetscFree(rtmp);CHKERRQ(ierr); 1301653a6975SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 1302653a6975SHong Zhang 1303653a6975SHong Zhang C->factor = FACTOR_CHOLESKY; 1304653a6975SHong Zhang C->assembled = PETSC_TRUE; 1305653a6975SHong Zhang C->preallocated = PETSC_TRUE; 1306653a6975SHong Zhang PetscLogFlops(b->mbs); 1307b00f7748SHong Zhang if (ndamp) { 1308b45a75daSHong Zhang PetscLogInfo(0,"MatCholeskyFactorNumerical_SeqSBAIJ_1_NaturalOrdering: number of damping tries %d damping value %g\n",ndamp,damping); 1309b00f7748SHong Zhang } 1310fb10cecfSBarry Smith if (nshift) { 1311fb10cecfSBarry Smith PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering diagonal shifted %d shifts\n",nshift); 1312fb10cecfSBarry Smith } 1313fb10cecfSBarry Smith 1314653a6975SHong Zhang PetscFunctionReturn(0); 1315653a6975SHong Zhang } 1316653a6975SHong Zhang 1317653a6975SHong Zhang #undef __FUNCT__ 13184a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ" 131915e8a5b3SHong Zhang int MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,MatFactorInfo *info) 132049b5e25fSSatish Balay { 13214445ddedSHong Zhang int ierr; 132249b5e25fSSatish Balay Mat C; 132349b5e25fSSatish Balay 132449b5e25fSSatish Balay PetscFunctionBegin; 132515e8a5b3SHong Zhang ierr = MatCholeskyFactorSymbolic(A,perm,info,&C);CHKERRQ(ierr); 1326a4ada70bSHong Zhang ierr = MatCholeskyFactorNumeric(A,&C);CHKERRQ(ierr); 13274445ddedSHong Zhang ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 132849b5e25fSSatish Balay PetscFunctionReturn(0); 132949b5e25fSSatish Balay } 133049b5e25fSSatish Balay 133149b5e25fSSatish Balay 1332