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