/*$Id: sbaijfact.c,v 1.61 2001/08/06 21:15:47 bsmith Exp $*/ /* Using Modified Sparse Row (MSR) storage. See page 85, "Iterative Methods ..." by Saad. */ /* Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP. */ #include "src/mat/impls/baij/seq/baij.h" #include "src/mat/impls/sbaij/seq/sbaij.h" #include "src/vec/vecimpl.h" #include "src/inline/ilu.h" #include "include/petscis.h" #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ" int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat A,IS perm,PetscReal f,Mat *B) { Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b; int *rip,ierr,i,mbs = a->mbs,*ai,*aj; int *jutmp,bs = a->bs,bs2=a->bs2; int m,nzi,realloc = 0; int *jl,*q,jumin,jmin,jmax,juptr,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd; PetscTruth perm_identity; PetscFunctionBegin; /* check whether perm is the identity mapping */ ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); if (!perm_identity) a->permute = PETSC_TRUE; ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); if (perm_identity){ /* without permutation */ ai = a->i; aj = a->j; } else { /* non-trivial permutation */ ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr); ai = a->inew; aj = a->jnew; } /* initialization */ /* Don't know how many column pointers are needed so estimate. Use Modified Sparse Row storage for u and ju, see Sasd pp.85 */ ierr = PetscMalloc((mbs+1)*sizeof(int),&iu);CHKERRQ(ierr); umax = (int)(f*ai[mbs] + 1); umax += mbs + 1; ierr = PetscMalloc(umax*sizeof(int),&ju);CHKERRQ(ierr); iu[0] = mbs+1; juptr = mbs; ierr = PetscMalloc(mbs*sizeof(int),&jl);CHKERRQ(ierr); ierr = PetscMalloc(mbs*sizeof(int),&q);CHKERRQ(ierr); for (i=0; i k){ qm = k; do { m = qm; qm = q[m]; } while(qm < vj); if (qm == vj) { printf(" error: duplicate entry in A\n"); break; } 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; /* printf(" iu[%d]=%d, umax=%d\n", k+1, iu[k+1],umax);*/ /* allocate more space to ju if needed */ if (iu[k+1] > umax) { printf("allocate more space, iu[%d]=%d > umax=%d\n",k+1, 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(int),&jutmp);CHKERRQ(ierr); ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(int));CHKERRQ(ierr); ierr = PetscFree(ju);CHKERRQ(ierr); ju = jutmp; realloc++; /* count how many times we realloc */ } /* save nonzero structure of k-th row in ju */ i=k; jumin = juptr + 1; juptr += nzk; for (j=jumin; jcomm,bs,bs*mbs,bs*mbs,0,PETSC_NULL,B);CHKERRQ(ierr); /* PetscLogObjectParent(*B,iperm); */ b = (Mat_SeqSBAIJ*)(*B)->data; ierr = PetscFree(b->imax);CHKERRQ(ierr); b->singlemalloc = PETSC_FALSE; /* the next line frees the default space generated by the Create() */ ierr = PetscFree(b->a);CHKERRQ(ierr); ierr = PetscFree(b->ilen);CHKERRQ(ierr); 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 MatCholeskyInfo */ 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 */ PetscLogObjectMemory(*B,(iu[mbs]-mbs)*(sizeof(int)+sizeof(MatScalar))); b->s_maxnz = b->s_nz = iu[mbs]; (*B)->factor = FACTOR_CHOLESKY; (*B)->info.factor_mallocs = realloc; (*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; } if (perm_identity){ switch (bs) { case 1: (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering; (*B)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=1\n"); break; case 2: (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; (*B)->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); break; case 3: (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; (*B)->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:sing special in-place natural ordering factor and solve BS=3\n"); break; case 4: (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; (*B)->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); break; case 5: (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; (*B)->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); break; case 6: (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; (*B)->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); break; case 7: (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; (*B)->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); break; default: (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering; (*B)->ops->solve = MatSolve_SeqSBAIJ_N_NaturalOrdering; PetscLogInfo(A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS>7\n"); break; } } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N" int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat A,Mat *B) { Mat C = *B; Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; IS perm = b->row; int *perm_ptr,ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; int *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; int bs=a->bs,bs2 = a->bs2; MatScalar *ba = b->a,*aa,*ap,*dk,*uik; MatScalar *u,*diag,*rtmp,*rtmp_ptr; MatScalar *work; int *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(int),&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(int),&a2anew);CHKERRQ(ierr); ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(int));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->factor = FACTOR_CHOLESKY; C->assembled = PETSC_TRUE; C->preallocated = PETSC_TRUE; PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */ PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering" int MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat A,Mat *B) { Mat C = *B; Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; int ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; int *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; int bs=a->bs,bs2 = a->bs2; MatScalar *ba = b->a,*aa,*ap,*dk,*uik; MatScalar *u,*diag,*rtmp,*rtmp_ptr; MatScalar *work; int *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(int),&il);CHKERRQ(ierr); jl = il + mbs; for (i=0; ii; aj = a->j; aa = a->a; /* for each row k */ for (k = 0; kfactor = FACTOR_CHOLESKY; C->assembled = PETSC_TRUE; C->preallocated = PETSC_TRUE; PetscLogFlops(1.3333*bs*bs2*b->mbs); /* 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" int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat A,Mat *B) { Mat C = *B; Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; IS perm = b->row; int *perm_ptr,ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; int *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; 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(int),&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(int),&a2anew);CHKERRQ(ierr); ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(int));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->factor = FACTOR_CHOLESKY; C->assembled = PETSC_TRUE; C->preallocated = PETSC_TRUE; PetscLogFlops(1.3333*8*b->mbs); /* 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" int MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat A,Mat *B) { Mat C = *B; Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; int ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; int *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; 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(int),&il);CHKERRQ(ierr); jl = il + mbs; for (i=0; ii; aj = a->j; aa = a->a; /* for each row k */ for (k = 0; kfactor = FACTOR_CHOLESKY; C->assembled = PETSC_TRUE; C->preallocated = PETSC_TRUE; PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */ PetscFunctionReturn(0); } /* Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP. Version for blocks are 1 by 1. */ #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1" int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat A,Mat *B) { Mat C = *B; Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; IS ip = b->row; int *rip,ierr,i,j,mbs = a->mbs,*bi = b->i,*bj = b->j; int *ai,*aj,*r; int k,jmin,jmax,*jl,*il,vj,nexti,ili; MatScalar *rtmp; MatScalar *ba = b->a,*aa,ak; MatScalar dk,uikdi; PetscFunctionBegin; 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; ierr = PetscMalloc(ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); ierr = PetscMemcpy(aa,a->a,ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); ierr = PetscMalloc(ai[mbs]*sizeof(int),&r);CHKERRQ(ierr); ierr= PetscMemcpy(r,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr); jmin = ai[0]; jmax = ai[mbs]; for (j=jmin; 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 = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr); jl = il + mbs; for (i=0; ipermute){ ierr = PetscFree(aa);CHKERRQ(ierr); } ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); C->factor = FACTOR_CHOLESKY; C->assembled = PETSC_TRUE; C->preallocated = PETSC_TRUE; PetscLogFlops(b->mbs); PetscFunctionReturn(0); } /* Version for when blocks are 1 by 1 Using natural ordering */ #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering" int MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat A,Mat *B) { Mat C = *B; Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; int ierr,i,j,mbs = a->mbs,*bi = b->i,*bj = b->j; int *ai,*aj; int k,jmin,jmax,*jl,*il,vj,nexti,ili; MatScalar *rtmp,*ba = b->a,*aa,dk,uikdi; 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(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr); jl = il + mbs; for (i=0; ii; aj = a->j; aa = a->a; /* for each row k */ for (k = 0; kfactor = FACTOR_CHOLESKY; C->assembled = PETSC_TRUE; C->preallocated = PETSC_TRUE; PetscLogFlops(b->mbs); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ" int MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,PetscReal f) { int ierr; Mat C; PetscFunctionBegin; ierr = MatCholeskyFactorSymbolic(A,perm,f,&C);CHKERRQ(ierr); ierr = MatCholeskyFactorNumeric(A,&C);CHKERRQ(ierr); ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); PetscFunctionReturn(0); }