/* Using Modified Sparse Row (MSR) storage. See page 85, "Iterative Methods ..." by Saad. */ /*$Id: sbaijfact.c,v 1.35 2000/10/30 19:58:55 hzhang Exp hzhang $*/ /* Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP. */ #include "sbaij.h" #include "src/mat/impls/baij/seq/baij.h" #include "src/vec/vecimpl.h" #include "src/inline/ilu.h" #include "include/petscis.h" #undef __FUNC__ #define __FUNC__ "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 ident; */ PetscFunctionBegin; PetscValidHeaderSpecific(perm,IS_COOKIE); if (A->M != A->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square"); /* check whether perm is the identity mapping */ /* ierr = ISView(perm, VIEWER_STDOUT_SELF);CHKERRA(ierr); ierr = ISIdentity(perm,&ident);CHKERRQ(ierr); printf("ident = %d\n", ident); */ ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); for (i=0; ipermute = PETSC_TRUE; /* printf("non-trivial perm\n"); */ break; } } if (!a->permute){ /* without permutation */ ai = a->i; aj = a->j; } else { /* non-trivial permutation */ ierr = MatReorderingSeqSBAIJ(A, perm);CHKERRA(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 */ iu = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(iu); umax = (int)(f*ai[mbs] + 1); umax += mbs + 1; ju = (int*)PetscMalloc(umax*sizeof(int));CHKPTRQ(ju); iu[0] = mbs+1; juptr = mbs; jl = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(jl); q = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(q); 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 */ jutmp = (int*)PetscMalloc(umax*sizeof(int));CHKPTRQ(jutmp); 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); /* PLogObjectParent(*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); b->a = (MatScalar*)PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*bs2);CHKPTRQ(b->a); b->j = ju; b->i = iu; b->diag = 0; b->ilen = 0; b->imax = 0; b->row = perm; ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); b->icol = perm; ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); b->solve_work = (Scalar*)PetscMalloc((bs*mbs+bs)*sizeof(Scalar));CHKPTRQ(b->solve_work); /* In b structure: Free imax, ilen, old a, old j. Allocate idnew, solve_work, new a, new j */ PLogObjectMemory(*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; } #ifdef TEMP for (k=0; ki[k+1] - b->i[k]; printf("\n b->i[%d]: %d, nzk: %d\n",k,b->i[k],nzk); jmin = b->i[k]; jmax = b->i[k+1]; for (j=jmin; jj[j]); } } #endif PetscFunctionReturn(0); } #undef __FUNC__ #define __FUNC__ "MatCholeskyFactorNumeric_SeqSBAIJ_N" int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat A,Mat *B) { Mat C = *B; Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; IS isrow = b->row,isicol = b->icol; int *r,*ic,ierr,i,j,n = a->mbs,*bi = b->i,*bj = b->j; int *ajtmpold,*ajtmp,nz,row,bslog,*ai=a->i,*aj=a->j,k,flg; int *diag_offset=b->diag,diag,bs=a->bs,bs2 = a->bs2,*v_pivots,*pj; MatScalar *ba = b->a,*aa = a->a,*pv,*v,*rtmp,*multiplier,*v_work,*pc,*w; PetscFunctionBegin; ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); rtmp = (MatScalar*)PetscMalloc(bs2*(n+1)*sizeof(MatScalar));CHKPTRQ(rtmp); ierr = PetscMemzero(rtmp,bs2*(n+1)*sizeof(MatScalar));CHKERRQ(ierr); /* generate work space needed by dense LU factorization */ v_work = (MatScalar*)PetscMalloc(bs*sizeof(int) + (bs+bs2)*sizeof(MatScalar));CHKPTRQ(v_work); multiplier = v_work + bs; v_pivots = (int*)(multiplier + bs2); /* flops in while loop */ bslog = 2*bs*bs2; for (i=0; ia */ pv = ba + bs2*bi[i]; pj = bj + bi[i]; nz = bi[i+1] - bi[i]; for (j=0; jfactor = FACTOR_LU; C->assembled = PETSC_TRUE; PLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */ PetscFunctionReturn(0); } /* Version for when blocks are 7 by 7 */ #undef __FUNC__ #define __FUNC__ "MatCholeskyFactorNumeric_SeqSBAIJ_7" int MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat A,Mat *B) { Mat C = *B; Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; IS isrow = b->row,isicol = b->icol; int *r,*ic,ierr,i,j,n = a->mbs,*bi = b->i,*bj = b->j; int *ajtmpold,*ajtmp,nz,row; int *diag_offset = b->diag,idx,*ai=a->i,*aj=a->j,*pj; MatScalar *pv,*v,*rtmp,*pc,*w,*x; MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4; MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,x16; MatScalar x17,x18,x19,x20,x21,x22,x23,x24,x25,p10,p11,p12,p13,p14; MatScalar p15,p16,p17,p18,p19,p20,p21,p22,p23,p24,p25,m10,m11,m12; MatScalar m13,m14,m15,m16,m17,m18,m19,m20,m21,m22,m23,m24,m25; MatScalar p26,p27,p28,p29,p30,p31,p32,p33,p34,p35,p36; MatScalar p37,p38,p39,p40,p41,p42,p43,p44,p45,p46,p47,p48,p49; MatScalar x26,x27,x28,x29,x30,x31,x32,x33,x34,x35,x36; MatScalar x37,x38,x39,x40,x41,x42,x43,x44,x45,x46,x47,x48,x49; MatScalar m26,m27,m28,m29,m30,m31,m32,m33,m34,m35,m36; MatScalar m37,m38,m39,m40,m41,m42,m43,m44,m45,m46,m47,m48,m49; MatScalar *ba = b->a,*aa = a->a; PetscFunctionBegin; ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); rtmp = (MatScalar*)PetscMalloc(49*(n+1)*sizeof(MatScalar));CHKPTRQ(rtmp); for (i=0; ia */ pv = ba + 49*bi[i]; pj = bj + bi[i]; nz = bi[i+1] - bi[i]; for (j=0; jfactor = FACTOR_LU; C->assembled = PETSC_TRUE; PLogFlops(1.3333*343*b->mbs); /* from inverting diagonal blocks */ PetscFunctionReturn(0); } /* Version for when blocks are 7 by 7 Using natural ordering */ #undef __FUNC__ #define __FUNC__ "MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering" int MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering(Mat A,Mat *B) { Mat C = *B; Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; int ierr,i,j,n = a->mbs,*bi = b->i,*bj = b->j; int *ajtmpold,*ajtmp,nz,row; int *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; MatScalar *pv,*v,*rtmp,*pc,*w,*x; MatScalar x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15; MatScalar x16,x17,x18,x19,x20,x21,x22,x23,x24,x25; MatScalar p1,p2,p3,p4,p5,p6,p7,p8,p9,p10,p11,p12,p13,p14,p15; MatScalar p16,p17,p18,p19,p20,p21,p22,p23,p24,p25; MatScalar m1,m2,m3,m4,m5,m6,m7,m8,m9,m10,m11,m12,m13,m14,m15; MatScalar m16,m17,m18,m19,m20,m21,m22,m23,m24,m25; MatScalar p26,p27,p28,p29,p30,p31,p32,p33,p34,p35,p36; MatScalar p37,p38,p39,p40,p41,p42,p43,p44,p45,p46,p47,p48,p49; MatScalar x26,x27,x28,x29,x30,x31,x32,x33,x34,x35,x36; MatScalar x37,x38,x39,x40,x41,x42,x43,x44,x45,x46,x47,x48,x49; MatScalar m26,m27,m28,m29,m30,m31,m32,m33,m34,m35,m36; MatScalar m37,m38,m39,m40,m41,m42,m43,m44,m45,m46,m47,m48,m49; MatScalar *ba = b->a,*aa = a->a; PetscFunctionBegin; rtmp = (MatScalar*)PetscMalloc(49*(n+1)*sizeof(MatScalar));CHKPTRQ(rtmp); for (i=0; ia */ pv = ba + 49*bi[i]; pj = bj + bi[i]; nz = bi[i+1] - bi[i]; for (j=0; jfactor = FACTOR_LU; C->assembled = PETSC_TRUE; PLogFlops(1.3333*343*b->mbs); /* from inverting diagonal blocks */ PetscFunctionReturn(0); } /* ------------------------------------------------------------*/ /* Version for when blocks are 6 by 6 */ #undef __FUNC__ #define __FUNC__ "MatCholeskyFactorNumeric_SeqSBAIJ_6" int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat A,Mat *B) { Mat C = *B; Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; IS isrow = b->row,isicol = b->icol; int *r,*ic,ierr,i,j,n = a->mbs,*bi = b->i,*bj = b->j; int *ajtmpold,*ajtmp,nz,row; int *diag_offset = b->diag,idx,*ai=a->i,*aj=a->j,*pj; MatScalar *pv,*v,*rtmp,*pc,*w,*x; MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4; MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,x16; MatScalar x17,x18,x19,x20,x21,x22,x23,x24,x25,p10,p11,p12,p13,p14; MatScalar p15,p16,p17,p18,p19,p20,p21,p22,p23,p24,p25,m10,m11,m12; MatScalar m13,m14,m15,m16,m17,m18,m19,m20,m21,m22,m23,m24,m25; MatScalar p26,p27,p28,p29,p30,p31,p32,p33,p34,p35,p36; MatScalar x26,x27,x28,x29,x30,x31,x32,x33,x34,x35,x36; MatScalar m26,m27,m28,m29,m30,m31,m32,m33,m34,m35,m36; MatScalar *ba = b->a,*aa = a->a; PetscFunctionBegin; ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); rtmp = (MatScalar*)PetscMalloc(36*(n+1)*sizeof(MatScalar));CHKPTRQ(rtmp); for (i=0; ia */ pv = ba + 36*bi[i]; pj = bj + bi[i]; nz = bi[i+1] - bi[i]; for (j=0; jfactor = FACTOR_LU; C->assembled = PETSC_TRUE; PLogFlops(1.3333*216*b->mbs); /* from inverting diagonal blocks */ PetscFunctionReturn(0); } /* 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_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; int ierr,i,j,n = a->mbs,*bi = b->i,*bj = b->j; int *ajtmpold,*ajtmp,nz,row; int *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; MatScalar *pv,*v,*rtmp,*pc,*w,*x; MatScalar x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15; MatScalar x16,x17,x18,x19,x20,x21,x22,x23,x24,x25; MatScalar p1,p2,p3,p4,p5,p6,p7,p8,p9,p10,p11,p12,p13,p14,p15; MatScalar p16,p17,p18,p19,p20,p21,p22,p23,p24,p25; MatScalar m1,m2,m3,m4,m5,m6,m7,m8,m9,m10,m11,m12,m13,m14,m15; MatScalar m16,m17,m18,m19,m20,m21,m22,m23,m24,m25; MatScalar p26,p27,p28,p29,p30,p31,p32,p33,p34,p35,p36; MatScalar x26,x27,x28,x29,x30,x31,x32,x33,x34,x35,x36; MatScalar m26,m27,m28,m29,m30,m31,m32,m33,m34,m35,m36; MatScalar *ba = b->a,*aa = a->a; PetscFunctionBegin; rtmp = (MatScalar*)PetscMalloc(36*(n+1)*sizeof(MatScalar));CHKPTRQ(rtmp); for (i=0; ia */ pv = ba + 36*bi[i]; pj = bj + bi[i]; nz = bi[i+1] - bi[i]; for (j=0; jfactor = FACTOR_LU; C->assembled = PETSC_TRUE; PLogFlops(1.3333*216*b->mbs); /* from inverting diagonal blocks */ PetscFunctionReturn(0); } /* Version for when blocks are 5 by 5 */ #undef __FUNC__ #define __FUNC__ "MatCholeskyFactorNumeric_SeqSBAIJ_5" int MatCholeskyFactorNumeric_SeqSBAIJ_5(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,*d,*rtmp,*rtmp_ptr; PetscFunctionBegin; /* initialization */ rtmp = (MatScalar*)PetscMalloc(25*mbs*sizeof(MatScalar));CHKPTRQ(rtmp); ierr = PetscMemzero(rtmp,25*mbs*sizeof(MatScalar));CHKERRQ(ierr); il = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(il); jl = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(jl); for (i=0; ipermute){ ai = a->i; aj = a->j; aa = a->a; } else { ai = a->inew; aj = a->jnew; aa = (MatScalar*)PetscMalloc(25*ai[mbs]*sizeof(MatScalar));CHKPTRQ(aa); ierr = PetscMemcpy(aa,a->a,25*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); a2anew = (int*)PetscMalloc(ai[mbs]*sizeof(int));CHKPTRQ(a2anew); 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*25; /* ptr to the beginning of j-th block of aa */ for (k=0; k<25; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */ for (k=0; k<5; k++){ /* j-th block of aa <- dk^T */ for (k1=0; k1<5; k1++) *ap++ = dk[k + 5*k1]; } } } } ierr = PetscFree(a2anew);CHKERRA(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; PLogFlops(1.3333*125*b->mbs); /* from inverting diagonal blocks */ PetscFunctionReturn(0); } #ifdef OLD /* Version for when blocks are 5 by 5 */ #undef __FUNC__ #define __FUNC__ "MatCholeskyFactorNumeric_SeqSBAIJ_5" int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat A,Mat *B) { Mat C = *B; Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; IS isrow = b->row,isicol = b->icol; int *r,*ic,ierr,i,j,n = a->mbs,*bi = b->i,*bj = b->j; int *ajtmpold,*ajtmp,nz,row; int *diag_offset = b->diag,idx,*ai=a->i,*aj=a->j,*pj; MatScalar *pv,*v,*rtmp,*pc,*w,*x; MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4; MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,x16; MatScalar x17,x18,x19,x20,x21,x22,x23,x24,x25,p10,p11,p12,p13,p14; MatScalar p15,p16,p17,p18,p19,p20,p21,p22,p23,p24,p25,m10,m11,m12; MatScalar m13,m14,m15,m16,m17,m18,m19,m20,m21,m22,m23,m24,m25; MatScalar *ba = b->a,*aa = a->a; PetscFunctionBegin; ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); rtmp = (MatScalar*)PetscMalloc(25*(n+1)*sizeof(MatScalar));CHKPTRQ(rtmp); for (i=0; ia */ pv = ba + 25*bi[i]; pj = bj + bi[i]; nz = bi[i+1] - bi[i]; for (j=0; jfactor = FACTOR_LU; C->assembled = PETSC_TRUE; PLogFlops(1.3333*125*b->mbs); /* from inverting diagonal blocks */ PetscFunctionReturn(0); } #endif /* Version for when blocks are 5 by 5 Using natural ordering */ #undef __FUNC__ #define __FUNC__ "MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering" int MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering(Mat A,Mat *B) { Mat C = *B; Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; int ierr,i,j,n = a->mbs,*bi = b->i,*bj = b->j; int *ajtmpold,*ajtmp,nz,row; int *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; MatScalar *pv,*v,*rtmp,*pc,*w,*x; MatScalar x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15; MatScalar x16,x17,x18,x19,x20,x21,x22,x23,x24,x25; MatScalar p1,p2,p3,p4,p5,p6,p7,p8,p9,p10,p11,p12,p13,p14,p15; MatScalar p16,p17,p18,p19,p20,p21,p22,p23,p24,p25; MatScalar m1,m2,m3,m4,m5,m6,m7,m8,m9,m10,m11,m12,m13,m14,m15; MatScalar m16,m17,m18,m19,m20,m21,m22,m23,m24,m25; MatScalar *ba = b->a,*aa = a->a; PetscFunctionBegin; rtmp = (MatScalar*)PetscMalloc(25*(n+1)*sizeof(MatScalar));CHKPTRQ(rtmp); for (i=0; ia */ pv = ba + 25*bi[i]; pj = bj + bi[i]; nz = bi[i+1] - bi[i]; for (j=0; jfactor = FACTOR_LU; C->assembled = PETSC_TRUE; PLogFlops(1.3333*125*b->mbs); /* from inverting diagonal blocks */ PetscFunctionReturn(0); } /* Version for when blocks are 4 by 4 Using natural ordering */ #undef __FUNC__ #define __FUNC__ "MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering" int MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering(Mat A,Mat *B) { Mat C = *B; Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; int ierr,i,j,n = a->mbs,*bi = b->i,*bj = b->j; int *ajtmpold,*ajtmp,nz,row; int *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; MatScalar *pv,*v,*rtmp,*pc,*w,*x; MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4; MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,x16; MatScalar p10,p11,p12,p13,p14,p15,p16,m10,m11,m12; MatScalar m13,m14,m15,m16; MatScalar *ba = b->a,*aa = a->a; PetscFunctionBegin; rtmp = (MatScalar*)PetscMalloc(16*(n+1)*sizeof(MatScalar));CHKPTRQ(rtmp); for (i=0; ia */ pv = ba + 16*bi[i]; pj = bj + bi[i]; nz = bi[i+1] - bi[i]; for (j=0; jfactor = FACTOR_LU; C->assembled = PETSC_TRUE; C->preallocated = PETSC_TRUE; PLogFlops(1.3333*64*b->mbs); /* from inverting diagonal blocks */ PetscFunctionReturn(0); } /* Version for when blocks are 4 by 4 */ #undef __FUNC__ #define __FUNC__ "MatCholeskyFactorNumeric_SeqSBAIJ_4" int MatCholeskyFactorNumeric_SeqSBAIJ_4(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 */ rtmp = (MatScalar*)PetscMalloc(16*mbs*sizeof(MatScalar));CHKPTRQ(rtmp); ierr = PetscMemzero(rtmp,16*mbs*sizeof(MatScalar));CHKERRQ(ierr); il = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(il); jl = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(jl); for (i=0; ipermute){ ai = a->i; aj = a->j; aa = a->a; } else { ai = a->inew; aj = a->jnew; aa = (MatScalar*)PetscMalloc(16*ai[mbs]*sizeof(MatScalar));CHKPTRQ(aa); ierr = PetscMemcpy(aa,a->a,16*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); a2anew = (int*)PetscMalloc(ai[mbs]*sizeof(int));CHKPTRQ(a2anew); 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*16; /* ptr to the beginning of j-th block of aa */ for (k=0; k<16; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */ for (k=0; k<4; k++){ /* j-th block of aa <- dk^T */ for (k1=0; k1<4; k1++) *ap++ = dk[k + 4*k1]; } } } } ierr = PetscFree(a2anew);CHKERRA(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; PLogFlops(1.3333*64*b->mbs); /* from inverting diagonal blocks */ PetscFunctionReturn(0); } /* Version for when blocks are 3 by 3 */ #undef __FUNC__ #define __FUNC__ "MatCholeskyFactorNumeric_SeqSBAIJ_3" int MatCholeskyFactorNumeric_SeqSBAIJ_3(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 */ rtmp = (MatScalar*)PetscMalloc(9*mbs*sizeof(MatScalar));CHKPTRQ(rtmp); ierr = PetscMemzero(rtmp,9*mbs*sizeof(MatScalar));CHKERRQ(ierr); il = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(il); jl = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(jl); for (i=0; ipermute){ ai = a->i; aj = a->j; aa = a->a; } else { ai = a->inew; aj = a->jnew; aa = (MatScalar*)PetscMalloc(9*ai[mbs]*sizeof(MatScalar));CHKPTRQ(aa); ierr = PetscMemcpy(aa,a->a,9*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); a2anew = (int*)PetscMalloc(ai[mbs]*sizeof(int));CHKPTRQ(a2anew); 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*9; /* ptr to the beginning of j-th block of aa */ for (k=0; k<9; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */ for (k=0; k<3; k++){ /* j-th block of aa <- dk^T */ for (k1=0; k1<3; k1++) *ap++ = dk[k + 3*k1]; } } } } ierr = PetscFree(a2anew);CHKERRA(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; PLogFlops(1.3333*27*b->mbs); /* from inverting diagonal blocks */ PetscFunctionReturn(0); } /* Version for when blocks are 3 by 3 Using natural ordering */ #undef __FUNC__ #define __FUNC__ "MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering" int MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering(Mat A,Mat *B) { Mat C = *B; Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; int ierr,i,j,n = a->mbs,*bi = b->i,*bj = b->j; int *ajtmpold,*ajtmp,nz,row; int *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; MatScalar *pv,*v,*rtmp,*pc,*w,*x; MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4; MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9; MatScalar *ba = b->a,*aa = a->a; PetscFunctionBegin; rtmp = (MatScalar*)PetscMalloc(9*(n+1)*sizeof(MatScalar));CHKPTRQ(rtmp); for (i=0; ia */ pv = ba + 9*bi[i]; pj = bj + bi[i]; nz = bi[i+1] - bi[i]; for (j=0; jfactor = FACTOR_CHOLESKY; C->assembled = PETSC_TRUE; C->preallocated = PETSC_TRUE; PLogFlops(1.3333*27*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 __FUNC__ #define __FUNC__ "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 */ rtmp = (MatScalar*)PetscMalloc(4*mbs*sizeof(MatScalar));CHKPTRQ(rtmp); ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr); il = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(il); jl = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(jl); for (i=0; ipermute){ ai = a->i; aj = a->j; aa = a->a; } else { ai = a->inew; aj = a->jnew; aa = (MatScalar*)PetscMalloc(4*ai[mbs]*sizeof(MatScalar));CHKPTRQ(aa); ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); a2anew = (int*)PetscMalloc(ai[mbs]*sizeof(int));CHKPTRQ(a2anew); 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);CHKERRA(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; PLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */ PetscFunctionReturn(0); } /* Version for when blocks are 2 by 2 Using natural ordering */ #undef __FUNC__ #define __FUNC__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering" int MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat A,Mat *B) { Mat C = *B; Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; int ierr,i,j,n = a->mbs,*bi = b->i,*bj = b->j; int *ajtmpold,*ajtmp,nz,row; int *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; MatScalar *pv,*v,*rtmp,*pc,*w,*x; MatScalar p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4; MatScalar *ba = b->a,*aa = a->a; PetscFunctionBegin; rtmp = (MatScalar*)PetscMalloc(4*(n+1)*sizeof(MatScalar));CHKPTRQ(rtmp); for (i=0; ia */ pv = ba + 4*bi[i]; pj = bj + bi[i]; nz = bi[i+1] - bi[i]; for (j=0; jfactor = FACTOR_LU; C->assembled = PETSC_TRUE; C->preallocated = PETSC_TRUE; PLogFlops(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 __FUNC__ #define __FUNC__ "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; MatScalar *rtmp; MatScalar *ba = b->a,*aa,ak; MatScalar dk,uikdi; int k,jmin,jmax,*jl,*il,vj,nexti,ili; 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; aa = (MatScalar*)PetscMalloc(ai[mbs]*sizeof(MatScalar));CHKPTRQ(aa); ierr = PetscMemcpy(aa,a->a,ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); r = (int*)PetscMalloc(ai[mbs]*sizeof(int));CHKPTRQ(r); 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 */ rtmp = (MatScalar*)PetscMalloc(mbs*sizeof(MatScalar));CHKPTRQ(rtmp); il = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(il); jl = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(jl); 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; PLogFlops(b->mbs); #ifdef TEMP printf("in factnum_1\n"); for (k=0; ki[k+1] - b->i[k]; printf("\n b->i[%d]: %d, nzk: %d, diag: %g\n",k,b->i[k],i,b->a[k]); jmin = b->i[k]; jmax = b->i[k+1]; for (j=jmin; jj[j],b->a[j]); } } #endif PetscFunctionReturn(0); } #undef __FUNC__ #define __FUNC__ "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); }