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 PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 21 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 22 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 23 24 /* Determine file from which we read the matrix A */ 25 PetscCall(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 PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd)); 30 PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 31 PetscCall(MatLoad(A,fd)); 32 PetscCall(VecCreate(PETSC_COMM_WORLD,&b)); 33 PetscCall(VecLoad(b,fd)); 34 PetscCall(PetscViewerDestroy(&fd)); 35 PetscCall(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 PetscCall(MatNorm(A,NORM_INFINITY,&Anorm)); 38 39 /* Create vectors */ 40 PetscCall(VecDuplicate(b,&x)); 41 PetscCall(VecDuplicate(x,&u)); /* save the true solution */ 42 43 /* Test Cholesky Factorization */ 44 PetscCall(MatGetOrdering(A,MATORDERINGNATURAL,&perm,&iperm)); 45 46 if (rank == 0) printf(" Clique Cholesky:\n"); 47 PetscCall(MatGetFactor(A,MATSOLVERCLIQUE,MAT_FACTOR_CHOLESKY,&F)); 48 49 info.fill = 5.0; 50 PetscCall(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 PetscCall(MatCholeskyFactorNumeric(F,A,&info)); 55 56 /* Test MatSolve() */ 57 if (testMatSolve && nfact == 2) { 58 PetscCall(MatSolve(F,b,x)); 59 60 /* Check the residual */ 61 PetscCall(MatMult(A,x,u)); 62 PetscCall(VecAXPY(u,-1.0,b)); 63 PetscCall(VecNorm(u,NORM_INFINITY,&norm)); 64 /* if (norm > tol) { */ 65 if (rank == 0) { 66 PetscCall(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 PetscCall(MatDestroy(&A)); 74 PetscCall(MatDestroy(&F)); 75 PetscCall(ISDestroy(&perm)); 76 PetscCall(ISDestroy(&iperm)); 77 PetscCall(VecDestroy(&x)); 78 PetscCall(VecDestroy(&b)); 79 PetscCall(VecDestroy(&u)); 80 PetscCall(PetscFinalize()); 81 return 0; 82 } 83