Actual source code: ex16.c
2: /* Usage: mpiexec ex16 [-help] [all PETSc options] */
4: static char help[] = "Solves a sequence of linear systems with different right-hand-side vectors.\n\
5: Input parameters include:\n\
6: -ntimes <ntimes> : number of linear systems to solve\n\
7: -view_exact_sol : write exact solution vector to stdout\n\
8: -m <mesh_x> : number of mesh points in x-direction\n\
9: -n <mesh_n> : number of mesh points in y-direction\n\n";
11: /*T
12: Concepts: KSP^repeatedly solving linear systems;
13: Concepts: KSP^Laplacian, 2d
14: Concepts: Laplacian, 2d
15: Processors: n
16: 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
30: int main(int argc,char **args)
31: {
32: Vec x,b,u; /* approx solution, RHS, exact solution */
33: Mat A; /* linear system matrix */
34: KSP ksp; /* linear solver context */
35: PetscReal norm; /* norm of solution error */
37: PetscInt ntimes,i,j,k,Ii,J,Istart,Iend;
38: PetscInt m = 8,n = 7,its;
39: PetscTruth flg = PETSC_FALSE;
40: PetscScalar v,one = 1.0,neg_one = -1.0,rhs;
42: PetscInitialize(&argc,&args,(char *)0,help);
43: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
44: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
46: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47: Compute the matrix for use in solving a series of
48: linear systems of the form, A x_i = b_i, for i=1,2,...
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.
55: */
56: MatCreate(PETSC_COMM_WORLD,&A);
57: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
58: MatSetFromOptions(A);
60: /*
61: Currently, all PETSc parallel matrix formats are partitioned by
62: contiguous chunks of rows across the processors. Determine which
63: rows of the matrix are locally owned.
64: */
65: MatGetOwnershipRange(A,&Istart,&Iend);
67: /*
68: Set matrix elements for the 2-D, five-point stencil in parallel.
69: - Each processor needs to insert only elements that it owns
70: locally (but any non-local elements will be sent to the
71: appropriate processor during matrix assembly).
72: - Always specify global rows and columns of matrix entries.
73: */
74: for (Ii=Istart; Ii<Iend; Ii++) {
75: v = -1.0; i = Ii/n; j = Ii - i*n;
76: if (i>0) {J = Ii - n; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
77: if (i<m-1) {J = Ii + n; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
78: if (j>0) {J = Ii - 1; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
79: if (j<n-1) {J = Ii + 1; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
80: v = 4.0; MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES);
81: }
83: /*
84: Assemble matrix, using the 2-step process:
85: MatAssemblyBegin(), MatAssemblyEnd()
86: Computations can be done while messages are in transition
87: by placing code between these two statements.
88: */
89: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
90: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
92: /*
93: Create parallel vectors.
94: - When using VecCreate(), VecSetSizes() and VecSetFromOptions(),
95: we specify only the vector's global
96: dimension; the parallel partitioning is determined at runtime.
97: - When solving a linear system, the vectors and matrices MUST
98: be partitioned accordingly. PETSc automatically generates
99: appropriately partitioned matrices and vectors when MatCreate()
100: and VecCreate() are used with the same communicator.
101: - Note: We form 1 vector from scratch and then duplicate as needed.
102: */
103: VecCreate(PETSC_COMM_WORLD,&u);
104: VecSetSizes(u,PETSC_DECIDE,m*n);
105: VecSetFromOptions(u);
106: VecDuplicate(u,&b);
107: VecDuplicate(b,&x);
109: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110: Create the linear solver and set various options
111: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113: /*
114: Create linear solver context
115: */
116: KSPCreate(PETSC_COMM_WORLD,&ksp);
118: /*
119: Set operators. Here the matrix that defines the linear system
120: also serves as the preconditioning matrix.
121: */
122: KSPSetOperators(ksp,A,A,SAME_PRECONDITIONER);
124: /*
125: Set runtime options, e.g.,
126: -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
127: These options will override those specified above as long as
128: KSPSetFromOptions() is called _after_ any other customization
129: routines.
130: */
131: KSPSetFromOptions(ksp);
133: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134: Solve several linear systems of the form A x_i = b_i
135: I.e., we retain the same matrix (A) for all systems, but
136: change the right-hand-side vector (b_i) at each step.
138: In this case, we simply call KSPSolve() multiple times. The
139: preconditioner setup operations (e.g., factorization for ILU)
140: be done during the first call to KSPSolve() only; such operations
141: will NOT be repeated for successive solves.
142: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144: ntimes = 2;
145: PetscOptionsGetInt(PETSC_NULL,"-ntimes",&ntimes,PETSC_NULL);
146: for (k=1; k<ntimes+1; k++) {
148: /*
149: Set exact solution; then compute right-hand-side vector. We use
150: an exact solution of a vector with all elements equal to 1.0*k.
151: */
152: rhs = one * (PetscReal)k;
153: VecSet(u,rhs);
154: MatMult(A,u,b);
156: /*
157: View the exact solution vector if desired
158: */
159: PetscOptionsGetTruth(PETSC_NULL,"-view_exact_sol",&flg,PETSC_NULL);
160: if (flg) {VecView(u,PETSC_VIEWER_STDOUT_WORLD);}
162: KSPSolve(ksp,b,x);
164: /*
165: Check the error
166: */
167: VecAXPY(x,neg_one,u);
168: VecNorm(x,NORM_2,&norm);
169: KSPGetIterationNumber(ksp,&its);
170: /*
171: Print convergence information. PetscPrintf() produces a single
172: print statement from all processes that share a communicator.
173: */
174: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A System %D: iterations %D\n",norm,k,its);
175: }
177: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178: Clean up
179: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
180: /*
181: Free work space. All PETSc objects should be destroyed when they
182: are no longer needed.
183: */
184: KSPDestroy(ksp);
185: VecDestroy(u); VecDestroy(x);
186: VecDestroy(b); MatDestroy(A);
188: /*
189: Always call PetscFinalize() before exiting a program. This routine
190: - finalizes the PETSc libraries as well as MPI
191: - provides summary and diagnostic information if certain runtime
192: options are chosen (e.g., -log_summary).
193: */
194: PetscFinalize();
195: return 0;
196: }