Actual source code: ex12.c
1: /*$Id: ex12.c,v 1.24 2001/08/07 21:30:54 bsmith Exp $*/
3: /* Program usage: mpirun -np <procs> ex12 [-help] [all PETSc options] */
5: static char help[] = "Solves a linear system in parallel with SLES.n
6: Input parameters include:n
7: -m <mesh_x> : number of mesh points in x-directionn
8: -n <mesh_n> : number of mesh points in y-directionnn";
10: /*T
11: Concepts: SLES^solving a system of linear equations
12: Concepts: SLES^Laplacian, 2d
13: Concepts: PC^registering preconditioners
14: Processors: n
15: T*/
17: /*
18: Demonstrates registering a new preconditioner (PC) type.
20: To register a PC type whose code is linked into the executable,
21: use PCRegister(). To register a PC type in a dynamic library use PCRegisterDynamic()
23: Also provide the prototype for your PCCreate_XXX() function. In
24: this example we use the PETSc implementation of the Jacobi method,
25: PCCreate_Jacobi() just as an example.
27: See the file src/sles/pc/impls/jacobi/jacobi.c for details on how to
28: write a new PC component.
30: See the manual page PCRegisterDynamic() for details on how to register a method.
31: */
33: /*
34: Include "petscsles.h" so that we can use SLES solvers. Note that this file
35: automatically includes:
36: petsc.h - base PETSc routines petscvec.h - vectors
37: petscsys.h - system routines petscmat.h - matrices
38: petscis.h - index sets petscksp.h - Krylov subspace methods
39: petscviewer.h - viewers petscpc.h - preconditioners
40: */
41: #include petscsles.h
43: EXTERN_C_BEGIN
44: extern int PCCreate_Jacobi(PC);
45: EXTERN_C_END
47: int main(int argc,char **args)
48: {
49: Vec x,b,u; /* approx solution, RHS, exact solution */
50: Mat A; /* linear system matrix */
51: SLES sles; /* linear solver context */
52: PetscReal norm; /* norm of solution error */
53: int i,j,I,J,Istart,Iend,ierr,m = 8,n = 7,its;
54: PetscScalar v,one = 1.0,neg_one = -1.0;
55: PC pc; /* preconditioner context */
57: PetscInitialize(&argc,&args,(char *)0,help);
58: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
59: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
61: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62: Compute the matrix and right-hand-side vector that define
63: the linear system, Ax = b.
64: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
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,&A);
72: MatSetFromOptions(A);
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(A,&Istart,&Iend);
81: /*
82: Set matrix elements for the 2-D, five-point stencil 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 rows 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(A,1,&I,1,&J,&v,INSERT_VALUES);}
91: if (i<m-1) {J = I + n; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
92: if (j>0) {J = I - 1; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
93: if (j<n-1) {J = I + 1; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
94: v = 4.0; MatSetValues(A,1,&I,1,&I,&v,INSERT_VALUES);
95: }
97: /*
98: Assemble matrix, using the 2-step process:
99: MatAssemblyBegin(), MatAssemblyEnd()
100: Computations can be done while messages are in transition
101: by placing code between these two statements.
102: */
103: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
104: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
106: /*
107: Create parallel vectors.
108: - When using VecCreate(), VecSetSizes() and VecSetFromOptions(),
109: we specify only the vector's global
110: dimension; the parallel partitioning is determined at runtime.
111: - When solving a linear system, the vectors and matrices MUST
112: be partitioned accordingly. PETSc automatically generates
113: appropriately partitioned matrices and vectors when MatCreate()
114: and VecCreate() are used with the same communicator.
115: - Note: We form 1 vector from scratch and then duplicate as needed.
116: */
117: VecCreate(PETSC_COMM_WORLD,&u);
118: VecSetSizes(u,PETSC_DECIDE,m*n);
119: VecSetFromOptions(u);
120: VecDuplicate(u,&b);
121: VecDuplicate(b,&x);
123: /*
124: Set exact solution; then compute right-hand-side vector.
125: Use an exact solution of a vector with all elements of 1.0;
126: */
127: VecSet(&one,u);
128: MatMult(A,u,b);
130: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131: Create the linear solver and set various options
132: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
134: /*
135: Create linear solver context
136: */
137: SLESCreate(PETSC_COMM_WORLD,&sles);
139: /*
140: Set operators. Here the matrix that defines the linear system
141: also serves as the preconditioning matrix.
142: */
143: SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN);
145: /*
146: First register a new PC type with the command PCRegister()
147: */
148: PCRegister("ourjacobi",0,"PCCreate_Jacobi",PCCreate_Jacobi);
149:
150: /*
151: Set the PC type to be the new method
152: */
153: SLESGetPC(sles,&pc);
154: PCSetType(pc,"ourjacobi");
156: /*
157: Set runtime options, e.g.,
158: -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
159: These options will override those specified above as long as
160: SLESSetFromOptions() is called _after_ any other customization
161: routines.
162: */
163: SLESSetFromOptions(sles);
165: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
166: Solve the linear system
167: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
169: SLESSolve(sles,b,x,&its);
171: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
172: Check solution and clean up
173: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
175: /*
176: Check the error
177: */
178: VecAXPY(&neg_one,u,x);
179: VecNorm(x,NORM_2,&norm);
181: /* Scale the norm */
182: /* norm *= sqrt(1.0/((m+1)*(n+1))); */
184: /*
185: Print convergence information. PetscPrintf() produces a single
186: print statement from all processes that share a communicator.
187: */
188: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A iterations %dn",norm,its);
190: /*
191: Free work space. All PETSc objects should be destroyed when they
192: are no longer needed.
193: */
194: SLESDestroy(sles);
195: VecDestroy(u); VecDestroy(x);
196: VecDestroy(b); MatDestroy(A);
198: /*
199: Always call PetscFinalize() before exiting a program. This routine
200: - finalizes the PETSc libraries as well as MPI
201: - provides summary and diagnostic information if certain runtime
202: options are chosen (e.g., -log_summary).
203: */
204: PetscFinalize();
205: return 0;
206: }