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:      petsc.h       - base PETSc routines   petscvec.h - vectors
 30:      petscsys.h    - system routines       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,I,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;

 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:   PetscOptionsHasName(PETSC_NULL,"-norandom",&flg);
 90:   if (flg) use_random = 0;
 91:   else     use_random = 1;
 92:   if (use_random) {
 93:     PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT_IMAGINARY,&rctx);
 94:   } else {
 95:     sigma2 = 10.0*PETSC_i;
 96:   }
 97:   h2 = 1.0/((n+1)*(n+1));
 98:   for (I=Istart; I<Iend; I++) {
 99:     v = -1.0; i = I/n; j = I - i*n;
100:     if (i>0) {
101:       J = I-n; MatSetValues(A,1,&I,1,&J,&v,ADD_VALUES);}
102:     if (i<n-1) {
103:       J = I+n; MatSetValues(A,1,&I,1,&J,&v,ADD_VALUES);}
104:     if (j>0) {
105:       J = I-1; MatSetValues(A,1,&I,1,&J,&v,ADD_VALUES);}
106:     if (j<n-1) {
107:       J = I+1; MatSetValues(A,1,&I,1,&J,&v,ADD_VALUES);}
108:     if (use_random) {PetscRandomGetValue(rctx,&sigma2);}
109:     v = 4.0 - sigma1*h2 + sigma2*h2;
110:     MatSetValues(A,1,&I,1,&I,&v,ADD_VALUES);
111:   }
112:   if (use_random) {PetscRandomDestroy(rctx);}

114:   /* 
115:      Assemble matrix, using the 2-step process:
116:        MatAssemblyBegin(), MatAssemblyEnd()
117:      Computations can be done while messages are in transition
118:      by placing code between these two statements.
119:   */
120:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
121:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

123:   /* 
124:      Create parallel vectors.
125:       - When using VecCreate(), VecSetSizes() and VecSetFromOptions(),
126:       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,&u);
131:   VecSetSizes(u,PETSC_DECIDE,dim);
132:   VecSetFromOptions(u);
133:   VecDuplicate(u,&b);
134:   VecDuplicate(b,&x);

136:   /* 
137:      Set exact solution; then compute right-hand-side vector.
138:   */
139: 
140:   if (use_random) {
141:     PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT,&rctx);
142:     VecSetRandom(u,rctx);
143:   } else {
144:     VecSet(u,pfive);
145:   }
146:   MatMult(A,u,b);

148:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
149:                 Create the linear solver and set various options
150:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

152:   /* 
153:      Create linear solver context
154:   */
155:   KSPCreate(PETSC_COMM_WORLD,&ksp);

157:   /* 
158:      Set operators. Here the matrix that defines the linear system
159:      also serves as the preconditioning matrix.
160:   */
161:   KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN);

163:   /* 
164:     Set runtime options, e.g.,
165:         -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
166:   */
167:   KSPSetFromOptions(ksp);

169:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
170:                       Solve the linear system
171:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

173:   KSPSolve(ksp,b,x);

175:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
176:                       Check solution and clean up
177:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

179:   /*
180:       Print the first 3 entries of x; this demonstrates extraction of the
181:       real and imaginary components of the complex vector, x.
182:   */
183:   PetscOptionsHasName(PETSC_NULL,"-print_x3",&flg);
184:   if (flg) {
185:     VecGetArray(x,&xa);
186:     PetscPrintf(PETSC_COMM_WORLD,"The first three entries of x are:\n");
187:     for (i=0; i<3; i++){
188:       PetscPrintf(PETSC_COMM_WORLD,"x[%D] = %g + %g i\n",i,PetscRealPart(xa[i]),PetscImaginaryPart(xa[i]));
189:   }
190:     VecRestoreArray(x,&xa);
191:   }

193:   /* 
194:      Check the error
195:   */
196:   VecAXPY(x,none,u);
197:   VecNorm(x,NORM_2,&norm);
198:   KSPGetIterationNumber(ksp,&its);
199:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A iterations %D\n",norm,its);

201:   /* 
202:      Free work space.  All PETSc objects should be destroyed when they
203:      are no longer needed.
204:   */
205:   KSPDestroy(ksp);
206:   if (use_random) {PetscRandomDestroy(rctx);}
207:   VecDestroy(u); VecDestroy(x);
208:   VecDestroy(b); MatDestroy(A);
209:   PetscFinalize();
210:   return 0;
211: }