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