Actual source code: ex9.c
2: static char help[] = "The solution of 2 different linear systems with different linear solvers.\n\
3: Also, this example illustrates the repeated\n\
4: solution of linear systems, while reusing matrix, vector, and solver data\n\
5: structures throughout the process. Note the various stages of event logging.\n\n";
7: /*T
8: Concepts: KSP^repeatedly solving linear systems;
9: Concepts: PetscLog^profiling multiple stages of code;
10: Concepts: PetscLog^user-defined event profiling;
11: Processors: n
12: T*/
14: /*
15: Include "petscksp.h" so that we can use KSP solvers. Note that this file
16: automatically includes:
17: petsc.h - base PETSc routines petscvec.h - vectors
18: petscsys.h - system routines petscmat.h - matrices
19: petscis.h - index sets petscksp.h - Krylov subspace methods
20: petscviewer.h - viewers petscpc.h - preconditioners
21: */
22: #include petscksp.h
24: /*
25: Declare user-defined routines
26: */
32: int main(int argc,char **args)
33: {
34: Vec x1,b1,x2,b2; /* solution and RHS vectors for systems #1 and #2 */
35: Vec u; /* exact solution vector */
36: Mat C1,C2; /* matrices for systems #1 and #2 */
37: KSP ksp1,ksp2; /* KSP contexts for systems #1 and #2 */
38: PetscInt ntimes = 3; /* number of times to solve the linear systems */
39: PetscEvent CHECK_ERROR; /* event number for error checking */
40: PetscInt ldim,low,high,iglobal,Istart,Iend,Istart2,Iend2;
41: PetscInt I,J,i,j,m = 3,n = 2,its,t;
43: int stages[3];
44: PetscTruth flg;
45: PetscScalar v;
46: PetscMPIInt rank,size;
48: PetscInitialize(&argc,&args,(char *)0,help);
49: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
50: PetscOptionsGetInt(PETSC_NULL,"-t",&ntimes,PETSC_NULL);
51: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
52: MPI_Comm_size(PETSC_COMM_WORLD,&size);
53: n = 2*size;
55: /*
56: Register various stages for profiling
57: */
58: PetscLogStageRegister(&stages[0],"Prelim setup");
59: PetscLogStageRegister(&stages[1],"Linear System 1");
60: PetscLogStageRegister(&stages[2],"Linear System 2");
62: /*
63: Register a user-defined event for profiling (error checking).
64: */
65: CHECK_ERROR = 0;
66: PetscLogEventRegister(&CHECK_ERROR,"Check Error",KSP_COOKIE);
68: /* - - - - - - - - - - - - Stage 0: - - - - - - - - - - - - - -
69: Preliminary Setup
70: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
72: PetscLogStagePush(stages[0]);
74: /*
75: Create data structures for first linear system.
76: - Create parallel matrix, specifying only its global dimensions.
77: When using MatCreate(), the matrix format can be specified at
78: runtime. Also, the parallel partitioning of the matrix is
79: determined by PETSc at runtime.
80: - Create parallel vectors.
81: - When using VecSetSizes(), we specify only the vector's global
82: dimension; the parallel partitioning is determined at runtime.
83: - Note: We form 1 vector from scratch and then duplicate as needed.
84: */
85: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,&C1);
86: MatSetFromOptions(C1);
87: MatGetOwnershipRange(C1,&Istart,&Iend);
88: VecCreate(PETSC_COMM_WORLD,&u);
89: VecSetSizes(u,PETSC_DECIDE,m*n);
90: VecSetFromOptions(u);
91: VecDuplicate(u,&b1);
92: VecDuplicate(u,&x1);
94: /*
95: Create first linear solver context.
96: Set runtime options (e.g., -pc_type <type>).
97: Note that the first linear system uses the default option
98: names, while the second linear systme uses a different
99: options prefix.
100: */
101: KSPCreate(PETSC_COMM_WORLD,&ksp1);
102: KSPSetFromOptions(ksp1);
104: /*
105: Set user-defined monitoring routine for first linear system.
106: */
107: PetscOptionsHasName(PETSC_NULL,"-my_ksp_monitor",&flg);
108: if (flg) {KSPSetMonitor(ksp1,MyKSPMonitor,PETSC_NULL,0);}
110: /*
111: Create data structures for second linear system.
112: */
113: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,&C2);
114: MatSetFromOptions(C2);
115: MatGetOwnershipRange(C2,&Istart2,&Iend2);
116: VecDuplicate(u,&b2);
117: VecDuplicate(u,&x2);
119: /*
120: Create second linear solver context
121: */
122: KSPCreate(PETSC_COMM_WORLD,&ksp2);
124: /*
125: Set different options prefix for second linear system.
126: Set runtime options (e.g., -s2_pc_type <type>)
127: */
128: KSPAppendOptionsPrefix(ksp2,"s2_");
129: KSPSetFromOptions(ksp2);
131: /*
132: Assemble exact solution vector in parallel. Note that each
133: processor needs to set only its local part of the vector.
134: */
135: VecGetLocalSize(u,&ldim);
136: VecGetOwnershipRange(u,&low,&high);
137: for (i=0; i<ldim; i++) {
138: iglobal = i + low;
139: v = (PetscScalar)(i + 100*rank);
140: VecSetValues(u,1,&iglobal,&v,ADD_VALUES);
141: }
142: VecAssemblyBegin(u);
143: VecAssemblyEnd(u);
145: /*
146: Log the number of flops for computing vector entries
147: */
148: PetscLogFlops(2*ldim);
150: /*
151: End curent profiling stage
152: */
153: PetscLogStagePop();
155: /* --------------------------------------------------------------
156: Linear solver loop:
157: Solve 2 different linear systems several times in succession
158: -------------------------------------------------------------- */
160: for (t=0; t<ntimes; t++) {
162: /* - - - - - - - - - - - - Stage 1: - - - - - - - - - - - - - -
163: Assemble and solve first linear system
164: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
166: /*
167: Begin profiling stage #1
168: */
169: PetscLogStagePush(stages[1]);
171: /*
172: Initialize all matrix entries to zero. MatZeroEntries() retains
173: the nonzero structure of the matrix for sparse formats.
174: */
175: if (t > 0) {MatZeroEntries(C1);}
177: /*
178: Set matrix entries in parallel. Also, log the number of flops
179: for computing matrix entries.
180: - Each processor needs to insert only elements that it owns
181: locally (but any non-local elements will be sent to the
182: appropriate processor during matrix assembly).
183: - Always specify global row and columns of matrix entries.
184: */
185: for (I=Istart; I<Iend; I++) {
186: v = -1.0; i = I/n; j = I - i*n;
187: if (i>0) {J = I - n; MatSetValues(C1,1,&I,1,&J,&v,ADD_VALUES);}
188: if (i<m-1) {J = I + n; MatSetValues(C1,1,&I,1,&J,&v,ADD_VALUES);}
189: if (j>0) {J = I - 1; MatSetValues(C1,1,&I,1,&J,&v,ADD_VALUES);}
190: if (j<n-1) {J = I + 1; MatSetValues(C1,1,&I,1,&J,&v,ADD_VALUES);}
191: v = 4.0; MatSetValues(C1,1,&I,1,&I,&v,ADD_VALUES);
192: }
193: for (I=Istart; I<Iend; I++) { /* Make matrix nonsymmetric */
194: v = -1.0*(t+0.5); i = I/n;
195: if (i>0) {J = I - n; MatSetValues(C1,1,&I,1,&J,&v,ADD_VALUES);}
196: }
197: PetscLogFlops(2*(Istart-Iend));
199: /*
200: Assemble matrix, using the 2-step process:
201: MatAssemblyBegin(), MatAssemblyEnd()
202: Computations can be done while messages are in transition
203: by placing code between these two statements.
204: */
205: MatAssemblyBegin(C1,MAT_FINAL_ASSEMBLY);
206: MatAssemblyEnd(C1,MAT_FINAL_ASSEMBLY);
208: /*
209: Indicate same nonzero structure of successive linear system matrices
210: */
211: MatSetOption(C1,MAT_NO_NEW_NONZERO_LOCATIONS);
213: /*
214: Compute right-hand-side vector
215: */
216: MatMult(C1,u,b1);
218: /*
219: Set operators. Here the matrix that defines the linear system
220: also serves as the preconditioning matrix.
221: - The flag SAME_NONZERO_PATTERN indicates that the
222: preconditioning matrix has identical nonzero structure
223: as during the last linear solve (although the values of
224: the entries have changed). Thus, we can save some
225: work in setting up the preconditioner (e.g., no need to
226: redo symbolic factorization for ILU/ICC preconditioners).
227: - If the nonzero structure of the matrix is different during
228: the second linear solve, then the flag DIFFERENT_NONZERO_PATTERN
229: must be used instead. If you are unsure whether the
230: matrix structure has changed or not, use the flag
231: DIFFERENT_NONZERO_PATTERN.
232: - Caution: If you specify SAME_NONZERO_PATTERN, PETSc
233: believes your assertion and does not check the structure
234: of the matrix. If you erroneously claim that the structure
235: is the same when it actually is not, the new preconditioner
236: will not function correctly. Thus, use this optimization
237: feature with caution!
238: */
239: KSPSetOperators(ksp1,C1,C1,SAME_NONZERO_PATTERN);
241: /*
242: Use the previous solution of linear system #1 as the initial
243: guess for the next solve of linear system #1. The user MUST
244: call KSPSetInitialGuessNonzero() in indicate use of an initial
245: guess vector; otherwise, an initial guess of zero is used.
246: */
247: if (t>0) {
248: KSPSetInitialGuessNonzero(ksp1,PETSC_TRUE);
249: }
251: /*
252: Solve the first linear system. Here we explicitly call
253: KSPSetUp() for more detailed performance monitoring of
254: certain preconditioners, such as ICC and ILU. This call
255: is optional, ase KSPSetUp() will automatically be called
256: within KSPSolve() if it hasn't been called already.
257: */
258: KSPSetUp(ksp1);
259: KSPSolve(ksp1,b1,x1);
260: KSPGetIterationNumber(ksp1,&its);
262: /*
263: Check error of solution to first linear system
264: */
265: CheckError(u,x1,b1,its,CHECK_ERROR);
267: /* - - - - - - - - - - - - Stage 2: - - - - - - - - - - - - - -
268: Assemble and solve second linear system
269: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
271: /*
272: Conclude profiling stage #1; begin profiling stage #2
273: */
274: PetscLogStagePop();
275: PetscLogStagePush(stages[2]);
277: /*
278: Initialize all matrix entries to zero
279: */
280: if (t > 0) {MatZeroEntries(C2);}
282: /*
283: Assemble matrix in parallel. Also, log the number of flops
284: for computing matrix entries.
285: - To illustrate the features of parallel matrix assembly, we
286: intentionally set the values differently from the way in
287: which the matrix is distributed across the processors. Each
288: entry that is not owned locally will be sent to the appropriate
289: processor during MatAssemblyBegin() and MatAssemblyEnd().
290: - For best efficiency the user should strive to set as many
291: entries locally as possible.
292: */
293: for (i=0; i<m; i++) {
294: for (j=2*rank; j<2*rank+2; j++) {
295: v = -1.0; I = j + n*i;
296: if (i>0) {J = I - n; MatSetValues(C2,1,&I,1,&J,&v,ADD_VALUES);}
297: if (i<m-1) {J = I + n; MatSetValues(C2,1,&I,1,&J,&v,ADD_VALUES);}
298: if (j>0) {J = I - 1; MatSetValues(C2,1,&I,1,&J,&v,ADD_VALUES);}
299: if (j<n-1) {J = I + 1; MatSetValues(C2,1,&I,1,&J,&v,ADD_VALUES);}
300: v = 6.0 + t*0.5; MatSetValues(C2,1,&I,1,&I,&v,ADD_VALUES);
301: }
302: }
303: for (I=Istart2; I<Iend2; I++) { /* Make matrix nonsymmetric */
304: v = -1.0*(t+0.5); i = I/n;
305: if (i>0) {J = I - n; MatSetValues(C2,1,&I,1,&J,&v,ADD_VALUES);}
306: }
307: MatAssemblyBegin(C2,MAT_FINAL_ASSEMBLY);
308: MatAssemblyEnd(C2,MAT_FINAL_ASSEMBLY);
309: PetscLogFlops(2*(Istart-Iend));
311: /*
312: Indicate same nonzero structure of successive linear system matrices
313: */
314: MatSetOption(C2,MAT_NO_NEW_NONZERO_LOCATIONS);
316: /*
317: Compute right-hand-side vector
318: */
319: MatMult(C2,u,b2);
321: /*
322: Set operators. Here the matrix that defines the linear system
323: also serves as the preconditioning matrix. Indicate same nonzero
324: structure of successive preconditioner matrices by setting flag
325: SAME_NONZERO_PATTERN.
326: */
327: KSPSetOperators(ksp2,C2,C2,SAME_NONZERO_PATTERN);
329: /*
330: Solve the second linear system
331: */
332: KSPSetUp(ksp2);
333: KSPSolve(ksp2,b2,x2);
334: KSPGetIterationNumber(ksp2,&its);
336: /*
337: Check error of solution to second linear system
338: */
339: CheckError(u,x2,b2,its,CHECK_ERROR);
341: /*
342: Conclude profiling stage #2
343: */
344: PetscLogStagePop();
345: }
346: /* --------------------------------------------------------------
347: End of linear solver loop
348: -------------------------------------------------------------- */
350: /*
351: Free work space. All PETSc objects should be destroyed when they
352: are no longer needed.
353: */
354: KSPDestroy(ksp1); KSPDestroy(ksp2);
355: VecDestroy(x1); VecDestroy(x2);
356: VecDestroy(b1); VecDestroy(b2);
357: MatDestroy(C1); MatDestroy(C2);
358: VecDestroy(u);
360: PetscFinalize();
361: return 0;
362: }
365: /* ------------------------------------------------------------- */
366: /*
367: CheckError - Checks the error of the solution.
369: Input Parameters:
370: u - exact solution
371: x - approximate solution
372: b - work vector
373: its - number of iterations for convergence
374: CHECK_ERROR - the event number for error checking
375: (for use with profiling)
377: Notes:
378: In order to profile this section of code separately from the
379: rest of the program, we register it as an "event" with
380: PetscLogEventRegister() in the main program. Then, we indicate
381: the start and end of this event by respectively calling
382: PetscLogEventBegin(CHECK_ERROR,u,x,b,0);
383: PetscLogEventEnd(CHECK_ERROR,u,x,b,0);
384: Here, we specify the objects most closely associated with
385: the event (the vectors u,x,b). Such information is optional;
386: we could instead just use 0 instead for all objects.
387: */
388: PetscErrorCode CheckError(Vec u,Vec x,Vec b,PetscInt its,PetscEvent CHECK_ERROR)
389: {
390: PetscScalar none = -1.0;
391: PetscReal norm;
394: PetscLogEventBegin(CHECK_ERROR,u,x,b,0);
396: /*
397: Compute error of the solution, using b as a work vector.
398: */
399: VecCopy(x,b);
400: VecAXPY(&none,u,b);
401: VecNorm(b,NORM_2,&norm);
402: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A, Iterations %D\n",norm,its);
403: PetscLogEventEnd(CHECK_ERROR,u,x,b,0);
404: return 0;
405: }
406: /* ------------------------------------------------------------- */
409: /*
410: MyKSPMonitor - This is a user-defined routine for monitoring
411: the KSP iterative solvers.
413: Input Parameters:
414: ksp - iterative context
415: n - iteration number
416: rnorm - 2-norm (preconditioned) residual value (may be estimated)
417: dummy - optional user-defined monitor context (unused here)
418: */
419: PetscErrorCode MyKSPMonitor(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
420: {
421: Vec x;
424: /*
425: Build the solution vector
426: */
427: KSPBuildSolution(ksp,PETSC_NULL,&x);
429: /*
430: Write the solution vector and residual norm to stdout.
431: - PetscPrintf() handles output for multiprocessor jobs
432: by printing from only one processor in the communicator.
433: - The parallel viewer PETSC_VIEWER_STDOUT_WORLD handles
434: data from multiple processors so that the output
435: is not jumbled.
436: */
437: PetscPrintf(PETSC_COMM_WORLD,"iteration %D solution vector:\n",n);
438: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
439: PetscPrintf(PETSC_COMM_WORLD,"iteration %D KSP Residual norm %14.12e \n",n,rnorm);
440: return 0;
441: }