Actual source code: ex2.c

petsc-master 2017-01-23
Report Typos and Errors

  2: static char help[] = "Solves a linear system in parallel with KSP.\n\
  3: Input parameters include:\n\
  4:   -random_exact_sol : use a random exact solution vector\n\
  5:   -view_exact_sol   : write exact solution vector to stdout\n\
  6:   -m <mesh_x>       : number of mesh points in x-direction\n\
  7:   -n <mesh_n>       : number of mesh points in y-direction\n\n";

  9: /*T
 10:    Concepts: KSP^basic parallel example;
 11:    Concepts: KSP^Laplacian, 2d
 12:    Concepts: Laplacian, 2d
 13:    Processors: n
 14: T*/

 16: /*
 17:   Include "petscksp.h" so that we can use KSP solvers.  Note that this file
 18:   automatically includes:
 19:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 20:      petscmat.h - matrices
 21:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 22:      petscviewer.h - viewers               petscpc.h  - preconditioners
 23: */
 24:  #include <petscksp.h>

 26: int main(int argc,char **args)
 27: {
 28:   Vec            x,b,u;    /* approx solution, RHS, exact solution */
 29:   Mat            A;        /* linear system matrix */
 30:   KSP            ksp;      /* linear solver context */
 31:   PetscRandom    rctx;     /* random number generator context */
 32:   PetscReal      norm;     /* norm of solution error */
 33:   PetscInt       i,j,Ii,J,Istart,Iend,m = 8,n = 7,its;
 35:   PetscBool      flg = PETSC_FALSE;
 36:   PetscScalar    v;
 37: #if defined(PETSC_USE_LOG)
 38:   PetscLogStage stage;
 39: #endif

 41:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 42:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 43:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 44:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 45:          Compute the matrix and right-hand-side vector that define
 46:          the linear system, Ax = b.
 47:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 48:   /*
 49:      Create parallel matrix, specifying only its global dimensions.
 50:      When using MatCreate(), the matrix format can be specified at
 51:      runtime. Also, the parallel partitioning of the matrix is
 52:      determined by PETSc at runtime.

 54:      Performance tuning note:  For problems of substantial size,
 55:      preallocation of matrix memory is crucial for attaining good
 56:      performance. See the matrix chapter of the users manual for details.
 57:   */
 58:   MatCreate(PETSC_COMM_WORLD,&A);
 59:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
 60:   MatSetFromOptions(A);
 61:   MatMPIAIJSetPreallocation(A,5,NULL,5,NULL);
 62:   MatSeqAIJSetPreallocation(A,5,NULL);
 63:   MatSeqSBAIJSetPreallocation(A,1,5,NULL);

 65:   /*
 66:      Currently, all PETSc parallel matrix formats are partitioned by
 67:      contiguous chunks of rows across the processors.  Determine which
 68:      rows of the matrix are locally owned.
 69:   */
 70:   MatGetOwnershipRange(A,&Istart,&Iend);

 72:   /*
 73:      Set matrix elements for the 2-D, five-point stencil in parallel.
 74:       - Each processor needs to insert only elements that it owns
 75:         locally (but any non-local elements will be sent to the
 76:         appropriate processor during matrix assembly).
 77:       - Always specify global rows and columns of matrix entries.

 79:      Note: this uses the less common natural ordering that orders first
 80:      all the unknowns for x = h then for x = 2h etc; Hence you see J = Ii +- n
 81:      instead of J = I +- m as you might expect. The more standard ordering
 82:      would first do all variables for y = h, then y = 2h etc.

 84:    */
 85:   PetscLogStageRegister("Assembly", &stage);
 86:   PetscLogStagePush(stage);
 87:   for (Ii=Istart; Ii<Iend; Ii++) {
 88:     v = -1.0; i = Ii/n; j = Ii - i*n;
 89:     if (i>0)   {J = Ii - n; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);}
 90:     if (i<m-1) {J = Ii + n; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);}
 91:     if (j>0)   {J = Ii - 1; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);}
 92:     if (j<n-1) {J = Ii + 1; MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);}
 93:     v = 4.0; MatSetValues(A,1,&Ii,1,&Ii,&v,ADD_VALUES);
 94:   }

 96:   /*
 97:      Assemble matrix, using the 2-step process:
 98:        MatAssemblyBegin(), MatAssemblyEnd()
 99:      Computations can be done while messages are in transition
100:      by placing code between these two statements.
101:   */
102:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
103:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
104:   PetscLogStagePop();

