Actual source code: ex23.c

petsc-3.8.4 2018-03-24
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>

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

 32:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 33:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);

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

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

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

 56:   VecGetOwnershipRange(x,&rstart,&rend);
 57:   VecGetLocalSize(x,&nlocal);

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

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

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

 75:   /*
 76:      Assemble matrix.

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


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

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

104:   /* Assemble the matrix */
105:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
106:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

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

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

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

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

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

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

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

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

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

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