Actual source code: ex11.c

petsc-3.7.7 2017-09-25
Report Typos and Errors
  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;

 51:   PetscInitialize(&argc,&args,(char*)0,help);
 52: #if !defined(PETSC_USE_COMPLEX)
 53:   SETERRQ(PETSC_COMM_WORLD,1,"This example requires complex numbers");
 54: #endif

 56:   PetscOptionsGetReal(NULL,NULL,"-sigma1",&sigma1,NULL);
 57:   PetscOptionsGetInt(NULL,NULL,"-n",&n,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);
 73:   MatSetUp(A);

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

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

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

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

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

143:   /*
144:      Set exact solution; then compute right-hand-side vector.
145:   */

147:   if (use_random) {
148:     PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
149:     PetscRandomSetFromOptions(rctx);
150:     VecSetRandom(u,rctx);
151:   } else {
152:     VecSet(u,pfive);
153:   }
154:   MatMult(A,u,b);

156:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157:                 Create the linear solver and set various options
158:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

160:   /*
161:      Create linear solver context
162:   */
163:   KSPCreate(PETSC_COMM_WORLD,&ksp);

165:   /*
166:      Set operators. Here the matrix that defines the linear system
167:      also serves as the preconditioning matrix.
168:   */
169:   KSPSetOperators(ksp,A,A);

171:   /*
172:     Set runtime options, e.g.,
173:         -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
174:   */
175:   KSPSetFromOptions(ksp);

177:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178:                       Solve the linear system
179:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

181:   KSPSolve(ksp,b,x);

183:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
184:                       Check solution and clean up
185:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

202:   /*
203:      Check the error
204:   */
205:   VecAXPY(x,none,u);
206:   VecNorm(x,NORM_2,&norm);
207:   KSPGetIterationNumber(ksp,&its);
208:   if (norm < 1.e-12) {
209:     PetscPrintf(PETSC_COMM_WORLD,"Norm of error < 1.e-12 iterations %D\n",its);
210:   } else {
211:     PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g iterations %D\n",(double)norm,its);
212:   }

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: }