Actual source code: ex1.c
petsc-3.10.5 2019-03-28
2: static char help[] = "Solves a tridiagonal linear system with KSP.\n\n";
4: /*T
5: Concepts: KSP^solving a system of linear equations
6: Processors: 1
7: T*/
11: /*
12: Include "petscksp.h" so that we can use KSP solvers. Note that this file
13: automatically includes:
14: petscsys.h - base PETSc routines petscvec.h - vectors
15: petscmat.h - matrices
16: petscis.h - index sets petscksp.h - Krylov subspace methods
17: petscviewer.h - viewers petscpc.h - preconditioners
19: Note: The corresponding parallel example is ex23.c
20: */
21: #include <petscksp.h>
23: int main(int argc,char **args)
24: {
25: Vec x, b, u; /* approx solution, RHS, exact solution */
26: Mat A; /* linear system matrix */
27: KSP ksp; /* linear solver context */
28: PC pc; /* preconditioner context */
29: PetscReal norm; /* norm of solution error */
31: PetscInt i,n = 10,col[3],its;
32: PetscMPIInt size;
33: PetscScalar one = 1.0,value[3];
34: PetscBool nonzeroguess = PETSC_FALSE,changepcside = PETSC_FALSE;
36: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
37: MPI_Comm_size(PETSC_COMM_WORLD,&size);
38: if (size != 1) SETERRQ(PETSC_COMM_WORLD,1,"This is a uniprocessor example only!");
39: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
40: PetscOptionsGetBool(NULL,NULL,"-nonzero_guess",&nonzeroguess,NULL);
43: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
44: Compute the matrix and right-hand-side vector that define
45: the linear system, Ax = b.
46: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
48: /*
49: Create vectors. Note that we form 1 vector from scratch and
50: then duplicate as needed.
51: */
52: VecCreate(PETSC_COMM_WORLD,&x);
53: PetscObjectSetName((PetscObject) x, "Solution");
54: VecSetSizes(x,PETSC_DECIDE,n);
55: VecSetFromOptions(x);
56: VecDuplicate(x,&b);
57: VecDuplicate(x,&u);
59: /*
60: Create matrix. When using MatCreate(), the matrix format can
61: be specified at runtime.
63: Performance tuning note: For problems of substantial size,
64: preallocation of matrix memory is crucial for attaining good
65: performance. See the matrix chapter of the users manual for details.
66: */
67: MatCreate(PETSC_COMM_WORLD,&A);
68: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
69: MatSetFromOptions(A);
70: MatSetUp(A);
72: /*
73: Assemble matrix
74: */
75: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
76: for (i=1; i<n-1; i++) {
77: col[0] = i-1; col[1] = i; col[2] = i+1;
78: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
79: }
80: i = n - 1; col[0] = n - 2; col[1] = n - 1;
81: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
82: i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
83: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
84: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
85: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
87: /*
88: Set exact solution; then compute right-hand-side vector.
89: */
90: VecSet(u,one);
91: MatMult(A,u,b);
93: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94: Create the linear solver and set various options
95: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96: /*
97: Create linear solver context
98: */
99: KSPCreate(PETSC_COMM_WORLD,&ksp);
101: /*
102: Set operators. Here the matrix that defines the linear system
103: also serves as the preconditioning matrix.
104: */
105: KSPSetOperators(ksp,A,A);
107: /*
108: Test if you can change the KSP side and type after they have been previously set
109: */
110: PetscOptionsGetBool(NULL,NULL,"-change_pc_side",&changepcside,NULL);
111: if (changepcside) {
112: KSPSetUp(ksp);
113: PetscOptionsInsertString(NULL,"-ksp_norm_type unpreconditioned -ksp_pc_side right");
114: KSPSetFromOptions(ksp);
115: KSPSetUp(ksp);
116: }
118: /*
119: Set linear solver defaults for this problem (optional).
120: - By extracting the KSP and PC contexts from the KSP context,
121: we can then directly call any KSP and PC routines to set
122: various options.
123: - The following four statements are optional; all of these
124: parameters could alternatively be specified at runtime via
125: KSPSetFromOptions();
126: */
127: KSPGetPC(ksp,&pc);
128: PCSetType(pc,PCJACOBI);
129: KSPSetTolerances(ksp,1.e-5,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
131: /*
132: Set runtime options, e.g.,
133: -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
134: These options will override those specified above as long as
135: KSPSetFromOptions() is called _after_ any other customization
136: routines.
137: */
138: KSPSetFromOptions(ksp);
140: if (nonzeroguess) {
141: PetscScalar p = .5;
142: VecSet(x,p);
143: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
144: }
146: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147: Solve the linear system
148: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149: /*
150: Solve linear system
151: */
152: KSPSolve(ksp,b,x);
154: /*
155: View solver info; we could instead use the option -ksp_view to
156: print this info to the screen at the conclusion of KSPSolve().
157: */
158: KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD);
160: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161: Check solution and clean up
162: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163: /*
164: Check the error
165: */
166: VecAXPY(x,-1.0,u);
167: VecNorm(x,NORM_2,&norm);
168: KSPGetIterationNumber(ksp,&its);
169: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %D\n",(double)norm,its);
171: /*
172: Free work space. All PETSc objects should be destroyed when they
173: are no longer needed.
174: */
175: VecDestroy(&x); VecDestroy(&u);
176: VecDestroy(&b); MatDestroy(&A);
177: KSPDestroy(&ksp);
179: /*
180: Always call PetscFinalize() before exiting a program. This routine
181: - finalizes the PETSc libraries as well as MPI
182: - provides summary and diagnostic information if certain runtime
183: options are chosen (e.g., -log_view).
184: */
185: PetscFinalize();
186: return ierr;
187: }
190: /*TEST
192: test:
193: args: -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always
195: test:
196: suffix: 2
197: args: -pc_type sor -pc_sor_symmetric -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always
199: test:
200: suffix: 2_aijcusparse
201: requires: veccuda
202: args: -pc_type sor -pc_sor_symmetric -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always -mat_type aijcusparse -vec_type cuda
204: test:
205: suffix: 3
206: args: -pc_type eisenstat -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always
208: test:
209: suffix: 3_aijcusparse
210: requires: veccuda
211: args: -pc_type eisenstat -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always -mat_type aijcusparse -vec_type cuda
213: test:
214: suffix: aijcusparse
215: requires: veccuda
216: args: -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always -mat_type aijcusparse -vec_type cuda
217: output_file: output/ex1_1_aijcusparse.out
219: test:
220: suffix: changepcside
221: args: -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always -change_pc_side
223: TEST*/