Actual source code: ex11.c
petsc-3.3-p7 2013-05-11
2: static char help[] = "Solves a linear system in parallel with KSP.\n\n";
4: /*T
5: Concepts: KSP^solving a Helmholtz equation
6: Concepts: complex numbers;
7: Concepts: Helmholtz equation
8: Processors: n
9: T*/
11: /*
12: Description: Solves a complex linear system in parallel with KSP.
14: The model problem:
15: Solve Helmholtz equation on the unit square: (0,1) x (0,1)
16: -delta u - sigma1*u + i*sigma2*u = f,
17: where delta = Laplace operator
18: Dirichlet b.c.'s on all sides
19: Use the 2-D, five-point finite difference stencil.
21: Compiling the code:
22: This code uses the complex numbers version of PETSc, so configure
23: must be run to enable this
24: */
26: /*
27: Include "petscksp.h" so that we can use KSP solvers. Note that this file
28: automatically includes:
29: petscsys.h - base PETSc routines petscvec.h - vectors
30: petscmat.h - matrices
31: petscis.h - index sets petscksp.h - Krylov subspace methods
32: petscviewer.h - viewers petscpc.h - preconditioners
33: */
34: #include <petscksp.h>
38: int main(int argc,char **args)
39: {
40: Vec x,b,u; /* approx solution, RHS, exact solution */
41: Mat A; /* linear system matrix */
42: KSP ksp; /* linear solver context */
43: PetscReal norm; /* norm of solution error */
44: PetscInt dim,i,j,Ii,J,Istart,Iend,n = 6,its,use_random;
46: PetscScalar v,none = -1.0,sigma2,pfive = 0.5,*xa;
47: PetscRandom rctx;
48: PetscReal h2,sigma1 = 100.0;
49: PetscBool flg = PETSC_FALSE;
50: PetscScalar a=1.0+PETSC_i;
52: PetscInitialize(&argc,&args,(char *)0,help);
53: #if !defined(PETSC_USE_COMPLEX)
54: SETERRQ(PETSC_COMM_WORLD,1,"This example requires complex numbers");
55: #endif
57: a=1.0+PETSC_i;
58: printf("%G+%Gi\n",PetscRealPart(a),PetscImaginaryPart(a));
60: PetscOptionsGetReal(PETSC_NULL,"-sigma1",&sigma1,PETSC_NULL);
61: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
62: dim = n*n;
64: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
65: Compute the matrix and right-hand-side vector that define
66: the linear system, Ax = b.
67: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
68: /*
69: Create parallel matrix, specifying only its global dimensions.
70: When using MatCreate(), the matrix format can be specified at
71: runtime. Also, the parallel partitioning of the matrix is
72: determined by PETSc at runtime.
73: */
74: MatCreate(PETSC_COMM_WORLD,&A);
75: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,dim,dim);
76: MatSetFromOptions(A);
78: /*
79: Currently, all PETSc parallel matrix formats are partitioned by
80: contiguous chunks of rows across the processors. Determine which
81: rows of the matrix are locally owned.
82: */
83: MatGetOwnershipRange(A,&Istart,&Iend);
85: /*
86: Set matrix elements in parallel.
87: - Each processor needs to insert only elements that it owns
88: locally (but any non-local elements will be sent to the
89: appropriate processor during matrix assembly).
90: - Always specify global rows and columns of matrix entries.
91: */
93: PetscOptionsGetBool(PETSC_NULL,"-norandom",&flg,PETSC_NULL);
94: if (flg) use_random = 0;
95: else use_random = 1;
96: if (use_random) {
97: PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
98: PetscRandomSetFromOptions(rctx);
99: PetscRandomSetInterval(rctx,0.0,PETSC_i);
100: } else {
101: sigma2 = 10.0*PETSC_i;
102: }
103: h2 = 1.0/((n+1)*(n+1));
104: for (Ii=Istart; Ii<Iend; Ii++) {
105: v = -1.0; i = Ii/n; j = Ii - i*n;
106: if (i>0) {
107: J = Ii-n; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);}
108: if (i<n-1) {
109: J = Ii+n; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);}
110: if (j>0) {
111: J = Ii-1; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);}
112: if (j<n-1) {
113: J = Ii+1; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);}
114: if (use_random) {PetscRandomGetValue(rctx,&sigma2);}
115: v = 4.0 - sigma1*h2 + sigma2*h2;
116: MatSetValues(A,1,&Ii,1,&Ii,&v,ADD_VALUES);
117: }
118: if (use_random) {PetscRandomDestroy(&rctx);}
120: /*
121: Assemble matrix, using the 2-step process:
122: MatAssemblyBegin(), MatAssemblyEnd()
123: Computations can be done while messages are in transition
124: by placing code between these two statements.
125: */
126: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
127: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
129: /*
130: Create parallel vectors.
131: - When using VecCreate(), VecSetSizes() and VecSetFromOptions(),
132: we specify only the vector's global
133: dimension; the parallel partitioning is determined at runtime.
134: - Note: We form 1 vector from scratch and then duplicate as needed.
135: */
136: VecCreate(PETSC_COMM_WORLD,&u);
137: VecSetSizes(u,PETSC_DECIDE,dim);
138: VecSetFromOptions(u);
139: VecDuplicate(u,&b);
140: VecDuplicate(b,&x);
142: /*
143: Set exact solution; then compute right-hand-side vector.
144: */
145:
146: if (use_random) {
147: PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
148: PetscRandomSetFromOptions(rctx);
149: VecSetRandom(u,rctx);
150: } else {
151: VecSet(u,pfive);
152: }
153: MatMult(A,u,b);
155: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
156: Create the linear solver and set various options
157: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
159: /*
160: Create linear solver context
161: */
162: KSPCreate(PETSC_COMM_WORLD,&ksp);
164: /*
165: Set operators. Here the matrix that defines the linear system
166: also serves as the preconditioning matrix.
167: */
168: KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN);
170: /*
171: Set runtime options, e.g.,
172: -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
173: */
174: KSPSetFromOptions(ksp);
176: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
177: Solve the linear system
178: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
180: KSPSolve(ksp,b,x);
182: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
183: Check solution and clean up
184: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
186: /*
187: Print the first 3 entries of x; this demonstrates extraction of the
188: real and imaginary components of the complex vector, x.
189: */
190: flg = PETSC_FALSE;
191: PetscOptionsGetBool(PETSC_NULL,"-print_x3",&flg,PETSC_NULL);
192: if (flg) {
193: VecGetArray(x,&xa);
194: PetscPrintf(PETSC_COMM_WORLD,"The first three entries of x are:\n");
195: for (i=0; i<3; i++){
196: PetscPrintf(PETSC_COMM_WORLD,"x[%D] = %G + %G i\n",i,PetscRealPart(xa[i]),PetscImaginaryPart(xa[i]));
197: }
198: VecRestoreArray(x,&xa);
199: }
201: /*
202: Check the error
203: */
204: VecAXPY(x,none,u);
205: VecNorm(x,NORM_2,&norm);
206: KSPGetIterationNumber(ksp,&its);
207: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %G iterations %D\n",norm,its);
209: /*
210: Free work space. All PETSc objects should be destroyed when they
211: are no longer needed.
212: */
213: KSPDestroy(&ksp);
214: if (use_random) {PetscRandomDestroy(&rctx);}
215: VecDestroy(&u); VecDestroy(&x);
216: VecDestroy(&b); MatDestroy(&A);
217: PetscFinalize();
218: return 0;
219: }