Actual source code: ex5.c
1: /*$Id: ex5.c,v 1.93 2001/09/11 16:33:29 bsmith Exp $*/
3: static char help[] = "Solves two linear systems in parallel with KSP. The code\n\
4: illustrates repeated solution of linear systems with the same preconditioner\n\
5: method but different matrices (having the same nonzero structure). The code\n\
6: also uses multiple profiling stages. Input arguments are\n\
7: -m <size> : problem size\n\
8: -mat_nonsym : use nonsymmetric matrix (default is symmetric)\n\n";
10: /*T
11: Concepts: KSP^repeatedly solving linear systems;
12: Concepts: PetscLog^profiling multiple stages of code;
13: Processors: n
14: T*/
16: /*
17: Include "petscksp.h" so that we can use KSP solvers. Note that this file
18: automatically includes:
19: petsc.h - base PETSc routines petscvec.h - vectors
20: petscsys.h - system routines petscmat.h - matrices
21: petscis.h - index sets petscksp.h - Krylov subspace methods
22: petscviewer.h - viewers petscpc.h - preconditioners
23: */
24: #include petscksp.h
28: int main(int argc,char **args)
29: {
30: KSP ksp; /* linear solver context */
31: Mat C; /* matrix */
32: Vec x,u,b; /* approx solution, RHS, exact solution */
33: PetscReal norm; /* norm of solution error */
34: PetscScalar v,none = -1.0;
35: int I,J,ldim,ierr,low,high,iglobal,Istart,Iend;
36: int i,j,m = 3,n = 2,rank,size,its;
37: PetscTruth mat_nonsymmetric;
38: int stages[2];
40: PetscInitialize(&argc,&args,(char *)0,help);
41: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
42: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
43: MPI_Comm_size(PETSC_COMM_WORLD,&size);
44: n = 2*size;
46: /*
47: Set flag if we are doing a nonsymmetric problem; the default is symmetric.
48: */
49: PetscOptionsHasName(PETSC_NULL,"-mat_nonsym",&mat_nonsymmetric);
51: /*
52: Register two stages for separate profiling of the two linear solves.
53: Use the runtime option -log_summary for a printout of performance
54: statistics at the program's conlusion.
55: */
56: PetscLogStageRegister(&stages[0],"Original Solve");
57: PetscLogStageRegister(&stages[1],"Second Solve");
59: /* -------------- Stage 0: Solve Original System ---------------------- */
60: /*
61: Indicate to PETSc profiling that we're beginning the first stage
62: */
63: PetscLogStagePush(stages[0]);
65: /*
66: Create parallel matrix, specifying only its global dimensions.
67: When using MatCreate(), the matrix format can be specified at
68: runtime. Also, the parallel partitioning of the matrix is
69: determined by PETSc at runtime.
70: */
71: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,&C);
72: MatSetFromOptions(C);
74: /*
75: Currently, all PETSc parallel matrix formats are partitioned by
76: contiguous chunks of rows across the processors. Determine which
77: rows of the matrix are locally owned.
78: */
79: MatGetOwnershipRange(C,&Istart,&Iend);
81: /*
82: Set matrix entries matrix in parallel.
83: - Each processor needs to insert only elements that it owns
84: locally (but any non-local elements will be sent to the
85: appropriate processor during matrix assembly).
86: - Always specify global row and columns of matrix entries.
87: */
88: for (I=Istart; I<Iend; I++) {
89: v = -1.0; i = I/n; j = I - i*n;
90: if (i>0) {J = I - n; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
91: if (i<m-1) {J = I + n; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
92: if (j>0) {J = I - 1; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
93: if (j<n-1) {J = I + 1; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
94: v = 4.0; MatSetValues(C,1,&I,1,&I,&v,ADD_VALUES);
95: }
97: /*
98: Make the matrix nonsymmetric if desired
99: */
100: if (mat_nonsymmetric) {
101: for (I=Istart; I<Iend; I++) {
102: v = -1.5; i = I/n;
103: if (i>1) {J = I-n-1; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
104: }
105: } else {
106: MatSetOption(C,MAT_SYMMETRIC);
107: MatSetOption(C,MAT_SYMMETRY_ETERNAL);
108: }
110: /*
111: Assemble matrix, using the 2-step process:
112: MatAssemblyBegin(), MatAssemblyEnd()
113: Computations can be done while messages are in transition
114: by placing code between these two statements.
115: */
116: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
117: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
119: /*
120: Create parallel vectors.
121: - When using VecSetSizes(), we specify only the vector's global
122: dimension; the parallel partitioning is determined at runtime.
123: - Note: We form 1 vector from scratch and then duplicate as needed.
124: */
125: VecCreate(PETSC_COMM_WORLD,&u);
126: VecSetSizes(u,PETSC_DECIDE,m*n);
127: VecSetFromOptions(u);
128: VecDuplicate(u,&b);
129: VecDuplicate(b,&x);
131: /*
132: Currently, all parallel PETSc vectors are partitioned by
133: contiguous chunks across the processors. Determine which
134: range of entries are locally owned.
135: */
136: VecGetOwnershipRange(x,&low,&high);
138: /*
139: Set elements within the exact solution vector in parallel.
140: - Each processor needs to insert only elements that it owns
141: locally (but any non-local entries will be sent to the
142: appropriate processor during vector assembly).
143: - Always specify global locations of vector entries.
144: */
145: VecGetLocalSize(x,&ldim);
146: for (i=0; i<ldim; i++) {
147: iglobal = i + low;
148: v = (PetscScalar)(i + 100*rank);
149: VecSetValues(u,1,&iglobal,&v,INSERT_VALUES);
150: }
152: /*
153: Assemble vector, using the 2-step process:
154: VecAssemblyBegin(), VecAssemblyEnd()
155: Computations can be done while messages are in transition,
156: by placing code between these two statements.
157: */
158: VecAssemblyBegin(u);
159: VecAssemblyEnd(u);
161: /*
162: Compute right-hand-side vector
163: */
164: MatMult(C,u,b);
165:
166: /*
167: Create linear solver context
168: */
169: KSPCreate(PETSC_COMM_WORLD,&ksp);
171: /*
172: Set operators. Here the matrix that defines the linear system
173: also serves as the preconditioning matrix.
174: */
175: KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
177: /*
178: Set runtime options (e.g., -ksp_type <type> -pc_type <type>)
179: */
181: KSPSetFromOptions(ksp);
183: /*
184: Solve linear system. Here we explicitly call KSPSetUp() for more
185: detailed performance monitoring of certain preconditioners, such
186: as ICC and ILU. This call is optional, as KSPSetUp() will
187: automatically be called within KSPSolve() if it hasn't been
188: called already.
189: */
190: KSPSetRhs(ksp,b);
191: KSPSetSolution(ksp,x);
192: KSPSetUp(ksp);
193: KSPSolve(ksp);
194:
195: /*
196: Check the error
197: */
198: VecAXPY(&none,u,x);
199: VecNorm(x,NORM_2,&norm);
200: KSPGetIterationNumber(ksp,&its);
201: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A, Iterations %d\n",norm,its);
203: /* -------------- Stage 1: Solve Second System ---------------------- */
204: /*
205: Solve another linear system with the same method. We reuse the KSP
206: context, matrix and vector data structures, and hence save the
207: overhead of creating new ones.
209: Indicate to PETSc profiling that we're concluding the first
210: stage with PetscLogStagePop(), and beginning the second stage with
211: PetscLogStagePush().
212: */
213: PetscLogStagePop();
214: PetscLogStagePush(stages[1]);
216: /*
217: Initialize all matrix entries to zero. MatZeroEntries() retains the
218: nonzero structure of the matrix for sparse formats.
219: */
220: MatZeroEntries(C);
222: /*
223: Assemble matrix again. Note that we retain the same matrix data
224: structure and the same nonzero pattern; we just change the values
225: of the matrix entries.
226: */
227: for (i=0; i<m; i++) {
228: for (j=2*rank; j<2*rank+2; j++) {
229: v = -1.0; I = j + n*i;
230: if (i>0) {J = I - n; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
231: if (i<m-1) {J = I + n; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
232: if (j>0) {J = I - 1; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
233: if (j<n-1) {J = I + 1; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
234: v = 6.0; MatSetValues(C,1,&I,1,&I,&v,ADD_VALUES);
235: }
236: }
237: if (mat_nonsymmetric) {
238: for (I=Istart; I<Iend; I++) {
239: v = -1.5; i = I/n;
240: if (i>1) {J = I-n-1; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
241: }
242: }
243: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
244: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
246: /*
247: Compute another right-hand-side vector
248: */
249: MatMult(C,u,b);
251: /*
252: Set operators. Here the matrix that defines the linear system
253: also serves as the preconditioning matrix.
254: - The flag SAME_NONZERO_PATTERN indicates that the
255: preconditioning matrix has identical nonzero structure
256: as during the last linear solve (although the values of
257: the entries have changed). Thus, we can save some
258: work in setting up the preconditioner (e.g., no need to
259: redo symbolic factorization for ILU/ICC preconditioners).
260: - If the nonzero structure of the matrix is different during
261: the second linear solve, then the flag DIFFERENT_NONZERO_PATTERN
262: must be used instead. If you are unsure whether the
263: matrix structure has changed or not, use the flag
264: DIFFERENT_NONZERO_PATTERN.
265: - Caution: If you specify SAME_NONZERO_PATTERN, PETSc
266: believes your assertion and does not check the structure
267: of the matrix. If you erroneously claim that the structure
268: is the same when it actually is not, the new preconditioner
269: will not function correctly. Thus, use this optimization
270: feature with caution!
271: */
272: KSPSetOperators(ksp,C,C,SAME_NONZERO_PATTERN);
274: /*
275: Solve linear system
276: */
277: KSPSetRhs(ksp,b);
278: KSPSetSolution(ksp,x);
279: KSPSetUp(ksp);
280: KSPSolve(ksp);
282: /*
283: Check the error
284: */
285: VecAXPY(&none,u,x);
286: VecNorm(x,NORM_2,&norm);
287: KSPGetIterationNumber(ksp,&its);
288: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A, Iterations %d\n",norm,its);
290: /*
291: Free work space. All PETSc objects should be destroyed when they
292: are no longer needed.
293: */
294: KSPDestroy(ksp);
295: VecDestroy(u);
296: VecDestroy(x);
297: VecDestroy(b);
298: MatDestroy(C);
300: /*
301: Indicate to PETSc profiling that we're concluding the second stage
302: */
303: PetscLogStagePop();
305: PetscFinalize();
306: return 0;
307: }