Actual source code: ex23.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: static char help[] = "Solves a tridiagonal linear system.\n\n";

  4: /*T
  5:    Concepts: KSP^basic parallel example;
  6:    Processors: n
  7: T*/

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

 17:   Note:  The corresponding uniprocessor example is ex1.c
 18: */
 19: #include <petscksp.h>

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

 34:   PetscInitialize(&argc,&args,(char*)0,help);
 35:   PetscOptionsGetInt(NULL,"-n",&n,NULL);

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

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

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

 58:   VecGetOwnershipRange(x,&rstart,&rend);
 59:   VecGetLocalSize(x,&nlocal);

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

 65:      Performance tuning note:  For problems of substantial size,
 66:      preallocation of matrix memory is crucial for attaining good
 67:      performance. See the matrix chapter of the users manual for details.

 69:      We pass in nlocal as the "local" size of the matrix to force it
 70:      to have the same parallel layout as the vector created above.
 71:   */
 72:   MatCreate(PETSC_COMM_WORLD,&A);
 73:   MatSetSizes(A,nlocal,nlocal,n,n);
 74:   MatSetFromOptions(A);
 75:   MatSetUp(A);

 77:   /*
 78:      Assemble matrix.

 80:      The linear system is distributed across the processors by
 81:      chunks of contiguous rows, which correspond to contiguous
 82:      sections of the mesh on which the problem is discretized.
 83:      For matrix assembly, each processor contributes entries for
 84:      the part that it owns locally.
 85:   */


 88:   if (!rstart) {
 89:     rstart = 1;
 90:     i      = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
 91:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 92:   }
 93:   if (rend == n) {
 94:     rend = n-1;
 95:     i    = n-1; col[0] = n-2; col[1] = n-1; value[0] = -1.0; value[1] = 2.0;
 96:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 97:   }

 99:   /* Set entries corresponding to the mesh interior */
100:   value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
101:   for (i=rstart; i<rend; i++) {
102:     col[0] = i-1; col[1] = i; col[2] = i+1;
103:     MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
104:   }

106:   /* Assemble the matrix */
107:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
108:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

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

116:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117:                 Create the linear solver and set various options
118:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119:   /*
120:      Create linear solver context
121:   */
122:   KSPCreate(PETSC_COMM_WORLD,&ksp);

124:   /*
125:      Set operators. Here the matrix that defines the linear system
126:      also serves as the preconditioning matrix.
127:   */
128:   KSPSetOperators(ksp,A,A);

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

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

152:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153:                       Solve the linear system
154:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155:   /*
156:      Solve linear system
157:   */
158:   KSPSolve(ksp,b,x);

160:   /*
161:      View solver info; we could instead use the option -ksp_view to
162:      print this info to the screen at the conclusion of KSPSolve().
163:   */
164:   KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD);

166:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167:                       Check solution and clean up
168:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
169:   /*
170:      Check the error
171:   */
172:   VecAXPY(x,neg_one,u);
173:   VecNorm(x,NORM_2,&norm);
174:   KSPGetIterationNumber(ksp,&its);
175:   if (norm > tol) {
176:     PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %D\n",(double)norm,its);
177:   }

179:   /*
180:      Free work space.  All PETSc objects should be destroyed when they
181:      are no longer needed.
182:   */
183:   VecDestroy(&x); VecDestroy(&u);
184:   VecDestroy(&b); MatDestroy(&A);
185:   KSPDestroy(&ksp);

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