xref: /petsc/src/mat/tests/ex168.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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   PetscMPIInt    rank,size;
12c4762a1bSJed Brown   PetscInt       m,n,nfact;
13c4762a1bSJed Brown   PetscReal      norm,tol=1.e-12,Anorm;
14c4762a1bSJed Brown   IS             perm,iperm;
15c4762a1bSJed Brown   MatFactorInfo  info;
16c4762a1bSJed Brown   PetscBool      flg,testMatSolve=PETSC_TRUE;
17c4762a1bSJed Brown   PetscViewer    fd;              /* viewer */
18c4762a1bSJed Brown   char           file[PETSC_MAX_PATH_LEN]; /* input file name */
19c4762a1bSJed Brown 
20*327415f7SBarry Smith   PetscFunctionBeginUser;
219566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
24c4762a1bSJed Brown 
25c4762a1bSJed Brown   /* Determine file from which we read the matrix A */
269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg));
2728b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option");
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   /* Load matrix A */
309566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
319566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
329566063dSJacob Faibussowitsch   PetscCall(MatLoad(A,fd));
339566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&b));
349566063dSJacob Faibussowitsch   PetscCall(VecLoad(b,fd));
359566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&fd));
369566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A,&m,&n));
3708401ef6SPierre Jolivet   PetscCheck(m == n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%d, %d)", m, n);
389566063dSJacob Faibussowitsch   PetscCall(MatNorm(A,NORM_INFINITY,&Anorm));
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   /* Create vectors */
419566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(b,&x));
429566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&u)); /* save the true solution */
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /* Test Cholesky Factorization */
459566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(A,MATORDERINGNATURAL,&perm,&iperm));
46c4762a1bSJed Brown 
47dd400576SPatrick Sanan   if (rank == 0) printf(" Clique Cholesky:\n");
489566063dSJacob Faibussowitsch   PetscCall(MatGetFactor(A,MATSOLVERCLIQUE,MAT_FACTOR_CHOLESKY,&F));
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   info.fill = 5.0;
519566063dSJacob Faibussowitsch   PetscCall(MatCholeskyFactorSymbolic(F,A,perm,&info));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   for (nfact = 0; nfact < 1; nfact++) {
54dd400576SPatrick Sanan     if (rank == 0) printf(" %d-the Cholesky numfactorization \n",nfact);
559566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorNumeric(F,A,&info));
56c4762a1bSJed Brown 
57c4762a1bSJed Brown     /* Test MatSolve() */
58c4762a1bSJed Brown     if (testMatSolve && nfact == 2) {
599566063dSJacob Faibussowitsch       PetscCall(MatSolve(F,b,x));
60c4762a1bSJed Brown 
61c4762a1bSJed Brown       /* Check the residual */
629566063dSJacob Faibussowitsch       PetscCall(MatMult(A,x,u));
639566063dSJacob Faibussowitsch       PetscCall(VecAXPY(u,-1.0,b));
649566063dSJacob Faibussowitsch       PetscCall(VecNorm(u,NORM_INFINITY,&norm));
65c4762a1bSJed Brown       /* if (norm > tol) { */
66dd400576SPatrick Sanan       if (rank == 0) {
679566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF,"MatSolve: rel residual %g/%g = %g, LU numfact %d\n",norm,Anorm,norm/Anorm,nfact));
68c4762a1bSJed Brown       }
69c4762a1bSJed Brown       /*} */
70c4762a1bSJed Brown     }
71c4762a1bSJed Brown   }
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /* Free data structures */
749566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
759566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&F));
769566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
779566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iperm));
789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
809566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
819566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
82b122ec5aSJacob Faibussowitsch   return 0;
83c4762a1bSJed Brown }
84