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