1be1d678aSKris Buschelman #define PETSCMAT_DLL 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" 1813f74950SBarry Smith PetscErrorCode MatGetInertia_SeqSBAIJ(Mat F,PetscInt *nneig,PetscInt *nzero,PetscInt *npos) 195f9f512dSHong Zhang { 20638f5ce0SDinesh Kaushik Mat_SeqSBAIJ *fact_ptr = (Mat_SeqSBAIJ*)F->data; 21dd6ea824SBarry Smith MatScalar *dd = fact_ptr->a; 22d0f46423SBarry Smith PetscInt mbs=fact_ptr->mbs,bs=F->rmap->bs,i,nneig_tmp,npos_tmp,*fi = fact_ptr->i; 235f9f512dSHong Zhang 245f9f512dSHong Zhang PetscFunctionBegin; 2577431f27SBarry Smith if (bs != 1) SETERRQ1(PETSC_ERR_SUP,"No support for bs: %D >1 yet",bs); 26eeeff2ecSHong Zhang nneig_tmp = 0; npos_tmp = 0; 27eeeff2ecSHong Zhang for (i=0; i<mbs; i++){ 28eeeff2ecSHong Zhang if (PetscRealPart(dd[*fi]) > 0.0){ 29eeeff2ecSHong Zhang npos_tmp++; 30eeeff2ecSHong Zhang } else if (PetscRealPart(dd[*fi]) < 0.0){ 31eeeff2ecSHong Zhang nneig_tmp++; 325f9f512dSHong Zhang } 33eeeff2ecSHong Zhang fi++; 343e0d88b5SBarry Smith } 35eeeff2ecSHong Zhang if (nneig) *nneig = nneig_tmp; 36eeeff2ecSHong Zhang if (npos) *npos = npos_tmp; 37eeeff2ecSHong Zhang if (nzero) *nzero = mbs - nneig_tmp - npos_tmp; 38eeeff2ecSHong Zhang 395f9f512dSHong Zhang PetscFunctionReturn(0); 405f9f512dSHong Zhang } 418dc9108eSHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */ 425f9f512dSHong Zhang 43db4efbfdSBarry Smith extern PetscErrorCode MatSeqSBAIJSetNumericFactorization(Mat,PetscTruth); 44db4efbfdSBarry Smith 455f9f512dSHong Zhang /* 465f9f512dSHong Zhang Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP. 4710c27e3dSHong Zhang Use Modified Sparse Row (MSR) storage for u and ju. See page 85, "Iterative Methods ..." by Saad. 485f9f512dSHong Zhang */ 4910c27e3dSHong Zhang #undef __FUNCT__ 5010c27e3dSHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ_MSR" 51*719d5645SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(Mat F,Mat A,IS perm,MatFactorInfo *info) 5210c27e3dSHong Zhang { 5310c27e3dSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b; 5410c27e3dSHong Zhang PetscErrorCode ierr; 5510c27e3dSHong Zhang PetscInt *rip,i,mbs = a->mbs,*ai,*aj; 56d0f46423SBarry Smith PetscInt *jutmp,bs = A->rmap->bs,bs2=a->bs2; 5710c27e3dSHong Zhang PetscInt m,reallocs = 0,prow; 5810c27e3dSHong Zhang PetscInt *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd; 5910c27e3dSHong Zhang PetscReal f = info->fill; 6010c27e3dSHong Zhang PetscTruth perm_identity; 6110c27e3dSHong Zhang 6210c27e3dSHong Zhang PetscFunctionBegin; 6310c27e3dSHong Zhang /* check whether perm is the identity mapping */ 6410c27e3dSHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 6510c27e3dSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 6610c27e3dSHong Zhang 6710c27e3dSHong Zhang if (perm_identity){ /* without permutation */ 6810c27e3dSHong Zhang a->permute = PETSC_FALSE; 6910c27e3dSHong Zhang ai = a->i; aj = a->j; 7010c27e3dSHong Zhang } else { /* non-trivial permutation */ 7110c27e3dSHong Zhang a->permute = PETSC_TRUE; 7210c27e3dSHong Zhang ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr); 7310c27e3dSHong Zhang ai = a->inew; aj = a->jnew; 7410c27e3dSHong Zhang } 7510c27e3dSHong Zhang 7610c27e3dSHong Zhang /* initialization */ 7710c27e3dSHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&iu);CHKERRQ(ierr); 7810c27e3dSHong Zhang umax = (PetscInt)(f*ai[mbs] + 1); umax += mbs + 1; 7910c27e3dSHong Zhang ierr = PetscMalloc(umax*sizeof(PetscInt),&ju);CHKERRQ(ierr); 8010c27e3dSHong Zhang iu[0] = mbs+1; 8110c27e3dSHong Zhang juidx = mbs + 1; /* index for ju */ 8210c27e3dSHong Zhang ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&jl);CHKERRQ(ierr); /* linked list for pivot row */ 8310c27e3dSHong Zhang q = jl + mbs; /* linked list for col index */ 8410c27e3dSHong Zhang for (i=0; i<mbs; i++){ 8510c27e3dSHong Zhang jl[i] = mbs; 8610c27e3dSHong Zhang q[i] = 0; 8710c27e3dSHong Zhang } 8810c27e3dSHong Zhang 8910c27e3dSHong Zhang /* for each row k */ 9010c27e3dSHong Zhang for (k=0; k<mbs; k++){ 9110c27e3dSHong Zhang for (i=0; i<mbs; i++) q[i] = 0; /* to be removed! */ 9210c27e3dSHong Zhang nzk = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */ 9310c27e3dSHong Zhang q[k] = mbs; 9410c27e3dSHong Zhang /* initialize nonzero structure of k-th row to row rip[k] of A */ 9510c27e3dSHong Zhang jmin = ai[rip[k]] +1; /* exclude diag[k] */ 9610c27e3dSHong Zhang jmax = ai[rip[k]+1]; 9710c27e3dSHong Zhang for (j=jmin; j<jmax; j++){ 9810c27e3dSHong Zhang vj = rip[aj[j]]; /* col. value */ 9910c27e3dSHong Zhang if(vj > k){ 10010c27e3dSHong Zhang qm = k; 10110c27e3dSHong Zhang do { 10210c27e3dSHong Zhang m = qm; qm = q[m]; 10310c27e3dSHong Zhang } while(qm < vj); 10410c27e3dSHong Zhang if (qm == vj) { 10510c27e3dSHong Zhang SETERRQ(PETSC_ERR_PLIB,"Duplicate entry in A\n"); 10610c27e3dSHong Zhang } 10710c27e3dSHong Zhang nzk++; 10810c27e3dSHong Zhang q[m] = vj; 10910c27e3dSHong Zhang q[vj] = qm; 11010c27e3dSHong Zhang } /* if(vj > k) */ 11110c27e3dSHong Zhang } /* for (j=jmin; j<jmax; j++) */ 11210c27e3dSHong Zhang 11310c27e3dSHong Zhang /* modify nonzero structure of k-th row by computing fill-in 11410c27e3dSHong Zhang for each row i to be merged in */ 11510c27e3dSHong Zhang prow = k; 11610c27e3dSHong Zhang prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */ 11710c27e3dSHong Zhang 11810c27e3dSHong Zhang while (prow < k){ 11910c27e3dSHong Zhang /* merge row prow into k-th row */ 12010c27e3dSHong Zhang jmin = iu[prow] + 1; jmax = iu[prow+1]; 12110c27e3dSHong Zhang qm = k; 12210c27e3dSHong Zhang for (j=jmin; j<jmax; j++){ 12310c27e3dSHong Zhang vj = ju[j]; 12410c27e3dSHong Zhang do { 12510c27e3dSHong Zhang m = qm; qm = q[m]; 12610c27e3dSHong Zhang } while (qm < vj); 12710c27e3dSHong Zhang if (qm != vj){ 12810c27e3dSHong Zhang nzk++; q[m] = vj; q[vj] = qm; qm = vj; 12910c27e3dSHong Zhang } 13010c27e3dSHong Zhang } 13110c27e3dSHong Zhang prow = jl[prow]; /* next pivot row */ 13210c27e3dSHong Zhang } 13310c27e3dSHong Zhang 13410c27e3dSHong Zhang /* add k to row list for first nonzero element in k-th row */ 13510c27e3dSHong Zhang if (nzk > 0){ 13610c27e3dSHong Zhang i = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */ 13710c27e3dSHong Zhang jl[k] = jl[i]; jl[i] = k; 13810c27e3dSHong Zhang } 13910c27e3dSHong Zhang iu[k+1] = iu[k] + nzk; 14010c27e3dSHong Zhang 14110c27e3dSHong Zhang /* allocate more space to ju if needed */ 14210c27e3dSHong Zhang if (iu[k+1] > umax) { 14310c27e3dSHong Zhang /* estimate how much additional space we will need */ 14410c27e3dSHong Zhang /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 14510c27e3dSHong Zhang /* just double the memory each time */ 14610c27e3dSHong Zhang maxadd = umax; 14710c27e3dSHong Zhang if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2; 14810c27e3dSHong Zhang umax += maxadd; 14910c27e3dSHong Zhang 15010c27e3dSHong Zhang /* allocate a longer ju */ 15110c27e3dSHong Zhang ierr = PetscMalloc(umax*sizeof(PetscInt),&jutmp);CHKERRQ(ierr); 15210c27e3dSHong Zhang ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(PetscInt));CHKERRQ(ierr); 15310c27e3dSHong Zhang ierr = PetscFree(ju);CHKERRQ(ierr); 15410c27e3dSHong Zhang ju = jutmp; 15510c27e3dSHong Zhang reallocs++; /* count how many times we realloc */ 15610c27e3dSHong Zhang } 15710c27e3dSHong Zhang 15810c27e3dSHong Zhang /* save nonzero structure of k-th row in ju */ 15910c27e3dSHong Zhang i=k; 16010c27e3dSHong Zhang while (nzk --) { 16110c27e3dSHong Zhang i = q[i]; 16210c27e3dSHong Zhang ju[juidx++] = i; 16310c27e3dSHong Zhang } 16410c27e3dSHong Zhang } 16510c27e3dSHong Zhang 1666cf91177SBarry Smith #if defined(PETSC_USE_INFO) 16710c27e3dSHong Zhang if (ai[mbs] != 0) { 16810c27e3dSHong Zhang PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]); 169ae15b995SBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr); 170ae15b995SBarry Smith ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 171ae15b995SBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr); 172ae15b995SBarry Smith ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 17310c27e3dSHong Zhang } else { 174ae15b995SBarry Smith ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 17510c27e3dSHong Zhang } 17663ba0a88SBarry Smith #endif 17710c27e3dSHong Zhang 17810c27e3dSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 17910c27e3dSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 18010c27e3dSHong Zhang 18110c27e3dSHong Zhang /* put together the new matrix */ 182*719d5645SBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(F,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 18310c27e3dSHong Zhang 184*719d5645SBarry Smith /* ierr = PetscLogObjectParent(B,iperm);CHKERRQ(ierr); */ 185*719d5645SBarry Smith b = (Mat_SeqSBAIJ*)(F)->data; 18610c27e3dSHong Zhang b->singlemalloc = PETSC_FALSE; 187e6b907acSBarry Smith b->free_a = PETSC_TRUE; 188e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 18910c27e3dSHong Zhang ierr = PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr); 19010c27e3dSHong Zhang b->j = ju; 19110c27e3dSHong Zhang b->i = iu; 19210c27e3dSHong Zhang b->diag = 0; 19310c27e3dSHong Zhang b->ilen = 0; 19410c27e3dSHong Zhang b->imax = 0; 19510c27e3dSHong Zhang b->row = perm; 19610c27e3dSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 19710c27e3dSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 19810c27e3dSHong Zhang b->icol = perm; 19910c27e3dSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 20010c27e3dSHong Zhang ierr = PetscMalloc((bs*mbs+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 20110c27e3dSHong Zhang /* In b structure: Free imax, ilen, old a, old j. 20210c27e3dSHong Zhang Allocate idnew, solve_work, new a, new j */ 203*719d5645SBarry Smith ierr = PetscLogObjectMemory(F,(iu[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 20410c27e3dSHong Zhang b->maxnz = b->nz = iu[mbs]; 20510c27e3dSHong Zhang 206*719d5645SBarry Smith (F)->info.factor_mallocs = reallocs; 207*719d5645SBarry Smith (F)->info.fill_ratio_given = f; 20810c27e3dSHong Zhang if (ai[mbs] != 0) { 209*719d5645SBarry Smith (F)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]); 21010c27e3dSHong Zhang } else { 211*719d5645SBarry Smith (F)->info.fill_ratio_needed = 0.0; 21210c27e3dSHong Zhang } 213*719d5645SBarry Smith ierr = MatSeqSBAIJSetNumericFactorization(F,perm_identity);CHKERRQ(ierr); 21410c27e3dSHong Zhang PetscFunctionReturn(0); 21510c27e3dSHong Zhang } 21610c27e3dSHong Zhang /* 21710c27e3dSHong Zhang Symbolic U^T*D*U factorization for SBAIJ format. 21810c27e3dSHong Zhang */ 21998a8e78dSHong Zhang #include "petscbt.h" 22098a8e78dSHong Zhang #include "src/mat/utils/freespace.h" 2214a2ae208SSatish Balay #undef __FUNCT__ 2224a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ" 223*719d5645SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat fact,Mat A,IS perm,MatFactorInfo *info) 22449b5e25fSSatish Balay { 22598a8e78dSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 22698a8e78dSHong Zhang Mat_SeqSBAIJ *b; 2276849ba73SBarry Smith PetscErrorCode ierr; 22858ebbce7SBarry Smith PetscTruth perm_identity,missing; 22998a8e78dSHong Zhang PetscReal fill = info->fill; 230d0f46423SBarry Smith PetscInt *rip,i,mbs=a->mbs,bs=A->rmap->bs,*ai,*aj,reallocs=0,prow,d; 23198a8e78dSHong Zhang PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 2327625dc9aSHong Zhang PetscInt nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr; 233a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 23498a8e78dSHong Zhang PetscBT lnkbt; 23549b5e25fSSatish Balay 23649b5e25fSSatish Balay PetscFunctionBegin; 23758ebbce7SBarry Smith ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr); 23858ebbce7SBarry Smith if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d); 23958ebbce7SBarry Smith 24010c27e3dSHong Zhang /* 24110c27e3dSHong Zhang This code originally uses Modified Sparse Row (MSR) storage 24210c27e3dSHong Zhang (see page 85, "Iterative Methods ..." by Saad) for the output matrix B - bad choise! 24310c27e3dSHong Zhang Then it is rewritten so the factor B takes seqsbaij format. However the associated 24410c27e3dSHong Zhang MatCholeskyFactorNumeric_() have not been modified for the cases of bs>1 or !perm_identity, 24510c27e3dSHong Zhang thus the original code in MSR format is still used for these cases. 24610c27e3dSHong Zhang The code below should replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR() whenever 24710c27e3dSHong Zhang MatCholeskyFactorNumeric_() is modified for using sbaij symbolic factor. 24810c27e3dSHong Zhang */ 249fff829cfSHong Zhang if (bs > 1){ 250*719d5645SBarry Smith ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(fact,A,perm,info);CHKERRQ(ierr); 251fff829cfSHong Zhang PetscFunctionReturn(0); 252fff829cfSHong Zhang } 25310c27e3dSHong Zhang 25498a8e78dSHong Zhang /* check whether perm is the identity mapping */ 25598a8e78dSHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 25698a8e78dSHong Zhang 257fff829cfSHong Zhang if (perm_identity){ 258fff829cfSHong Zhang a->permute = PETSC_FALSE; 259fff829cfSHong Zhang ai = a->i; aj = a->j; 260fff829cfSHong Zhang } else { 261d33dd558SHong Zhang SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format"); 262d33dd558SHong Zhang /* There are bugs for reordeing. Needs further work. 263d33dd558SHong Zhang MatReordering for sbaij cannot be efficient. User should use aij formt! */ 264fff829cfSHong Zhang a->permute = PETSC_TRUE; 265fff829cfSHong Zhang ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr); 266fff829cfSHong Zhang ai = a->inew; aj = a->jnew; 267fff829cfSHong Zhang } 268fff829cfSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 269fff829cfSHong Zhang 270fff829cfSHong Zhang /* initialization */ 2717625dc9aSHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 27298a8e78dSHong Zhang ui[0] = 0; 27398a8e78dSHong Zhang 27498a8e78dSHong Zhang /* jl: linked list for storing indices of the pivot rows 2757625dc9aSHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */ 2767625dc9aSHong Zhang ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr); 2777625dc9aSHong Zhang il = jl + mbs; 2787625dc9aSHong Zhang cols = il + mbs; 2797625dc9aSHong Zhang ui_ptr = (PetscInt**)(cols + mbs); 2807625dc9aSHong Zhang 2817625dc9aSHong Zhang for (i=0; i<mbs; i++){ 2827625dc9aSHong Zhang jl[i] = mbs; il[i] = 0; 283fff829cfSHong Zhang } 284fff829cfSHong Zhang 28598a8e78dSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 2867625dc9aSHong Zhang nlnk = mbs + 1; 2877625dc9aSHong Zhang ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 288fff829cfSHong Zhang 2897625dc9aSHong Zhang /* initial FreeSpace size is fill*(ai[mbs]+1) */ 290a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr); 29198a8e78dSHong Zhang current_space = free_space; 29298a8e78dSHong Zhang 2937625dc9aSHong Zhang for (k=0; k<mbs; k++){ /* for each active row k */ 29498a8e78dSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 29598a8e78dSHong Zhang nzk = 0; 29698a8e78dSHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 2977625dc9aSHong Zhang for (j=0; j<ncols; j++){ 2987625dc9aSHong Zhang i = *(aj + ai[rip[k]] + j); 2997625dc9aSHong Zhang cols[j] = rip[i]; 3007625dc9aSHong Zhang } 3017625dc9aSHong Zhang ierr = PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 30298a8e78dSHong Zhang nzk += nlnk; 30398a8e78dSHong Zhang 30498a8e78dSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 30598a8e78dSHong Zhang prow = jl[k]; /* 1st pivot row */ 306fff829cfSHong Zhang 307fff829cfSHong Zhang while (prow < k){ 308fff829cfSHong Zhang nextprow = jl[prow]; 30998a8e78dSHong Zhang /* merge prow into k-th row */ 3107625dc9aSHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */ 31198a8e78dSHong Zhang jmax = ui[prow+1]; 31298a8e78dSHong Zhang ncols = jmax-jmin; 3137625dc9aSHong Zhang uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */ 3145a8e39fbSHong Zhang ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 31598a8e78dSHong Zhang nzk += nlnk; 316fff829cfSHong Zhang 31798a8e78dSHong Zhang /* update il and jl for prow */ 318fff829cfSHong Zhang if (jmin < jmax){ 319fff829cfSHong Zhang il[prow] = jmin; 3207625dc9aSHong Zhang j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 321fff829cfSHong Zhang } 322fff829cfSHong Zhang prow = nextprow; 323fff829cfSHong Zhang } 324fff829cfSHong Zhang 32598a8e78dSHong Zhang /* if free space is not available, make more free space */ 32698a8e78dSHong Zhang if (current_space->local_remaining<nzk) { 3277625dc9aSHong Zhang i = mbs - k + 1; /* num of unfactored rows */ 32898a8e78dSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 329a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 33098a8e78dSHong Zhang reallocs++; 33198a8e78dSHong Zhang } 33298a8e78dSHong Zhang 33398a8e78dSHong Zhang /* copy data into free space, then initialize lnk */ 3347625dc9aSHong Zhang ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 33598a8e78dSHong Zhang 33698a8e78dSHong Zhang /* add the k-th row into il and jl */ 33798a8e78dSHong Zhang if (nzk-1 > 0){ 3387625dc9aSHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */ 339fff829cfSHong Zhang jl[k] = jl[i]; jl[i] = k; 34098a8e78dSHong Zhang il[k] = ui[k] + 1; 341fff829cfSHong Zhang } 3427625dc9aSHong Zhang ui_ptr[k] = current_space->array; 34398a8e78dSHong Zhang current_space->array += nzk; 34498a8e78dSHong Zhang current_space->local_used += nzk; 34598a8e78dSHong Zhang current_space->local_remaining -= nzk; 346fff829cfSHong Zhang 34798a8e78dSHong Zhang ui[k+1] = ui[k] + nzk; 348fff829cfSHong Zhang } 349fff829cfSHong Zhang 3506cf91177SBarry Smith #if defined(PETSC_USE_INFO) 3517625dc9aSHong Zhang if (ai[mbs] != 0) { 3527625dc9aSHong Zhang PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 353ae15b995SBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 354ae15b995SBarry Smith ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 355ae15b995SBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 356fff829cfSHong Zhang } else { 357ae15b995SBarry Smith ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 358fff829cfSHong Zhang } 35963ba0a88SBarry Smith #endif 360fff829cfSHong Zhang 361fff829cfSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 362fff829cfSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 363fff829cfSHong Zhang 36498a8e78dSHong Zhang /* destroy list of free space and other temporary array(s) */ 3657625dc9aSHong Zhang ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 366a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 36798a8e78dSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 368fff829cfSHong Zhang 36998a8e78dSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 370*719d5645SBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(fact,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 37198a8e78dSHong Zhang 372*719d5645SBarry Smith b = (Mat_SeqSBAIJ*)(fact)->data; 373fff829cfSHong Zhang b->singlemalloc = PETSC_FALSE; 374e6b907acSBarry Smith b->free_a = PETSC_TRUE; 375e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 3767625dc9aSHong Zhang ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 37798a8e78dSHong Zhang b->j = uj; 37898a8e78dSHong Zhang b->i = ui; 379fff829cfSHong Zhang b->diag = 0; 380fff829cfSHong Zhang b->ilen = 0; 381fff829cfSHong Zhang b->imax = 0; 382fff829cfSHong Zhang b->row = perm; 383fff829cfSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 384fff829cfSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 385fff829cfSHong Zhang b->icol = perm; 386fff829cfSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 3877625dc9aSHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 388*719d5645SBarry Smith ierr = PetscLogObjectMemory(fact,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 3897625dc9aSHong Zhang b->maxnz = b->nz = ui[mbs]; 390fff829cfSHong Zhang 391*719d5645SBarry Smith (fact)->info.factor_mallocs = reallocs; 392*719d5645SBarry Smith (fact)->info.fill_ratio_given = fill; 3937625dc9aSHong Zhang if (ai[mbs] != 0) { 394*719d5645SBarry Smith (fact)->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 395fff829cfSHong Zhang } else { 396*719d5645SBarry Smith (fact)->info.fill_ratio_needed = 0.0; 397fff829cfSHong Zhang } 398*719d5645SBarry Smith ierr = MatSeqSBAIJSetNumericFactorization(fact,perm_identity);CHKERRQ(ierr); 399064503c5SHong Zhang PetscFunctionReturn(0); 400064503c5SHong Zhang } 4014a2ae208SSatish Balay #undef __FUNCT__ 4024a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N" 403*719d5645SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat C,Mat A,MatFactorInfo *info) 40449b5e25fSSatish Balay { 4054c16a6a6SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 4064c16a6a6SHong Zhang IS perm = b->row; 4076849ba73SBarry Smith PetscErrorCode ierr; 40813f74950SBarry Smith PetscInt *perm_ptr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 40913f74950SBarry Smith PetscInt *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 410d0f46423SBarry Smith PetscInt bs=A->rmap->bs,bs2 = a->bs2,bslog = 0; 4114c16a6a6SHong Zhang MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 4124c16a6a6SHong Zhang MatScalar *u,*diag,*rtmp,*rtmp_ptr; 41328de702eSHong Zhang MatScalar *work; 41413f74950SBarry Smith PetscInt *pivots; 4154c16a6a6SHong Zhang 4164c16a6a6SHong Zhang PetscFunctionBegin; 4174c16a6a6SHong Zhang /* initialization */ 41882502324SSatish Balay ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 4194c16a6a6SHong Zhang ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr); 42013f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); 42128de702eSHong Zhang jl = il + mbs; 4224c16a6a6SHong Zhang for (i=0; i<mbs; i++) { 4234c16a6a6SHong Zhang jl[i] = mbs; il[0] = 0; 4244c16a6a6SHong Zhang } 425b0a32e0cSBarry Smith ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr); 42628de702eSHong Zhang uik = dk + bs2; 42728de702eSHong Zhang work = uik + bs2; 42813f74950SBarry Smith ierr = PetscMalloc(bs*sizeof(PetscInt),&pivots);CHKERRQ(ierr); 4294c16a6a6SHong Zhang 4304c16a6a6SHong Zhang ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr); 4314c16a6a6SHong Zhang 4324c16a6a6SHong Zhang /* check permutation */ 4334c16a6a6SHong Zhang if (!a->permute){ 4344c16a6a6SHong Zhang ai = a->i; aj = a->j; aa = a->a; 4354c16a6a6SHong Zhang } else { 4364c16a6a6SHong Zhang ai = a->inew; aj = a->jnew; 43782502324SSatish Balay ierr = PetscMalloc(bs2*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); 4384c16a6a6SHong Zhang ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); 43913f74950SBarry Smith ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr); 44013f74950SBarry Smith ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr); 4414c16a6a6SHong Zhang 442187a9f4bSHong Zhang /* flops in while loop */ 443187a9f4bSHong Zhang bslog = 2*bs*bs2; 444187a9f4bSHong Zhang 4454c16a6a6SHong Zhang for (i=0; i<mbs; i++){ 4464c16a6a6SHong Zhang jmin = ai[i]; jmax = ai[i+1]; 4474c16a6a6SHong Zhang for (j=jmin; j<jmax; j++){ 4484c16a6a6SHong Zhang while (a2anew[j] != j){ 4494c16a6a6SHong Zhang k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k; 4504c16a6a6SHong Zhang for (k1=0; k1<bs2; k1++){ 4514c16a6a6SHong Zhang dk[k1] = aa[k*bs2+k1]; 4524c16a6a6SHong Zhang aa[k*bs2+k1] = aa[j*bs2+k1]; 4534c16a6a6SHong Zhang aa[j*bs2+k1] = dk[k1]; 4544c16a6a6SHong Zhang } 4554c16a6a6SHong Zhang } 4564c16a6a6SHong Zhang /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */ 4574c16a6a6SHong Zhang if (i > aj[j]){ 4584c16a6a6SHong Zhang /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */ 4594c16a6a6SHong Zhang ap = aa + j*bs2; /* ptr to the beginning of j-th block of aa */ 4604c16a6a6SHong Zhang for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */ 4614c16a6a6SHong Zhang for (k=0; k<bs; k++){ /* j-th block of aa <- dk^T */ 4624c16a6a6SHong Zhang for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1]; 4634c16a6a6SHong Zhang } 4644c16a6a6SHong Zhang } 4654c16a6a6SHong Zhang } 4664c16a6a6SHong Zhang } 467323b833fSBarry Smith ierr = PetscFree(a2anew);CHKERRQ(ierr); 4684c16a6a6SHong Zhang } 4694c16a6a6SHong Zhang 4704c16a6a6SHong Zhang /* for each row k */ 4714c16a6a6SHong Zhang for (k = 0; k<mbs; k++){ 4724c16a6a6SHong Zhang 4734c16a6a6SHong Zhang /*initialize k-th row with elements nonzero in row perm(k) of A */ 4744c16a6a6SHong Zhang jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1]; 475057f5ba7SHong Zhang 4764c16a6a6SHong Zhang ap = aa + jmin*bs2; 4774c16a6a6SHong Zhang for (j = jmin; j < jmax; j++){ 4784c16a6a6SHong Zhang vj = perm_ptr[aj[j]]; /* block col. index */ 4794c16a6a6SHong Zhang rtmp_ptr = rtmp + vj*bs2; 4804c16a6a6SHong Zhang for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++; 4814c16a6a6SHong Zhang } 4824c16a6a6SHong Zhang 4834c16a6a6SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 4844c16a6a6SHong Zhang ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 4854c16a6a6SHong Zhang i = jl[k]; /* first row to be added to k_th row */ 4864c16a6a6SHong Zhang 487057f5ba7SHong Zhang while (i < k){ 4884c16a6a6SHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 4894c16a6a6SHong Zhang 4904c16a6a6SHong Zhang /* compute multiplier */ 4914c16a6a6SHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 4924c16a6a6SHong Zhang 4934c16a6a6SHong Zhang /* uik = -inv(Di)*U_bar(i,k) */ 4944c16a6a6SHong Zhang diag = ba + i*bs2; 4954c16a6a6SHong Zhang u = ba + ili*bs2; 4964c16a6a6SHong Zhang ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 4974c16a6a6SHong Zhang Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u); 4984c16a6a6SHong Zhang 4994c16a6a6SHong Zhang /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 5004c16a6a6SHong Zhang Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u); 501187a9f4bSHong Zhang ierr = PetscLogFlops(bslog*2);CHKERRQ(ierr); 5024c16a6a6SHong Zhang 5034c16a6a6SHong Zhang /* update -U(i,k) */ 5044c16a6a6SHong Zhang ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 5054c16a6a6SHong Zhang 5064c16a6a6SHong Zhang /* add multiple of row i to k-th row ... */ 5074c16a6a6SHong Zhang jmin = ili + 1; jmax = bi[i+1]; 5084c16a6a6SHong Zhang if (jmin < jmax){ 5094c16a6a6SHong Zhang for (j=jmin; j<jmax; j++) { 5104c16a6a6SHong Zhang /* rtmp += -U(i,k)^T * U_bar(i,j) */ 5114c16a6a6SHong Zhang rtmp_ptr = rtmp + bj[j]*bs2; 5124c16a6a6SHong Zhang u = ba + j*bs2; 5134c16a6a6SHong Zhang Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u); 5144c16a6a6SHong Zhang } 515187a9f4bSHong Zhang ierr = PetscLogFlops(bslog*(jmax-jmin));CHKERRQ(ierr); 5164c16a6a6SHong Zhang 5174c16a6a6SHong Zhang /* ... add i to row list for next nonzero entry */ 5184c16a6a6SHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 5194c16a6a6SHong Zhang j = bj[jmin]; 5204c16a6a6SHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 5214c16a6a6SHong Zhang } 5224c16a6a6SHong Zhang i = nexti; 5234c16a6a6SHong Zhang } 5244c16a6a6SHong Zhang 5254c16a6a6SHong Zhang /* save nonzero entries in k-th row of U ... */ 5264c16a6a6SHong Zhang 5274c16a6a6SHong Zhang /* invert diagonal block */ 5284c16a6a6SHong Zhang diag = ba+k*bs2; 5294c16a6a6SHong Zhang ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr); 530d230e6fdSBarry Smith ierr = Kernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr); 5314c16a6a6SHong Zhang 5324c16a6a6SHong Zhang jmin = bi[k]; jmax = bi[k+1]; 5334c16a6a6SHong Zhang if (jmin < jmax) { 5344c16a6a6SHong Zhang for (j=jmin; j<jmax; j++){ 5354c16a6a6SHong Zhang vj = bj[j]; /* block col. index of U */ 5364c16a6a6SHong Zhang u = ba + j*bs2; 5374c16a6a6SHong Zhang rtmp_ptr = rtmp + vj*bs2; 5384c16a6a6SHong Zhang for (k1=0; k1<bs2; k1++){ 5394c16a6a6SHong Zhang *u++ = *rtmp_ptr; 5404c16a6a6SHong Zhang *rtmp_ptr++ = 0.0; 5414c16a6a6SHong Zhang } 5424c16a6a6SHong Zhang } 5434c16a6a6SHong Zhang 5444c16a6a6SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 5454c16a6a6SHong Zhang il[k] = jmin; 5464c16a6a6SHong Zhang i = bj[jmin]; 5474c16a6a6SHong Zhang jl[k] = jl[i]; jl[i] = k; 5484c16a6a6SHong Zhang } 5494c16a6a6SHong Zhang } 5504c16a6a6SHong Zhang 5514c16a6a6SHong Zhang ierr = PetscFree(rtmp);CHKERRQ(ierr); 5524c16a6a6SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 5534c16a6a6SHong Zhang ierr = PetscFree(dk);CHKERRQ(ierr); 5544c16a6a6SHong Zhang ierr = PetscFree(pivots);CHKERRQ(ierr); 5554c16a6a6SHong Zhang if (a->permute){ 5564c16a6a6SHong Zhang ierr = PetscFree(aa);CHKERRQ(ierr); 5574c16a6a6SHong Zhang } 5584c16a6a6SHong Zhang 5594c16a6a6SHong Zhang ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr); 560db4efbfdSBarry Smith C->ops->solve = MatSolve_SeqSBAIJ_N; 561db4efbfdSBarry Smith C->ops->solvetranspose = MatSolve_SeqSBAIJ_N; 562db4efbfdSBarry Smith C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_N; 563db4efbfdSBarry Smith C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_N; 564db4efbfdSBarry Smith 5654c16a6a6SHong Zhang C->assembled = PETSC_TRUE; 5664c16a6a6SHong Zhang C->preallocated = PETSC_TRUE; 567efee365bSSatish Balay ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 5684c16a6a6SHong Zhang PetscFunctionReturn(0); 5694c16a6a6SHong Zhang } 570d4edadd4SHong Zhang 5714a2ae208SSatish Balay #undef __FUNCT__ 5724a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering" 573*719d5645SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat C,Mat A,MatFactorInfo *info) 574671cb588SHong Zhang { 575671cb588SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 576dfbe8321SBarry Smith PetscErrorCode ierr; 57713f74950SBarry Smith PetscInt i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 57813f74950SBarry Smith PetscInt *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 579d0f46423SBarry Smith PetscInt bs=A->rmap->bs,bs2 = a->bs2,bslog; 580671cb588SHong Zhang MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 581671cb588SHong Zhang MatScalar *u,*diag,*rtmp,*rtmp_ptr; 58228de702eSHong Zhang MatScalar *work; 58313f74950SBarry Smith PetscInt *pivots; 584671cb588SHong Zhang 585671cb588SHong Zhang PetscFunctionBegin; 58682502324SSatish Balay ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 587671cb588SHong Zhang ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr); 58813f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); 58928de702eSHong Zhang jl = il + mbs; 590671cb588SHong Zhang for (i=0; i<mbs; i++) { 591671cb588SHong Zhang jl[i] = mbs; il[0] = 0; 592671cb588SHong Zhang } 593b0a32e0cSBarry Smith ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr); 59428de702eSHong Zhang uik = dk + bs2; 59528de702eSHong Zhang work = uik + bs2; 59613f74950SBarry Smith ierr = PetscMalloc(bs*sizeof(PetscInt),&pivots);CHKERRQ(ierr); 597671cb588SHong Zhang 598671cb588SHong Zhang ai = a->i; aj = a->j; aa = a->a; 599671cb588SHong Zhang 600187a9f4bSHong Zhang /* flops in while loop */ 601187a9f4bSHong Zhang bslog = 2*bs*bs2; 602187a9f4bSHong Zhang 603671cb588SHong Zhang /* for each row k */ 604671cb588SHong Zhang for (k = 0; k<mbs; k++){ 605671cb588SHong Zhang 606671cb588SHong Zhang /*initialize k-th row with elements nonzero in row k of A */ 607671cb588SHong Zhang jmin = ai[k]; jmax = ai[k+1]; 608671cb588SHong Zhang ap = aa + jmin*bs2; 609671cb588SHong Zhang for (j = jmin; j < jmax; j++){ 610671cb588SHong Zhang vj = aj[j]; /* block col. index */ 611671cb588SHong Zhang rtmp_ptr = rtmp + vj*bs2; 612671cb588SHong Zhang for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++; 613671cb588SHong Zhang } 614671cb588SHong Zhang 615671cb588SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 616671cb588SHong Zhang ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 617671cb588SHong Zhang i = jl[k]; /* first row to be added to k_th row */ 618671cb588SHong Zhang 619057f5ba7SHong Zhang while (i < k){ 620671cb588SHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 621671cb588SHong Zhang 622671cb588SHong Zhang /* compute multiplier */ 623671cb588SHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 624671cb588SHong Zhang 625671cb588SHong Zhang /* uik = -inv(Di)*U_bar(i,k) */ 626671cb588SHong Zhang diag = ba + i*bs2; 627671cb588SHong Zhang u = ba + ili*bs2; 628671cb588SHong Zhang ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 629671cb588SHong Zhang Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u); 630671cb588SHong Zhang 631671cb588SHong Zhang /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 632671cb588SHong Zhang Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u); 633187a9f4bSHong Zhang ierr = PetscLogFlops(bslog*2);CHKERRQ(ierr); 634671cb588SHong Zhang 635671cb588SHong Zhang /* update -U(i,k) */ 636671cb588SHong Zhang ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 637671cb588SHong Zhang 638671cb588SHong Zhang /* add multiple of row i to k-th row ... */ 639671cb588SHong Zhang jmin = ili + 1; jmax = bi[i+1]; 640671cb588SHong Zhang if (jmin < jmax){ 641671cb588SHong Zhang for (j=jmin; j<jmax; j++) { 642671cb588SHong Zhang /* rtmp += -U(i,k)^T * U_bar(i,j) */ 643671cb588SHong Zhang rtmp_ptr = rtmp + bj[j]*bs2; 644671cb588SHong Zhang u = ba + j*bs2; 645671cb588SHong Zhang Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u); 646671cb588SHong Zhang } 647187a9f4bSHong Zhang ierr = PetscLogFlops(bslog*(jmax-jmin));CHKERRQ(ierr); 648671cb588SHong Zhang 649671cb588SHong Zhang /* ... add i to row list for next nonzero entry */ 650671cb588SHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 651671cb588SHong Zhang j = bj[jmin]; 652671cb588SHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 653671cb588SHong Zhang } 654671cb588SHong Zhang i = nexti; 655671cb588SHong Zhang } 656671cb588SHong Zhang 657671cb588SHong Zhang /* save nonzero entries in k-th row of U ... */ 658671cb588SHong Zhang 659671cb588SHong Zhang /* invert diagonal block */ 660671cb588SHong Zhang diag = ba+k*bs2; 661671cb588SHong Zhang ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr); 662d230e6fdSBarry Smith ierr = Kernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr); 663671cb588SHong Zhang 664671cb588SHong Zhang jmin = bi[k]; jmax = bi[k+1]; 665671cb588SHong Zhang if (jmin < jmax) { 666671cb588SHong Zhang for (j=jmin; j<jmax; j++){ 667671cb588SHong Zhang vj = bj[j]; /* block col. index of U */ 668671cb588SHong Zhang u = ba + j*bs2; 669671cb588SHong Zhang rtmp_ptr = rtmp + vj*bs2; 670671cb588SHong Zhang for (k1=0; k1<bs2; k1++){ 671671cb588SHong Zhang *u++ = *rtmp_ptr; 672671cb588SHong Zhang *rtmp_ptr++ = 0.0; 673671cb588SHong Zhang } 674671cb588SHong Zhang } 675671cb588SHong Zhang 676671cb588SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 677671cb588SHong Zhang il[k] = jmin; 678671cb588SHong Zhang i = bj[jmin]; 679671cb588SHong Zhang jl[k] = jl[i]; jl[i] = k; 680671cb588SHong Zhang } 681671cb588SHong Zhang } 682671cb588SHong Zhang 683671cb588SHong Zhang ierr = PetscFree(rtmp);CHKERRQ(ierr); 684671cb588SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 685671cb588SHong Zhang ierr = PetscFree(dk);CHKERRQ(ierr); 686671cb588SHong Zhang ierr = PetscFree(pivots);CHKERRQ(ierr); 687671cb588SHong Zhang 688db4efbfdSBarry Smith C->ops->solve = MatSolve_SeqSBAIJ_N_NaturalOrdering; 689db4efbfdSBarry Smith C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_NaturalOrdering; 690db4efbfdSBarry Smith C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering; 691db4efbfdSBarry Smith C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering; 692671cb588SHong Zhang C->assembled = PETSC_TRUE; 693671cb588SHong Zhang C->preallocated = PETSC_TRUE; 694efee365bSSatish Balay ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 695671cb588SHong Zhang PetscFunctionReturn(0); 696671cb588SHong Zhang } 697671cb588SHong Zhang 69849b5e25fSSatish Balay /* 699fcf159c0SHong Zhang Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP. 700cc0c071aSHong Zhang Version for blocks 2 by 2. 70149b5e25fSSatish Balay */ 7024a2ae208SSatish Balay #undef __FUNCT__ 7034a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2" 704*719d5645SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat C,Mat A,MatFactorInfo *info) 70549b5e25fSSatish Balay { 70691602c66SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 707cc0c071aSHong Zhang IS perm = b->row; 7086849ba73SBarry Smith PetscErrorCode ierr; 70913f74950SBarry Smith PetscInt *perm_ptr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 71013f74950SBarry Smith PetscInt *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 711a1723e09SHong Zhang MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 712cc0c071aSHong Zhang MatScalar *u,*diag,*rtmp,*rtmp_ptr; 71362bba022SBarry Smith PetscReal shift = info->shiftinblocks; 71449b5e25fSSatish Balay 71549b5e25fSSatish Balay PetscFunctionBegin; 71691602c66SHong Zhang /* initialization */ 71791602c66SHong Zhang /* il and jl record the first nonzero element in each row of the accessing 71891602c66SHong Zhang window U(0:k, k:mbs-1). 71991602c66SHong Zhang jl: list of rows to be added to uneliminated rows 72091602c66SHong Zhang i>= k: jl(i) is the first row to be added to row i 72191602c66SHong Zhang i< k: jl(i) is the row following row i in some list of rows 72291602c66SHong Zhang jl(i) = mbs indicates the end of a list 72391602c66SHong Zhang il(i): points to the first nonzero element in columns k,...,mbs-1 of 72491602c66SHong Zhang row i of U */ 72582502324SSatish Balay ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 726cc0c071aSHong Zhang ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr); 72713f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); 72828de702eSHong Zhang jl = il + mbs; 72991602c66SHong Zhang for (i=0; i<mbs; i++) { 7303845f261SHong Zhang jl[i] = mbs; il[0] = 0; 73191602c66SHong Zhang } 73282502324SSatish Balay ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr); 73328de702eSHong Zhang uik = dk + 4; 734cc0c071aSHong Zhang ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr); 735cc0c071aSHong Zhang 736cc0c071aSHong Zhang /* check permutation */ 737cc0c071aSHong Zhang if (!a->permute){ 738cc0c071aSHong Zhang ai = a->i; aj = a->j; aa = a->a; 739cc0c071aSHong Zhang } else { 740cc0c071aSHong Zhang ai = a->inew; aj = a->jnew; 74182502324SSatish Balay ierr = PetscMalloc(4*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); 742cc0c071aSHong Zhang ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); 74313f74950SBarry Smith ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr); 74413f74950SBarry Smith ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr); 745cc0c071aSHong Zhang 746cc0c071aSHong Zhang for (i=0; i<mbs; i++){ 747cc0c071aSHong Zhang jmin = ai[i]; jmax = ai[i+1]; 748cc0c071aSHong Zhang for (j=jmin; j<jmax; j++){ 749cc0c071aSHong Zhang while (a2anew[j] != j){ 750cc0c071aSHong Zhang k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k; 751cc0c071aSHong Zhang for (k1=0; k1<4; k1++){ 752cc0c071aSHong Zhang dk[k1] = aa[k*4+k1]; 753cc0c071aSHong Zhang aa[k*4+k1] = aa[j*4+k1]; 754cc0c071aSHong Zhang aa[j*4+k1] = dk[k1]; 755cc0c071aSHong Zhang } 756cc0c071aSHong Zhang } 757cc0c071aSHong Zhang /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */ 758cc0c071aSHong Zhang if (i > aj[j]){ 759a1723e09SHong Zhang /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */ 760cc0c071aSHong Zhang ap = aa + j*4; /* ptr to the beginning of the block */ 761cc0c071aSHong Zhang dk[1] = ap[1]; /* swap ap[1] and ap[2] */ 762cc0c071aSHong Zhang ap[1] = ap[2]; 763cc0c071aSHong Zhang ap[2] = dk[1]; 764cc0c071aSHong Zhang } 765cc0c071aSHong Zhang } 766cc0c071aSHong Zhang } 767ac355199SBarry Smith ierr = PetscFree(a2anew);CHKERRQ(ierr); 768cc0c071aSHong Zhang } 7693845f261SHong Zhang 77091602c66SHong Zhang /* for each row k */ 77191602c66SHong Zhang for (k = 0; k<mbs; k++){ 77291602c66SHong Zhang 77391602c66SHong Zhang /*initialize k-th row with elements nonzero in row perm(k) of A */ 774cc0c071aSHong Zhang jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1]; 775cc0c071aSHong Zhang ap = aa + jmin*4; 77691602c66SHong Zhang for (j = jmin; j < jmax; j++){ 777cc0c071aSHong Zhang vj = perm_ptr[aj[j]]; /* block col. index */ 778cc0c071aSHong Zhang rtmp_ptr = rtmp + vj*4; 779cc0c071aSHong Zhang for (i=0; i<4; i++) *rtmp_ptr++ = *ap++; 78091602c66SHong Zhang } 78191602c66SHong Zhang 78291602c66SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 783cc0c071aSHong Zhang ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr); 78491602c66SHong Zhang i = jl[k]; /* first row to be added to k_th row */ 78591602c66SHong Zhang 786057f5ba7SHong Zhang while (i < k){ 78791602c66SHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 78891602c66SHong Zhang 7893845f261SHong Zhang /* compute multiplier */ 79091602c66SHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 7913845f261SHong Zhang 7923845f261SHong Zhang /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */ 793cc0c071aSHong Zhang diag = ba + i*4; 794cc0c071aSHong Zhang u = ba + ili*4; 795cc0c071aSHong Zhang uik[0] = -(diag[0]*u[0] + diag[2]*u[1]); 796cc0c071aSHong Zhang uik[1] = -(diag[1]*u[0] + diag[3]*u[1]); 797cc0c071aSHong Zhang uik[2] = -(diag[0]*u[2] + diag[2]*u[3]); 798cc0c071aSHong Zhang uik[3] = -(diag[1]*u[2] + diag[3]*u[3]); 7993845f261SHong Zhang 8003845f261SHong Zhang /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */ 801cc0c071aSHong Zhang dk[0] += uik[0]*u[0] + uik[1]*u[1]; 802cc0c071aSHong Zhang dk[1] += uik[2]*u[0] + uik[3]*u[1]; 803cc0c071aSHong Zhang dk[2] += uik[0]*u[2] + uik[1]*u[3]; 804cc0c071aSHong Zhang dk[3] += uik[2]*u[2] + uik[3]*u[3]; 8053845f261SHong Zhang 806187a9f4bSHong Zhang ierr = PetscLogFlops(16*2);CHKERRQ(ierr); 807187a9f4bSHong Zhang 8083845f261SHong Zhang /* update -U(i,k): ba[ili] = uik */ 809cc0c071aSHong Zhang ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr); 81091602c66SHong Zhang 81191602c66SHong Zhang /* add multiple of row i to k-th row ... */ 81291602c66SHong Zhang jmin = ili + 1; jmax = bi[i+1]; 81391602c66SHong Zhang if (jmin < jmax){ 8143845f261SHong Zhang for (j=jmin; j<jmax; j++) { 8153845f261SHong Zhang /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */ 816cc0c071aSHong Zhang rtmp_ptr = rtmp + bj[j]*4; 817cc0c071aSHong Zhang u = ba + j*4; 818cc0c071aSHong Zhang rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1]; 819cc0c071aSHong Zhang rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1]; 820cc0c071aSHong Zhang rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3]; 821cc0c071aSHong Zhang rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3]; 8223845f261SHong Zhang } 823187a9f4bSHong Zhang ierr = PetscLogFlops(16*(jmax-jmin));CHKERRQ(ierr); 8243845f261SHong Zhang 82591602c66SHong Zhang /* ... add i to row list for next nonzero entry */ 82691602c66SHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 82791602c66SHong Zhang j = bj[jmin]; 82891602c66SHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 82991602c66SHong Zhang } 830a1723e09SHong Zhang i = nexti; 83191602c66SHong Zhang } 832cc0c071aSHong Zhang 83391602c66SHong Zhang /* save nonzero entries in k-th row of U ... */ 8343845f261SHong Zhang 835cc0c071aSHong Zhang /* invert diagonal block */ 836cc0c071aSHong Zhang diag = ba+k*4; 837cc0c071aSHong Zhang ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr); 83862bba022SBarry Smith ierr = Kernel_A_gets_inverse_A_2(diag,shift);CHKERRQ(ierr); 8393845f261SHong Zhang 84091602c66SHong Zhang jmin = bi[k]; jmax = bi[k+1]; 84191602c66SHong Zhang if (jmin < jmax) { 84291602c66SHong Zhang for (j=jmin; j<jmax; j++){ 843cc0c071aSHong Zhang vj = bj[j]; /* block col. index of U */ 844cc0c071aSHong Zhang u = ba + j*4; 845cc0c071aSHong Zhang rtmp_ptr = rtmp + vj*4; 846cc0c071aSHong Zhang for (k1=0; k1<4; k1++){ 847cc0c071aSHong Zhang *u++ = *rtmp_ptr; 848cc0c071aSHong Zhang *rtmp_ptr++ = 0.0; 849cc0c071aSHong Zhang } 850cc0c071aSHong Zhang } 8513845f261SHong Zhang 85291602c66SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 85391602c66SHong Zhang il[k] = jmin; 85491602c66SHong Zhang i = bj[jmin]; 85591602c66SHong Zhang jl[k] = jl[i]; jl[i] = k; 85691602c66SHong Zhang } 85791602c66SHong Zhang } 8583845f261SHong Zhang 85949b5e25fSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 86091602c66SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 8613845f261SHong Zhang ierr = PetscFree(dk);CHKERRQ(ierr); 86291602c66SHong Zhang if (a->permute) { 86391602c66SHong Zhang ierr = PetscFree(aa);CHKERRQ(ierr); 86491602c66SHong Zhang } 865cc0c071aSHong Zhang ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr); 866db4efbfdSBarry Smith C->ops->solve = MatSolve_SeqSBAIJ_2; 867db4efbfdSBarry Smith C->ops->solvetranspose = MatSolve_SeqSBAIJ_2; 86849b5e25fSSatish Balay C->assembled = PETSC_TRUE; 8695c0bcdfcSHong Zhang C->preallocated = PETSC_TRUE; 870efee365bSSatish Balay ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 87149b5e25fSSatish Balay PetscFunctionReturn(0); 87249b5e25fSSatish Balay } 87391602c66SHong Zhang 87449b5e25fSSatish Balay /* 87549b5e25fSSatish Balay Version for when blocks are 2 by 2 Using natural ordering 87649b5e25fSSatish Balay */ 8774a2ae208SSatish Balay #undef __FUNCT__ 8784a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering" 879*719d5645SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat C,Mat A,MatFactorInfo *info) 88049b5e25fSSatish Balay { 881ab27746eSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 882dfbe8321SBarry Smith PetscErrorCode ierr; 88313f74950SBarry Smith PetscInt i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 88413f74950SBarry Smith PetscInt *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 885ab27746eSHong Zhang MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 886ab27746eSHong Zhang MatScalar *u,*diag,*rtmp,*rtmp_ptr; 88762bba022SBarry Smith PetscReal shift = info->shiftinblocks; 88849b5e25fSSatish Balay 88949b5e25fSSatish Balay PetscFunctionBegin; 890ab27746eSHong Zhang /* initialization */ 891ab27746eSHong Zhang /* il and jl record the first nonzero element in each row of the accessing 892ab27746eSHong Zhang window U(0:k, k:mbs-1). 893ab27746eSHong Zhang jl: list of rows to be added to uneliminated rows 894ab27746eSHong Zhang i>= k: jl(i) is the first row to be added to row i 895ab27746eSHong Zhang i< k: jl(i) is the row following row i in some list of rows 896ab27746eSHong Zhang jl(i) = mbs indicates the end of a list 897ab27746eSHong Zhang il(i): points to the first nonzero element in columns k,...,mbs-1 of 898ab27746eSHong Zhang row i of U */ 89982502324SSatish Balay ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 900ab27746eSHong Zhang ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr); 90113f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); 90228de702eSHong Zhang jl = il + mbs; 903ab27746eSHong Zhang for (i=0; i<mbs; i++) { 904ab27746eSHong Zhang jl[i] = mbs; il[0] = 0; 90549b5e25fSSatish Balay } 90682502324SSatish Balay ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr); 90728de702eSHong Zhang uik = dk + 4; 908ab27746eSHong Zhang 909ab27746eSHong Zhang ai = a->i; aj = a->j; aa = a->a; 910ab27746eSHong Zhang 911ab27746eSHong Zhang /* for each row k */ 912ab27746eSHong Zhang for (k = 0; k<mbs; k++){ 913ab27746eSHong Zhang 914ab27746eSHong Zhang /*initialize k-th row with elements nonzero in row k of A */ 915ab27746eSHong Zhang jmin = ai[k]; jmax = ai[k+1]; 916ab27746eSHong Zhang ap = aa + jmin*4; 917ab27746eSHong Zhang for (j = jmin; j < jmax; j++){ 918ab27746eSHong Zhang vj = aj[j]; /* block col. index */ 919ab27746eSHong Zhang rtmp_ptr = rtmp + vj*4; 920ab27746eSHong Zhang for (i=0; i<4; i++) *rtmp_ptr++ = *ap++; 92149b5e25fSSatish Balay } 922ab27746eSHong Zhang 923ab27746eSHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 924ab27746eSHong Zhang ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr); 925ab27746eSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 926ab27746eSHong Zhang 927057f5ba7SHong Zhang while (i < k){ 928ab27746eSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 929ab27746eSHong Zhang 930ab27746eSHong Zhang /* compute multiplier */ 931ab27746eSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 932ab27746eSHong Zhang 933ab27746eSHong Zhang /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */ 934ab27746eSHong Zhang diag = ba + i*4; 935ab27746eSHong Zhang u = ba + ili*4; 936ab27746eSHong Zhang uik[0] = -(diag[0]*u[0] + diag[2]*u[1]); 937ab27746eSHong Zhang uik[1] = -(diag[1]*u[0] + diag[3]*u[1]); 938ab27746eSHong Zhang uik[2] = -(diag[0]*u[2] + diag[2]*u[3]); 939ab27746eSHong Zhang uik[3] = -(diag[1]*u[2] + diag[3]*u[3]); 940ab27746eSHong Zhang 941ab27746eSHong Zhang /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */ 942ab27746eSHong Zhang dk[0] += uik[0]*u[0] + uik[1]*u[1]; 943ab27746eSHong Zhang dk[1] += uik[2]*u[0] + uik[3]*u[1]; 944ab27746eSHong Zhang dk[2] += uik[0]*u[2] + uik[1]*u[3]; 945ab27746eSHong Zhang dk[3] += uik[2]*u[2] + uik[3]*u[3]; 946ab27746eSHong Zhang 947187a9f4bSHong Zhang ierr = PetscLogFlops(16*2);CHKERRQ(ierr); 948187a9f4bSHong Zhang 949ab27746eSHong Zhang /* update -U(i,k): ba[ili] = uik */ 950ab27746eSHong Zhang ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr); 951ab27746eSHong Zhang 952ab27746eSHong Zhang /* add multiple of row i to k-th row ... */ 953ab27746eSHong Zhang jmin = ili + 1; jmax = bi[i+1]; 954ab27746eSHong Zhang if (jmin < jmax){ 955ab27746eSHong Zhang for (j=jmin; j<jmax; j++) { 956ab27746eSHong Zhang /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */ 957ab27746eSHong Zhang rtmp_ptr = rtmp + bj[j]*4; 958ab27746eSHong Zhang u = ba + j*4; 959ab27746eSHong Zhang rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1]; 960ab27746eSHong Zhang rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1]; 961ab27746eSHong Zhang rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3]; 962ab27746eSHong Zhang rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3]; 96349b5e25fSSatish Balay } 964187a9f4bSHong Zhang ierr = PetscLogFlops(16*(jmax-jmin));CHKERRQ(ierr); 965ab27746eSHong Zhang 966ab27746eSHong Zhang /* ... add i to row list for next nonzero entry */ 967ab27746eSHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 968ab27746eSHong Zhang j = bj[jmin]; 969ab27746eSHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 97049b5e25fSSatish Balay } 971ab27746eSHong Zhang i = nexti; 97249b5e25fSSatish Balay } 973ab27746eSHong Zhang 974ab27746eSHong Zhang /* save nonzero entries in k-th row of U ... */ 975ab27746eSHong Zhang 97649b5e25fSSatish Balay /* invert diagonal block */ 977ab27746eSHong Zhang diag = ba+k*4; 978ab27746eSHong Zhang ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr); 97962bba022SBarry Smith ierr = Kernel_A_gets_inverse_A_2(diag,shift);CHKERRQ(ierr); 980ab27746eSHong Zhang 981ab27746eSHong Zhang jmin = bi[k]; jmax = bi[k+1]; 982ab27746eSHong Zhang if (jmin < jmax) { 983ab27746eSHong Zhang for (j=jmin; j<jmax; j++){ 984ab27746eSHong Zhang vj = bj[j]; /* block col. index of U */ 985ab27746eSHong Zhang u = ba + j*4; 986ab27746eSHong Zhang rtmp_ptr = rtmp + vj*4; 987ab27746eSHong Zhang for (k1=0; k1<4; k1++){ 988ab27746eSHong Zhang *u++ = *rtmp_ptr; 989ab27746eSHong Zhang *rtmp_ptr++ = 0.0; 990ab27746eSHong Zhang } 991ab27746eSHong Zhang } 992ab27746eSHong Zhang 993ab27746eSHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 994ab27746eSHong Zhang il[k] = jmin; 995ab27746eSHong Zhang i = bj[jmin]; 996ab27746eSHong Zhang jl[k] = jl[i]; jl[i] = k; 997ab27746eSHong Zhang } 99849b5e25fSSatish Balay } 99949b5e25fSSatish Balay 100049b5e25fSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 1001ab27746eSHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 1002ab27746eSHong Zhang ierr = PetscFree(dk);CHKERRQ(ierr); 1003ab27746eSHong Zhang 1004db4efbfdSBarry Smith C->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 1005db4efbfdSBarry Smith C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; 1006db4efbfdSBarry Smith C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering; 1007db4efbfdSBarry Smith C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering; 100849b5e25fSSatish Balay C->assembled = PETSC_TRUE; 10095c0bcdfcSHong Zhang C->preallocated = PETSC_TRUE; 1010efee365bSSatish Balay ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 101149b5e25fSSatish Balay PetscFunctionReturn(0); 101249b5e25fSSatish Balay } 101349b5e25fSSatish Balay 101449b5e25fSSatish Balay /* 101598a8e78dSHong Zhang Numeric U^T*D*U factorization for SBAIJ format. 101691602c66SHong Zhang Version for blocks are 1 by 1. 101749b5e25fSSatish Balay */ 10184a2ae208SSatish Balay #undef __FUNCT__ 10194a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1" 1020*719d5645SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat C,Mat A,MatFactorInfo *info) 102149b5e25fSSatish Balay { 102249b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data; 102349b5e25fSSatish Balay IS ip=b->row; 10246849ba73SBarry Smith PetscErrorCode ierr; 1025997a0949SHong Zhang PetscInt *rip,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j,*bcol; 1026fff829cfSHong Zhang PetscInt *ai,*aj,*a2anew; 1027997a0949SHong Zhang PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 1028997a0949SHong Zhang MatScalar *rtmp,*ba=b->a,*bval,*aa,dk,uikdi; 10293cea4cbeSHong Zhang PetscReal zeropivot,rs,shiftnz; 1030fbf22428SSatish Balay PetscReal shiftpd; 10313cea4cbeSHong Zhang ChShift_Ctx sctx; 10323cea4cbeSHong Zhang PetscInt newshift; 103349b5e25fSSatish Balay 103449b5e25fSSatish Balay PetscFunctionBegin; 10353cea4cbeSHong Zhang /* initialization */ 10363cea4cbeSHong Zhang shiftnz = info->shiftnz; 10373cea4cbeSHong Zhang shiftpd = info->shiftpd; 10383cea4cbeSHong Zhang zeropivot = info->zeropivot; 10393cea4cbeSHong Zhang 104049b5e25fSSatish Balay ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 1041cb718733SHong Zhang if (!a->permute){ 10422d007305SHong Zhang ai = a->i; aj = a->j; aa = a->a; 10432d007305SHong Zhang } else { 10442d007305SHong Zhang ai = a->inew; aj = a->jnew; 1045fff829cfSHong Zhang nz = ai[mbs]; 1046fff829cfSHong Zhang ierr = PetscMalloc(nz*sizeof(MatScalar),&aa);CHKERRQ(ierr); 1047fff829cfSHong Zhang a2anew = a->a2anew; 1048997a0949SHong Zhang bval = a->a; 1049fff829cfSHong Zhang for (j=0; j<nz; j++){ 1050997a0949SHong Zhang aa[a2anew[j]] = *(bval++); 10512d007305SHong Zhang } 10522d007305SHong Zhang } 10532d007305SHong Zhang 105491602c66SHong Zhang /* initialization */ 105549b5e25fSSatish Balay /* il and jl record the first nonzero element in each row of the accessing 105649b5e25fSSatish Balay window U(0:k, k:mbs-1). 105749b5e25fSSatish Balay jl: list of rows to be added to uneliminated rows 105849b5e25fSSatish Balay i>= k: jl(i) is the first row to be added to row i 105949b5e25fSSatish Balay i< k: jl(i) is the row following row i in some list of rows 106049b5e25fSSatish Balay jl(i) = mbs indicates the end of a list 106149b5e25fSSatish Balay il(i): points to the first nonzero element in columns k,...,mbs-1 of 106249b5e25fSSatish Balay row i of U */ 1063997a0949SHong Zhang nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); 1064997a0949SHong Zhang ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 106528de702eSHong Zhang jl = il + mbs; 1066997a0949SHong Zhang rtmp = (MatScalar*)(jl + mbs); 1067997a0949SHong Zhang 10683cea4cbeSHong Zhang sctx.shift_amount = 0; 10693cea4cbeSHong Zhang sctx.nshift = 0; 1070997a0949SHong Zhang do { 10713cea4cbeSHong Zhang sctx.chshift = PETSC_FALSE; 107249b5e25fSSatish Balay for (i=0; i<mbs; i++) { 107349b5e25fSSatish Balay rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 107449b5e25fSSatish Balay } 1075997a0949SHong Zhang 107649b5e25fSSatish Balay for (k = 0; k<mbs; k++){ 1077997a0949SHong Zhang /*initialize k-th row by the perm[k]-th row of A */ 107849b5e25fSSatish Balay jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 10797701de02SHong Zhang bval = ba + bi[k]; 108049b5e25fSSatish Balay for (j = jmin; j < jmax; j++){ 1081997a0949SHong Zhang col = rip[aj[j]]; 1082997a0949SHong Zhang rtmp[col] = aa[j]; 10837701de02SHong Zhang *bval++ = 0.0; /* for in-place factorization */ 108449b5e25fSSatish Balay } 10853cea4cbeSHong Zhang 10863cea4cbeSHong Zhang /* shift the diagonal of the matrix */ 10873cea4cbeSHong Zhang if (sctx.nshift) rtmp[k] += sctx.shift_amount; 108849b5e25fSSatish Balay 108991602c66SHong Zhang /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 109049b5e25fSSatish Balay dk = rtmp[k]; 109149b5e25fSSatish Balay i = jl[k]; /* first row to be added to k_th row */ 109249b5e25fSSatish Balay 1093057f5ba7SHong Zhang while (i < k){ 109449b5e25fSSatish Balay nexti = jl[i]; /* next row to be added to k_th row */ 109549b5e25fSSatish Balay 1096fff829cfSHong Zhang /* compute multiplier, update diag(k) and U(i,k) */ 109749b5e25fSSatish Balay ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 1098997a0949SHong Zhang uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 109949b5e25fSSatish Balay dk += uikdi*ba[ili]; 1100658e7b3eSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 110149b5e25fSSatish Balay 1102997a0949SHong Zhang /* add multiple of row i to k-th row */ 110349b5e25fSSatish Balay jmin = ili + 1; jmax = bi[i+1]; 110449b5e25fSSatish Balay if (jmin < jmax){ 110549b5e25fSSatish Balay for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 1106187a9f4bSHong Zhang ierr = PetscLogFlops(2*(jmax-jmin));CHKERRQ(ierr); 1107187a9f4bSHong Zhang 1108fff829cfSHong Zhang /* update il and jl for row i */ 1109fff829cfSHong Zhang il[i] = jmin; 1110fff829cfSHong Zhang j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 111149b5e25fSSatish Balay } 1112ab27746eSHong Zhang i = nexti; 111349b5e25fSSatish Balay } 111449b5e25fSSatish Balay 11153cea4cbeSHong Zhang /* shift the diagonals when zero pivot is detected */ 11163cea4cbeSHong Zhang /* compute rs=sum of abs(off-diagonal) */ 11173cea4cbeSHong Zhang rs = 0.0; 1118997a0949SHong Zhang jmin = bi[k]+1; 1119997a0949SHong Zhang nz = bi[k+1] - jmin; 1120997a0949SHong Zhang if (nz){ 1121997a0949SHong Zhang bcol = bj + jmin; 1122997a0949SHong Zhang while (nz--){ 112389f845c8SHong Zhang rs += PetscAbsScalar(rtmp[*bcol]); 112489f845c8SHong Zhang bcol++; 1125997a0949SHong Zhang } 1126997a0949SHong Zhang } 11273cea4cbeSHong Zhang 11283cea4cbeSHong Zhang sctx.rs = rs; 11293cea4cbeSHong Zhang sctx.pv = dk; 113045938b79SHong Zhang ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr); 113145938b79SHong Zhang if (newshift == 1) break; /* sctx.shift_amount is updated */ 113249b5e25fSSatish Balay 1133997a0949SHong Zhang /* copy data into U(k,:) */ 113498a8e78dSHong Zhang ba[bi[k]] = 1.0/dk; /* U(k,k) */ 1135fff829cfSHong Zhang jmin = bi[k]+1; jmax = bi[k+1]; 113649b5e25fSSatish Balay if (jmin < jmax) { 113749b5e25fSSatish Balay for (j=jmin; j<jmax; j++){ 1138997a0949SHong Zhang col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 113949b5e25fSSatish Balay } 114098a8e78dSHong Zhang /* add the k-th row into il and jl */ 114149b5e25fSSatish Balay il[k] = jmin; 114298a8e78dSHong Zhang i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 114349b5e25fSSatish Balay } 114449b5e25fSSatish Balay } 11453cea4cbeSHong Zhang } while (sctx.chshift); 114649b5e25fSSatish Balay ierr = PetscFree(il);CHKERRQ(ierr); 114798a8e78dSHong Zhang if (a->permute){ierr = PetscFree(aa);CHKERRQ(ierr);} 114849b5e25fSSatish Balay 114949b5e25fSSatish Balay ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 1150db4efbfdSBarry Smith C->ops->solve = MatSolve_SeqSBAIJ_1; 1151db4efbfdSBarry Smith C->ops->solves = MatSolves_SeqSBAIJ_1; 1152db4efbfdSBarry Smith C->ops->solvetranspose = MatSolve_SeqSBAIJ_1; 1153db4efbfdSBarry Smith C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1; 1154db4efbfdSBarry Smith C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1; 115549b5e25fSSatish Balay C->assembled = PETSC_TRUE; 11565c0bcdfcSHong Zhang C->preallocated = PETSC_TRUE; 1157d0f46423SBarry Smith ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); 11583cea4cbeSHong Zhang if (sctx.nshift){ 11593cea4cbeSHong Zhang if (shiftnz) { 11601e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 11613cea4cbeSHong Zhang } else if (shiftpd) { 11621e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 1163997a0949SHong Zhang } 1164997a0949SHong Zhang } 116549b5e25fSSatish Balay PetscFunctionReturn(0); 116649b5e25fSSatish Balay } 116749b5e25fSSatish Balay 1168671cb588SHong Zhang /* 1169671cb588SHong Zhang Version for when blocks are 1 by 1 Using natural ordering 1170671cb588SHong Zhang */ 11714a2ae208SSatish Balay #undef __FUNCT__ 11724a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering" 1173*719d5645SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat C,Mat A,MatFactorInfo *info) 1174671cb588SHong Zhang { 1175671cb588SHong Zhang Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data; 1176dfbe8321SBarry Smith PetscErrorCode ierr; 117713f74950SBarry Smith PetscInt i,j,mbs = a->mbs; 117813f74950SBarry Smith PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 11793cea4cbeSHong Zhang PetscInt k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz; 1180653a6975SHong Zhang MatScalar *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval; 11813cea4cbeSHong Zhang PetscReal zeropivot,rs,shiftnz; 1182fbf22428SSatish Balay PetscReal shiftpd; 11833cea4cbeSHong Zhang ChShift_Ctx sctx; 11843cea4cbeSHong Zhang PetscInt newshift; 1185653a6975SHong Zhang 1186653a6975SHong Zhang PetscFunctionBegin; 1187653a6975SHong Zhang /* initialization */ 11883cea4cbeSHong Zhang shiftnz = info->shiftnz; 11893cea4cbeSHong Zhang shiftpd = info->shiftpd; 11903cea4cbeSHong Zhang zeropivot = info->zeropivot; 11913cea4cbeSHong Zhang 1192653a6975SHong Zhang /* il and jl record the first nonzero element in each row of the accessing 1193653a6975SHong Zhang window U(0:k, k:mbs-1). 1194653a6975SHong Zhang jl: list of rows to be added to uneliminated rows 1195653a6975SHong Zhang i>= k: jl(i) is the first row to be added to row i 1196653a6975SHong Zhang i< k: jl(i) is the row following row i in some list of rows 1197653a6975SHong Zhang jl(i) = mbs indicates the end of a list 1198653a6975SHong Zhang il(i): points to the first nonzero element in U(i,k:mbs-1) 1199653a6975SHong Zhang */ 1200653a6975SHong Zhang ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 120113f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); 1202653a6975SHong Zhang jl = il + mbs; 1203b00f7748SHong Zhang 12043cea4cbeSHong Zhang sctx.shift_amount = 0; 12053cea4cbeSHong Zhang sctx.nshift = 0; 1206b00f7748SHong Zhang do { 12073cea4cbeSHong Zhang sctx.chshift = PETSC_FALSE; 1208653a6975SHong Zhang for (i=0; i<mbs; i++) { 1209653a6975SHong Zhang rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 1210653a6975SHong Zhang } 1211653a6975SHong Zhang 1212997a0949SHong Zhang for (k = 0; k<mbs; k++){ 1213653a6975SHong Zhang /*initialize k-th row with elements nonzero in row perm(k) of A */ 1214653a6975SHong Zhang nz = ai[k+1] - ai[k]; 1215653a6975SHong Zhang acol = aj + ai[k]; 1216653a6975SHong Zhang aval = aa + ai[k]; 1217653a6975SHong Zhang bval = ba + bi[k]; 1218653a6975SHong Zhang while (nz -- ){ 1219653a6975SHong Zhang rtmp[*acol++] = *aval++; 1220653a6975SHong Zhang *bval++ = 0.0; /* for in-place factorization */ 1221653a6975SHong Zhang } 12223cea4cbeSHong Zhang 12233cea4cbeSHong Zhang /* shift the diagonal of the matrix */ 12243cea4cbeSHong Zhang if (sctx.nshift) rtmp[k] += sctx.shift_amount; 1225653a6975SHong Zhang 1226653a6975SHong Zhang /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 1227653a6975SHong Zhang dk = rtmp[k]; 1228653a6975SHong Zhang i = jl[k]; /* first row to be added to k_th row */ 1229653a6975SHong Zhang 1230653a6975SHong Zhang while (i < k){ 1231653a6975SHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 1232653a6975SHong Zhang /* compute multiplier, update D(k) and U(i,k) */ 1233653a6975SHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 1234653a6975SHong Zhang uikdi = - ba[ili]*ba[bi[i]]; 1235653a6975SHong Zhang dk += uikdi*ba[ili]; 1236653a6975SHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 1237653a6975SHong Zhang 1238653a6975SHong Zhang /* add multiple of row i to k-th row ... */ 1239653a6975SHong Zhang jmin = ili + 1; 1240653a6975SHong Zhang nz = bi[i+1] - jmin; 1241653a6975SHong Zhang if (nz > 0){ 1242653a6975SHong Zhang bcol = bj + jmin; 1243653a6975SHong Zhang bval = ba + jmin; 1244187a9f4bSHong Zhang ierr = PetscLogFlops(2*nz);CHKERRQ(ierr); 1245653a6975SHong Zhang while (nz --) rtmp[*bcol++] += uikdi*(*bval++); 1246187a9f4bSHong Zhang 1247997a0949SHong Zhang /* update il and jl for i-th row */ 1248997a0949SHong Zhang il[i] = jmin; 1249997a0949SHong Zhang j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 1250653a6975SHong Zhang } 1251653a6975SHong Zhang i = nexti; 1252653a6975SHong Zhang } 1253653a6975SHong Zhang 12543cea4cbeSHong Zhang /* shift the diagonals when zero pivot is detected */ 12553cea4cbeSHong Zhang /* compute rs=sum of abs(off-diagonal) */ 12563cea4cbeSHong Zhang rs = 0.0; 125721c26570Svictorle jmin = bi[k]+1; 125821c26570Svictorle nz = bi[k+1] - jmin; 125921c26570Svictorle if (nz){ 126021c26570Svictorle bcol = bj + jmin; 126121c26570Svictorle while (nz--){ 126289f845c8SHong Zhang rs += PetscAbsScalar(rtmp[*bcol]); 126389f845c8SHong Zhang bcol++; 126421c26570Svictorle } 126521c26570Svictorle } 12663cea4cbeSHong Zhang 12673cea4cbeSHong Zhang sctx.rs = rs; 12683cea4cbeSHong Zhang sctx.pv = dk; 126945938b79SHong Zhang ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr); 127045938b79SHong Zhang if (newshift == 1) break; /* sctx.shift_amount is updated */ 1271653a6975SHong Zhang 1272997a0949SHong Zhang /* copy data into U(k,:) */ 1273653a6975SHong Zhang ba[bi[k]] = 1.0/dk; 1274653a6975SHong Zhang jmin = bi[k]+1; 1275653a6975SHong Zhang nz = bi[k+1] - jmin; 1276653a6975SHong Zhang if (nz){ 1277653a6975SHong Zhang bcol = bj + jmin; 1278653a6975SHong Zhang bval = ba + jmin; 1279653a6975SHong Zhang while (nz--){ 1280653a6975SHong Zhang *bval++ = rtmp[*bcol]; 1281653a6975SHong Zhang rtmp[*bcol++] = 0.0; 1282653a6975SHong Zhang } 1283997a0949SHong Zhang /* add k-th row into il and jl */ 1284653a6975SHong Zhang il[k] = jmin; 1285997a0949SHong Zhang i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 1286653a6975SHong Zhang } 1287b00f7748SHong Zhang } /* end of for (k = 0; k<mbs; k++) */ 12883cea4cbeSHong Zhang } while (sctx.chshift); 1289653a6975SHong Zhang ierr = PetscFree(rtmp);CHKERRQ(ierr); 1290653a6975SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 1291653a6975SHong Zhang 1292db4efbfdSBarry Smith C->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1293db4efbfdSBarry Smith C->ops->solves = MatSolves_SeqSBAIJ_1; 1294db4efbfdSBarry Smith C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1295db4efbfdSBarry Smith C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering; 1296db4efbfdSBarry Smith C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering; 1297db4efbfdSBarry Smith 1298653a6975SHong Zhang C->assembled = PETSC_TRUE; 1299653a6975SHong Zhang C->preallocated = PETSC_TRUE; 1300d0f46423SBarry Smith ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); 13013cea4cbeSHong Zhang if (sctx.nshift){ 13023cea4cbeSHong Zhang if (shiftnz) { 13031e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 13043cea4cbeSHong Zhang } else if (shiftpd) { 13051e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 1306b00f7748SHong Zhang } 1307fb10cecfSBarry Smith } 1308653a6975SHong Zhang PetscFunctionReturn(0); 1309653a6975SHong Zhang } 1310653a6975SHong Zhang 1311653a6975SHong Zhang #undef __FUNCT__ 13124a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ" 1313dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,MatFactorInfo *info) 131449b5e25fSSatish Balay { 1315dfbe8321SBarry Smith PetscErrorCode ierr; 131649b5e25fSSatish Balay Mat C; 131749b5e25fSSatish Balay 131849b5e25fSSatish Balay PetscFunctionBegin; 1319*719d5645SBarry Smith ierr = MatGetFactor(A,"petsc",MAT_FACTOR_CHOLESKY,&C);CHKERRQ(ierr); 1320*719d5645SBarry Smith ierr = MatCholeskyFactorSymbolic(C,A,perm,info);CHKERRQ(ierr); 1321*719d5645SBarry Smith ierr = MatCholeskyFactorNumeric(C,A,info);CHKERRQ(ierr); 1322db4efbfdSBarry Smith A->ops->solve = C->ops->solve; 1323db4efbfdSBarry Smith A->ops->solvetranspose = C->ops->solvetranspose; 13244445ddedSHong Zhang ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 132549b5e25fSSatish Balay PetscFunctionReturn(0); 132649b5e25fSSatish Balay } 132749b5e25fSSatish Balay 132849b5e25fSSatish Balay 1329