Actual source code: ex23.c

  2: /* Program usage:  mpirun ex23 [-help] [all PETSc options] */

  4: static char help[] = "Solves a tridiagonal linear system.\n\n";

  6: /*T
  7:    Concepts: KSP^basic parallel example;
  8:    Processors: n
  9: T*/

 11: /* 
 12:   Include "petscksp.h" so that we can use KSP solvers.  Note that this file
 13:   automatically includes:
 14:      petsc.h       - base PETSc routines   petscvec.h - vectors
 15:      petscsys.h    - system routines       petscmat.h - matrices
 16:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 17:      petscviewer.h - viewers               petscpc.h  - preconditioners

 19:   Note:  The corresponding uniprocessor example is ex1.c
 20: */
 21:  #include petscksp.h

 25: int main(int argc,char **args)
 26: {
 27:   Vec            x, b, u;      /* approx solution, RHS, exact solution */
 28:   Mat            A;            /* linear system matrix */
 29:   KSP            ksp;         /* linear solver context */
 30:   PC             pc;           /* preconditioner context */
 31:   PetscReal      norm;         /* norm of solution error */
 33:   PetscInt       i,n = 10,col[3],its,rstart,rend,nlocal;
 34:   PetscScalar    neg_one = -1.0,one = 1.0,value[3];

 36:   PetscInitialize(&argc,&args,(char *)0,help);
 37:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);

 39:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 40:          Compute the matrix and right-hand-side vector that define
 41:          the linear system, Ax = b.
 42:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 44:   /* 
 45:      Create vectors.  Note that we form 1 vector from scratch and
 46:      then duplicate as needed. For this simple case let PETSc decide how
 47:      many elements of the vector are stored on each processor. The second
 48:      argument to VecSetSizes() below causes PETSc to decide.
 49:   */
 50:   VecCreate(PETSC_COMM_WORLD,&x);
 51:   VecSetSizes(x,PETSC_DECIDE,n);
 52:   VecSetFromOptions(x);
 53:   VecDuplicate(x,&b);
 54:   VecDuplicate(x,&u);

 56:   /* Identify the starting and ending mesh points on each
 57:      processor for the interior part of the mesh. We let PETSc decide
 58:      above. */

 60:   VecGetOwnershipRange(x,&rstart,&rend);
 61:   VecGetLocalSize(x,&nlocal);

 63:   /* 
 64:      Create matrix.  When using MatCreate(), the matrix format can
 65:      be specified at runtime.

 67:      Performance tuning note:  For problems of substantial size,
 68:      preallocation of matrix memory is crucial for attaining good 
 69:      performance.  Since preallocation is not possible via the generic
 70:      matrix creation routine MatCreate(), we recommend for practical 
 71:      problems instead to use the creation routine for a particular matrix
 72:      format, e.g.,
 73:          MatCreateMPIAIJ() - sequential AIJ (compressed sparse row)
 74:          MatCreateMPIBAIJ() - block AIJ
 75:      See the matrix chapter of the users manual for details.

 77:      We pass in nlocal as the "local" size of the matrix to force it
 78:      to have the same parallel layout as the vector created above.
 79:   */
 80:   MatCreate(PETSC_COMM_WORLD,nlocal,nlocal,n,n,&A);
 81:   MatSetFromOptions(A);

 83:   /* 
 84:      Assemble matrix.  

 86:      The linear system is distributed across the processors by 
 87:      chunks of contiguous rows, which correspond to contiguous
 88:      sections of the mesh on which the problem is discretized.  
 89:      For matrix assembly, each processor contributes entries for
 90:      the part that it owns locally.
 91:   */


 94:   if (!rstart) {
 95:     rstart = 1;
 96:     i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
 97:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 98:   }
 99:   if (rend == n) {
100:     rend = n-1;
101:     i = n-1; col[0] = n-2; col[1] = n-1; value[0] = -1.0; value[1] = 2.0;
102:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
103:   }

105:   /* Set entries corresponding to the mesh interior */
106:   value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
107:   for (i=rstart; i<rend; i++) {
108:     col[0] = i-1; col[1] = i; col[2] = i+1;
109:     MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
110:   }

112:   /* Assemble the matrix */
113:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
114:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

116:   /* 
117:      Set exact solution; then compute right-hand-side vector.
118:   */
119:   VecSet(&one,u);
120:   MatMult(A,u,b);

122:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
123:                 Create the linear solver and set various options
124:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
125:   /* 
126:      Create linear solver context
127:   */
128:   KSPCreate(PETSC_COMM_WORLD,&ksp);

130:   /* 
131:      Set operators. Here the matrix that defines the linear system
132:      also serves as the preconditioning matrix.
133:   */
134:   KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN);

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

149:   /* 
150:     Set runtime options, e.g.,
151:         -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
152:     These options will override those specified above as long as
153:     KSPSetFromOptions() is called _after_ any other customization
154:     routines.
155:   */
156:   KSPSetFromOptions(ksp);
157: 
158:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
159:                       Solve the linear system
160:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
161:   /* 
162:      Solve linear system
163:   */
164:   KSPSolve(ksp,b,x);

166:   /* 
167:      View solver info; we could instead use the option -ksp_view to
168:      print this info to the screen at the conclusion of KSPSolve().
169:   */
170:   KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD);

172:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
173:                       Check solution and clean up
174:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
175:   /* 
176:      Check the error
177:   */
178:   VecAXPY(&neg_one,u,x);
179:   VecNorm(x,NORM_2,&norm);
180:   KSPGetIterationNumber(ksp,&its);
181:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A, Iterations %D\n",norm,its);
182:   /* 
183:      Free work space.  All PETSc objects should be destroyed when they
184:      are no longer needed.
185:   */
186:   VecDestroy(x); VecDestroy(u);
187:   VecDestroy(b); MatDestroy(A);
188:   KSPDestroy(ksp);

190:   /*
191:      Always call PetscFinalize() before exiting a program.  This routine
192:        - finalizes the PETSc libraries as well as MPI
193:        - provides summary and diagnostic information if certain runtime
194:          options are chosen (e.g., -log_summary).
195:   */
196:   PetscFinalize();
197:   return 0;
198: }