106:   /* A is symmetric. Set symmetric flag to enable ICC/Cholesky preconditioner */
107:   MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);

109:   /*
110:      Create parallel vectors.
111:       - We form 1 vector from scratch and then duplicate as needed.
112:       - When using VecCreate(), VecSetSizes and VecSetFromOptions()
113:         in this example, we specify only the
114:         vector's global dimension; the parallel partitioning is determined
115:         at runtime.
116:       - When solving a linear system, the vectors and matrices MUST
117:         be partitioned accordingly.  PETSc automatically generates
118:         appropriately partitioned matrices and vectors when MatCreate()
119:         and VecCreate() are used with the same communicator.
120:       - The user can alternatively specify the local vector and matrix
121:         dimensions when more sophisticated partitioning is needed
122:         (replacing the PETSC_DECIDE argument in the VecSetSizes() statement
123:         below).
124:   */
125:   VecCreate(PETSC_COMM_WORLD,&u);
126:   VecSetSizes(u,PETSC_DECIDE,m*n);
127:   VecSetFromOptions(u);
128:   VecDuplicate(u,&b);
129:   VecDuplicate(b,&x);

131:   /*
132:      Set exact solution; then compute right-hand-side vector.
133:      By default we use an exact solution of a vector with all
134:      elements of 1.0;  Alternatively, using the runtime option
135:      -random_sol forms a solution vector with random components.
136:   */
137:   PetscOptionsGetBool(NULL,NULL,"-random_exact_sol",&flg,NULL);
138:   if (flg) {
139:     PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
140:     PetscRandomSetFromOptions(rctx);
141:     VecSetRandom(u,rctx);
142:     PetscRandomDestroy(&rctx);
143:   } else {
144:     VecSet(u,1.0);
145:   }
146:   MatMult(A,u,b);

148:   /*
149:      View the exact solution vector if desired
150:   */
151:   flg  = PETSC_FALSE;
152:   PetscOptionsGetBool(NULL,NULL,"-view_exact_sol",&flg,NULL);
153:   if (flg) {VecView(u,PETSC_VIEWER_STDOUT_WORLD);}

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

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

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

170:   /*
171:      Set linear solver defaults for this problem (optional).
172:      - By extracting the KSP and PC contexts from the KSP context,
173:        we can then directly call any KSP and PC routines to set
174:        various options.
175:      - The following two statements are optional; all of these
176:        parameters could alternatively be specified at runtime via
177:        KSPSetFromOptions().  All of these defaults can be
178:        overridden at runtime, as indicated below.
179:   */
180:   KSPSetTolerances(ksp,1.e-2/((m+1)*(n+1)),1.e-50,PETSC_DEFAULT,
181:                           PETSC_DEFAULT);

183:   /*
184:     Set runtime options, e.g.,
185:         -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
186:     These options will override those specified above as long as
187:     KSPSetFromOptions() is called _after_ any other customization
188:     routines.
189:   */
190:   KSPSetFromOptions(ksp);

192:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
193:                       Solve the linear system
194:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

196:   KSPSolve(ksp,b,x);

198:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
199:                       Check solution and clean up
200:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

202:   /*
203:      Check the error
204:   */
205:   VecAXPY(x,-1.0,u);
206:   VecNorm(x,NORM_2,&norm);
207:   KSPGetIterationNumber(ksp,&its);

209:   /*
210:      Print convergence information.  PetscPrintf() produces a single
211:      print statement from all processes that share a communicator.
212:      An alternative is PetscFPrintf(), which prints to a file.
213:   */
214:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g iterations %D\n",(double)norm,its);

216:   /*
217:      Free work space.  All PETSc objects should be destroyed when they
218:      are no longer needed.
219:   */
220:   KSPDestroy(&ksp);
221:   VecDestroy(&u);  VecDestroy(&x);
222:   VecDestroy(&b);  MatDestroy(&A);

224:   /*
225:      Always call PetscFinalize() before exiting a program.  This routine
226:        - finalizes the PETSc libraries as well as MPI
227:        - provides summary and diagnostic information if certain runtime
228:          options are chosen (e.g., -log_view).
229:   */
230:   PetscFinalize();
231:   return ierr;
232: }