xref: /petsc/src/ksp/ksp/tutorials/ex75.c (revision ff58222e5b8d6da69c593a20a2ec03b6a9c84ee1)
1 #include <petsc.h>
2 
3 static char help[] = "Solves a series of linear systems using KSPHPDDM.\n\n";
4 
5 int main(int argc, char **args) {
6   Vec x, b; /* computed solution and RHS */
7   Mat A;    /* linear system matrix */
8   KSP ksp;  /* linear solver context */
9 #if defined(PETSC_HAVE_HPDDM)
10   Mat U; /* deflation space */
11 #endif
12   PetscInt    i, j, nmat = 10;
13   PetscViewer viewer;
14   char        dir[PETSC_MAX_PATH_LEN], name[256];
15   PetscBool   flg, reset = PETSC_FALSE;
16 
17   PetscFunctionBeginUser;
18   PetscCall(PetscInitialize(&argc, &args, NULL, help));
19   PetscCall(PetscStrcpy(dir, "."));
20   PetscCall(PetscOptionsGetString(NULL, NULL, "-load_dir", dir, sizeof(dir), NULL));
21   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nmat", &nmat, NULL));
22   PetscCall(PetscOptionsGetBool(NULL, NULL, "-reset", &reset, NULL));
23   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
24   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
25   PetscCall(KSPSetOperators(ksp, A, A));
26   for (i = 0; i < nmat; i++) {
27     j = i + 400;
28     PetscCall(PetscSNPrintf(name, sizeof(name), "%s/A_%" PetscInt_FMT ".dat", dir, j));
29     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, name, FILE_MODE_READ, &viewer));
30     PetscCall(MatLoad(A, viewer));
31     PetscCall(PetscViewerDestroy(&viewer));
32     if (i == 0) PetscCall(MatCreateVecs(A, &x, &b));
33     PetscCall(PetscSNPrintf(name, sizeof(name), "%s/rhs_%" PetscInt_FMT ".dat", dir, j));
34     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, name, FILE_MODE_READ, &viewer));
35     PetscCall(VecLoad(b, viewer));
36     PetscCall(PetscViewerDestroy(&viewer));
37     PetscCall(KSPSetFromOptions(ksp));
38     PetscCall(KSPSolve(ksp, b, x));
39     PetscCall(PetscObjectTypeCompare((PetscObject)ksp, KSPHPDDM, &flg));
40 #if defined(PETSC_HAVE_HPDDM)
41     if (flg && reset) {
42       PetscCall(KSPHPDDMGetDeflationMat(ksp, &U));
43       PetscCall(KSPReset(ksp));
44       PetscCall(KSPSetOperators(ksp, A, A));
45       PetscCall(KSPSetFromOptions(ksp));
46       PetscCall(KSPSetUp(ksp));
47       if (U) {
48         PetscCall(KSPHPDDMSetDeflationMat(ksp, U));
49         PetscCall(MatDestroy(&U));
50       }
51     }
52 #endif
53   }
54   PetscCall(VecDestroy(&x));
55   PetscCall(VecDestroy(&b));
56   PetscCall(MatDestroy(&A));
57   PetscCall(KSPDestroy(&ksp));
58   PetscCall(PetscFinalize());
59   return 0;
60 }
61 
62 /*TEST
63 
64    test:
65       suffix: 1
66       nsize: 1
67       requires: hpddm datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)
68       args: -nmat 1 -pc_type none -ksp_converged_reason -ksp_type {{gmres hpddm}shared output} -ksp_max_it 1000 -ksp_gmres_restart 1000 -ksp_rtol 1e-10 -ksp_hpddm_type {{gmres bgmres}shared output} -options_left no -load_dir ${DATAFILESPATH}/matrices/hpddm/GCRODR
69 
70    test:
71       requires: hpddm datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)
72       suffix: 1_icc
73       nsize: 1
74       args: -nmat 1 -pc_type icc -ksp_converged_reason -ksp_type {{gmres hpddm}shared output} -ksp_max_it 1000 -ksp_gmres_restart 1000 -ksp_rtol 1e-10 -load_dir ${DATAFILESPATH}/matrices/hpddm/GCRODR
75 
76    testset:
77       requires: hpddm datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)
78       args: -nmat 3 -pc_type none -ksp_converged_reason -ksp_type hpddm -ksp_max_it 1000 -ksp_gmres_restart 40 -ksp_rtol 1e-10 -ksp_hpddm_type {{gcrodr bgcrodr}shared output} -ksp_hpddm_recycle 20 -load_dir ${DATAFILESPATH}/matrices/hpddm/GCRODR
79       test:
80         nsize: 1
81         suffix: 2_seq
82         output_file: output/ex75_2.out
83       test:
84         nsize: 2
85         suffix: 2_par
86         output_file: output/ex75_2.out
87 
88    testset:
89       requires: hpddm datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)
90       nsize: 1
91       args: -nmat 3 -pc_type icc -ksp_converged_reason -ksp_type hpddm -ksp_max_it 1000 -ksp_gmres_restart 40 -ksp_rtol 1e-10 -ksp_hpddm_type {{gcrodr bgcrodr}shared output} -ksp_hpddm_recycle 20 -reset {{false true}shared output} -load_dir ${DATAFILESPATH}/matrices/hpddm/GCRODR
92       test:
93         suffix: 2_icc
94         args:
95       test:
96         suffix: 2_icc_atol
97         args: -ksp_atol 1e-12
98 
99    test:
100       requires: hpddm datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) slepc defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
101       nsize: 2
102       suffix: symmetric
103       args: -nmat 3 -pc_type jacobi -ksp_converged_reason -ksp_type hpddm -ksp_max_it 1000 -ksp_gmres_restart 40 -ksp_atol 1e-11 -ksp_hpddm_type bgcrodr -ksp_hpddm_recycle 20 -reset {{false true}shared output} -load_dir ${DATAFILESPATH}/matrices/hpddm/GCRODR -ksp_hpddm_recycle_symmetric true
104 
105 TEST*/
106