#define PETSCMAT_DLL #include "src/mat/impls/baij/seq/baij.h" #include "src/mat/impls/sbaij/seq/sbaij.h" #include "src/inline/ilu.h" #include "include/petscis.h" #if !defined(PETSC_USE_COMPLEX) /* 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->i; PetscFunctionBegin; if (bs != 1) SETERRQ1(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); } #endif /* !defined(PETSC_USE_COMPLEX) */ extern PetscErrorCode MatSeqSBAIJSetNumericFactorization(Mat,PetscTruth); /* 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 A,IS perm,MatFactorInfo *info,Mat *B) { Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b; PetscErrorCode ierr; PetscInt *rip,i,mbs = a->mbs,*ai,*aj; PetscInt *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; PetscTruth 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 = PetscMalloc((mbs+1)*sizeof(PetscInt),&iu);CHKERRQ(ierr); umax = (PetscInt)(f*ai[mbs] + 1); umax += mbs + 1; ierr = PetscMalloc(umax*sizeof(PetscInt),&ju);CHKERRQ(ierr); iu[0] = mbs+1; juidx = mbs + 1; /* index for ju */ ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&jl);CHKERRQ(ierr); /* linked list for pivot row */ q = jl + mbs; /* linked list for col index */ for (i=0; i k){ qm = k; do { m = qm; qm = q[m]; } while(qm < vj); if (qm == vj) { SETERRQ(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 = PetscMalloc(umax*sizeof(PetscInt),&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,f,af);CHKERRQ(ierr); ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",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 = PetscFree(jl);CHKERRQ(ierr); /* put together the new matrix */ ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); /* ierr = PetscLogObjectParent(*B,iperm);CHKERRQ(ierr); */ b = (Mat_SeqSBAIJ*)(*B)->data; b->singlemalloc = PETSC_FALSE; b->free_a = PETSC_TRUE; b->free_ij = PETSC_TRUE; ierr = PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*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 = PetscMalloc((bs*mbs+bs)*sizeof(PetscScalar),&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(*B,(iu[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); b->maxnz = b->nz = iu[mbs]; (*B)->factor = MAT_FACTOR_CHOLESKY; (*B)->info.factor_mallocs = reallocs; (*B)->info.fill_ratio_given = f; if (ai[mbs] != 0) { (*B)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]); } else { (*B)->info.fill_ratio_needed = 0.0; } ierr = MatSeqSBAIJSetNumericFactorization(*B,perm_identity);CHKERRQ(ierr); PetscFunctionReturn(0); } /* Symbolic U^T*D*U factorization for SBAIJ format. */ #include "petscbt.h" #include "src/mat/utils/freespace.h" #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ" PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) { Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; Mat_SeqSBAIJ *b; PetscErrorCode ierr; PetscTruth perm_identity,missing; PetscReal fill = info->fill; PetscInt *rip,i,mbs=a->mbs,bs=A->rmap->bs,*ai,*aj,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=PETSC_NULL,current_space=PETSC_NULL; PetscBT lnkbt; PetscFunctionBegin; ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr); if (missing) SETERRQ1(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(A,perm,info,fact);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_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format"); /* There are bugs for reordeing. Needs further work. MatReordering for sbaij cannot be efficient. User should use aij formt! */ a->permute = PETSC_TRUE; ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr); ai = a->inew; aj = a->jnew; } ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); /* initialization */ ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&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 = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr); il = jl + mbs; cols = il + mbs; ui_ptr = (PetscInt**)(cols + mbs); 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; } #if defined(PETSC_USE_INFO) if (ai[mbs] != 0) { PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); } else { ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); } #endif ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); ierr = PetscFree(jl);CHKERRQ(ierr); /* destroy list of free space and other temporary array(s) */ ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&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,PETSC_NULL);CHKERRQ(ierr); b = (Mat_SeqSBAIJ*)(*fact)->data; b->singlemalloc = PETSC_FALSE; b->free_a = PETSC_TRUE; b->free_ij = PETSC_TRUE; ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&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 = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); ierr = PetscLogObjectMemory(*fact,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); b->maxnz = b->nz = ui[mbs]; (*fact)->factor = MAT_FACTOR_CHOLESKY; (*fact)->info.factor_mallocs = reallocs; (*fact)->info.fill_ratio_given = fill; if (ai[mbs] != 0) { (*fact)->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); } else { (*fact)->info.fill_ratio_needed = 0.0; } ierr = MatSeqSBAIJSetNumericFactorization(*fact,perm_identity);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N" PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat A,MatFactorInfo *info,Mat *B) { Mat C = *B; Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; IS perm = b->row; PetscErrorCode ierr; PetscInt *perm_ptr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; PetscInt *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; PetscInt bs=A->rmap->bs,bs2 = a->bs2,bslog = 0; MatScalar *ba = b->a,*aa,*ap,*dk,*uik; MatScalar *u,*diag,*rtmp,*rtmp_ptr; MatScalar *work; PetscInt *pivots; PetscFunctionBegin; /* initialization */ ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr); ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); jl = il + mbs; for (i=0; ipermute){ ai = a->i; aj = a->j; aa = a->a; } else { ai = a->inew; aj = a->jnew; ierr = PetscMalloc(bs2*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr); ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr); /* flops in while loop */ bslog = 2*bs*bs2; 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; C->ops->solvetranspose = MatSolve_SeqSBAIJ_N; C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_N; C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_N; 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 A,MatFactorInfo *info,Mat *B) { Mat C = *B; 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,bslog; MatScalar *ba = b->a,*aa,*ap,*dk,*uik; MatScalar *u,*diag,*rtmp,*rtmp_ptr; MatScalar *work; PetscInt *pivots; PetscFunctionBegin; ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr); ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); jl = il + mbs; for (i=0; ii; aj = a->j; aa = a->a; /* flops in while loop */ bslog = 2*bs*bs2; /* for each row k */ for (k = 0; kops->solve = MatSolve_SeqSBAIJ_N_NaturalOrdering; C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_NaturalOrdering; C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering; C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering; 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 A,MatFactorInfo *info,Mat *B) { Mat C = *B; Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; IS perm = b->row; PetscErrorCode ierr; PetscInt *perm_ptr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; PetscInt *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; MatScalar *ba = b->a,*aa,*ap,*dk,*uik; MatScalar *u,*diag,*rtmp,*rtmp_ptr; PetscReal shift = info->shiftinblocks; 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 = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr); ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); jl = il + mbs; for (i=0; ipermute){ ai = a->i; aj = a->j; aa = a->a; } else { ai = a->inew; aj = a->jnew; ierr = PetscMalloc(4*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&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; C->ops->solvetranspose = MatSolve_SeqSBAIJ_2; 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 A,MatFactorInfo *info,Mat *B) { Mat C = *B; 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,*uik; MatScalar *u,*diag,*rtmp,*rtmp_ptr; PetscReal shift = info->shiftinblocks; 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 = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr); ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); jl = il + mbs; for (i=0; ii; aj = a->j; aa = a->a; /* for each row k */ for (k = 0; kops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering; C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering; 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" PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat A,MatFactorInfo *info,Mat *B) { Mat C = *B; Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data; IS ip=b->row; PetscErrorCode ierr; PetscInt *rip,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j,*bcol; PetscInt *ai,*aj,*a2anew; PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; MatScalar *rtmp,*ba=b->a,*bval,*aa,dk,uikdi; PetscReal zeropivot,rs,shiftnz; PetscReal shiftpd; ChShift_Ctx sctx; PetscInt newshift; PetscFunctionBegin; /* initialization */ shiftnz = info->shiftnz; shiftpd = info->shiftpd; zeropivot = info->zeropivot; 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 = PetscMalloc(nz*sizeof(MatScalar),&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 */ nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); jl = il + mbs; rtmp = (MatScalar*)(jl + mbs); sctx.shift_amount = 0; sctx.nshift = 0; do { sctx.chshift = PETSC_FALSE; for (i=0; ipermute){ierr = PetscFree(aa);CHKERRQ(ierr);} ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); C->ops->solve = MatSolve_SeqSBAIJ_1; C->ops->solves = MatSolves_SeqSBAIJ_1; C->ops->solvetranspose = MatSolve_SeqSBAIJ_1; C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1; C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1; C->assembled = PETSC_TRUE; C->preallocated = PETSC_TRUE; ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); if (sctx.nshift){ if (shiftnz) { ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); } else if (shiftpd) { ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); } } PetscFunctionReturn(0); } /* Version for when blocks are 1 by 1 Using natural ordering */ #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering" PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B) { Mat C = *B; 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 zeropivot,rs,shiftnz; PetscReal shiftpd; ChShift_Ctx sctx; PetscInt newshift; PetscFunctionBegin; /* initialization */ shiftnz = info->shiftnz; shiftpd = info->shiftpd; zeropivot = info->zeropivot; /* 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 = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); jl = il + mbs; sctx.shift_amount = 0; sctx.nshift = 0; do { sctx.chshift = PETSC_FALSE; for (i=0; i 0){ bcol = bj + jmin; bval = ba + jmin; ierr = PetscLogFlops(2*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 = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr); if (newshift == 1) break; /* sctx.shift_amount is updated */ /* 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; C->ops->solves = MatSolves_SeqSBAIJ_1; C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering; C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering; C->assembled = PETSC_TRUE; C->preallocated = PETSC_TRUE; ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); if (sctx.nshift){ if (shiftnz) { ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); } else if (shiftpd) { ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); } } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ" PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,MatFactorInfo *info) { PetscErrorCode ierr; Mat C; PetscFunctionBegin; ierr = MatCholeskyFactorSymbolic(A,perm,info,&C);CHKERRQ(ierr); ierr = MatCholeskyFactorNumeric(A,info,&C);CHKERRQ(ierr); A->ops->solve = C->ops->solve; A->ops->solvetranspose = C->ops->solvetranspose; ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); PetscFunctionReturn(0); }