#include <../src/mat/impls/baij/seq/baij.h> #include <../src/mat/impls/sbaij/seq/sbaij.h> #include #include /* input: F -- numeric factor output: nneg, nzero, npos: matrix inertia */ #undef __FUNCT__ #define __FUNCT__ "MatGetInertia_SeqSBAIJ" PetscErrorCode MatGetInertia_SeqSBAIJ(Mat F,PetscInt *nneig,PetscInt *nzero,PetscInt *npos) { Mat_SeqSBAIJ *fact_ptr = (Mat_SeqSBAIJ*)F->data; MatScalar *dd = fact_ptr->a; PetscInt mbs =fact_ptr->mbs,bs=F->rmap->bs,i,nneig_tmp,npos_tmp,*fi = fact_ptr->diag; PetscFunctionBegin; if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for bs: %D >1 yet",bs); nneig_tmp = 0; npos_tmp = 0; for (i=0; i 0.0) npos_tmp++; else if (PetscRealPart(dd[*fi]) < 0.0) nneig_tmp++; fi++; } if (nneig) *nneig = nneig_tmp; if (npos) *npos = npos_tmp; if (nzero) *nzero = mbs - nneig_tmp - npos_tmp; PetscFunctionReturn(0); } /* Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP. Use Modified Sparse Row (MSR) storage for u and ju. See page 85, "Iterative Methods ..." by Saad. */ #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ_MSR" PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(Mat F,Mat A,IS perm,const MatFactorInfo *info) { Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b; PetscErrorCode ierr; const PetscInt *rip,*ai,*aj; PetscInt i,mbs = a->mbs,*jutmp,bs = A->rmap->bs,bs2=a->bs2; PetscInt m,reallocs = 0,prow; PetscInt *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd; PetscReal f = info->fill; PetscBool perm_identity; PetscFunctionBegin; /* check whether perm is the identity mapping */ ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); if (perm_identity) { /* without permutation */ a->permute = PETSC_FALSE; ai = a->i; aj = a->j; } else { /* non-trivial permutation */ a->permute = PETSC_TRUE; ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr); ai = a->inew; aj = a->jnew; } /* initialization */ ierr = PetscMalloc1(mbs+1,&iu);CHKERRQ(ierr); umax = (PetscInt)(f*ai[mbs] + 1); umax += mbs + 1; ierr = PetscMalloc1(umax,&ju);CHKERRQ(ierr); iu[0] = mbs+1; juidx = mbs + 1; /* index for ju */ /* jl linked list for pivot row -- linked list for col index */ ierr = PetscMalloc2(mbs,&jl,mbs,&q);CHKERRQ(ierr); for (i=0; i k) { qm = k; do { m = qm; qm = q[m]; } while (qm < vj); if (qm == vj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Duplicate entry in A\n"); nzk++; q[m] = vj; q[vj] = qm; } /* if (vj > k) */ } /* for (j=jmin; j 0) { i = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */ jl[k] = jl[i]; jl[i] = k; } iu[k+1] = iu[k] + nzk; /* allocate more space to ju if needed */ if (iu[k+1] > umax) { /* estimate how much additional space we will need */ /* use the strategy suggested by David Hysom */ /* just double the memory each time */ maxadd = umax; if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2; umax += maxadd; /* allocate a longer ju */ ierr = PetscMalloc1(umax,&jutmp);CHKERRQ(ierr); ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(PetscInt));CHKERRQ(ierr); ierr = PetscFree(ju);CHKERRQ(ierr); ju = jutmp; reallocs++; /* count how many times we realloc */ } /* save nonzero structure of k-th row in ju */ i=k; while (nzk--) { i = q[i]; ju[juidx++] = i; } } #if defined(PETSC_USE_INFO) if (ai[mbs] != 0) { PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]); ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);CHKERRQ(ierr); ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr); ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g);\n",(double)af);CHKERRQ(ierr); ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); } else { ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); } #endif ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); ierr = PetscFree2(jl,q);CHKERRQ(ierr); /* put together the new matrix */ ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(F,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); /* ierr = PetscLogObjectParent((PetscObject)B,(PetscObject)iperm);CHKERRQ(ierr); */ b = (Mat_SeqSBAIJ*)(F)->data; b->singlemalloc = PETSC_FALSE; b->free_a = PETSC_TRUE; b->free_ij = PETSC_TRUE; ierr = PetscMalloc1((iu[mbs]+1)*bs2,&b->a);CHKERRQ(ierr); b->j = ju; b->i = iu; b->diag = 0; b->ilen = 0; b->imax = 0; b->row = perm; b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); b->icol = perm; ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); ierr = PetscMalloc1(bs*mbs+bs,&b->solve_work);CHKERRQ(ierr); /* In b structure: Free imax, ilen, old a, old j. Allocate idnew, solve_work, new a, new j */ ierr = PetscLogObjectMemory((PetscObject)F,(iu[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); b->maxnz = b->nz = iu[mbs]; (F)->info.factor_mallocs = reallocs; (F)->info.fill_ratio_given = f; if (ai[mbs] != 0) { (F)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]); } else { (F)->info.fill_ratio_needed = 0.0; } ierr = MatSeqSBAIJSetNumericFactorization_inplace(F,perm_identity);CHKERRQ(ierr); PetscFunctionReturn(0); } /* Symbolic U^T*D*U factorization for SBAIJ format. See MatICCFactorSymbolic_SeqAIJ() for description of its data structure. */ #include #include <../src/mat/utils/freespace.h> #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ" PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info) { Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; Mat_SeqSBAIJ *b; PetscErrorCode ierr; PetscBool perm_identity,missing; PetscReal fill = info->fill; const PetscInt *rip,*ai=a->i,*aj=a->j; PetscInt i,mbs=a->mbs,bs=A->rmap->bs,reallocs=0,prow; PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; PetscInt nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr,*udiag; PetscFreeSpaceList free_space=NULL,current_space=NULL; PetscBT lnkbt; PetscFunctionBegin; if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n); ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr); if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i); if (bs > 1) { ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(fact,A,perm,info);CHKERRQ(ierr); PetscFunctionReturn(0); } /* check whether perm is the identity mapping */ ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); if (!perm_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format"); a->permute = PETSC_FALSE; ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); /* initialization */ ierr = PetscMalloc1(mbs+1,&ui);CHKERRQ(ierr); ierr = PetscMalloc1(mbs+1,&udiag);CHKERRQ(ierr); ui[0] = 0; /* jl: linked list for storing indices of the pivot rows il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */ ierr = PetscMalloc4(mbs,&ui_ptr,mbs,&il,mbs,&jl,mbs,&cols);CHKERRQ(ierr); for (i=0; ilocal_remainingarray,lnkbt);CHKERRQ(ierr); /* add the k-th row into il and jl */ if (nzk > 1) { i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */ jl[k] = jl[i]; jl[i] = k; il[k] = ui[k] + 1; } ui_ptr[k] = current_space->array; current_space->array += nzk; current_space->local_used += nzk; current_space->local_remaining -= nzk; ui[k+1] = ui[k] + nzk; } ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr); /* destroy list of free space and other temporary array(s) */ ierr = PetscMalloc1(ui[mbs]+1,&uj);CHKERRQ(ierr); ierr = PetscFreeSpaceContiguous_Cholesky(&free_space,uj,mbs,ui,udiag);CHKERRQ(ierr); /* store matrix factor */ ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); /* put together the new matrix in MATSEQSBAIJ format */ ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); b = (Mat_SeqSBAIJ*)fact->data; b->singlemalloc = PETSC_FALSE; b->free_a = PETSC_TRUE; b->free_ij = PETSC_TRUE; ierr = PetscMalloc1(ui[mbs]+1,&b->a);CHKERRQ(ierr); b->j = uj; b->i = ui; b->diag = udiag; b->free_diag = PETSC_TRUE; b->ilen = 0; b->imax = 0; b->row = perm; b->icol = perm; ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ ierr = PetscMalloc1(mbs+1,&b->solve_work);CHKERRQ(ierr); ierr = PetscLogObjectMemory((PetscObject)fact,ui[mbs]*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); b->maxnz = b->nz = ui[mbs]; fact->info.factor_mallocs = reallocs; fact->info.fill_ratio_given = fill; if (ai[mbs] != 0) { fact->info.fill_ratio_needed = ((PetscReal)ui[mbs])/ai[mbs]; } else { fact->info.fill_ratio_needed = 0.0; } #if defined(PETSC_USE_INFO) if (ai[mbs] != 0) { PetscReal af = fact->info.fill_ratio_needed; ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);CHKERRQ(ierr); ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr); ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);CHKERRQ(ierr); } else { ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); } #endif fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ_inplace" PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(Mat fact,Mat A,IS perm,const MatFactorInfo *info) { Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; Mat_SeqSBAIJ *b; PetscErrorCode ierr; PetscBool perm_identity,missing; PetscReal fill = info->fill; const PetscInt *rip,*ai,*aj; PetscInt i,mbs=a->mbs,bs=A->rmap->bs,reallocs=0,prow,d; PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; PetscInt nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr; PetscFreeSpaceList free_space=NULL,current_space=NULL; PetscBT lnkbt; PetscFunctionBegin; ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr); if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d); /* This code originally uses Modified Sparse Row (MSR) storage (see page 85, "Iterative Methods ..." by Saad) for the output matrix B - bad choise! Then it is rewritten so the factor B takes seqsbaij format. However the associated MatCholeskyFactorNumeric_() have not been modified for the cases of bs>1 or !perm_identity, thus the original code in MSR format is still used for these cases. The code below should replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR() whenever MatCholeskyFactorNumeric_() is modified for using sbaij symbolic factor. */ if (bs > 1) { ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(fact,A,perm,info);CHKERRQ(ierr); PetscFunctionReturn(0); } /* check whether perm is the identity mapping */ ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); if (perm_identity) { a->permute = PETSC_FALSE; ai = a->i; aj = a->j; } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format"); ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); /* initialization */ ierr = PetscMalloc1(mbs+1,&ui);CHKERRQ(ierr); ui[0] = 0; /* jl: linked list for storing indices of the pivot rows il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */ ierr = PetscMalloc4(mbs,&ui_ptr,mbs,&il,mbs,&jl,mbs,&cols);CHKERRQ(ierr); for (i=0; ilocal_remainingarray,lnkbt);CHKERRQ(ierr); /* add the k-th row into il and jl */ if (nzk-1 > 0) { i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */ jl[k] = jl[i]; jl[i] = k; il[k] = ui[k] + 1; } ui_ptr[k] = current_space->array; current_space->array += nzk; current_space->local_used += nzk; current_space->local_remaining -= nzk; ui[k+1] = ui[k] + nzk; } ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr); /* destroy list of free space and other temporary array(s) */ ierr = PetscMalloc1(ui[mbs]+1,&uj);CHKERRQ(ierr); ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); /* put together the new matrix in MATSEQSBAIJ format */ ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); b = (Mat_SeqSBAIJ*)fact->data; b->singlemalloc = PETSC_FALSE; b->free_a = PETSC_TRUE; b->free_ij = PETSC_TRUE; ierr = PetscMalloc1(ui[mbs]+1,&b->a);CHKERRQ(ierr); b->j = uj; b->i = ui; b->diag = 0; b->ilen = 0; b->imax = 0; b->row = perm; b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); b->icol = perm; ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); ierr = PetscMalloc1(mbs+1,&b->solve_work);CHKERRQ(ierr); ierr = PetscLogObjectMemory((PetscObject)fact,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); b->maxnz = b->nz = ui[mbs]; fact->info.factor_mallocs = reallocs; fact->info.fill_ratio_given = fill; if (ai[mbs] != 0) { fact->info.fill_ratio_needed = ((PetscReal)ui[mbs])/ai[mbs]; } else { fact->info.fill_ratio_needed = 0.0; } #if defined(PETSC_USE_INFO) if (ai[mbs] != 0) { PetscReal af = fact->info.fill_ratio_needed; ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);CHKERRQ(ierr); ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr); ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);CHKERRQ(ierr); } else { ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); } #endif ierr = MatSeqSBAIJSetNumericFactorization_inplace(fact,perm_identity);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N" PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat C,Mat A,const MatFactorInfo *info) { Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data; IS perm = b->row; PetscErrorCode ierr; const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j; PetscInt i,j; PetscInt *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; PetscInt bs =A->rmap->bs,bs2 = a->bs2; MatScalar *ba = b->a,*aa,*ap,*dk,*uik; MatScalar *u,*diag,*rtmp,*rtmp_ptr; MatScalar *work; PetscInt *pivots; PetscFunctionBegin; /* initialization */ ierr = PetscCalloc1(bs2*mbs,&rtmp);CHKERRQ(ierr); ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr); for (i=0; ipermute) { ai = a->i; aj = a->j; aa = a->a; } else { ai = a->inew; aj = a->jnew; ierr = PetscMalloc1(bs2*ai[mbs],&aa);CHKERRQ(ierr); ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); ierr = PetscMalloc1(ai[mbs],&a2anew);CHKERRQ(ierr); ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr); for (i=0; i aj[j]) { /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */ ap = aa + j*bs2; /* ptr to the beginning of j-th block of aa */ for (k=0; kpermute) { ierr = PetscFree(aa);CHKERRQ(ierr); } ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr); C->ops->solve = MatSolve_SeqSBAIJ_N_inplace; C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_inplace; C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_N_inplace; C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_N_inplace; C->assembled = PETSC_TRUE; C->preallocated = PETSC_TRUE; ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering" PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info) { Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data; PetscErrorCode ierr; PetscInt i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; PetscInt *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; PetscInt bs =A->rmap->bs,bs2 = a->bs2; MatScalar *ba = b->a,*aa,*ap,*dk,*uik; MatScalar *u,*diag,*rtmp,*rtmp_ptr; MatScalar *work; PetscInt *pivots; PetscFunctionBegin; ierr = PetscCalloc1(bs2*mbs,&rtmp);CHKERRQ(ierr); ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr); for (i=0; ii; aj = a->j; aa = a->a; /* for each row k */ for (k = 0; kops->solve = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace; C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace; C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace; C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace; C->assembled = PETSC_TRUE; C->preallocated = PETSC_TRUE; ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ PetscFunctionReturn(0); } /* Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP. Version for blocks 2 by 2. */ #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2" PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat C,Mat A,const MatFactorInfo *info) { Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data; IS perm = b->row; PetscErrorCode ierr; const PetscInt *ai,*aj,*perm_ptr; PetscInt i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; PetscInt *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; MatScalar *ba = b->a,*aa,*ap; MatScalar *u,*diag,*rtmp,*rtmp_ptr,dk[4],uik[4]; PetscReal shift = info->shiftamount; PetscFunctionBegin; /* initialization */ /* il and jl record the first nonzero element in each row of the accessing window U(0:k, k:mbs-1). jl: list of rows to be added to uneliminated rows i>= k: jl(i) is the first row to be added to row i i< k: jl(i) is the row following row i in some list of rows jl(i) = mbs indicates the end of a list il(i): points to the first nonzero element in columns k,...,mbs-1 of row i of U */ ierr = PetscCalloc1(4*mbs,&rtmp);CHKERRQ(ierr); ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr); for (i=0; ipermute) { ai = a->i; aj = a->j; aa = a->a; } else { ai = a->inew; aj = a->jnew; ierr = PetscMalloc1(4*ai[mbs],&aa);CHKERRQ(ierr); ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); ierr = PetscMalloc1(ai[mbs],&a2anew);CHKERRQ(ierr); ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr); for (i=0; i aj[j]) { /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */ ap = aa + j*4; /* ptr to the beginning of the block */ dk[1] = ap[1]; /* swap ap[1] and ap[2] */ ap[1] = ap[2]; ap[2] = dk[1]; } } } ierr = PetscFree(a2anew);CHKERRQ(ierr); } /* for each row k */ for (k = 0; kpermute) { ierr = PetscFree(aa);CHKERRQ(ierr); } ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr); C->ops->solve = MatSolve_SeqSBAIJ_2_inplace; C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_inplace; C->assembled = PETSC_TRUE; C->preallocated = PETSC_TRUE; ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ PetscFunctionReturn(0); } /* Version for when blocks are 2 by 2 Using natural ordering */ #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering" PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info) { Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data; PetscErrorCode ierr; PetscInt i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; PetscInt *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; MatScalar *ba = b->a,*aa,*ap,dk[8],uik[8]; MatScalar *u,*diag,*rtmp,*rtmp_ptr; PetscReal shift = info->shiftamount; PetscFunctionBegin; /* initialization */ /* il and jl record the first nonzero element in each row of the accessing window U(0:k, k:mbs-1). jl: list of rows to be added to uneliminated rows i>= k: jl(i) is the first row to be added to row i i< k: jl(i) is the row following row i in some list of rows jl(i) = mbs indicates the end of a list il(i): points to the first nonzero element in columns k,...,mbs-1 of row i of U */ ierr = PetscCalloc1(4*mbs,&rtmp);CHKERRQ(ierr); ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr); for (i=0; ii; aj = a->j; aa = a->a; /* for each row k */ for (k = 0; kops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace; C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace; C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace; C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace; C->assembled = PETSC_TRUE; C->preallocated = PETSC_TRUE; ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ PetscFunctionReturn(0); } /* Numeric U^T*D*U factorization for SBAIJ format. Version for blocks are 1 by 1. */ #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace" PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace(Mat C,Mat A,const MatFactorInfo *info) { Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ*)C->data; IS ip=b->row; PetscErrorCode ierr; const PetscInt *ai,*aj,*rip; PetscInt *a2anew,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j,*bcol; PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; MatScalar *rtmp,*ba=b->a,*bval,*aa,dk,uikdi; PetscReal rs; FactorShiftCtx sctx; PetscFunctionBegin; /* MatPivotSetUp(): initialize shift context sctx */ ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr); ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); if (!a->permute) { ai = a->i; aj = a->j; aa = a->a; } else { ai = a->inew; aj = a->jnew; nz = ai[mbs]; ierr = PetscMalloc1(nz,&aa);CHKERRQ(ierr); a2anew = a->a2anew; bval = a->a; for (j=0; j= k: jl(i) is the first row to be added to row i i< k: jl(i) is the row following row i in some list of rows jl(i) = mbs indicates the end of a list il(i): points to the first nonzero element in columns k,...,mbs-1 of row i of U */ ierr = PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&jl);CHKERRQ(ierr); do { sctx.newshift = PETSC_FALSE; for (i=0; ipermute) {ierr = PetscFree(aa);CHKERRQ(ierr);} ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); C->ops->solve = MatSolve_SeqSBAIJ_1_inplace; C->ops->solves = MatSolves_SeqSBAIJ_1_inplace; C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace; C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_inplace; C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_inplace; C->assembled = PETSC_TRUE; C->preallocated = PETSC_TRUE; ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); if (sctx.nshift) { if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) { ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr); } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr); } } PetscFunctionReturn(0); } /* Version for when blocks are 1 by 1 Using natural ordering under new datastructure Modified from MatCholeskyFactorNumeric_SeqAIJ() */ #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering" PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info) { Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)B->data; PetscErrorCode ierr; PetscInt i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bdiag=b->diag,*bjtmp; PetscInt *ai=a->i,*aj=a->j,*ajtmp; PetscInt k,jmin,jmax,*c2r,*il,col,nexti,ili,nz; MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; FactorShiftCtx sctx; PetscReal rs; MatScalar d,*v; PetscFunctionBegin; ierr = PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&c2r);CHKERRQ(ierr); /* MatPivotSetUp(): initialize shift context sctx */ ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr); if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */ sctx.shift_top = info->zeropivot; ierr = PetscMemzero(rtmp,mbs*sizeof(MatScalar));CHKERRQ(ierr); for (i=0; idiag[i]]; rtmp[i] += -PetscRealPart(d); /* diagonal entry */ ajtmp = aj + ai[i] + 1; /* exclude diagonal */ v = aa + ai[i] + 1; nz = ai[i+1] - ai[i] - 1; for (j=0; j sctx.shift_top) sctx.shift_top = PetscRealPart(rtmp[i]); } sctx.shift_top *= 1.1; sctx.nshift_max = 5; sctx.shift_lo = 0.; sctx.shift_hi = 1.; } /* allocate working arrays c2r: linked list, keep track of pivot rows for a given column. c2r[col]: head of the list for a given col il: for active k row, il[i] gives the index of the 1st nonzero entry in U[i,k:n-1] in bj and ba arrays */ do { sctx.newshift = PETSC_FALSE; for (i=0; iops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; B->ops->solves = MatSolves_SeqSBAIJ_1; B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; B->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering; B->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering; B->assembled = PETSC_TRUE; B->preallocated = PETSC_TRUE; ierr = PetscLogFlops(B->rmap->n);CHKERRQ(ierr); /* MatPivotView() */ if (sctx.nshift) { if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %g, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,(double)sctx.shift_amount,(double)sctx.shift_fraction,(double)sctx.shift_top);CHKERRQ(ierr); } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) { ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr); } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) { ierr = PetscInfo2(A,"number of shift_inblocks applied %D, each shift_amount %g\n",sctx.nshift,(double)info->shiftamount);CHKERRQ(ierr); } } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace" PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info) { Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ*)C->data; PetscErrorCode ierr; PetscInt i,j,mbs = a->mbs; PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; PetscInt k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz; MatScalar *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval; PetscReal rs; FactorShiftCtx sctx; PetscFunctionBegin; /* MatPivotSetUp(): initialize shift context sctx */ ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr); /* initialization */ /* il and jl record the first nonzero element in each row of the accessing window U(0:k, k:mbs-1). jl: list of rows to be added to uneliminated rows i>= k: jl(i) is the first row to be added to row i i< k: jl(i) is the row following row i in some list of rows jl(i) = mbs indicates the end of a list il(i): points to the first nonzero element in U(i,k:mbs-1) */ ierr = PetscMalloc1(mbs,&rtmp);CHKERRQ(ierr); ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr); do { sctx.newshift = PETSC_FALSE; for (i=0; i 0) { bcol = bj + jmin; bval = ba + jmin; ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr); while (nz--) rtmp[*bcol++] += uikdi*(*bval++); /* update il and jl for i-th row */ il[i] = jmin; j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; } i = nexti; } /* shift the diagonals when zero pivot is detected */ /* compute rs=sum of abs(off-diagonal) */ rs = 0.0; jmin = bi[k]+1; nz = bi[k+1] - jmin; if (nz) { bcol = bj + jmin; while (nz--) { rs += PetscAbsScalar(rtmp[*bcol]); bcol++; } } sctx.rs = rs; sctx.pv = dk; ierr = MatPivotCheck(C,A,info,&sctx,k);CHKERRQ(ierr); if (sctx.newshift) break; /* sctx.shift_amount is updated */ dk = sctx.pv; /* copy data into U(k,:) */ ba[bi[k]] = 1.0/dk; jmin = bi[k]+1; nz = bi[k+1] - jmin; if (nz) { bcol = bj + jmin; bval = ba + jmin; while (nz--) { *bval++ = rtmp[*bcol]; rtmp[*bcol++] = 0.0; } /* add k-th row into il and jl */ il[k] = jmin; i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; } } /* end of for (k = 0; kops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; C->ops->solves = MatSolves_SeqSBAIJ_1_inplace; C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace; C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace; C->assembled = PETSC_TRUE; C->preallocated = PETSC_TRUE; ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); if (sctx.nshift) { if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) { ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr); } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr); } } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ" PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,const MatFactorInfo *info) { PetscErrorCode ierr; Mat C; PetscFunctionBegin; ierr = MatGetFactor(A,"petsc",MAT_FACTOR_CHOLESKY,&C);CHKERRQ(ierr); ierr = MatCholeskyFactorSymbolic(C,A,perm,info);CHKERRQ(ierr); ierr = MatCholeskyFactorNumeric(C,A,info);CHKERRQ(ierr); A->ops->solve = C->ops->solve; A->ops->solvetranspose = C->ops->solvetranspose; ierr = MatHeaderMerge(A,&C);CHKERRQ(ierr); PetscFunctionReturn(0); }