Actual source code: ex2.c
petsc-3.10.5 2019-03-28
2: static char help[] = "Solves a linear system in parallel with KSP.\n\
3: Input parameters include:\n\
4: -random_exact_sol : use a random exact solution vector\n\
5: -view_exact_sol : write exact solution vector to stdout\n\
6: -m <mesh_x> : number of mesh points in x-direction\n\
7: -n <mesh_n> : number of mesh points in y-direction\n\n";
9: /*T
10: Concepts: KSP^basic parallel example;
11: Concepts: KSP^Laplacian, 2d
12: Concepts: Laplacian, 2d
13: Processors: n
14: T*/
18: /*
19: Include "petscksp.h" so that we can use KSP solvers. Note that this file
20: automatically includes:
21: petscsys.h - base PETSc routines petscvec.h - vectors
22: petscmat.h - matrices
23: petscis.h - index sets petscksp.h - Krylov subspace methods
24: petscviewer.h - viewers petscpc.h - preconditioners
25: */
26: #include <petscksp.h>
28: int main(int argc,char **args)
29: {
30: Vec x,b,u; /* approx solution, RHS, exact solution */
31: Mat A; /* linear system matrix */
32: KSP ksp; /* linear solver context */
33: PetscRandom rctx; /* random number generator context */
34: PetscReal norm; /* norm of solution error */
35: PetscInt i,j,Ii,J,Istart,Iend,m = 8,n = 7,its;
37: PetscBool flg = PETSC_FALSE;
38: PetscScalar v;
39: #if defined(PETSC_USE_LOG)
40: PetscLogStage stage;
41: #endif
43: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
44: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
45: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
46: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47: Compute the matrix and right-hand-side vector that define
48: the linear system, Ax = b.
49: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50: /*
51: Create parallel matrix, specifying only its global dimensions.
52: When using MatCreate(), the matrix format can be specified at
53: runtime. Also, the parallel partitioning of the matrix is
54: determined by PETSc at runtime.
56: Performance tuning note: For problems of substantial size,
57: preallocation of matrix memory is crucial for attaining good
58: performance. See the matrix chapter of the users manual for details.
59: */
60: MatCreate(PETSC_COMM_WORLD,&A);
61: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
62: MatSetFromOptions(A);
63: MatMPIAIJSetPreallocation(A,5,NULL,5,NULL);
64: MatSeqAIJSetPreallocation(A,5,NULL);
65: MatSeqSBAIJSetPreallocation(A,1,5,NULL);
66: MatMPISBAIJSetPreallocation(A,1,5,NULL,5,NULL);
67: MatMPISELLSetPreallocation(A,5,NULL,5,NULL);
68: MatSeqSELLSetPreallocation(A,5,NULL);
70: /*
71: Currently, all PETSc parallel matrix formats are partitioned by
72: contiguous chunks of rows across the processors. Determine which
73: rows of the matrix are locally owned.
74: */
75: MatGetOwnershipRange(A,&Istart,&Iend);
77: /*
78: Set matrix elements for the 2-D, five-point stencil in parallel.
79: - Each processor needs to insert only elements that it owns
80: locally (but any non-local elements will be sent to the
81: appropriate processor during matrix assembly).
82: - Always specify global rows and columns of matrix entries.
84: Note: this uses the less common natural ordering that orders first
85: all the unknowns for x = h then for x = 2h etc; Hence you see J = Ii +- n
86: instead of J = I +- m as you might expect. The more standard ordering
87: would first do all variables for y = h, then y = 2h etc.
89: */
90: PetscLogStageRegister("Assembly", &stage);
91: PetscLogStagePush(stage);
92: for (Ii=Istart; Ii<Iend; Ii++) {
93: v = -1.0; i = Ii/n; j = Ii - i*n;
94: if (i>0) {J = Ii - n; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);}
95: if (i<m-1) {J = Ii + n; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);}
96: if (j>0) {J = Ii - 1; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);}
97: if (j<n-1) {J = Ii + 1; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);}
98: v = 4.0; MatSetValues(A,1,&Ii,1,&Ii,&v,ADD_VALUES);
99: }
101: /*
102: Assemble matrix, using the 2-step process:
103: MatAssemblyBegin(), MatAssemblyEnd()
104: Computations can be done while messages are in transition
105: by placing code between these two statements.
106: */
107: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
108: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
109: PetscLogStagePop();
111: /* A is symmetric. Set symmetric flag to enable ICC/Cholesky preconditioner */
112: MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);
114: /*
115: Create parallel vectors.
116: - We form 1 vector from scratch and then duplicate as needed.
117: - When using VecCreate(), VecSetSizes and VecSetFromOptions()
118: in this example, we specify only the
119: vector's global dimension; the parallel partitioning is determined
120: at runtime.
121: - When solving a linear system, the vectors and matrices MUST
122: be partitioned accordingly. PETSc automatically generates
123: appropriately partitioned matrices and vectors when MatCreate()
124: and VecCreate() are used with the same communicator.
125: - The user can alternatively specify the local vector and matrix
126: dimensions when more sophisticated partitioning is needed
127: (replacing the PETSC_DECIDE argument in the VecSetSizes() statement
128: below).
129: */
130: VecCreate(PETSC_COMM_WORLD,&u);
131: VecSetSizes(u,PETSC_DECIDE,m*n);
132: VecSetFromOptions(u);
133: VecDuplicate(u,&b);
134: VecDuplicate(b,&x);
136: /*
137: Set exact solution; then compute right-hand-side vector.
138: By default we use an exact solution of a vector with all
139: elements of 1.0; Alternatively, using the runtime option
140: -random_sol forms a solution vector with random components.
141: */
142: PetscOptionsGetBool(NULL,NULL,"-random_exact_sol",&flg,NULL);
143: if (flg) {
144: PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
145: PetscRandomSetFromOptions(rctx);
146: VecSetRandom(u,rctx);
147: PetscRandomDestroy(&rctx);
148: } else {
149: VecSet(u,1.0);
150: }
151: MatMult(A,u,b);
153: /*
154: View the exact solution vector if desired
155: */
156: flg = PETSC_FALSE;
157: PetscOptionsGetBool(NULL,NULL,"-view_exact_sol",&flg,NULL);
158: if (flg) {VecView(u,PETSC_VIEWER_STDOUT_WORLD);}
160: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161: Create the linear solver and set various options
162: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
164: /*
165: Create linear solver context
166: */
167: KSPCreate(PETSC_COMM_WORLD,&ksp);
169: /*
170: Set operators. Here the matrix that defines the linear system
171: also serves as the preconditioning matrix.
172: */
173: KSPSetOperators(ksp,A,A);
175: /*
176: Set linear solver defaults for this problem (optional).
177: - By extracting the KSP and PC contexts from the KSP context,
178: we can then directly call any KSP and PC routines to set
179: various options.
180: - The following two statements are optional; all of these
181: parameters could alternatively be specified at runtime via
182: KSPSetFromOptions(). All of these defaults can be
183: overridden at runtime, as indicated below.
184: */
185: KSPSetTolerances(ksp,1.e-2/((m+1)*(n+1)),1.e-50,PETSC_DEFAULT,
186: PETSC_DEFAULT);
188: /*
189: Set runtime options, e.g.,
190: -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
191: These options will override those specified above as long as
192: KSPSetFromOptions() is called _after_ any other customization
193: routines.
194: */
195: KSPSetFromOptions(ksp);
197: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198: Solve the linear system
199: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
201: KSPSolve(ksp,b,x);
203: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204: Check solution and clean up
205: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
207: /*
208: Check the error
209: */
210: VecAXPY(x,-1.0,u);
211: VecNorm(x,NORM_2,&norm);
212: KSPGetIterationNumber(ksp,&its);
214: /*
215: Print convergence information. PetscPrintf() produces a single
216: print statement from all processes that share a communicator.
217: An alternative is PetscFPrintf(), which prints to a file.
218: */
219: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g iterations %D\n",(double)norm,its);
221: /*
222: Free work space. All PETSc objects should be destroyed when they
223: are no longer needed.
224: */
225: KSPDestroy(&ksp);
226: VecDestroy(&u); VecDestroy(&x);
227: VecDestroy(&b); MatDestroy(&A);
229: /*
230: Always call PetscFinalize() before exiting a program. This routine
231: - finalizes the PETSc libraries as well as MPI
232: - provides summary and diagnostic information if certain runtime
233: options are chosen (e.g., -log_view).
234: */
235: PetscFinalize();
236: return ierr;
237: }
240: /*TEST
242: build:
243: requires: !complex !single
245: test:
246: suffix: chebyest_1
247: args: -m 80 -n 80 -ksp_pc_side right -pc_type ksp -ksp_ksp_type chebyshev -ksp_ksp_max_it 5 -ksp_ksp_chebyshev_esteig 0.9,0,0,1.1 -ksp_monitor_short
249: test:
250: suffix: chebyest_2
251: args: -m 80 -n 80 -ksp_pc_side right -pc_type ksp -ksp_ksp_type chebyshev -ksp_ksp_max_it 5 -ksp_ksp_chebyshev_esteig 0.9,0,0,1.1 -ksp_esteig_ksp_type cg -ksp_monitor_short
253: test:
254: args: -ksp_monitor_short -m 5 -n 5 -ksp_gmres_cgs_refinement_type refine_always
256: test:
257: suffix: 2
258: nsize: 2
259: args: -ksp_monitor_short -m 5 -n 5 -ksp_gmres_cgs_refinement_type refine_always
261: test:
262: suffix: 3
263: args: -pc_type sor -pc_sor_symmetric -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always
265: test:
266: suffix: 4
267: args: -pc_type eisenstat -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always
269: test:
270: suffix: 5
271: nsize: 2
272: args: -ksp_monitor_short -m 5 -n 5 -mat_view draw -ksp_gmres_cgs_refinement_type refine_always -nox
273: output_file: output/ex2_2.out
275: test:
276: suffix: bjacobi
277: nsize: 4
278: args: -pc_type bjacobi -pc_bjacobi_blocks 1 -ksp_monitor_short -sub_pc_type jacobi -sub_ksp_type gmres
280: test:
281: suffix: bjacobi_2
282: nsize: 4
283: args: -pc_type bjacobi -pc_bjacobi_blocks 2 -ksp_monitor_short -sub_pc_type jacobi -sub_ksp_type gmres -ksp_view
285: test:
286: suffix: bjacobi_3
287: nsize: 4
288: args: -pc_type bjacobi -pc_bjacobi_blocks 4 -ksp_monitor_short -sub_pc_type jacobi -sub_ksp_type gmres
290: test:
291: suffix: fbcgs
292: args: -ksp_type fbcgs -pc_type ilu
294: test:
295: suffix: fbcgs_2
296: nsize: 3
297: args: -ksp_type fbcgsr -pc_type bjacobi
299: test:
300: suffix: groppcg
301: args: -ksp_monitor_short -ksp_type groppcg -m 9 -n 9
303: test:
304: suffix: mkl_pardiso_cholesky
305: requires: mkl_pardiso
306: args: -ksp_type preonly -pc_type cholesky -mat_type sbaij -pc_factor_mat_solver_type mkl_pardiso
308: test:
309: suffix: mkl_pardiso_lu
310: requires: mkl_pardiso
311: args: -ksp_type preonly -pc_type lu -pc_factor_mat_solver_type mkl_pardiso
313: test:
314: suffix: pipebcgs
315: args: -ksp_monitor_short -ksp_type pipebcgs -m 9 -n 9
317: test:
318: suffix: pipecg
319: args: -ksp_monitor_short -ksp_type pipecg -m 9 -n 9
321: test:
322: suffix: pipecgrr
323: args: -ksp_monitor_short -ksp_type pipecgrr -m 9 -n 9
325: test:
326: suffix: pipecr
327: args: -ksp_monitor_short -ksp_type pipecr -m 9 -n 9
329: test:
330: suffix: pipelcg
331: args: -ksp_monitor_short -ksp_type pipelcg -m 9 -n 9 -pc_type none -ksp_pipelcg_pipel 2 -ksp_pipelcg_lmax 2
332: filter: grep -v "sqrt breakdown in iteration"
334: test:
335: suffix: sell
336: args: -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always -m 9 -n 9 -mat_type sell
338: test:
339: requires: mumps
340: suffix: sell_mumps
341: args: -ksp_type preonly -m 9 -n 12 -mat_type sell -pc_type lu -pc_factor_mat_solver_type mumps -pc_factor_mat_ordering_type natural
343: test:
344: suffix: telescope
345: nsize: 4
346: args: -m 100 -n 100 -ksp_converged_reason -pc_type telescope -pc_telescope_reduction_factor 4 -telescope_pc_type bjacobi
348: test:
349: suffix: umfpack
350: requires: suitesparse
351: args: -ksp_type preonly -pc_type lu -pc_factor_mat_solver_type umfpack
353: test:
354: suffix: pc_symmetric
355: args: -m 10 -n 9 -ksp_converged_reason -ksp_type gmres -ksp_pc_side symmetric -pc_type cholesky
356: TEST*/