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