Actual source code: ex23.c


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

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

 12:   Note:  The corresponding uniprocessor example is ex1.c
 13: */
 14: #include <petscksp.h>

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

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

 29:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 30:          Compute the matrix and right-hand-side vector that define
 31:          the linear system, Ax = b.
 32:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

 46:   /* Identify the starting and ending mesh points on each
 47:      processor for the interior part of the mesh. We let PETSc decide
 48:      above. */

 50:   VecGetOwnershipRange(x,&rstart,&rend);
 51:   VecGetLocalSize(x,&nlocal);

 53:   /*
 54:      Create matrix.  When using MatCreate(), the matrix format can
 55:      be specified at runtime.

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

 61:      We pass in nlocal as the "local" size of the matrix to force it
 62:      to have the same parallel layout as the vector created above.
 63:   */
 64:   MatCreate(PETSC_COMM_WORLD,&A);
 65:   MatSetSizes(A,nlocal,nlocal,n,n);
 66:   MatSetFromOptions(A);
 67:   MatSetUp(A);

 69:   /*
 70:      Assemble matrix.

 72:      The linear system is distributed across the processors by
 73:      chunks of contiguous rows, which correspond to contiguous
 74:      sections of the mesh on which the problem is discretized.
 75:      For matrix assembly, each processor contributes entries for
 76:      the part that it owns locally.
 77:   */

 79:   if (!rstart) {
 80:     rstart = 1;
 81:     i      = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
 82:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 83:   }
 84:   if (rend == n) {
 85:     rend = n-1;
 86:     i    = n-1; col[0] = n-2; col[1] = n-1; value[0] = -1.0; value[1] = 2.0;
 87:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 88:   }

 90:   /* Set entries corresponding to the mesh interior */
 91:   value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
 92:   for (i=rstart; i<rend; i++) {
 93:     col[0] = i-1; col[1] = i; col[2] = i+1;
 94:     MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 95:   }

 97:   /* Assemble the matrix */
 98:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 99:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

101:   /*
102:      Set exact solution; then compute right-hand-side vector.
103:   */
104:   VecSet(u,one);
105:   MatMult(A,u,b);

107:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108:                 Create the linear solver and set various options
109:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110:   /*
111:      Create linear solver context
112:   */
113:   KSPCreate(PETSC_COMM_WORLD,&ksp);

115:   /*
116:      Set operators. Here the matrix that defines the linear system
117:      also serves as the preconditioning matrix.
118:   */
119:   KSPSetOperators(ksp,A,A);

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

134:   /*
135:     Set runtime options, e.g.,
136:         -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
137:     These options will override those specified above as long as
138:     KSPSetFromOptions() is called _after_ any other customization
139:     routines.
140:   */
141:   KSPSetFromOptions(ksp);

143:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144:                       Solve the linear system
145:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
146:   /*
147:      Solve linear system
148:   */
149:   KSPSolve(ksp,b,x);

151:   /*
152:      View solver info; we could instead use the option -ksp_view to
153:      print this info to the screen at the conclusion of KSPSolve().
154:   */
155:   KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD);

157:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158:                       Check solution and clean up
159:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
160:   /*
161:      Check the error
162:   */
163:   VecAXPY(x,-1.0,u);
164:   VecNorm(x,NORM_2,&norm);
165:   KSPGetIterationNumber(ksp,&its);
166:   if (norm > tol) {
167:     PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %D\n",(double)norm,its);
168:   }

170:   /*
171:      Free work space.  All PETSc objects should be destroyed when they
172:      are no longer needed.
173:   */
174:   VecDestroy(&x)); PetscCall(VecDestroy(&u);
175:   VecDestroy(&b)); PetscCall(MatDestroy(&A);
176:   KSPDestroy(&ksp);

178:   /*
179:      Always call PetscFinalize() before exiting a program.  This routine
180:        - finalizes the PETSc libraries as well as MPI
181:        - provides summary and diagnostic information if certain runtime
182:          options are chosen (e.g., -log_view).
183:   */
184:   PetscFinalize();
185:   return 0;
186: }

188: /*TEST

190:    build:
191:       requires: !complex !single

193:    test:
194:       args: -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always

196:    test:
197:       suffix: 2
198:       nsize: 3
199:       args: -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always

201:    test:
202:       suffix: 3
203:       nsize: 2
204:       args: -ksp_monitor_short -ksp_rtol 1e-6 -ksp_type pipefgmres

206: TEST*/