Actual source code: ex11.c

petsc-3.4.5 2014-06-29
  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(NULL,"-sigma1",&sigma1,NULL);
 61:   PetscOptionsGetInt(NULL,"-n",&n,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);
 77:   MatSetUp(A);

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

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

 94:   PetscOptionsGetBool(NULL,"-norandom",&flg,NULL);
 95:   if (flg) use_random = 0;
 96:   else use_random = 1;
 97:   if (use_random) {
 98:     PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
 99:     PetscRandomSetFromOptions(rctx);
100:     PetscRandomSetInterval(rctx,0.0,PETSC_i);
101:   } else {
102:     sigma2 = 10.0*PETSC_i;
103:   }
104:   h2 = 1.0/((n+1)*(n+1));
105:   for (Ii=Istart; Ii<Iend; Ii++) {
106:     v = -1.0; i = Ii/n; j = Ii - i*n;
107:     if (i>0) {
108:       J = Ii-n; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);
109:     }
110:     if (i<n-1) {
111:       J = Ii+n; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);
112:     }
113:     if (j>0) {
114:       J = Ii-1; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);
115:     }
116:     if (j<n-1) {
117:       J = Ii+1; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);
118:     }
119:     if (use_random) {PetscRandomGetValue(rctx,&sigma2);}
120:     v    = 4.0 - sigma1*h2 + sigma2*h2;
121:     MatSetValues(A,1,&Ii,1,&Ii,&v,ADD_VALUES);
122:   }
123:   if (use_random) {PetscRandomDestroy(&rctx);}

125:   /*
126:      Assemble matrix, using the 2-step process:
127:        MatAssemblyBegin(), MatAssemblyEnd()
128:      Computations can be done while messages are in transition
129:      by placing code between these two statements.
130:   */
131:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
132:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

134:   /*
135:      Create parallel vectors.
136:       - When using VecCreate(), VecSetSizes() and VecSetFromOptions(),
137:       we specify only the vector's global
138:         dimension; the parallel partitioning is determined at runtime.
139:       - Note: We form 1 vector from scratch and then duplicate as needed.
140:   */
141:   VecCreate(PETSC_COMM_WORLD,&u);
142:   VecSetSizes(u,PETSC_DECIDE,dim);
143:   VecSetFromOptions(u);
144:   VecDuplicate(u,&b);
145:   VecDuplicate(b,&x);

147:   /*
148:      Set exact solution; then compute right-hand-side vector.
149:   */

151:   if (use_random) {
152:     PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
153:     PetscRandomSetFromOptions(rctx);
154:     VecSetRandom(u,rctx);
155:   } else {
156:     VecSet(u,pfive);
157:   }
158:   MatMult(A,u,b);

160:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161:                 Create the linear solver and set various options
162:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

164:   /*
165:      Create linear solver context
166:   */
167:   KSPCreate(PETSC_COMM_WORLD,&ksp);

169:   /*
170:      Set operators. Here the matrix that defines the linear system
171:      also serves as the preconditioning matrix.
172:   */
173:   KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN);

175:   /*
176:     Set runtime options, e.g.,
177:         -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
178:   */
179:   KSPSetFromOptions(ksp);

181:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
182:                       Solve the linear system
183:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

185:   KSPSolve(ksp,b,x);

187:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
188:                       Check solution and clean up
189:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

191:   /*
192:       Print the first 3 entries of x; this demonstrates extraction of the
193:       real and imaginary components of the complex vector, x.
194:   */
195:   flg  = PETSC_FALSE;
196:   PetscOptionsGetBool(NULL,"-print_x3",&flg,NULL);
197:   if (flg) {
198:     VecGetArray(x,&xa);
199:     PetscPrintf(PETSC_COMM_WORLD,"The first three entries of x are:\n");
200:     for (i=0; i<3; i++) {
201:       PetscPrintf(PETSC_COMM_WORLD,"x[%D] = %G + %G i\n",i,PetscRealPart(xa[i]),PetscImaginaryPart(xa[i]));
202:     }
203:     VecRestoreArray(x,&xa);
204:   }

206:   /*
207:      Check the error
208:   */
209:   VecAXPY(x,none,u);
210:   VecNorm(x,NORM_2,&norm);
211:   KSPGetIterationNumber(ksp,&its);
212:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %G iterations %D\n",norm,its);

214:   /*
215:      Free work space.  All PETSc objects should be destroyed when they
216:      are no longer needed.
217:   */
218:   KSPDestroy(&ksp);
219:   if (use_random) {PetscRandomDestroy(&rctx);}
220:   VecDestroy(&u); VecDestroy(&x);
221:   VecDestroy(&b); MatDestroy(&A);
222:   PetscFinalize();
223:   return 0;
224: }