Actual source code: ex23.c
1: /*$Id: ex23.c,v 1.11 2001/08/07 21:30:54 bsmith Exp $*/
3: /* Program usage: mpirun ex23 [-help] [all PETSc options] */
5: static char help[] = "Solves a tridiagonal linear system.nn";
7: /*T
8: Concepts: SLES^basic parallel example;
9: Processors: n
10: T*/
12: /*
13: Include "petscsles.h" so that we can use SLES solvers. Note that this file
14: automatically includes:
15: petsc.h - base PETSc routines petscvec.h - vectors
16: petscsys.h - system routines petscmat.h - matrices
17: petscis.h - index sets petscksp.h - Krylov subspace methods
18: petscviewer.h - viewers petscpc.h - preconditioners
20: Note: The corresponding uniprocessor example is ex1.c
21: */
22: #include petscsles.h
24: int main(int argc,char **args)
25: {
26: Vec x, b, u; /* approx solution, RHS, exact solution */
27: Mat A; /* linear system matrix */
28: SLES sles; /* linear solver context */
29: PC pc; /* preconditioner context */
30: KSP ksp; /* Krylov subspace method context */
31: PetscReal norm; /* norm of solution error */
32: int ierr,i,n = 10,col[3],its,rstart,rend,nlocal;
33: PetscScalar neg_one = -1.0,one = 1.0,value[3];
35: PetscInitialize(&argc,&args,(char *)0,help);
36: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
38: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
39: Compute the matrix and right-hand-side vector that define
40: the linear system, Ax = b.
41: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
43: /*
44: Create vectors. Note that we form 1 vector from scratch and
45: then duplicate as needed. For this simple case let PETSc decide how
46: many elements of the vector are stored on each processor. The second
47: argument to VecSetSizes() below causes PETSc to decide.
48: */
49: VecCreate(PETSC_COMM_WORLD,&x);
50: VecSetSizes(x,PETSC_DECIDE,n);
51: VecSetFromOptions(x);
52: VecDuplicate(x,&b);
53: VecDuplicate(x,&u);
55: /* Identify the starting and ending mesh points on each
56: processor for the interior part of the mesh. We let PETSc decide
57: above. */
59: VecGetOwnershipRange(x,&rstart,&rend);
60: VecGetLocalSize(x,&nlocal);
62: /*
63: Create matrix. When using MatCreate(), the matrix format can
64: be specified at runtime.
66: Performance tuning note: For problems of substantial size,
67: preallocation of matrix memory is crucial for attaining good
68: performance. Since preallocation is not possible via the generic
69: matrix creation routine MatCreate(), we recommend for practical
70: problems instead to use the creation routine for a particular matrix
71: format, e.g.,
72: MatCreateMPIAIJ() - sequential AIJ (compressed sparse row)
73: MatCreateMPIBAIJ() - block AIJ
74: See the matrix chapter of the users manual for details.
76: We pass in nlocal as the "local" size of the matrix to force it
77: to have the same parallel layout as the vector created above.
78: */
79: MatCreate(PETSC_COMM_WORLD,nlocal,nlocal,n,n,&A);
80: MatSetFromOptions(A);
82: /*
83: Assemble matrix.
85: The linear system is distributed across the processors by
86: chunks of contiguous rows, which correspond to contiguous
87: sections of the mesh on which the problem is discretized.
88: For matrix assembly, each processor contributes entries for
89: the part that it owns locally.
90: */
93: if (!rstart) {
94: rstart = 1;
95: i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
96: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
97: }
98: if (rend == n) {
99: rend = n-1;
100: i = n-1; col[0] = n-2; col[1] = n-1; value[0] = -1.0; value[1] = 2.0;
101: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
102: }
104: /* Set entries corresponding to the mesh interior */
105: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
106: for (i=rstart; i<rend; i++) {
107: col[0] = i-1; col[1] = i; col[2] = i+1;
108: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
109: }
111: /* Assemble the matrix */
112: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
113: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
115: /*
116: Set exact solution; then compute right-hand-side vector.
117: */
118: VecSet(&one,u);
119: MatMult(A,u,b);
121: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
122: Create the linear solver and set various options
123: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
124: /*
125: Create linear solver context
126: */
127: SLESCreate(PETSC_COMM_WORLD,&sles);
129: /*
130: Set operators. Here the matrix that defines the linear system
131: also serves as the preconditioning matrix.
132: */
133: SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN);
135: /*
136: Set linear solver defaults for this problem (optional).
137: - By extracting the KSP and PC contexts from the SLES context,
138: we can then directly call any KSP and PC routines to set
139: various options.
140: - The following four statements are optional; all of these
141: parameters could alternatively be specified at runtime via
142: SLESSetFromOptions();
143: */
144: SLESGetKSP(sles,&ksp);
145: SLESGetPC(sles,&pc);
146: PCSetType(pc,PCJACOBI);
147: KSPSetTolerances(ksp,1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
149: /*
150: Set runtime options, e.g.,
151: -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
152: These options will override those specified above as long as
153: SLESSetFromOptions() is called _after_ any other customization
154: routines.
155: */
156: SLESSetFromOptions(sles);
157:
158: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
159: Solve the linear system
160: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
161: /*
162: Solve linear system
163: */
164: SLESSolve(sles,b,x,&its);
166: /*
167: View solver info; we could instead use the option -sles_view to
168: print this info to the screen at the conclusion of SLESSolve().
169: */
170: SLESView(sles,PETSC_VIEWER_STDOUT_WORLD);
172: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
173: Check solution and clean up
174: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
175: /*
176: Check the error
177: */
178: VecAXPY(&neg_one,u,x);
179: ierr = VecNorm(x,NORM_2,&norm);
180: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A, Iterations %dn",norm,its);
181: /*
182: Free work space. All PETSc objects should be destroyed when they
183: are no longer needed.
184: */
185: VecDestroy(x); VecDestroy(u);
186: VecDestroy(b); MatDestroy(A);
187: SLESDestroy(sles);
189: /*
190: Always call PetscFinalize() before exiting a program. This routine
191: - finalizes the PETSc libraries as well as MPI
192: - provides summary and diagnostic information if certain runtime
193: options are chosen (e.g., -log_summary).
194: */
195: PetscFinalize();
196: return 0;
197: }