Actual source code: ex1.c
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);
38: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
39: Compute the matrix and right-hand-side vector that define
40: the linear system, Ax = b.
41: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
43: /*
44: Create vectors. Note that we form 1 vector from scratch and
45: then duplicate as needed.
46: */
47: VecCreate(PETSC_COMM_WORLD,&x);
48: PetscObjectSetName((PetscObject) x, "Solution");
49: VecSetSizes(x,PETSC_DECIDE,n);
50: VecSetFromOptions(x);
51: VecDuplicate(x,&b);
52: VecDuplicate(x,&u);
54: /*
55: Create matrix. When using MatCreate(), the matrix format can
56: be specified at runtime.
58: Performance tuning note: For problems of substantial size,
59: preallocation of matrix memory is crucial for attaining good
60: performance. See the matrix chapter of the users manual for details.
61: */
62: MatCreate(PETSC_COMM_WORLD,&A);
63: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
64: MatSetFromOptions(A);
65: MatSetUp(A);
67: /*
68: Assemble matrix
69: */
70: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
71: for (i=1; i<n-1; i++) {
72: col[0] = i-1; col[1] = i; col[2] = i+1;
73: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
74: }
75: i = n - 1; col[0] = n - 2; col[1] = n - 1;
76: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
77: i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
78: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
79: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
80: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
82: /*
83: Set exact solution; then compute right-hand-side vector.
84: */
85: VecSet(u,1.0);
86: MatMult(A,u,b);
88: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
89: Create the linear solver and set various options
90: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
91: KSPCreate(PETSC_COMM_WORLD,&ksp);
93: /*
94: Set operators. Here the matrix that defines the linear system
95: also serves as the matrix that defines the preconditioner.
96: */
97: KSPSetOperators(ksp,A,A);
99: /*
100: Set linear solver defaults for this problem (optional).
101: - By extracting the KSP and PC contexts from the KSP context,
102: we can then directly call any KSP and PC routines to set
103: various options.
104: - The following four statements are optional; all of these
105: parameters could alternatively be specified at runtime via
106: KSPSetFromOptions();
107: */
108: KSPGetPC(ksp,&pc);
109: PCSetType(pc,PCJACOBI);
110: KSPSetTolerances(ksp,1.e-5,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
112: /*
113: Set runtime options, e.g.,
114: -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
115: These options will override those specified above as long as
116: KSPSetFromOptions() is called _after_ any other customization
117: routines.
118: */
119: KSPSetFromOptions(ksp);
121: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
122: Solve the linear system
123: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
124: KSPSolve(ksp,b,x);
126: /*
127: View solver info; we could instead use the option -ksp_view to
128: print this info to the screen at the conclusion of KSPSolve().
129: */
130: KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD);
132: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133: Check the solution and clean up
134: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
135: VecAXPY(x,-1.0,u);
136: VecNorm(x,NORM_2,&norm);
137: KSPGetIterationNumber(ksp,&its);
138: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %D\n",(double)norm,its);
140: /*
141: Free work space. All PETSc objects should be destroyed when they
142: are no longer needed.
143: */
144: VecDestroy(&x); VecDestroy(&u);
145: VecDestroy(&b); MatDestroy(&A);
146: KSPDestroy(&ksp);
148: /*
149: Always call PetscFinalize() before exiting a program. This routine
150: - finalizes the PETSc libraries as well as MPI
151: - provides summary and diagnostic information if certain runtime
152: options are chosen (e.g., -log_view).
153: */
154: PetscFinalize();
155: return ierr;
156: }
158: /*TEST
160: test:
161: args: -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always
163: test:
164: suffix: 2
165: args: -pc_type sor -pc_sor_symmetric -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always
167: test:
168: suffix: 2_aijcusparse
169: requires: cuda
170: args: -pc_type sor -pc_sor_symmetric -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always -mat_type aijcusparse -vec_type cuda
172: test:
173: suffix: 3
174: args: -pc_type eisenstat -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always
176: test:
177: suffix: 3_aijcusparse
178: requires: cuda
179: args: -pc_type eisenstat -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always -mat_type aijcusparse -vec_type cuda
181: test:
182: suffix: aijcusparse
183: requires: cuda
184: args: -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always -mat_type aijcusparse -vec_type cuda
185: output_file: output/ex1_1_aijcusparse.out
187: TEST*/