static char help[] = "Tests cholesky, icc factorization and solve on sequential aij, baij and sbaij matrices. \n"; #include "petscmat.h" #undef __FUNCT__ #define __FUNCT__ "main" int main(int argc,char **args) { Vec x,y,b; Mat A; /* linear system matrix */ Mat sA,sC; /* symmetric part of the matrices */ PetscInt n,mbs=16,bs=1,nz=3,prob=1,i,j,col[3],block, row,Ii,J,n1,lvl; PetscErrorCode ierr; PetscMPIInt size; PetscReal norm2,tol=1.e-10,err[10]; PetscScalar neg_one = -1.0,four=4.0,value[3]; IS perm,cperm; PetscRandom rdm; PetscInt reorder=0,displ=0; MatFactorInfo factinfo; PetscTruth equal; PetscTruth TestAIJ=PETSC_FALSE,TestBAIJ=PETSC_TRUE; PetscInt TestShift=0; PetscInitialize(&argc,&args,(char *)0,help); ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); if (size != 1) SETERRQ(1,"This is a uniprocessor example only!"); ierr = PetscOptionsGetInt(PETSC_NULL,"-bs",&bs,PETSC_NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(PETSC_NULL,"-mbs",&mbs,PETSC_NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(PETSC_NULL,"-reorder",&reorder,PETSC_NULL);CHKERRQ(ierr); ierr = PetscOptionsGetTruth(PETSC_NULL,"-testaij",&TestAIJ,PETSC_NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(PETSC_NULL,"-testShift",&TestShift,PETSC_NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(PETSC_NULL,"-displ",&displ,PETSC_NULL);CHKERRQ(ierr); n = mbs*bs; if (TestAIJ){ /* A is in aij format */ ierr = MatCreateSeqAIJ(PETSC_COMM_WORLD,n,n,nz,PETSC_NULL,&A);CHKERRQ(ierr); TestBAIJ = PETSC_FALSE; } else { /* A is in baij format */ ierr=MatCreateSeqBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,PETSC_NULL,&A);CHKERRQ(ierr); TestAIJ = PETSC_FALSE; } /* Assemble matrix */ if (bs == 1){ ierr = PetscOptionsGetInt(PETSC_NULL,"-test_problem",&prob,PETSC_NULL);CHKERRQ(ierr); if (prob == 1){ /* tridiagonal matrix */ value[0] = -1.0; value[1] = 2.0; value[2] = -1.0; for (i=1; i0) { J = Ii - n1; ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); } if (i0) { J = Ii - 1; ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); } if (j 1 */ for (block=0; block0){ierr = PetscPrintf(PETSC_COMM_SELF,"AIJ: \n");} i = 0; for (lvl=-1; lvl<10; lvl++){ if (lvl==-1) { /* Cholesky factor */ factinfo.fill = 5.0; ierr = MatCholeskyFactorSymbolic(A,perm,&factinfo,&sC);CHKERRQ(ierr); } else { /* incomplete Cholesky factor */ factinfo.fill = 5.0; factinfo.levels = lvl; ierr = MatICCFactorSymbolic(A,perm,&factinfo,&sC);CHKERRQ(ierr); } ierr = MatCholeskyFactorNumeric(A,&factinfo,&sC);CHKERRQ(ierr); ierr = MatMult(A,x,b);CHKERRQ(ierr); ierr = MatSolve(sC,b,y);CHKERRQ(ierr); ierr = MatDestroy(sC);CHKERRQ(ierr); /* Check the error */ ierr = VecAXPY(y,neg_one,x);CHKERRQ(ierr); ierr = VecNorm(y,NORM_2,&norm2);CHKERRQ(ierr); if (displ>0){ierr = PetscPrintf(PETSC_COMM_SELF," lvl: %d, error: %G\n", lvl,norm2);} err[i++] = norm2; } } /* Test baij matrix A */ if (TestBAIJ){ if (displ>0){ierr = PetscPrintf(PETSC_COMM_SELF,"BAIJ: \n");} i = 0; for (lvl=-1; lvl<10; lvl++){ if (lvl==-1) { /* Cholesky factor */ factinfo.fill = 5.0; ierr = MatCholeskyFactorSymbolic(A,perm,&factinfo,&sC);CHKERRQ(ierr); } else { /* incomplete Cholesky factor */ factinfo.fill = 5.0; factinfo.levels = lvl; ierr = MatICCFactorSymbolic(A,perm,&factinfo,&sC);CHKERRQ(ierr); } ierr = MatCholeskyFactorNumeric(A,&factinfo,&sC);CHKERRQ(ierr); ierr = MatMult(A,x,b);CHKERRQ(ierr); ierr = MatSolve(sC,b,y);CHKERRQ(ierr); ierr = MatDestroy(sC);CHKERRQ(ierr); /* Check the error */ ierr = VecAXPY(y,neg_one,x);CHKERRQ(ierr); ierr = VecNorm(y,NORM_2,&norm2);CHKERRQ(ierr); if (displ>0){ierr = PetscPrintf(PETSC_COMM_SELF," lvl: %d, error: %G\n", lvl,norm2);} err[i++] = norm2; } } /* Test sbaij matrix sA */ if (displ>0){ierr = PetscPrintf(PETSC_COMM_SELF,"SBAIJ: \n");} i = 0; for (lvl=-1; lvl<10; lvl++){ if (lvl==-1) { /* Cholesky factor */ factinfo.fill = 5.0; ierr = MatCholeskyFactorSymbolic(sA,perm,&factinfo,&sC);CHKERRQ(ierr); } else { /* incomplete Cholesky factor */ factinfo.fill = 5.0; factinfo.levels = lvl; ierr = MatICCFactorSymbolic(sA,perm,&factinfo,&sC);CHKERRQ(ierr); } ierr = MatCholeskyFactorNumeric(sA,&factinfo,&sC);CHKERRQ(ierr); if (lvl==0 && bs==1){ /* Test inplace ICC(0) for sbaij sA */ Mat B; ierr = MatDuplicate(sA,MAT_COPY_VALUES,&B);CHKERRQ(ierr); ierr = MatICCFactor(B,perm,&factinfo);CHKERRQ(ierr); ierr = MatEqual(sC,B,&equal);CHKERRQ(ierr); if (!equal){ SETERRQ(PETSC_ERR_USER,"in-place Cholesky factor != out-place Cholesky factor"); } ierr = MatDestroy(B);CHKERRQ(ierr); } ierr = MatMult(sA,x,b);CHKERRQ(ierr); ierr = MatSolve(sC,b,y);CHKERRQ(ierr); /* Test MatSolves() */ if (bs == 1) { Vecs xx,bb; ierr = VecsCreateSeq(PETSC_COMM_SELF,n,4,&xx);CHKERRQ(ierr); ierr = VecsDuplicate(xx,&bb);CHKERRQ(ierr); ierr = MatSolves(sC,bb,xx);CHKERRQ(ierr); ierr = VecsDestroy(xx);CHKERRQ(ierr); ierr = VecsDestroy(bb);CHKERRQ(ierr); } ierr = MatDestroy(sC);CHKERRQ(ierr); /* Check the error */ ierr = VecAXPY(y,neg_one,x);CHKERRQ(ierr); ierr = VecNorm(y,NORM_2,&norm2);CHKERRQ(ierr); if (displ>0){ierr = PetscPrintf(PETSC_COMM_SELF," lvl: %d, error: %G\n", lvl,norm2); } err[i] -= norm2; if (err[i] > tol) SETERRQ2(PETSC_ERR_USER," level: %d, err: %G\n", lvl,err[i]); } ierr = ISDestroy(perm);CHKERRQ(ierr); ierr = MatDestroy(A);CHKERRQ(ierr); ierr = MatDestroy(sA);CHKERRQ(ierr); ierr = VecDestroy(x);CHKERRQ(ierr); ierr = VecDestroy(y);CHKERRQ(ierr); ierr = VecDestroy(b);CHKERRQ(ierr); ierr = PetscRandomDestroy(rdm);CHKERRQ(ierr); ierr = PetscFinalize();CHKERRQ(ierr); return 0; }