Actual source code: ex11.c
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: PetscTruth flg = PETSC_FALSE;
51: PetscInitialize(&argc,&args,(char *)0,help);
52: #if !defined(PETSC_USE_COMPLEX)
53: SETERRQ(1,"This example requires complex numbers");
54: #endif
56: PetscOptionsGetReal(PETSC_NULL,"-sigma1",&sigma1,PETSC_NULL);
57: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
58: dim = n*n;
60: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61: Compute the matrix and right-hand-side vector that define
62: the linear system, Ax = b.
63: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
64: /*
65: Create parallel matrix, specifying only its global dimensions.
66: When using MatCreate(), the matrix format can be specified at
67: runtime. Also, the parallel partitioning of the matrix is
68: determined by PETSc at runtime.
69: */
70: MatCreate(PETSC_COMM_WORLD,&A);
71: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,dim,dim);
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 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: */
89: PetscOptionsGetTruth(PETSC_NULL,"-norandom",&flg,PETSC_NULL);
90: if (flg) use_random = 0;
91: else use_random = 1;
92: if (use_random) {
93: PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
94: PetscRandomSetFromOptions(rctx);
95: PetscRandomSetInterval(rctx,0.0,PETSC_i);
96: } else {
97: sigma2 = 10.0*PETSC_i;
98: }
99: h2 = 1.0/((n+1)*(n+1));
100: for (Ii=Istart; Ii<Iend; Ii++) {
101: v = -1.0; i = Ii/n; j = Ii - i*n;
102: if (i>0) {
103: J = Ii-n; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);}
104: if (i<n-1) {
105: J = Ii+n; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);}
106: if (j>0) {
107: J = Ii-1; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);}
108: if (j<n-1) {
109: J = Ii+1; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);}
110: if (use_random) {PetscRandomGetValue(rctx,&sigma2);}
111: v = 4.0 - sigma1*h2 + sigma2*h2;
112: MatSetValues(A,1,&Ii,1,&Ii,&v,ADD_VALUES);
113: }
114: if (use_random) {PetscRandomDestroy(rctx);}
116: /*
117: Assemble matrix, using the 2-step process:
118: MatAssemblyBegin(), MatAssemblyEnd()
119: Computations can be done while messages are in transition
120: by placing code between these two statements.
121: */
122: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
123: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
125: /*
126: Create parallel vectors.
127: - When using VecCreate(), VecSetSizes() and VecSetFromOptions(),
128: we specify only the vector's global
129: dimension; the parallel partitioning is determined at runtime.
130: - Note: We form 1 vector from scratch and then duplicate as needed.
131: */
132: VecCreate(PETSC_COMM_WORLD,&u);
133: VecSetSizes(u,PETSC_DECIDE,dim);
134: VecSetFromOptions(u);
135: VecDuplicate(u,&b);
136: VecDuplicate(b,&x);
138: /*
139: Set exact solution; then compute right-hand-side vector.
140: */
141:
142: if (use_random) {
143: PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
144: PetscRandomSetFromOptions(rctx);
145: VecSetRandom(u,rctx);
146: } else {
147: VecSet(u,pfive);
148: }
149: MatMult(A,u,b);
151: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152: Create the linear solver and set various options
153: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155: /*
156: Create linear solver context
157: */
158: KSPCreate(PETSC_COMM_WORLD,&ksp);
160: /*
161: Set operators. Here the matrix that defines the linear system
162: also serves as the preconditioning matrix.
163: */
164: KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN);
166: /*
167: Set runtime options, e.g.,
168: -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
169: */
170: KSPSetFromOptions(ksp);
172: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
173: Solve the linear system
174: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
176: KSPSolve(ksp,b,x);
178: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
179: Check solution and clean up
180: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
182: /*
183: Print the first 3 entries of x; this demonstrates extraction of the
184: real and imaginary components of the complex vector, x.
185: */
186: flg = PETSC_FALSE;
187: PetscOptionsGetTruth(PETSC_NULL,"-print_x3",&flg,PETSC_NULL);
188: if (flg) {
189: VecGetArray(x,&xa);
190: PetscPrintf(PETSC_COMM_WORLD,"The first three entries of x are:\n");
191: for (i=0; i<3; i++){
192: PetscPrintf(PETSC_COMM_WORLD,"x[%D] = %G + %G i\n",i,PetscRealPart(xa[i]),PetscImaginaryPart(xa[i]));
193: }
194: VecRestoreArray(x,&xa);
195: }
197: /*
198: Check the error
199: */
200: VecAXPY(x,none,u);
201: VecNorm(x,NORM_2,&norm);
202: KSPGetIterationNumber(ksp,&its);
203: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A iterations %D\n",norm,its);
205: /*
206: Free work space. All PETSc objects should be destroyed when they
207: are no longer needed.
208: */
209: KSPDestroy(ksp);
210: if (use_random) {PetscRandomDestroy(rctx);}
211: VecDestroy(u); VecDestroy(x);
212: VecDestroy(b); MatDestroy(A);
213: PetscFinalize();
214: return 0;
215: }