1be1d678aSKris Buschelman #define PETSCMAT_DLL 249b5e25fSSatish Balay 37c4f633dSBarry Smith #include "../src/mat/impls/baij/seq/baij.h" 47c4f633dSBarry Smith #include "../src/mat/impls/sbaij/seq/sbaij.h" 57c4f633dSBarry Smith #include "../src/inline/ilu.h" 67c4f633dSBarry Smith #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" 510481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(Mat F,Mat A,IS perm,const MatFactorInfo *info) 5210c27e3dSHong Zhang { 5310c27e3dSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b; 5410c27e3dSHong Zhang PetscErrorCode ierr; 555d0c19d7SBarry Smith const PetscInt *rip,*ai,*aj; 565d0c19d7SBarry Smith PetscInt i,mbs = a->mbs,*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 */ 182719d5645SBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(F,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 18310c27e3dSHong Zhang 184719d5645SBarry Smith /* ierr = PetscLogObjectParent(B,iperm);CHKERRQ(ierr); */ 185719d5645SBarry 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 */ 203719d5645SBarry Smith ierr = PetscLogObjectMemory(F,(iu[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 20410c27e3dSHong Zhang b->maxnz = b->nz = iu[mbs]; 20510c27e3dSHong Zhang 206719d5645SBarry Smith (F)->info.factor_mallocs = reallocs; 207719d5645SBarry Smith (F)->info.fill_ratio_given = f; 20810c27e3dSHong Zhang if (ai[mbs] != 0) { 209719d5645SBarry Smith (F)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]); 21010c27e3dSHong Zhang } else { 211719d5645SBarry Smith (F)->info.fill_ratio_needed = 0.0; 21210c27e3dSHong Zhang } 213719d5645SBarry 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" 2207c4f633dSBarry Smith #include "../src/mat/utils/freespace.h" 2214a2ae208SSatish Balay #undef __FUNCT__ 2224a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ" 2230481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat fact,Mat A,IS perm,const 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; 2305d0c19d7SBarry Smith const PetscInt *rip,*ai,*aj; 2315d0c19d7SBarry Smith PetscInt i,mbs=a->mbs,bs=A->rmap->bs,reallocs=0,prow,d; 23298a8e78dSHong Zhang PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 2337625dc9aSHong Zhang PetscInt nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr; 234a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 23598a8e78dSHong Zhang PetscBT lnkbt; 23649b5e25fSSatish Balay 23749b5e25fSSatish Balay PetscFunctionBegin; 23858ebbce7SBarry Smith ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr); 23958ebbce7SBarry Smith if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d); 24058ebbce7SBarry Smith 24110c27e3dSHong Zhang /* 24210c27e3dSHong Zhang This code originally uses Modified Sparse Row (MSR) storage 24310c27e3dSHong Zhang (see page 85, "Iterative Methods ..." by Saad) for the output matrix B - bad choise! 24410c27e3dSHong Zhang Then it is rewritten so the factor B takes seqsbaij format. However the associated 24510c27e3dSHong Zhang MatCholeskyFactorNumeric_() have not been modified for the cases of bs>1 or !perm_identity, 24610c27e3dSHong Zhang thus the original code in MSR format is still used for these cases. 24710c27e3dSHong Zhang The code below should replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR() whenever 24810c27e3dSHong Zhang MatCholeskyFactorNumeric_() is modified for using sbaij symbolic factor. 24910c27e3dSHong Zhang */ 250fff829cfSHong Zhang if (bs > 1){ 251719d5645SBarry Smith ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(fact,A,perm,info);CHKERRQ(ierr); 252fff829cfSHong Zhang PetscFunctionReturn(0); 253fff829cfSHong Zhang } 25410c27e3dSHong Zhang 25598a8e78dSHong Zhang /* check whether perm is the identity mapping */ 25698a8e78dSHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 25798a8e78dSHong Zhang 258fff829cfSHong Zhang if (perm_identity){ 259fff829cfSHong Zhang a->permute = PETSC_FALSE; 260fff829cfSHong Zhang ai = a->i; aj = a->j; 261fff829cfSHong Zhang } else { 262d33dd558SHong Zhang SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format"); 263d33dd558SHong Zhang /* There are bugs for reordeing. Needs further work. 264d33dd558SHong Zhang MatReordering for sbaij cannot be efficient. User should use aij formt! */ 265fff829cfSHong Zhang a->permute = PETSC_TRUE; 266fff829cfSHong Zhang ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr); 267fff829cfSHong Zhang ai = a->inew; aj = a->jnew; 268fff829cfSHong Zhang } 269fff829cfSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 270fff829cfSHong Zhang 271fff829cfSHong Zhang /* initialization */ 2727625dc9aSHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 27398a8e78dSHong Zhang ui[0] = 0; 27498a8e78dSHong Zhang 27598a8e78dSHong Zhang /* jl: linked list for storing indices of the pivot rows 2767625dc9aSHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */ 2777625dc9aSHong Zhang ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr); 2787625dc9aSHong Zhang il = jl + mbs; 2797625dc9aSHong Zhang cols = il + mbs; 2807625dc9aSHong Zhang ui_ptr = (PetscInt**)(cols + mbs); 2817625dc9aSHong Zhang 2827625dc9aSHong Zhang for (i=0; i<mbs; i++){ 2837625dc9aSHong Zhang jl[i] = mbs; il[i] = 0; 284fff829cfSHong Zhang } 285fff829cfSHong Zhang 28698a8e78dSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 2877625dc9aSHong Zhang nlnk = mbs + 1; 2887625dc9aSHong Zhang ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 289fff829cfSHong Zhang 2907625dc9aSHong Zhang /* initial FreeSpace size is fill*(ai[mbs]+1) */ 291a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr); 29298a8e78dSHong Zhang current_space = free_space; 29398a8e78dSHong Zhang 2947625dc9aSHong Zhang for (k=0; k<mbs; k++){ /* for each active row k */ 29598a8e78dSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 29698a8e78dSHong Zhang nzk = 0; 29798a8e78dSHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 2987625dc9aSHong Zhang for (j=0; j<ncols; j++){ 2997625dc9aSHong Zhang i = *(aj + ai[rip[k]] + j); 3007625dc9aSHong Zhang cols[j] = rip[i]; 3017625dc9aSHong Zhang } 3027625dc9aSHong Zhang ierr = PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 30398a8e78dSHong Zhang nzk += nlnk; 30498a8e78dSHong Zhang 30598a8e78dSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 30698a8e78dSHong Zhang prow = jl[k]; /* 1st pivot row */ 307fff829cfSHong Zhang 308fff829cfSHong Zhang while (prow < k){ 309fff829cfSHong Zhang nextprow = jl[prow]; 31098a8e78dSHong Zhang /* merge prow into k-th row */ 3117625dc9aSHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */ 31298a8e78dSHong Zhang jmax = ui[prow+1]; 31398a8e78dSHong Zhang ncols = jmax-jmin; 3147625dc9aSHong Zhang uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */ 3155a8e39fbSHong Zhang ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 31698a8e78dSHong Zhang nzk += nlnk; 317fff829cfSHong Zhang 31898a8e78dSHong Zhang /* update il and jl for prow */ 319fff829cfSHong Zhang if (jmin < jmax){ 320fff829cfSHong Zhang il[prow] = jmin; 3217625dc9aSHong Zhang j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 322fff829cfSHong Zhang } 323fff829cfSHong Zhang prow = nextprow; 324fff829cfSHong Zhang } 325fff829cfSHong Zhang 32698a8e78dSHong Zhang /* if free space is not available, make more free space */ 32798a8e78dSHong Zhang if (current_space->local_remaining<nzk) { 3287625dc9aSHong Zhang i = mbs - k + 1; /* num of unfactored rows */ 32998a8e78dSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 330a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 33198a8e78dSHong Zhang reallocs++; 33298a8e78dSHong Zhang } 33398a8e78dSHong Zhang 33498a8e78dSHong Zhang /* copy data into free space, then initialize lnk */ 3357625dc9aSHong Zhang ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 33698a8e78dSHong Zhang 33798a8e78dSHong Zhang /* add the k-th row into il and jl */ 33898a8e78dSHong Zhang if (nzk-1 > 0){ 3397625dc9aSHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */ 340fff829cfSHong Zhang jl[k] = jl[i]; jl[i] = k; 34198a8e78dSHong Zhang il[k] = ui[k] + 1; 342fff829cfSHong Zhang } 3437625dc9aSHong Zhang ui_ptr[k] = current_space->array; 34498a8e78dSHong Zhang current_space->array += nzk; 34598a8e78dSHong Zhang current_space->local_used += nzk; 34698a8e78dSHong Zhang current_space->local_remaining -= nzk; 347fff829cfSHong Zhang 34898a8e78dSHong Zhang ui[k+1] = ui[k] + nzk; 349fff829cfSHong Zhang } 350fff829cfSHong Zhang 3516cf91177SBarry Smith #if defined(PETSC_USE_INFO) 3527625dc9aSHong Zhang if (ai[mbs] != 0) { 3537625dc9aSHong Zhang PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 354ae15b995SBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 355ae15b995SBarry Smith ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 356ae15b995SBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 357fff829cfSHong Zhang } else { 358ae15b995SBarry Smith ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 359fff829cfSHong Zhang } 36063ba0a88SBarry Smith #endif 361fff829cfSHong Zhang 362fff829cfSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 363fff829cfSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 364fff829cfSHong Zhang 36598a8e78dSHong Zhang /* destroy list of free space and other temporary array(s) */ 3667625dc9aSHong Zhang ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 367a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 36898a8e78dSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 369fff829cfSHong Zhang 37098a8e78dSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 371719d5645SBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(fact,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 37298a8e78dSHong Zhang 373719d5645SBarry Smith b = (Mat_SeqSBAIJ*)(fact)->data; 374fff829cfSHong Zhang b->singlemalloc = PETSC_FALSE; 375e6b907acSBarry Smith b->free_a = PETSC_TRUE; 376e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 3777625dc9aSHong Zhang ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 37898a8e78dSHong Zhang b->j = uj; 37998a8e78dSHong Zhang b->i = ui; 380fff829cfSHong Zhang b->diag = 0; 381fff829cfSHong Zhang b->ilen = 0; 382fff829cfSHong Zhang b->imax = 0; 383fff829cfSHong Zhang b->row = perm; 384fff829cfSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 385fff829cfSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 386fff829cfSHong Zhang b->icol = perm; 387fff829cfSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 3887625dc9aSHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 389719d5645SBarry Smith ierr = PetscLogObjectMemory(fact,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 3907625dc9aSHong Zhang b->maxnz = b->nz = ui[mbs]; 391fff829cfSHong Zhang 392719d5645SBarry Smith (fact)->info.factor_mallocs = reallocs; 393719d5645SBarry Smith (fact)->info.fill_ratio_given = fill; 3947625dc9aSHong Zhang if (ai[mbs] != 0) { 395719d5645SBarry Smith (fact)->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 396fff829cfSHong Zhang } else { 397719d5645SBarry Smith (fact)->info.fill_ratio_needed = 0.0; 398fff829cfSHong Zhang } 399719d5645SBarry Smith ierr = MatSeqSBAIJSetNumericFactorization(fact,perm_identity);CHKERRQ(ierr); 400064503c5SHong Zhang PetscFunctionReturn(0); 401064503c5SHong Zhang } 4024a2ae208SSatish Balay #undef __FUNCT__ 4034a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N" 4040481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat C,Mat A,const MatFactorInfo *info) 40549b5e25fSSatish Balay { 4064c16a6a6SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 4074c16a6a6SHong Zhang IS perm = b->row; 4086849ba73SBarry Smith PetscErrorCode ierr; 4095d0c19d7SBarry Smith const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j; 4105d0c19d7SBarry Smith PetscInt i,j; 4115d0c19d7SBarry Smith PetscInt *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 412d0f46423SBarry Smith PetscInt bs=A->rmap->bs,bs2 = a->bs2,bslog = 0; 4134c16a6a6SHong Zhang MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 4144c16a6a6SHong Zhang MatScalar *u,*diag,*rtmp,*rtmp_ptr; 41528de702eSHong Zhang MatScalar *work; 41613f74950SBarry Smith PetscInt *pivots; 4174c16a6a6SHong Zhang 4184c16a6a6SHong Zhang PetscFunctionBegin; 4194c16a6a6SHong Zhang /* initialization */ 42082502324SSatish Balay ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 4214c16a6a6SHong Zhang ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr); 42213f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); 42328de702eSHong Zhang jl = il + mbs; 4244c16a6a6SHong Zhang for (i=0; i<mbs; i++) { 4254c16a6a6SHong Zhang jl[i] = mbs; il[0] = 0; 4264c16a6a6SHong Zhang } 427b0a32e0cSBarry Smith ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr); 42828de702eSHong Zhang uik = dk + bs2; 42928de702eSHong Zhang work = uik + bs2; 43013f74950SBarry Smith ierr = PetscMalloc(bs*sizeof(PetscInt),&pivots);CHKERRQ(ierr); 4314c16a6a6SHong Zhang 4324c16a6a6SHong Zhang ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr); 4334c16a6a6SHong Zhang 4344c16a6a6SHong Zhang /* check permutation */ 4354c16a6a6SHong Zhang if (!a->permute){ 4364c16a6a6SHong Zhang ai = a->i; aj = a->j; aa = a->a; 4374c16a6a6SHong Zhang } else { 4384c16a6a6SHong Zhang ai = a->inew; aj = a->jnew; 43982502324SSatish Balay ierr = PetscMalloc(bs2*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); 4404c16a6a6SHong Zhang ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); 44113f74950SBarry Smith ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr); 44213f74950SBarry Smith ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr); 4434c16a6a6SHong Zhang 444187a9f4bSHong Zhang /* flops in while loop */ 445187a9f4bSHong Zhang bslog = 2*bs*bs2; 446187a9f4bSHong Zhang 4474c16a6a6SHong Zhang for (i=0; i<mbs; i++){ 4484c16a6a6SHong Zhang jmin = ai[i]; jmax = ai[i+1]; 4494c16a6a6SHong Zhang for (j=jmin; j<jmax; j++){ 4504c16a6a6SHong Zhang while (a2anew[j] != j){ 4514c16a6a6SHong Zhang k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k; 4524c16a6a6SHong Zhang for (k1=0; k1<bs2; k1++){ 4534c16a6a6SHong Zhang dk[k1] = aa[k*bs2+k1]; 4544c16a6a6SHong Zhang aa[k*bs2+k1] = aa[j*bs2+k1]; 4554c16a6a6SHong Zhang aa[j*bs2+k1] = dk[k1]; 4564c16a6a6SHong Zhang } 4574c16a6a6SHong Zhang } 4584c16a6a6SHong Zhang /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */ 4594c16a6a6SHong Zhang if (i > aj[j]){ 4604c16a6a6SHong Zhang /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */ 4614c16a6a6SHong Zhang ap = aa + j*bs2; /* ptr to the beginning of j-th block of aa */ 4624c16a6a6SHong Zhang for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */ 4634c16a6a6SHong Zhang for (k=0; k<bs; k++){ /* j-th block of aa <- dk^T */ 4644c16a6a6SHong Zhang for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1]; 4654c16a6a6SHong Zhang } 4664c16a6a6SHong Zhang } 4674c16a6a6SHong Zhang } 4684c16a6a6SHong Zhang } 469323b833fSBarry Smith ierr = PetscFree(a2anew);CHKERRQ(ierr); 4704c16a6a6SHong Zhang } 4714c16a6a6SHong Zhang 4724c16a6a6SHong Zhang /* for each row k */ 4734c16a6a6SHong Zhang for (k = 0; k<mbs; k++){ 4744c16a6a6SHong Zhang 4754c16a6a6SHong Zhang /*initialize k-th row with elements nonzero in row perm(k) of A */ 4764c16a6a6SHong Zhang jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1]; 477057f5ba7SHong Zhang 4784c16a6a6SHong Zhang ap = aa + jmin*bs2; 4794c16a6a6SHong Zhang for (j = jmin; j < jmax; j++){ 4804c16a6a6SHong Zhang vj = perm_ptr[aj[j]]; /* block col. index */ 4814c16a6a6SHong Zhang rtmp_ptr = rtmp + vj*bs2; 4824c16a6a6SHong Zhang for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++; 4834c16a6a6SHong Zhang } 4844c16a6a6SHong Zhang 4854c16a6a6SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 4864c16a6a6SHong Zhang ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 4874c16a6a6SHong Zhang i = jl[k]; /* first row to be added to k_th row */ 4884c16a6a6SHong Zhang 489057f5ba7SHong Zhang while (i < k){ 4904c16a6a6SHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 4914c16a6a6SHong Zhang 4924c16a6a6SHong Zhang /* compute multiplier */ 4934c16a6a6SHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 4944c16a6a6SHong Zhang 4954c16a6a6SHong Zhang /* uik = -inv(Di)*U_bar(i,k) */ 4964c16a6a6SHong Zhang diag = ba + i*bs2; 4974c16a6a6SHong Zhang u = ba + ili*bs2; 4984c16a6a6SHong Zhang ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 4994c16a6a6SHong Zhang Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u); 5004c16a6a6SHong Zhang 5014c16a6a6SHong Zhang /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 5024c16a6a6SHong Zhang Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u); 503*dc0b31edSSatish Balay ierr = PetscLogFlops(bslog*2.0);CHKERRQ(ierr); 5044c16a6a6SHong Zhang 5054c16a6a6SHong Zhang /* update -U(i,k) */ 5064c16a6a6SHong Zhang ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 5074c16a6a6SHong Zhang 5084c16a6a6SHong Zhang /* add multiple of row i to k-th row ... */ 5094c16a6a6SHong Zhang jmin = ili + 1; jmax = bi[i+1]; 5104c16a6a6SHong Zhang if (jmin < jmax){ 5114c16a6a6SHong Zhang for (j=jmin; j<jmax; j++) { 5124c16a6a6SHong Zhang /* rtmp += -U(i,k)^T * U_bar(i,j) */ 5134c16a6a6SHong Zhang rtmp_ptr = rtmp + bj[j]*bs2; 5144c16a6a6SHong Zhang u = ba + j*bs2; 5154c16a6a6SHong Zhang Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u); 5164c16a6a6SHong Zhang } 517187a9f4bSHong Zhang ierr = PetscLogFlops(bslog*(jmax-jmin));CHKERRQ(ierr); 5184c16a6a6SHong Zhang 5194c16a6a6SHong Zhang /* ... add i to row list for next nonzero entry */ 5204c16a6a6SHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 5214c16a6a6SHong Zhang j = bj[jmin]; 5224c16a6a6SHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 5234c16a6a6SHong Zhang } 5244c16a6a6SHong Zhang i = nexti; 5254c16a6a6SHong Zhang } 5264c16a6a6SHong Zhang 5274c16a6a6SHong Zhang /* save nonzero entries in k-th row of U ... */ 5284c16a6a6SHong Zhang 5294c16a6a6SHong Zhang /* invert diagonal block */ 5304c16a6a6SHong Zhang diag = ba+k*bs2; 5314c16a6a6SHong Zhang ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr); 532d230e6fdSBarry Smith ierr = Kernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr); 5334c16a6a6SHong Zhang 5344c16a6a6SHong Zhang jmin = bi[k]; jmax = bi[k+1]; 5354c16a6a6SHong Zhang if (jmin < jmax) { 5364c16a6a6SHong Zhang for (j=jmin; j<jmax; j++){ 5374c16a6a6SHong Zhang vj = bj[j]; /* block col. index of U */ 5384c16a6a6SHong Zhang u = ba + j*bs2; 5394c16a6a6SHong Zhang rtmp_ptr = rtmp + vj*bs2; 5404c16a6a6SHong Zhang for (k1=0; k1<bs2; k1++){ 5414c16a6a6SHong Zhang *u++ = *rtmp_ptr; 5424c16a6a6SHong Zhang *rtmp_ptr++ = 0.0; 5434c16a6a6SHong Zhang } 5444c16a6a6SHong Zhang } 5454c16a6a6SHong Zhang 5464c16a6a6SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 5474c16a6a6SHong Zhang il[k] = jmin; 5484c16a6a6SHong Zhang i = bj[jmin]; 5494c16a6a6SHong Zhang jl[k] = jl[i]; jl[i] = k; 5504c16a6a6SHong Zhang } 5514c16a6a6SHong Zhang } 5524c16a6a6SHong Zhang 5534c16a6a6SHong Zhang ierr = PetscFree(rtmp);CHKERRQ(ierr); 5544c16a6a6SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 5554c16a6a6SHong Zhang ierr = PetscFree(dk);CHKERRQ(ierr); 5564c16a6a6SHong Zhang ierr = PetscFree(pivots);CHKERRQ(ierr); 5574c16a6a6SHong Zhang if (a->permute){ 5584c16a6a6SHong Zhang ierr = PetscFree(aa);CHKERRQ(ierr); 5594c16a6a6SHong Zhang } 5604c16a6a6SHong Zhang 5614c16a6a6SHong Zhang ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr); 562db4efbfdSBarry Smith C->ops->solve = MatSolve_SeqSBAIJ_N; 563db4efbfdSBarry Smith C->ops->solvetranspose = MatSolve_SeqSBAIJ_N; 564db4efbfdSBarry Smith C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_N; 565db4efbfdSBarry Smith C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_N; 566db4efbfdSBarry Smith 5674c16a6a6SHong Zhang C->assembled = PETSC_TRUE; 5684c16a6a6SHong Zhang C->preallocated = PETSC_TRUE; 569efee365bSSatish Balay ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 5704c16a6a6SHong Zhang PetscFunctionReturn(0); 5714c16a6a6SHong Zhang } 572d4edadd4SHong Zhang 5734a2ae208SSatish Balay #undef __FUNCT__ 5744a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering" 5750481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info) 576671cb588SHong Zhang { 577671cb588SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 578dfbe8321SBarry Smith PetscErrorCode ierr; 57913f74950SBarry Smith PetscInt i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 58013f74950SBarry Smith PetscInt *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 581d0f46423SBarry Smith PetscInt bs=A->rmap->bs,bs2 = a->bs2,bslog; 582671cb588SHong Zhang MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 583671cb588SHong Zhang MatScalar *u,*diag,*rtmp,*rtmp_ptr; 58428de702eSHong Zhang MatScalar *work; 58513f74950SBarry Smith PetscInt *pivots; 586671cb588SHong Zhang 587671cb588SHong Zhang PetscFunctionBegin; 58882502324SSatish Balay ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 589671cb588SHong Zhang ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr); 59013f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); 59128de702eSHong Zhang jl = il + mbs; 592671cb588SHong Zhang for (i=0; i<mbs; i++) { 593671cb588SHong Zhang jl[i] = mbs; il[0] = 0; 594671cb588SHong Zhang } 595b0a32e0cSBarry Smith ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr); 59628de702eSHong Zhang uik = dk + bs2; 59728de702eSHong Zhang work = uik + bs2; 59813f74950SBarry Smith ierr = PetscMalloc(bs*sizeof(PetscInt),&pivots);CHKERRQ(ierr); 599671cb588SHong Zhang 600671cb588SHong Zhang ai = a->i; aj = a->j; aa = a->a; 601671cb588SHong Zhang 602187a9f4bSHong Zhang /* flops in while loop */ 603187a9f4bSHong Zhang bslog = 2*bs*bs2; 604187a9f4bSHong Zhang 605671cb588SHong Zhang /* for each row k */ 606671cb588SHong Zhang for (k = 0; k<mbs; k++){ 607671cb588SHong Zhang 608671cb588SHong Zhang /*initialize k-th row with elements nonzero in row k of A */ 609671cb588SHong Zhang jmin = ai[k]; jmax = ai[k+1]; 610671cb588SHong Zhang ap = aa + jmin*bs2; 611671cb588SHong Zhang for (j = jmin; j < jmax; j++){ 612671cb588SHong Zhang vj = aj[j]; /* block col. index */ 613671cb588SHong Zhang rtmp_ptr = rtmp + vj*bs2; 614671cb588SHong Zhang for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++; 615671cb588SHong Zhang } 616671cb588SHong Zhang 617671cb588SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 618671cb588SHong Zhang ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 619671cb588SHong Zhang i = jl[k]; /* first row to be added to k_th row */ 620671cb588SHong Zhang 621057f5ba7SHong Zhang while (i < k){ 622671cb588SHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 623671cb588SHong Zhang 624671cb588SHong Zhang /* compute multiplier */ 625671cb588SHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 626671cb588SHong Zhang 627671cb588SHong Zhang /* uik = -inv(Di)*U_bar(i,k) */ 628671cb588SHong Zhang diag = ba + i*bs2; 629671cb588SHong Zhang u = ba + ili*bs2; 630671cb588SHong Zhang ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 631671cb588SHong Zhang Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u); 632671cb588SHong Zhang 633671cb588SHong Zhang /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 634671cb588SHong Zhang Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u); 635*dc0b31edSSatish Balay ierr = PetscLogFlops(bslog*2.0);CHKERRQ(ierr); 636671cb588SHong Zhang 637671cb588SHong Zhang /* update -U(i,k) */ 638671cb588SHong Zhang ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 639671cb588SHong Zhang 640671cb588SHong Zhang /* add multiple of row i to k-th row ... */ 641671cb588SHong Zhang jmin = ili + 1; jmax = bi[i+1]; 642671cb588SHong Zhang if (jmin < jmax){ 643671cb588SHong Zhang for (j=jmin; j<jmax; j++) { 644671cb588SHong Zhang /* rtmp += -U(i,k)^T * U_bar(i,j) */ 645671cb588SHong Zhang rtmp_ptr = rtmp + bj[j]*bs2; 646671cb588SHong Zhang u = ba + j*bs2; 647671cb588SHong Zhang Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u); 648671cb588SHong Zhang } 649187a9f4bSHong Zhang ierr = PetscLogFlops(bslog*(jmax-jmin));CHKERRQ(ierr); 650671cb588SHong Zhang 651671cb588SHong Zhang /* ... add i to row list for next nonzero entry */ 652671cb588SHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 653671cb588SHong Zhang j = bj[jmin]; 654671cb588SHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 655671cb588SHong Zhang } 656671cb588SHong Zhang i = nexti; 657671cb588SHong Zhang } 658671cb588SHong Zhang 659671cb588SHong Zhang /* save nonzero entries in k-th row of U ... */ 660671cb588SHong Zhang 661671cb588SHong Zhang /* invert diagonal block */ 662671cb588SHong Zhang diag = ba+k*bs2; 663671cb588SHong Zhang ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr); 664d230e6fdSBarry Smith ierr = Kernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr); 665671cb588SHong Zhang 666671cb588SHong Zhang jmin = bi[k]; jmax = bi[k+1]; 667671cb588SHong Zhang if (jmin < jmax) { 668671cb588SHong Zhang for (j=jmin; j<jmax; j++){ 669671cb588SHong Zhang vj = bj[j]; /* block col. index of U */ 670671cb588SHong Zhang u = ba + j*bs2; 671671cb588SHong Zhang rtmp_ptr = rtmp + vj*bs2; 672671cb588SHong Zhang for (k1=0; k1<bs2; k1++){ 673671cb588SHong Zhang *u++ = *rtmp_ptr; 674671cb588SHong Zhang *rtmp_ptr++ = 0.0; 675671cb588SHong Zhang } 676671cb588SHong Zhang } 677671cb588SHong Zhang 678671cb588SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 679671cb588SHong Zhang il[k] = jmin; 680671cb588SHong Zhang i = bj[jmin]; 681671cb588SHong Zhang jl[k] = jl[i]; jl[i] = k; 682671cb588SHong Zhang } 683671cb588SHong Zhang } 684671cb588SHong Zhang 685671cb588SHong Zhang ierr = PetscFree(rtmp);CHKERRQ(ierr); 686671cb588SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 687671cb588SHong Zhang ierr = PetscFree(dk);CHKERRQ(ierr); 688671cb588SHong Zhang ierr = PetscFree(pivots);CHKERRQ(ierr); 689671cb588SHong Zhang 690db4efbfdSBarry Smith C->ops->solve = MatSolve_SeqSBAIJ_N_NaturalOrdering; 691db4efbfdSBarry Smith C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_NaturalOrdering; 692db4efbfdSBarry Smith C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering; 693db4efbfdSBarry Smith C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering; 694671cb588SHong Zhang C->assembled = PETSC_TRUE; 695671cb588SHong Zhang C->preallocated = PETSC_TRUE; 696efee365bSSatish Balay ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 697671cb588SHong Zhang PetscFunctionReturn(0); 698671cb588SHong Zhang } 699671cb588SHong Zhang 70049b5e25fSSatish Balay /* 701fcf159c0SHong Zhang Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP. 702cc0c071aSHong Zhang Version for blocks 2 by 2. 70349b5e25fSSatish Balay */ 7044a2ae208SSatish Balay #undef __FUNCT__ 7054a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2" 7060481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat C,Mat A,const MatFactorInfo *info) 70749b5e25fSSatish Balay { 70891602c66SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 709cc0c071aSHong Zhang IS perm = b->row; 7106849ba73SBarry Smith PetscErrorCode ierr; 7115d0c19d7SBarry Smith const PetscInt *ai,*aj,*perm_ptr; 7125d0c19d7SBarry Smith PetscInt i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 7135d0c19d7SBarry Smith PetscInt *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 714a1723e09SHong Zhang MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 715cc0c071aSHong Zhang MatScalar *u,*diag,*rtmp,*rtmp_ptr; 71662bba022SBarry Smith PetscReal shift = info->shiftinblocks; 71749b5e25fSSatish Balay 71849b5e25fSSatish Balay PetscFunctionBegin; 71991602c66SHong Zhang /* initialization */ 72091602c66SHong Zhang /* il and jl record the first nonzero element in each row of the accessing 72191602c66SHong Zhang window U(0:k, k:mbs-1). 72291602c66SHong Zhang jl: list of rows to be added to uneliminated rows 72391602c66SHong Zhang i>= k: jl(i) is the first row to be added to row i 72491602c66SHong Zhang i< k: jl(i) is the row following row i in some list of rows 72591602c66SHong Zhang jl(i) = mbs indicates the end of a list 72691602c66SHong Zhang il(i): points to the first nonzero element in columns k,...,mbs-1 of 72791602c66SHong Zhang row i of U */ 72882502324SSatish Balay ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 729cc0c071aSHong Zhang ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr); 73013f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); 73128de702eSHong Zhang jl = il + mbs; 73291602c66SHong Zhang for (i=0; i<mbs; i++) { 7333845f261SHong Zhang jl[i] = mbs; il[0] = 0; 73491602c66SHong Zhang } 73582502324SSatish Balay ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr); 73628de702eSHong Zhang uik = dk + 4; 737cc0c071aSHong Zhang ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr); 738cc0c071aSHong Zhang 739cc0c071aSHong Zhang /* check permutation */ 740cc0c071aSHong Zhang if (!a->permute){ 741cc0c071aSHong Zhang ai = a->i; aj = a->j; aa = a->a; 742cc0c071aSHong Zhang } else { 743cc0c071aSHong Zhang ai = a->inew; aj = a->jnew; 74482502324SSatish Balay ierr = PetscMalloc(4*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); 745cc0c071aSHong Zhang ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); 74613f74950SBarry Smith ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr); 74713f74950SBarry Smith ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr); 748cc0c071aSHong Zhang 749cc0c071aSHong Zhang for (i=0; i<mbs; i++){ 750cc0c071aSHong Zhang jmin = ai[i]; jmax = ai[i+1]; 751cc0c071aSHong Zhang for (j=jmin; j<jmax; j++){ 752cc0c071aSHong Zhang while (a2anew[j] != j){ 753cc0c071aSHong Zhang k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k; 754cc0c071aSHong Zhang for (k1=0; k1<4; k1++){ 755cc0c071aSHong Zhang dk[k1] = aa[k*4+k1]; 756cc0c071aSHong Zhang aa[k*4+k1] = aa[j*4+k1]; 757cc0c071aSHong Zhang aa[j*4+k1] = dk[k1]; 758cc0c071aSHong Zhang } 759cc0c071aSHong Zhang } 760cc0c071aSHong Zhang /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */ 761cc0c071aSHong Zhang if (i > aj[j]){ 762a1723e09SHong Zhang /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */ 763cc0c071aSHong Zhang ap = aa + j*4; /* ptr to the beginning of the block */ 764cc0c071aSHong Zhang dk[1] = ap[1]; /* swap ap[1] and ap[2] */ 765cc0c071aSHong Zhang ap[1] = ap[2]; 766cc0c071aSHong Zhang ap[2] = dk[1]; 767cc0c071aSHong Zhang } 768cc0c071aSHong Zhang } 769cc0c071aSHong Zhang } 770ac355199SBarry Smith ierr = PetscFree(a2anew);CHKERRQ(ierr); 771cc0c071aSHong Zhang } 7723845f261SHong Zhang 77391602c66SHong Zhang /* for each row k */ 77491602c66SHong Zhang for (k = 0; k<mbs; k++){ 77591602c66SHong Zhang 77691602c66SHong Zhang /*initialize k-th row with elements nonzero in row perm(k) of A */ 777cc0c071aSHong Zhang jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1]; 778cc0c071aSHong Zhang ap = aa + jmin*4; 77991602c66SHong Zhang for (j = jmin; j < jmax; j++){ 780cc0c071aSHong Zhang vj = perm_ptr[aj[j]]; /* block col. index */ 781cc0c071aSHong Zhang rtmp_ptr = rtmp + vj*4; 782cc0c071aSHong Zhang for (i=0; i<4; i++) *rtmp_ptr++ = *ap++; 78391602c66SHong Zhang } 78491602c66SHong Zhang 78591602c66SHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 786cc0c071aSHong Zhang ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr); 78791602c66SHong Zhang i = jl[k]; /* first row to be added to k_th row */ 78891602c66SHong Zhang 789057f5ba7SHong Zhang while (i < k){ 79091602c66SHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 79191602c66SHong Zhang 7923845f261SHong Zhang /* compute multiplier */ 79391602c66SHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 7943845f261SHong Zhang 7953845f261SHong Zhang /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */ 796cc0c071aSHong Zhang diag = ba + i*4; 797cc0c071aSHong Zhang u = ba + ili*4; 798cc0c071aSHong Zhang uik[0] = -(diag[0]*u[0] + diag[2]*u[1]); 799cc0c071aSHong Zhang uik[1] = -(diag[1]*u[0] + diag[3]*u[1]); 800cc0c071aSHong Zhang uik[2] = -(diag[0]*u[2] + diag[2]*u[3]); 801cc0c071aSHong Zhang uik[3] = -(diag[1]*u[2] + diag[3]*u[3]); 8023845f261SHong Zhang 8033845f261SHong Zhang /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */ 804cc0c071aSHong Zhang dk[0] += uik[0]*u[0] + uik[1]*u[1]; 805cc0c071aSHong Zhang dk[1] += uik[2]*u[0] + uik[3]*u[1]; 806cc0c071aSHong Zhang dk[2] += uik[0]*u[2] + uik[1]*u[3]; 807cc0c071aSHong Zhang dk[3] += uik[2]*u[2] + uik[3]*u[3]; 8083845f261SHong Zhang 809*dc0b31edSSatish Balay ierr = PetscLogFlops(16.0*2.0);CHKERRQ(ierr); 810187a9f4bSHong Zhang 8113845f261SHong Zhang /* update -U(i,k): ba[ili] = uik */ 812cc0c071aSHong Zhang ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr); 81391602c66SHong Zhang 81491602c66SHong Zhang /* add multiple of row i to k-th row ... */ 81591602c66SHong Zhang jmin = ili + 1; jmax = bi[i+1]; 81691602c66SHong Zhang if (jmin < jmax){ 8173845f261SHong Zhang for (j=jmin; j<jmax; j++) { 8183845f261SHong Zhang /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */ 819cc0c071aSHong Zhang rtmp_ptr = rtmp + bj[j]*4; 820cc0c071aSHong Zhang u = ba + j*4; 821cc0c071aSHong Zhang rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1]; 822cc0c071aSHong Zhang rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1]; 823cc0c071aSHong Zhang rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3]; 824cc0c071aSHong Zhang rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3]; 8253845f261SHong Zhang } 826*dc0b31edSSatish Balay ierr = PetscLogFlops(16.0*(jmax-jmin));CHKERRQ(ierr); 8273845f261SHong Zhang 82891602c66SHong Zhang /* ... add i to row list for next nonzero entry */ 82991602c66SHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 83091602c66SHong Zhang j = bj[jmin]; 83191602c66SHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 83291602c66SHong Zhang } 833a1723e09SHong Zhang i = nexti; 83491602c66SHong Zhang } 835cc0c071aSHong Zhang 83691602c66SHong Zhang /* save nonzero entries in k-th row of U ... */ 8373845f261SHong Zhang 838cc0c071aSHong Zhang /* invert diagonal block */ 839cc0c071aSHong Zhang diag = ba+k*4; 840cc0c071aSHong Zhang ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr); 84162bba022SBarry Smith ierr = Kernel_A_gets_inverse_A_2(diag,shift);CHKERRQ(ierr); 8423845f261SHong Zhang 84391602c66SHong Zhang jmin = bi[k]; jmax = bi[k+1]; 84491602c66SHong Zhang if (jmin < jmax) { 84591602c66SHong Zhang for (j=jmin; j<jmax; j++){ 846cc0c071aSHong Zhang vj = bj[j]; /* block col. index of U */ 847cc0c071aSHong Zhang u = ba + j*4; 848cc0c071aSHong Zhang rtmp_ptr = rtmp + vj*4; 849cc0c071aSHong Zhang for (k1=0; k1<4; k1++){ 850cc0c071aSHong Zhang *u++ = *rtmp_ptr; 851cc0c071aSHong Zhang *rtmp_ptr++ = 0.0; 852cc0c071aSHong Zhang } 853cc0c071aSHong Zhang } 8543845f261SHong Zhang 85591602c66SHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 85691602c66SHong Zhang il[k] = jmin; 85791602c66SHong Zhang i = bj[jmin]; 85891602c66SHong Zhang jl[k] = jl[i]; jl[i] = k; 85991602c66SHong Zhang } 86091602c66SHong Zhang } 8613845f261SHong Zhang 86249b5e25fSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 86391602c66SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 8643845f261SHong Zhang ierr = PetscFree(dk);CHKERRQ(ierr); 86591602c66SHong Zhang if (a->permute) { 86691602c66SHong Zhang ierr = PetscFree(aa);CHKERRQ(ierr); 86791602c66SHong Zhang } 868cc0c071aSHong Zhang ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr); 869db4efbfdSBarry Smith C->ops->solve = MatSolve_SeqSBAIJ_2; 870db4efbfdSBarry Smith C->ops->solvetranspose = MatSolve_SeqSBAIJ_2; 87149b5e25fSSatish Balay C->assembled = PETSC_TRUE; 8725c0bcdfcSHong Zhang C->preallocated = PETSC_TRUE; 873efee365bSSatish Balay ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 87449b5e25fSSatish Balay PetscFunctionReturn(0); 87549b5e25fSSatish Balay } 87691602c66SHong Zhang 87749b5e25fSSatish Balay /* 87849b5e25fSSatish Balay Version for when blocks are 2 by 2 Using natural ordering 87949b5e25fSSatish Balay */ 8804a2ae208SSatish Balay #undef __FUNCT__ 8814a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering" 8820481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info) 88349b5e25fSSatish Balay { 884ab27746eSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 885dfbe8321SBarry Smith PetscErrorCode ierr; 88613f74950SBarry Smith PetscInt i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 88713f74950SBarry Smith PetscInt *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 888ab27746eSHong Zhang MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 889ab27746eSHong Zhang MatScalar *u,*diag,*rtmp,*rtmp_ptr; 89062bba022SBarry Smith PetscReal shift = info->shiftinblocks; 89149b5e25fSSatish Balay 89249b5e25fSSatish Balay PetscFunctionBegin; 893ab27746eSHong Zhang /* initialization */ 894ab27746eSHong Zhang /* il and jl record the first nonzero element in each row of the accessing 895ab27746eSHong Zhang window U(0:k, k:mbs-1). 896ab27746eSHong Zhang jl: list of rows to be added to uneliminated rows 897ab27746eSHong Zhang i>= k: jl(i) is the first row to be added to row i 898ab27746eSHong Zhang i< k: jl(i) is the row following row i in some list of rows 899ab27746eSHong Zhang jl(i) = mbs indicates the end of a list 900ab27746eSHong Zhang il(i): points to the first nonzero element in columns k,...,mbs-1 of 901ab27746eSHong Zhang row i of U */ 90282502324SSatish Balay ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 903ab27746eSHong Zhang ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr); 90413f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); 90528de702eSHong Zhang jl = il + mbs; 906ab27746eSHong Zhang for (i=0; i<mbs; i++) { 907ab27746eSHong Zhang jl[i] = mbs; il[0] = 0; 90849b5e25fSSatish Balay } 90982502324SSatish Balay ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr); 91028de702eSHong Zhang uik = dk + 4; 911ab27746eSHong Zhang 912ab27746eSHong Zhang ai = a->i; aj = a->j; aa = a->a; 913ab27746eSHong Zhang 914ab27746eSHong Zhang /* for each row k */ 915ab27746eSHong Zhang for (k = 0; k<mbs; k++){ 916ab27746eSHong Zhang 917ab27746eSHong Zhang /*initialize k-th row with elements nonzero in row k of A */ 918ab27746eSHong Zhang jmin = ai[k]; jmax = ai[k+1]; 919ab27746eSHong Zhang ap = aa + jmin*4; 920ab27746eSHong Zhang for (j = jmin; j < jmax; j++){ 921ab27746eSHong Zhang vj = aj[j]; /* block col. index */ 922ab27746eSHong Zhang rtmp_ptr = rtmp + vj*4; 923ab27746eSHong Zhang for (i=0; i<4; i++) *rtmp_ptr++ = *ap++; 92449b5e25fSSatish Balay } 925ab27746eSHong Zhang 926ab27746eSHong Zhang /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 927ab27746eSHong Zhang ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr); 928ab27746eSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 929ab27746eSHong Zhang 930057f5ba7SHong Zhang while (i < k){ 931ab27746eSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 932ab27746eSHong Zhang 933ab27746eSHong Zhang /* compute multiplier */ 934ab27746eSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 935ab27746eSHong Zhang 936ab27746eSHong Zhang /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */ 937ab27746eSHong Zhang diag = ba + i*4; 938ab27746eSHong Zhang u = ba + ili*4; 939ab27746eSHong Zhang uik[0] = -(diag[0]*u[0] + diag[2]*u[1]); 940ab27746eSHong Zhang uik[1] = -(diag[1]*u[0] + diag[3]*u[1]); 941ab27746eSHong Zhang uik[2] = -(diag[0]*u[2] + diag[2]*u[3]); 942ab27746eSHong Zhang uik[3] = -(diag[1]*u[2] + diag[3]*u[3]); 943ab27746eSHong Zhang 944ab27746eSHong Zhang /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */ 945ab27746eSHong Zhang dk[0] += uik[0]*u[0] + uik[1]*u[1]; 946ab27746eSHong Zhang dk[1] += uik[2]*u[0] + uik[3]*u[1]; 947ab27746eSHong Zhang dk[2] += uik[0]*u[2] + uik[1]*u[3]; 948ab27746eSHong Zhang dk[3] += uik[2]*u[2] + uik[3]*u[3]; 949ab27746eSHong Zhang 950*dc0b31edSSatish Balay ierr = PetscLogFlops(16.0*2.0);CHKERRQ(ierr); 951187a9f4bSHong Zhang 952ab27746eSHong Zhang /* update -U(i,k): ba[ili] = uik */ 953ab27746eSHong Zhang ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr); 954ab27746eSHong Zhang 955ab27746eSHong Zhang /* add multiple of row i to k-th row ... */ 956ab27746eSHong Zhang jmin = ili + 1; jmax = bi[i+1]; 957ab27746eSHong Zhang if (jmin < jmax){ 958ab27746eSHong Zhang for (j=jmin; j<jmax; j++) { 959ab27746eSHong Zhang /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */ 960ab27746eSHong Zhang rtmp_ptr = rtmp + bj[j]*4; 961ab27746eSHong Zhang u = ba + j*4; 962ab27746eSHong Zhang rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1]; 963ab27746eSHong Zhang rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1]; 964ab27746eSHong Zhang rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3]; 965ab27746eSHong Zhang rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3]; 96649b5e25fSSatish Balay } 967*dc0b31edSSatish Balay ierr = PetscLogFlops(16.0*(jmax-jmin));CHKERRQ(ierr); 968ab27746eSHong Zhang 969ab27746eSHong Zhang /* ... add i to row list for next nonzero entry */ 970ab27746eSHong Zhang il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 971ab27746eSHong Zhang j = bj[jmin]; 972ab27746eSHong Zhang jl[i] = jl[j]; jl[j] = i; /* update jl */ 97349b5e25fSSatish Balay } 974ab27746eSHong Zhang i = nexti; 97549b5e25fSSatish Balay } 976ab27746eSHong Zhang 977ab27746eSHong Zhang /* save nonzero entries in k-th row of U ... */ 978ab27746eSHong Zhang 97949b5e25fSSatish Balay /* invert diagonal block */ 980ab27746eSHong Zhang diag = ba+k*4; 981ab27746eSHong Zhang ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr); 98262bba022SBarry Smith ierr = Kernel_A_gets_inverse_A_2(diag,shift);CHKERRQ(ierr); 983ab27746eSHong Zhang 984ab27746eSHong Zhang jmin = bi[k]; jmax = bi[k+1]; 985ab27746eSHong Zhang if (jmin < jmax) { 986ab27746eSHong Zhang for (j=jmin; j<jmax; j++){ 987ab27746eSHong Zhang vj = bj[j]; /* block col. index of U */ 988ab27746eSHong Zhang u = ba + j*4; 989ab27746eSHong Zhang rtmp_ptr = rtmp + vj*4; 990ab27746eSHong Zhang for (k1=0; k1<4; k1++){ 991ab27746eSHong Zhang *u++ = *rtmp_ptr; 992ab27746eSHong Zhang *rtmp_ptr++ = 0.0; 993ab27746eSHong Zhang } 994ab27746eSHong Zhang } 995ab27746eSHong Zhang 996ab27746eSHong Zhang /* ... add k to row list for first nonzero entry in k-th row */ 997ab27746eSHong Zhang il[k] = jmin; 998ab27746eSHong Zhang i = bj[jmin]; 999ab27746eSHong Zhang jl[k] = jl[i]; jl[i] = k; 1000ab27746eSHong Zhang } 100149b5e25fSSatish Balay } 100249b5e25fSSatish Balay 100349b5e25fSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 1004ab27746eSHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 1005ab27746eSHong Zhang ierr = PetscFree(dk);CHKERRQ(ierr); 1006ab27746eSHong Zhang 1007db4efbfdSBarry Smith C->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 1008db4efbfdSBarry Smith C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; 1009db4efbfdSBarry Smith C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering; 1010db4efbfdSBarry Smith C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering; 101149b5e25fSSatish Balay C->assembled = PETSC_TRUE; 10125c0bcdfcSHong Zhang C->preallocated = PETSC_TRUE; 1013efee365bSSatish Balay ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 101449b5e25fSSatish Balay PetscFunctionReturn(0); 101549b5e25fSSatish Balay } 101649b5e25fSSatish Balay 101749b5e25fSSatish Balay /* 101898a8e78dSHong Zhang Numeric U^T*D*U factorization for SBAIJ format. 101991602c66SHong Zhang Version for blocks are 1 by 1. 102049b5e25fSSatish Balay */ 10214a2ae208SSatish Balay #undef __FUNCT__ 10224a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1" 10230481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat C,Mat A,const MatFactorInfo *info) 102449b5e25fSSatish Balay { 102549b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data; 102649b5e25fSSatish Balay IS ip=b->row; 10276849ba73SBarry Smith PetscErrorCode ierr; 10285d0c19d7SBarry Smith const PetscInt *ai,*aj,*rip; 10295d0c19d7SBarry Smith PetscInt *a2anew,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j,*bcol; 1030997a0949SHong Zhang PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 1031997a0949SHong Zhang MatScalar *rtmp,*ba=b->a,*bval,*aa,dk,uikdi; 10323cea4cbeSHong Zhang PetscReal zeropivot,rs,shiftnz; 1033fbf22428SSatish Balay PetscReal shiftpd; 10343cea4cbeSHong Zhang ChShift_Ctx sctx; 10353cea4cbeSHong Zhang PetscInt newshift; 103649b5e25fSSatish Balay 103749b5e25fSSatish Balay PetscFunctionBegin; 10383cea4cbeSHong Zhang /* initialization */ 10393cea4cbeSHong Zhang shiftnz = info->shiftnz; 10403cea4cbeSHong Zhang shiftpd = info->shiftpd; 10413cea4cbeSHong Zhang zeropivot = info->zeropivot; 10423cea4cbeSHong Zhang 104349b5e25fSSatish Balay ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 1044cb718733SHong Zhang if (!a->permute){ 10452d007305SHong Zhang ai = a->i; aj = a->j; aa = a->a; 10462d007305SHong Zhang } else { 10472d007305SHong Zhang ai = a->inew; aj = a->jnew; 1048fff829cfSHong Zhang nz = ai[mbs]; 1049fff829cfSHong Zhang ierr = PetscMalloc(nz*sizeof(MatScalar),&aa);CHKERRQ(ierr); 1050fff829cfSHong Zhang a2anew = a->a2anew; 1051997a0949SHong Zhang bval = a->a; 1052fff829cfSHong Zhang for (j=0; j<nz; j++){ 1053997a0949SHong Zhang aa[a2anew[j]] = *(bval++); 10542d007305SHong Zhang } 10552d007305SHong Zhang } 10562d007305SHong Zhang 105791602c66SHong Zhang /* initialization */ 105849b5e25fSSatish Balay /* il and jl record the first nonzero element in each row of the accessing 105949b5e25fSSatish Balay window U(0:k, k:mbs-1). 106049b5e25fSSatish Balay jl: list of rows to be added to uneliminated rows 106149b5e25fSSatish Balay i>= k: jl(i) is the first row to be added to row i 106249b5e25fSSatish Balay i< k: jl(i) is the row following row i in some list of rows 106349b5e25fSSatish Balay jl(i) = mbs indicates the end of a list 106449b5e25fSSatish Balay il(i): points to the first nonzero element in columns k,...,mbs-1 of 106549b5e25fSSatish Balay row i of U */ 1066997a0949SHong Zhang nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); 1067997a0949SHong Zhang ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 106828de702eSHong Zhang jl = il + mbs; 1069997a0949SHong Zhang rtmp = (MatScalar*)(jl + mbs); 1070997a0949SHong Zhang 10713cea4cbeSHong Zhang sctx.shift_amount = 0; 10723cea4cbeSHong Zhang sctx.nshift = 0; 1073997a0949SHong Zhang do { 10743cea4cbeSHong Zhang sctx.chshift = PETSC_FALSE; 107549b5e25fSSatish Balay for (i=0; i<mbs; i++) { 107649b5e25fSSatish Balay rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 107749b5e25fSSatish Balay } 1078997a0949SHong Zhang 107949b5e25fSSatish Balay for (k = 0; k<mbs; k++){ 1080997a0949SHong Zhang /*initialize k-th row by the perm[k]-th row of A */ 108149b5e25fSSatish Balay jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 10827701de02SHong Zhang bval = ba + bi[k]; 108349b5e25fSSatish Balay for (j = jmin; j < jmax; j++){ 1084997a0949SHong Zhang col = rip[aj[j]]; 1085997a0949SHong Zhang rtmp[col] = aa[j]; 10867701de02SHong Zhang *bval++ = 0.0; /* for in-place factorization */ 108749b5e25fSSatish Balay } 10883cea4cbeSHong Zhang 10893cea4cbeSHong Zhang /* shift the diagonal of the matrix */ 10903cea4cbeSHong Zhang if (sctx.nshift) rtmp[k] += sctx.shift_amount; 109149b5e25fSSatish Balay 109291602c66SHong Zhang /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 109349b5e25fSSatish Balay dk = rtmp[k]; 109449b5e25fSSatish Balay i = jl[k]; /* first row to be added to k_th row */ 109549b5e25fSSatish Balay 1096057f5ba7SHong Zhang while (i < k){ 109749b5e25fSSatish Balay nexti = jl[i]; /* next row to be added to k_th row */ 109849b5e25fSSatish Balay 1099fff829cfSHong Zhang /* compute multiplier, update diag(k) and U(i,k) */ 110049b5e25fSSatish Balay ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 1101997a0949SHong Zhang uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 110249b5e25fSSatish Balay dk += uikdi*ba[ili]; 1103658e7b3eSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 110449b5e25fSSatish Balay 1105997a0949SHong Zhang /* add multiple of row i to k-th row */ 110649b5e25fSSatish Balay jmin = ili + 1; jmax = bi[i+1]; 110749b5e25fSSatish Balay if (jmin < jmax){ 110849b5e25fSSatish Balay for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 1109*dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*(jmax-jmin));CHKERRQ(ierr); 1110187a9f4bSHong Zhang 1111fff829cfSHong Zhang /* update il and jl for row i */ 1112fff829cfSHong Zhang il[i] = jmin; 1113fff829cfSHong Zhang j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 111449b5e25fSSatish Balay } 1115ab27746eSHong Zhang i = nexti; 111649b5e25fSSatish Balay } 111749b5e25fSSatish Balay 11183cea4cbeSHong Zhang /* shift the diagonals when zero pivot is detected */ 11193cea4cbeSHong Zhang /* compute rs=sum of abs(off-diagonal) */ 11203cea4cbeSHong Zhang rs = 0.0; 1121997a0949SHong Zhang jmin = bi[k]+1; 1122997a0949SHong Zhang nz = bi[k+1] - jmin; 1123997a0949SHong Zhang if (nz){ 1124997a0949SHong Zhang bcol = bj + jmin; 1125997a0949SHong Zhang while (nz--){ 112689f845c8SHong Zhang rs += PetscAbsScalar(rtmp[*bcol]); 112789f845c8SHong Zhang bcol++; 1128997a0949SHong Zhang } 1129997a0949SHong Zhang } 11303cea4cbeSHong Zhang 11313cea4cbeSHong Zhang sctx.rs = rs; 11323cea4cbeSHong Zhang sctx.pv = dk; 113345938b79SHong Zhang ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr); 113445938b79SHong Zhang if (newshift == 1) break; /* sctx.shift_amount is updated */ 113549b5e25fSSatish Balay 1136997a0949SHong Zhang /* copy data into U(k,:) */ 113798a8e78dSHong Zhang ba[bi[k]] = 1.0/dk; /* U(k,k) */ 1138fff829cfSHong Zhang jmin = bi[k]+1; jmax = bi[k+1]; 113949b5e25fSSatish Balay if (jmin < jmax) { 114049b5e25fSSatish Balay for (j=jmin; j<jmax; j++){ 1141997a0949SHong Zhang col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 114249b5e25fSSatish Balay } 114398a8e78dSHong Zhang /* add the k-th row into il and jl */ 114449b5e25fSSatish Balay il[k] = jmin; 114598a8e78dSHong Zhang i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 114649b5e25fSSatish Balay } 114749b5e25fSSatish Balay } 11483cea4cbeSHong Zhang } while (sctx.chshift); 114949b5e25fSSatish Balay ierr = PetscFree(il);CHKERRQ(ierr); 115098a8e78dSHong Zhang if (a->permute){ierr = PetscFree(aa);CHKERRQ(ierr);} 115149b5e25fSSatish Balay 115249b5e25fSSatish Balay ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 1153db4efbfdSBarry Smith C->ops->solve = MatSolve_SeqSBAIJ_1; 1154db4efbfdSBarry Smith C->ops->solves = MatSolves_SeqSBAIJ_1; 1155db4efbfdSBarry Smith C->ops->solvetranspose = MatSolve_SeqSBAIJ_1; 1156db4efbfdSBarry Smith C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1; 1157db4efbfdSBarry Smith C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1; 115849b5e25fSSatish Balay C->assembled = PETSC_TRUE; 11595c0bcdfcSHong Zhang C->preallocated = PETSC_TRUE; 1160d0f46423SBarry Smith ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); 11613cea4cbeSHong Zhang if (sctx.nshift){ 11623cea4cbeSHong Zhang if (shiftnz) { 11631e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 11643cea4cbeSHong Zhang } else if (shiftpd) { 11651e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 1166997a0949SHong Zhang } 1167997a0949SHong Zhang } 116849b5e25fSSatish Balay PetscFunctionReturn(0); 116949b5e25fSSatish Balay } 117049b5e25fSSatish Balay 1171671cb588SHong Zhang /* 1172671cb588SHong Zhang Version for when blocks are 1 by 1 Using natural ordering 1173671cb588SHong Zhang */ 11744a2ae208SSatish Balay #undef __FUNCT__ 11754a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering" 11760481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info) 1177671cb588SHong Zhang { 1178671cb588SHong Zhang Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data; 1179dfbe8321SBarry Smith PetscErrorCode ierr; 118013f74950SBarry Smith PetscInt i,j,mbs = a->mbs; 118113f74950SBarry Smith PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 11823cea4cbeSHong Zhang PetscInt k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz; 1183653a6975SHong Zhang MatScalar *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval; 11843cea4cbeSHong Zhang PetscReal zeropivot,rs,shiftnz; 1185fbf22428SSatish Balay PetscReal shiftpd; 11863cea4cbeSHong Zhang ChShift_Ctx sctx; 11873cea4cbeSHong Zhang PetscInt newshift; 1188653a6975SHong Zhang 1189653a6975SHong Zhang PetscFunctionBegin; 1190653a6975SHong Zhang /* initialization */ 11913cea4cbeSHong Zhang shiftnz = info->shiftnz; 11923cea4cbeSHong Zhang shiftpd = info->shiftpd; 11933cea4cbeSHong Zhang zeropivot = info->zeropivot; 11943cea4cbeSHong Zhang 1195653a6975SHong Zhang /* il and jl record the first nonzero element in each row of the accessing 1196653a6975SHong Zhang window U(0:k, k:mbs-1). 1197653a6975SHong Zhang jl: list of rows to be added to uneliminated rows 1198653a6975SHong Zhang i>= k: jl(i) is the first row to be added to row i 1199653a6975SHong Zhang i< k: jl(i) is the row following row i in some list of rows 1200653a6975SHong Zhang jl(i) = mbs indicates the end of a list 1201653a6975SHong Zhang il(i): points to the first nonzero element in U(i,k:mbs-1) 1202653a6975SHong Zhang */ 1203653a6975SHong Zhang ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 120413f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); 1205653a6975SHong Zhang jl = il + mbs; 1206b00f7748SHong Zhang 12073cea4cbeSHong Zhang sctx.shift_amount = 0; 12083cea4cbeSHong Zhang sctx.nshift = 0; 1209b00f7748SHong Zhang do { 12103cea4cbeSHong Zhang sctx.chshift = PETSC_FALSE; 1211653a6975SHong Zhang for (i=0; i<mbs; i++) { 1212653a6975SHong Zhang rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 1213653a6975SHong Zhang } 1214653a6975SHong Zhang 1215997a0949SHong Zhang for (k = 0; k<mbs; k++){ 1216653a6975SHong Zhang /*initialize k-th row with elements nonzero in row perm(k) of A */ 1217653a6975SHong Zhang nz = ai[k+1] - ai[k]; 1218653a6975SHong Zhang acol = aj + ai[k]; 1219653a6975SHong Zhang aval = aa + ai[k]; 1220653a6975SHong Zhang bval = ba + bi[k]; 1221653a6975SHong Zhang while (nz -- ){ 1222653a6975SHong Zhang rtmp[*acol++] = *aval++; 1223653a6975SHong Zhang *bval++ = 0.0; /* for in-place factorization */ 1224653a6975SHong Zhang } 12253cea4cbeSHong Zhang 12263cea4cbeSHong Zhang /* shift the diagonal of the matrix */ 12273cea4cbeSHong Zhang if (sctx.nshift) rtmp[k] += sctx.shift_amount; 1228653a6975SHong Zhang 1229653a6975SHong Zhang /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 1230653a6975SHong Zhang dk = rtmp[k]; 1231653a6975SHong Zhang i = jl[k]; /* first row to be added to k_th row */ 1232653a6975SHong Zhang 1233653a6975SHong Zhang while (i < k){ 1234653a6975SHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 1235653a6975SHong Zhang /* compute multiplier, update D(k) and U(i,k) */ 1236653a6975SHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 1237653a6975SHong Zhang uikdi = - ba[ili]*ba[bi[i]]; 1238653a6975SHong Zhang dk += uikdi*ba[ili]; 1239653a6975SHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 1240653a6975SHong Zhang 1241653a6975SHong Zhang /* add multiple of row i to k-th row ... */ 1242653a6975SHong Zhang jmin = ili + 1; 1243653a6975SHong Zhang nz = bi[i+1] - jmin; 1244653a6975SHong Zhang if (nz > 0){ 1245653a6975SHong Zhang bcol = bj + jmin; 1246653a6975SHong Zhang bval = ba + jmin; 1247*dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr); 1248653a6975SHong Zhang while (nz --) rtmp[*bcol++] += uikdi*(*bval++); 1249187a9f4bSHong Zhang 1250997a0949SHong Zhang /* update il and jl for i-th row */ 1251997a0949SHong Zhang il[i] = jmin; 1252997a0949SHong Zhang j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 1253653a6975SHong Zhang } 1254653a6975SHong Zhang i = nexti; 1255653a6975SHong Zhang } 1256653a6975SHong Zhang 12573cea4cbeSHong Zhang /* shift the diagonals when zero pivot is detected */ 12583cea4cbeSHong Zhang /* compute rs=sum of abs(off-diagonal) */ 12593cea4cbeSHong Zhang rs = 0.0; 126021c26570Svictorle jmin = bi[k]+1; 126121c26570Svictorle nz = bi[k+1] - jmin; 126221c26570Svictorle if (nz){ 126321c26570Svictorle bcol = bj + jmin; 126421c26570Svictorle while (nz--){ 126589f845c8SHong Zhang rs += PetscAbsScalar(rtmp[*bcol]); 126689f845c8SHong Zhang bcol++; 126721c26570Svictorle } 126821c26570Svictorle } 12693cea4cbeSHong Zhang 12703cea4cbeSHong Zhang sctx.rs = rs; 12713cea4cbeSHong Zhang sctx.pv = dk; 127245938b79SHong Zhang ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr); 127345938b79SHong Zhang if (newshift == 1) break; /* sctx.shift_amount is updated */ 1274653a6975SHong Zhang 1275997a0949SHong Zhang /* copy data into U(k,:) */ 1276653a6975SHong Zhang ba[bi[k]] = 1.0/dk; 1277653a6975SHong Zhang jmin = bi[k]+1; 1278653a6975SHong Zhang nz = bi[k+1] - jmin; 1279653a6975SHong Zhang if (nz){ 1280653a6975SHong Zhang bcol = bj + jmin; 1281653a6975SHong Zhang bval = ba + jmin; 1282653a6975SHong Zhang while (nz--){ 1283653a6975SHong Zhang *bval++ = rtmp[*bcol]; 1284653a6975SHong Zhang rtmp[*bcol++] = 0.0; 1285653a6975SHong Zhang } 1286997a0949SHong Zhang /* add k-th row into il and jl */ 1287653a6975SHong Zhang il[k] = jmin; 1288997a0949SHong Zhang i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 1289653a6975SHong Zhang } 1290b00f7748SHong Zhang } /* end of for (k = 0; k<mbs; k++) */ 12913cea4cbeSHong Zhang } while (sctx.chshift); 1292653a6975SHong Zhang ierr = PetscFree(rtmp);CHKERRQ(ierr); 1293653a6975SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 1294653a6975SHong Zhang 1295db4efbfdSBarry Smith C->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1296db4efbfdSBarry Smith C->ops->solves = MatSolves_SeqSBAIJ_1; 1297db4efbfdSBarry Smith C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1298db4efbfdSBarry Smith C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering; 1299db4efbfdSBarry Smith C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering; 1300db4efbfdSBarry Smith 1301653a6975SHong Zhang C->assembled = PETSC_TRUE; 1302653a6975SHong Zhang C->preallocated = PETSC_TRUE; 1303d0f46423SBarry Smith ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); 13043cea4cbeSHong Zhang if (sctx.nshift){ 13053cea4cbeSHong Zhang if (shiftnz) { 13061e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 13073cea4cbeSHong Zhang } else if (shiftpd) { 13081e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 1309b00f7748SHong Zhang } 1310fb10cecfSBarry Smith } 1311653a6975SHong Zhang PetscFunctionReturn(0); 1312653a6975SHong Zhang } 1313653a6975SHong Zhang 1314653a6975SHong Zhang #undef __FUNCT__ 13154a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ" 13160481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,const MatFactorInfo *info) 131749b5e25fSSatish Balay { 1318dfbe8321SBarry Smith PetscErrorCode ierr; 131949b5e25fSSatish Balay Mat C; 132049b5e25fSSatish Balay 132149b5e25fSSatish Balay PetscFunctionBegin; 1322719d5645SBarry Smith ierr = MatGetFactor(A,"petsc",MAT_FACTOR_CHOLESKY,&C);CHKERRQ(ierr); 1323719d5645SBarry Smith ierr = MatCholeskyFactorSymbolic(C,A,perm,info);CHKERRQ(ierr); 1324719d5645SBarry Smith ierr = MatCholeskyFactorNumeric(C,A,info);CHKERRQ(ierr); 1325db4efbfdSBarry Smith A->ops->solve = C->ops->solve; 1326db4efbfdSBarry Smith A->ops->solvetranspose = C->ops->solvetranspose; 13274445ddedSHong Zhang ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 132849b5e25fSSatish Balay PetscFunctionReturn(0); 132949b5e25fSSatish Balay } 133049b5e25fSSatish Balay 133149b5e25fSSatish Balay 1332