Actual source code: ex2.c

  1: /*$Id: ex2.c,v 1.26 2001/08/07 21:30:50 bsmith Exp $*/

  3: static char help[] = "Demonstrates running several independent tasks in PETSc.nn";

  5: /*T
  6:    Concepts: SLES^solving linear equations
  7:    Processors: n

  9:    Comments: Demonstrates how to use PetscSetCommWorld() to tell a subset of
 10:              processors (in this case each individual processor) to run 
 11:              as if it was all the processors that PETSc is using. ADVANCED
 12:              example, not for beginning PETSc users.

 14: T*/

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

 26: EXTERN int slesex(int,char**);

 28: int main(int argc,char **argv)
 29: {
 30:     MPI_Init(&argc,&argv);
 31:     slesex(argc,argv);
 32:     MPI_Finalize();
 33:     return 0;
 34: }

 36: int slesex(int argc,char **args)
 37: {
 38:   Vec         x,b,u;      /* approx solution, RHS, exact solution */
 39:   Mat         A;            /* linear system matrix */
 40:   SLES        sles;         /* linear solver context */
 41:   PC          pc;           /* preconditioner context */
 42:   KSP         ksp;          /* Krylov subspace method context */
 43:   PetscReal   norm;         /* norm of solution error */
 44:   int         i,j,I,J,Istart,Iend,ierr,m = 8,n = 7,its;
 45:   PetscScalar v,one = 1.0,none = -1.0;

 47:   PetscSetCommWorld(PETSC_COMM_SELF);
 48:   PetscInitialize(&argc,&args,(char *)0,help);

 50:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
 51:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);

 53:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 54:          Compute the matrix and right-hand-side vector that define
 55:          the linear system, Ax = b.
 56:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 57:   /* 
 58:      Create parallel matrix, specifying only its global dimensions.
 59:      When using MatCreate(), the matrix format can be specified at
 60:      runtime. Also, the parallel partitioning of the matrix is
 61:      determined by PETSc at runtime.
 62:   */
 63:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,&A);
 64:   MatSetFromOptions(A);

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

 73:   /* 
 74:      Set matrix elements for the 2-D, five-point stencil in parallel.
 75:       - Each processor needs to insert only elements that it owns
 76:         locally (but any non-local elements will be sent to the
 77:         appropriate processor during matrix assembly). 
 78:       - Always specify global row and columns of matrix entries.
 79:    */
 80:   for (I=Istart; I<Iend; I++) {
 81:     v = -1.0; i = I/n; j = I - i*n;
 82:     if (i>0)   {J = I - n; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
 83:     if (i<m-1) {J = I + n; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
 84:     if (j>0)   {J = I - 1; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
 85:     if (j<n-1) {J = I + 1; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
 86:     v = 4.0; MatSetValues(A,1,&I,1,&I,&v,INSERT_VALUES);
 87:   }

 89:   /* 
 90:      Assemble matrix, using the 2-step process:
 91:        MatAssemblyBegin(), MatAssemblyEnd()
 92:      Computations can be done while messages are in transition,
 93:      by placing code between these two statements.
 94:   */
 95:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 96:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 98:   /* 
 99:      Create parallel vectors.
100:       - When using VecCreate() VecSetSizes() and VecSetFromOptions(),
101:         we specify only the vector's global
102:         dimension; the parallel partitioning is determined at runtime. 
103:       - Note: We form 1 vector from scratch and then duplicate as needed.
104:   */
105:   VecCreate(PETSC_COMM_WORLD,&u);
106:   VecSetSizes(u,PETSC_DECIDE,m*n);
107:   VecSetFromOptions(u);
108:   VecDuplicate(u,&b);
109:   VecDuplicate(b,&x);

111:   /* 
112:      Set exact solution; then compute right-hand-side vector.
113:   */
114:   VecSet(&one,u);
115:   MatMult(A,u,b);

117:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
118:                 Create the linear solver and set various options
119:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

121:   /* 
122:      Create linear solver context
123:   */
124:   SLESCreate(PETSC_COMM_WORLD,&sles);

126:   /* 
127:      Set operators. Here the matrix that defines the linear system
128:      also serves as the preconditioning matrix.
129:   */
130:   SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN);

132:   /* 
133:      Set linear solver defaults for this problem (optional).
134:      - By extracting the KSP and PC contexts from the SLES context,
135:        we can then directly directly call any KSP and PC routines
136:        to set various options.
137:      - The following four statements are optional; all of these
138:        parameters could alternatively be specified at runtime via
139:        SLESSetFromOptions();
140:   */
141:   SLESGetKSP(sles,&ksp);
142:   SLESGetPC(sles,&pc);
143:   PCSetType(pc,PCJACOBI);
144:   KSPSetTolerances(ksp,1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);

146:   /* 
147:     Set runtime options, e.g.,
148:         -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
149:     These options will override those specified above as long as
150:     SLESSetFromOptions() is called _after_ any other customization
151:     routines.
152:   */
153:   SLESSetFromOptions(sles);

155:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
156:                       Solve the linear system
157:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

161:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
162:                       Check solution and clean up
163:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

165:   /* 
166:      Check the error
167:   */
168:   VecAXPY(&none,u,x);
169:   VecNorm(x,NORM_2,&norm);
170:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A iterations %dn",norm,its);

172:   /* 
173:      Free work space.  All PETSc objects should be destroyed when they
174:      are no longer needed.
175:   */
176:   SLESDestroy(sles);
177:   VecDestroy(u);  VecDestroy(x);
178:   VecDestroy(b);  MatDestroy(A);
179:   PetscFinalize();
180:   return 0;
181: }