Actual source code: ex11.c

  1: /*$Id: ex11.c,v 1.37 2001/08/07 21:30:54 bsmith 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: #undef __FUNCT__
 43: int main(int argc,char **args)
 44: {
 45:   Vec         x,b,u;      /* approx solution, RHS, exact solution */
 46:   Mat         A;            /* linear system matrix */
 47:   SLES        sles;         /* linear solver context */
 48:   PetscReal   norm;         /* norm of solution error */
 49:   int         dim,i,j,I,J,Istart,Iend,ierr,n = 6,its,use_random;
 50:   PetscScalar v,none = -1.0,sigma2,pfive = 0.5,*xa;
 51:   PetscRandom rctx;
 52:   PetscReal   h2,sigma1 = 100.0;
 53:   PetscTruth  flg;

 55:   PetscInitialize(&argc,&args,(char *)0,help);
 56: #if !defined(PETSC_USE_COMPLEX)
 57:   SETERRQ(1,"This example requires complex numbers");
 58: #endif

 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,PETSC_DECIDE,PETSC_DECIDE,dim,dim,&A);
 75:   MatSetFromOptions(A);

 77:   /* 
 78:      Currently, all PETSc parallel matrix formats are partitioned by
 79:      contiguous chunks of rows across the processors.  Determine which
 80:      rows of the matrix are locally owned. 
 81:   */
 82:   MatGetOwnershipRange(A,&Istart,&Iend);

 84:   /* 
 85:      Set matrix elements in parallel.
 86:       - Each processor needs to insert only elements that it owns
 87:         locally (but any non-local elements will be sent to the
 88:         appropriate processor during matrix assembly). 
 89:       - Always specify global rows and columns of matrix entries.
 90:   */

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

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

126:   /* 
127:      Create parallel vectors.
128:       - When using VecCreate(), VecSetSizes() and VecSetFromOptions(),
129:       we specify only the vector's global
130:         dimension; the parallel partitioning is determined at runtime. 
131:       - Note: We form 1 vector from scratch and then duplicate as needed.
132:   */
133:   VecCreate(PETSC_COMM_WORLD,&u);
134:   VecSetSizes(u,PETSC_DECIDE,dim);
135:   VecSetFromOptions(u);
136:   VecDuplicate(u,&b);
137:   VecDuplicate(b,&x);

139:   /* 
140:      Set exact solution; then compute right-hand-side vector.
141:   */
142: 
143:   if (use_random) {
144:     PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT,&rctx);
145:     VecSetRandom(rctx,u);
146:   } else {
147:     VecSet(&pfive,u);
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:   SLESCreate(PETSC_COMM_WORLD,&sles);

160:   /* 
161:      Set operators. Here the matrix that defines the linear system
162:      also serves as the preconditioning matrix.
163:   */
164:   SLESSetOperators(sles,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:   SLESSetFromOptions(sles);

172:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
173:                       Solve the linear system
174:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

176:   SLESSolve(sles,b,x,&its);

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:   PetscOptionsHasName(PETSC_NULL,"-print_x3",&flg);
187:   if (flg) {
188:     VecGetArray(x,&xa);
189:     PetscPrintf(PETSC_COMM_WORLD,"The first three entries of x are:n");
190:     for (i=0; i<3; i++){
191:       PetscPrintf(PETSC_COMM_WORLD,"x[%d] = %g + %g in",i,PetscRealPart(xa[i]),PetscImaginaryPart(xa[i]));
192:   }
193:     VecRestoreArray(x,&xa);
194:   }

196:   /* 
197:      Check the error
198:   */
199:   VecAXPY(&none,u,x);
200:   VecNorm(x,NORM_2,&norm);
201:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A iterations %dn",norm,its);

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