Actual source code: ex1.c
petsc-3.14.5 2021-03-03
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*/
9: /*
10: Include "petscksp.h" so that we can use KSP solvers. Note that this file
11: automatically includes:
12: petscsys.h - base PETSc routines petscvec.h - vectors
13: petscmat.h - matrices petscpc.h - preconditioners
14: petscis.h - index sets
15: petscviewer.h - viewers
17: Note: The corresponding parallel example is ex23.c
18: */
19: #include <petscksp.h>
21: int main(int argc,char **args)
22: {
23: Vec x, b, u; /* approx solution, RHS, exact solution */
24: Mat A; /* linear system matrix */
25: KSP ksp; /* linear solver context */
26: PC pc; /* preconditioner context */
27: PetscReal norm; /* norm of solution error */
29: PetscInt i,n = 10,col[3],its;
30: PetscMPIInt size;
31: PetscScalar value[3];
33: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
34: MPI_Comm_size(PETSC_COMM_WORLD,&size);
35: if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
36: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
39: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
40: Compute the matrix and right-hand-side vector that define
41: the linear system, Ax = b.
42: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
44: /*
45: Create vectors. Note that we form 1 vector from scratch and
46: then duplicate as needed.
47: */
48: VecCreate(PETSC_COMM_WORLD,&x);
49: PetscObjectSetName((PetscObject) x, "Solution");
50: VecSetSizes(x,PETSC_DECIDE,n);
51: VecSetFromOptions(x);
52: VecDuplicate(x,&b);
53: VecDuplicate(x,&u);
55: /*
56: Create matrix. When using MatCreate(), the matrix format can
57: be specified at runtime.
59: Performance tuning note: For problems of substantial size,
60: preallocation of matrix memory is crucial for attaining good
61: performance. See the matrix chapter of the users manual for details.
62: */
63: MatCreate(PETSC_COMM_WORLD,&A);
64: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
65: MatSetFromOptions(A);
66: MatSetUp(A);
68: /*
69: Assemble matrix
70: */
71: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
72: for (i=1; i<n-1; i++) {
73: col[0] = i-1; col[1] = i; col[2] = i+1;
74: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
75: }
76: i = n - 1; col[0] = n - 2; col[1] = n - 1;
77: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
78: i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
79: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
80: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
81: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
83: /*
84: Set exact solution; then compute right-hand-side vector.
85: */
86: VecSet(u,1.0);
87: MatMult(A,u,b);
89: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
90: Create the linear solver and set various options
91: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
92: KSPCreate(PETSC_COMM_WORLD,&ksp);
94: /*
95: Set operators. Here the matrix that defines the linear system
96: also serves as the matrix that defines the preconditioner.
97: */
98: KSPSetOperators(ksp,A,A);
100: /*
101: Set linear solver defaults for this problem (optional).
102: - By extracting the KSP and PC contexts from the KSP context,
103: we can then directly call any KSP and PC routines to set
104: various options.
105: - The following four statements are optional; all of these
106: parameters could alternatively be specified at runtime via
107: KSPSetFromOptions();
108: */
109: KSPGetPC(ksp,&pc);
110: PCSetType(pc,PCJACOBI);
111: KSPSetTolerances(ksp,1.e-5,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
113: /*
114: Set runtime options, e.g.,
115: -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
116: These options will override those specified above as long as
117: KSPSetFromOptions() is called _after_ any other customization
118: routines.
119: */
120: KSPSetFromOptions(ksp);
122: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123: Solve the linear system
124: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
125: KSPSolve(ksp,b,x);
127: /*
128: View solver info; we could instead use the option -ksp_view to
129: print this info to the screen at the conclusion of KSPSolve().
130: */
131: KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD);
133: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134: Check the solution and clean up
135: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136: VecAXPY(x,-1.0,u);
137: VecNorm(x,NORM_2,&norm);
138: KSPGetIterationNumber(ksp,&its);
139: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %D\n",(double)norm,its);
141: /*
142: Free work space. All PETSc objects should be destroyed when they
143: are no longer needed.
144: */
145: VecDestroy(&x); VecDestroy(&u);
146: VecDestroy(&b); MatDestroy(&A);
147: KSPDestroy(&ksp);
149: /*
150: Always call PetscFinalize() before exiting a program. This routine
151: - finalizes the PETSc libraries as well as MPI
152: - provides summary and diagnostic information if certain runtime
153: options are chosen (e.g., -log_view).
154: */
155: PetscFinalize();
156: return ierr;
157: }
159: /*TEST
161: test:
162: args: -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always
164: test:
165: suffix: 2
166: args: -pc_type sor -pc_sor_symmetric -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always
168: test:
169: suffix: 2_aijcusparse
170: requires: cuda
171: args: -pc_type sor -pc_sor_symmetric -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always -mat_type aijcusparse -vec_type cuda
173: test:
174: suffix: 3
175: args: -pc_type eisenstat -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always
177: test:
178: suffix: 3_aijcusparse
179: requires: cuda
180: args: -pc_type eisenstat -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always -mat_type aijcusparse -vec_type cuda
182: test:
183: suffix: aijcusparse
184: requires: cuda
185: args: -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always -mat_type aijcusparse -vec_type cuda
186: output_file: output/ex1_1_aijcusparse.out
188: TEST*/