1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests external Clique direct solvers. Simplified from ex130.c\n\ 3*c4762a1bSJed Brown Example: mpiexec -n <np> ./ex168 -f <matrix binary file> \n\n"; 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown #include <petscmat.h> 6*c4762a1bSJed Brown 7*c4762a1bSJed Brown int main(int argc,char **args) 8*c4762a1bSJed Brown { 9*c4762a1bSJed Brown Mat A,F; 10*c4762a1bSJed Brown Vec u,x,b; 11*c4762a1bSJed Brown PetscErrorCode ierr; 12*c4762a1bSJed Brown PetscMPIInt rank,size; 13*c4762a1bSJed Brown PetscInt m,n,nfact; 14*c4762a1bSJed Brown PetscReal norm,tol=1.e-12,Anorm; 15*c4762a1bSJed Brown IS perm,iperm; 16*c4762a1bSJed Brown MatFactorInfo info; 17*c4762a1bSJed Brown PetscBool flg,testMatSolve=PETSC_TRUE; 18*c4762a1bSJed Brown PetscViewer fd; /* viewer */ 19*c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; /* input file name */ 20*c4762a1bSJed Brown 21*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 22*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRQ(ierr); 23*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRQ(ierr); 24*c4762a1bSJed Brown 25*c4762a1bSJed Brown /* Determine file from which we read the matrix A */ 26*c4762a1bSJed Brown ierr = PetscOptionsGetString(NULL,NULL,"-f",file,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 27*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 28*c4762a1bSJed Brown 29*c4762a1bSJed Brown /* Load matrix A */ 30*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 32*c4762a1bSJed Brown ierr = MatLoad(A,fd);CHKERRQ(ierr); 33*c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&b);CHKERRQ(ierr); 34*c4762a1bSJed Brown ierr = VecLoad(b,fd);CHKERRQ(ierr); 35*c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 36*c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 37*c4762a1bSJed Brown if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%d, %d)", m, n); 38*c4762a1bSJed Brown ierr = MatNorm(A,NORM_INFINITY,&Anorm);CHKERRQ(ierr); 39*c4762a1bSJed Brown 40*c4762a1bSJed Brown /* Create vectors */ 41*c4762a1bSJed Brown ierr = VecDuplicate(b,&x);CHKERRQ(ierr); 42*c4762a1bSJed Brown ierr = VecDuplicate(x,&u);CHKERRQ(ierr); /* save the true solution */ 43*c4762a1bSJed Brown 44*c4762a1bSJed Brown /* Test Cholesky Factorization */ 45*c4762a1bSJed Brown ierr = MatGetOrdering(A,MATORDERINGNATURAL,&perm,&iperm);CHKERRQ(ierr); 46*c4762a1bSJed Brown 47*c4762a1bSJed Brown if (!rank) printf(" Clique Cholesky:\n"); 48*c4762a1bSJed Brown ierr = MatGetFactor(A,MATSOLVERCLIQUE,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr); 49*c4762a1bSJed Brown 50*c4762a1bSJed Brown info.fill = 5.0; 51*c4762a1bSJed Brown ierr = MatCholeskyFactorSymbolic(F,A,perm,&info);CHKERRQ(ierr); 52*c4762a1bSJed Brown 53*c4762a1bSJed Brown for (nfact = 0; nfact < 1; nfact++) { 54*c4762a1bSJed Brown if (!rank) printf(" %d-the Cholesky numfactorization \n",nfact); 55*c4762a1bSJed Brown ierr = MatCholeskyFactorNumeric(F,A,&info);CHKERRQ(ierr); 56*c4762a1bSJed Brown 57*c4762a1bSJed Brown /* Test MatSolve() */ 58*c4762a1bSJed Brown if (testMatSolve && nfact == 2) { 59*c4762a1bSJed Brown ierr = MatSolve(F,b,x);CHKERRQ(ierr); 60*c4762a1bSJed Brown 61*c4762a1bSJed Brown /* Check the residual */ 62*c4762a1bSJed Brown ierr = MatMult(A,x,u);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = VecAXPY(u,-1.0,b);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = VecNorm(u,NORM_INFINITY,&norm);CHKERRQ(ierr); 65*c4762a1bSJed Brown /* if (norm > tol) { */ 66*c4762a1bSJed Brown if (!rank) { 67*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve: rel residual %g/%g = %g, LU numfact %d\n",norm,Anorm,norm/Anorm,nfact);CHKERRQ(ierr); 68*c4762a1bSJed Brown } 69*c4762a1bSJed Brown /*} */ 70*c4762a1bSJed Brown } 71*c4762a1bSJed Brown } 72*c4762a1bSJed Brown 73*c4762a1bSJed Brown /* Free data structures */ 74*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = MatDestroy(&F);CHKERRQ(ierr); 76*c4762a1bSJed Brown ierr = ISDestroy(&perm);CHKERRQ(ierr); 77*c4762a1bSJed Brown ierr = ISDestroy(&iperm);CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = VecDestroy(&b);CHKERRQ(ierr); 80*c4762a1bSJed Brown ierr = VecDestroy(&u);CHKERRQ(ierr); 81*c4762a1bSJed Brown ierr = PetscFinalize(); 82*c4762a1bSJed Brown return ierr; 83*c4762a1bSJed Brown } 84