/*$Id: sbaijfact10.c,v 1.2 2001/02/14 16:02:11 bsmith Exp bsmith $*/ #include "sbaij.h" #include "src/inline/ilu.h" /* Version for when blocks are 6 by 6 Using natural ordering */ #undef __FUNC__ #define __FUNC__ "MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering" int MatCholeskyFactorNumeric_SeqSBAIJ_6_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,*d,*w,*wp,u0,u1,u2,u3,u4,u5,u6,u7,u8,u9,u10,u11,u12; MatScalar u13,u14,u15,u16,u17,u18,u19,u20,u21,u22,u23,u24,u25,u26,u27; MatScalar u28,u29,u30,u31,u32,u33,u34,u35; PetscFunctionBegin; /* initialization */ ierr = PetscMalloc(36*mbs*sizeof(MatScalar),&w);CHKERRQ(ierr); ierr = PetscMemzero(w,36*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*216*b->mbs); /* from inverting diagonal blocks */ PetscFunctionReturn(0); }