xref: /petsc/src/ksp/ksp/tests/ex1.c (revision 8564c97f3fe574910f676ffb31bf76fa44548916)
1 static char help[] = "Tests solving linear system on 0 by 0 matrix, and KSPLSQR convergence test handling.\n\n";
2 
3 #include <petscksp.h>
4 
5 static PetscErrorCode GetConvergenceTestName(PetscErrorCode (*converged)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), char name[], size_t n)
6 {
7   PetscFunctionBegin;
8   if (converged == KSPConvergedDefault) {
9     PetscCall(PetscStrncpy(name, "default", n));
10   } else if (converged == KSPConvergedSkip) {
11     PetscCall(PetscStrncpy(name, "skip", n));
12   } else if (converged == KSPLSQRConvergedDefault) {
13     PetscCall(PetscStrncpy(name, "lsqr", n));
14   } else {
15     PetscCall(PetscStrncpy(name, "other", n));
16   }
17   PetscFunctionReturn(PETSC_SUCCESS);
18 }
19 
20 int main(int argc, char **args)
21 {
22   Mat       C;
23   PetscInt  N = 0;
24   Vec       u, b, x;
25   KSP       ksp;
26   PetscReal norm;
27   PetscBool flg = PETSC_FALSE;
28 
29   PetscFunctionBeginUser;
30   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
31 
32   /* create stiffness matrix */
33   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
34   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, N, N));
35   PetscCall(MatSetFromOptions(C));
36   PetscCall(MatSetUp(C));
37   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
38   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
39 
40   /* create right-hand side and solution */
41   PetscCall(VecCreate(PETSC_COMM_WORLD, &u));
42   PetscCall(VecSetSizes(u, PETSC_DECIDE, N));
43   PetscCall(VecSetFromOptions(u));
44   PetscCall(VecDuplicate(u, &b));
45   PetscCall(VecDuplicate(u, &x));
46   PetscCall(VecSet(u, 0.0));
47   PetscCall(VecSet(b, 0.0));
48 
49   PetscCall(VecAssemblyBegin(b));
50   PetscCall(VecAssemblyEnd(b));
51 
52   /* solve linear system */
53   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
54   PetscCall(KSPSetOperators(ksp, C, C));
55   PetscCall(KSPSetFromOptions(ksp));
56   PetscCall(KSPSolve(ksp, b, u));
57 
58   /* test proper handling of convergence test by KSPLSQR */
59   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_lsqr", &flg, NULL));
60   if (flg) {
61     char     *type;
62     char      convtestname[16];
63     PetscBool islsqr;
64     PetscErrorCode (*converged)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
65     PetscErrorCode (*converged1)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
66     PetscErrorCode (*destroy)(void *), (*destroy1)(void *);
67     void *ctx, *ctx1;
68 
69     {
70       const char *typeP;
71       PetscCall(KSPGetType(ksp, &typeP));
72       PetscCall(PetscStrallocpy(typeP, &type));
73     }
74     PetscCall(PetscStrcmp(type, KSPLSQR, &islsqr));
75     PetscCall(KSPGetConvergenceTest(ksp, &converged, &ctx, &destroy));
76     PetscCall(GetConvergenceTestName(converged, convtestname, 16));
77     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "convergence test: %s\n", convtestname));
78     PetscCall(KSPSetType(ksp, KSPLSQR));
79     PetscCall(KSPGetConvergenceTest(ksp, &converged1, &ctx1, &destroy1));
80     PetscCheck(converged1 == KSPLSQRConvergedDefault, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "convergence test should be KSPLSQRConvergedDefault");
81     PetscCheck(destroy1 == KSPConvergedDefaultDestroy, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "convergence test destroy function should be KSPConvergedDefaultDestroy");
82     if (islsqr) {
83       PetscCheck(converged1 == converged, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "convergence test should be kept");
84       PetscCheck(destroy1 == destroy, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "convergence test destroy function should be kept");
85       PetscCheck(ctx1 == ctx, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "convergence test context should be kept");
86     }
87     PetscCall(GetConvergenceTestName(converged1, convtestname, 16));
88     PetscCall(KSPViewFromOptions(ksp, NULL, "-ksp1_view"));
89     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "convergence test: %s\n", convtestname));
90     PetscCall(KSPSetType(ksp, type));
91     PetscCall(KSPGetConvergenceTest(ksp, &converged1, &ctx1, &destroy1));
92     PetscCheck(converged1 == converged, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "convergence test not reverted properly");
93     PetscCheck(destroy1 == destroy, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "convergence test destroy function not reverted properly");
94     PetscCheck(ctx1 == ctx, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "convergence test context not reverted properly");
95     PetscCall(GetConvergenceTestName(converged1, convtestname, 16));
96     PetscCall(KSPViewFromOptions(ksp, NULL, "-ksp2_view"));
97     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "convergence test: %s\n", convtestname));
98     PetscCall(PetscFree(type));
99   }
100 
101   PetscCall(MatMult(C, u, x));
102   PetscCall(VecAXPY(x, -1.0, b));
103   PetscCall(VecNorm(x, NORM_2, &norm));
104 
105   PetscCall(KSPDestroy(&ksp));
106   PetscCall(VecDestroy(&u));
107   PetscCall(VecDestroy(&x));
108   PetscCall(VecDestroy(&b));
109   PetscCall(MatDestroy(&C));
110   PetscCall(PetscFinalize());
111   return 0;
112 }
113 
114 /*TEST
115 
116     test:
117       args: -pc_type jacobi -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always
118 
119     test:
120       suffix: 2
121       nsize: 2
122       args: -pc_type jacobi -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always
123 
124     test:
125       suffix: 3
126       args: -pc_type sor -pc_sor_symmetric -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always
127 
128     test:
129       suffix: 5
130       args: -pc_type eisenstat -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always
131 
132     testset:
133       args: -test_lsqr -ksp{,1,2}_view -pc_type jacobi
134       filter: grep -E "(^  type:|preconditioning|norm type|convergence test:)"
135       test:
136         suffix: lsqr_0
137         args: -ksp_convergence_test {{default skip}separate output}
138       test:
139         suffix: lsqr_1
140         args: -ksp_type cg -ksp_convergence_test {{default skip}separate output}
141       test:
142         suffix: lsqr_2
143         args: -ksp_type lsqr
144 
145 TEST*/
146