Actual source code: ex11.c
1: /*$Id: ex11.c,v 1.34 2001/03/23 23:23:55 balay Exp $*/
3: static char help[] = "Solves a linear system in parallel with SLES.nn";
5: /*T
6: Concepts: SLES^solving a Helmholtz equation
7: Concepts: complex numbers;
8: Concepts: Helmholtz equation
9: Processors: n
10: T*/
12: /*
13: Description: Solves a complex linear system in parallel with SLES.
15: The model problem:
16: Solve Helmholtz equation on the unit square: (0,1) x (0,1)
17: -delta u - sigma1*u + i*sigma2*u = f,
18: where delta = Laplace operator
19: Dirichlet b.c.'s on all sides
20: Use the 2-D, five-point finite difference stencil.
22: Compiling the code:
23: This code uses the complex numbers version of PETSc, so one of the
24: following values of BOPT must be used for compiling the PETSc libraries
25: and this example:
26: BOPT=g_complex - debugging version
27: BOPT=O_complex - optimized version
28: BOPT=Opg_complex - profiling version
29: */
31: /*
32: Include "petscsles.h" so that we can use SLES solvers. Note that this file
33: automatically includes:
34: petsc.h - base PETSc routines petscvec.h - vectors
35: petscsys.h - system routines petscmat.h - matrices
36: petscis.h - index sets petscksp.h - Krylov subspace methods
37: petscviewer.h - viewers petscpc.h - preconditioners
38: */
39: #include petscsles.h
41: int main(int argc,char **args)
42: {
43: Vec x,b,u; /* approx solution, RHS, exact solution */
44: Mat A; /* linear system matrix */
45: SLES sles; /* linear solver context */
46: double norm; /* norm of solution error */
47: int dim,i,j,I,J,Istart,Iend,ierr,n = 6,its,use_random;
48: Scalar v,none = -1.0,sigma2,pfive = 0.5,*xa;
49: PetscRandom rctx;
50: double h2,sigma1 = 100.0;
51: PetscTruth flg;
53: PetscInitialize(&argc,&args,(char *)0,help);
54: #if !defined(PETSC_USE_COMPLEX)
55: SETERRQ(1,"This example requires complex numbers");
56: #endif
58: PetscOptionsGetDouble(PETSC_NULL,"-sigma1",&sigma1,PETSC_NULL);
59: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
60: dim = n*n;
62: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63: Compute the matrix and right-hand-side vector that define
64: the linear system, Ax = b.
65: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
66: /*
67: Create parallel matrix, specifying only its global dimensions.
68: When using MatCreate(), the matrix format can be specified at
69: runtime. Also, the parallel partitioning of the matrix is
70: determined by PETSc at runtime.
71: */
72: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,dim,dim,&A);
73: MatSetFromOptions(A);
75: /*
76: Currently, all PETSc parallel matrix formats are partitioned by
77: contiguous chunks of rows across the processors. Determine which
78: rows of the matrix are locally owned.
79: */
80: MatGetOwnershipRange(A,&Istart,&Iend);
82: /*
83: Set matrix elements in parallel.
84: - Each processor needs to insert only elements that it owns
85: locally (but any non-local elements will be sent to the
86: appropriate processor during matrix assembly).
87: - Always specify global rows and columns of matrix entries.
88: */
90: PetscOptionsHasName(PETSC_NULL,"-norandom",&flg);
91: if (flg) use_random = 0;
92: else use_random = 1;
93: if (use_random) {
94: PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT_IMAGINARY,&rctx);
95: } else {
96: sigma2 = 10.0*PETSC_i;
97: }
98: h2 = 1.0/((n+1)*(n+1));
99: for (I=Istart; I<Iend; I++) {
100: v = -1.0; i = I/n; j = I - i*n;
101: if (i>0) {
102: J = I-n; MatSetValues(A,1,&I,1,&J,&v,ADD_VALUES);}
103: if (i<n-1) {
104: J = I+n; MatSetValues(A,1,&I,1,&J,&v,ADD_VALUES);}
105: if (j>0) {
106: J = I-1; MatSetValues(A,1,&I,1,&J,&v,ADD_VALUES);}
107: if (j<n-1) {
108: J = I+1; MatSetValues(A,1,&I,1,&J,&v,ADD_VALUES);}
109: if (use_random) {PetscRandomGetValue(rctx,&sigma2);}
110: v = 4.0 - sigma1*h2 + sigma2*h2;
111: MatSetValues(A,1,&I,1,&I,&v,ADD_VALUES);
112: }
113: if (use_random) {PetscRandomDestroy(rctx);}
115: /*
116: Assemble matrix, using the 2-step process:
117: MatAssemblyBegin(), MatAssemblyEnd()
118: Computations can be done while messages are in transition
119: by placing code between these two statements.
120: */
121: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
122: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
124: /*
125: Create parallel vectors.
126: - When using VecCreate() and VecSetFromOptions(), we specify only the vector's global
127: dimension; the parallel partitioning is determined at runtime.
128: - Note: We form 1 vector from scratch and then duplicate as needed.
129: */
130: VecCreate(PETSC_COMM_WORLD,PETSC_DECIDE,dim,&u);
131: VecSetFromOptions(u);
132: VecDuplicate(u,&b);
133: VecDuplicate(b,&x);
135: /*
136: Set exact solution; then compute right-hand-side vector.
137: */
138:
139: if (use_random) {
140: PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT,&rctx);
141: VecSetRandom(rctx,u);
142: } else {
143: VecSet(&pfive,u);
144: }
145: MatMult(A,u,b);
147: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148: Create the linear solver and set various options
149: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
151: /*
152: Create linear solver context
153: */
154: SLESCreate(PETSC_COMM_WORLD,&sles);
156: /*
157: Set operators. Here the matrix that defines the linear system
158: also serves as the preconditioning matrix.
159: */
160: SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN);
162: /*
163: Set runtime options, e.g.,
164: -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
165: */
166: SLESSetFromOptions(sles);
168: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
169: Solve the linear system
170: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
172: SLESSolve(sles,b,x,&its);
174: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175: Check solution and clean up
176: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
178: /*
179: Print the first 3 entries of x; this demonstrates extraction of the
180: real and imaginary components of the complex vector, x.
181: */
182: PetscOptionsHasName(PETSC_NULL,"-print_x3",&flg);
183: if (flg) {
184: VecGetArray(x,&xa);
185: PetscPrintf(PETSC_COMM_WORLD,"The first three entries of x are:n");
186: for (i=0; i<3; i++){
187: PetscPrintf(PETSC_COMM_WORLD,"x[%d] = %g + %g in",i,PetscRealPart(xa[i]),PetscImaginaryPart(xa[i]));
188: }
189: VecRestoreArray(x,&xa);
190: }
192: /*
193: Check the error
194: */
195: VecAXPY(&none,u,x);
196: VecNorm(x,NORM_2,&norm);
197: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A iterations %dn",norm,its);
199: /*
200: Free work space. All PETSc objects should be destroyed when they
201: are no longer needed.
202: */
203: SLESDestroy(sles);
204: if (use_random) {PetscRandomDestroy(rctx);}
205: VecDestroy(u); VecDestroy(x);
206: VecDestroy(b); MatDestroy(A);
207: PetscFinalize();
208: return 0;
209: }