Actual source code: ex1.c


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

  4: /*T
  5:    Concepts: KSP^solving a system of linear equations
  6:    Processors: 1
  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              petscpc.h  - preconditioners
 14:      petscis.h     - index sets
 15:      petscviewer.h - viewers

 17:   Note:  The corresponding parallel example is ex23.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;         /* norm of solution error */
 28:   PetscInt       i,n = 10,col[3],its;
 29:   PetscMPIInt    size;
 30:   PetscScalar    value[3];

 32:   PetscInitialize(&argc,&args,(char*)0,help);
 33:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 35:   PetscOptionsGetInt(NULL,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.
 45:   */
 46:   VecCreate(PETSC_COMM_WORLD,&x);
 47:   PetscObjectSetName((PetscObject) x, "Solution");
 48:   VecSetSizes(x,PETSC_DECIDE,n);
 49:   VecSetFromOptions(x);
 50:   VecDuplicate(x,&b);
 51:   VecDuplicate(x,&u);

 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.
 60:   */
 61:   MatCreate(PETSC_COMM_WORLD,&A);
 62:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 63:   MatSetFromOptions(A);
 64:   MatSetUp(A);

 66:   /*
 67:      Assemble matrix
 68:   */
 69:   value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
 70:   for (i=1; i<n-1; i++) {
 71:     col[0] = i-1; col[1] = i; col[2] = i+1;
 72:     MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 73:   }
 74:   i    = n - 1; col[0] = n - 2; col[1] = n - 1;
 75:   MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 76:   i    = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
 77:   MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 78:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 79:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 81:   /*
 82:      Set exact solution; then compute right-hand-side vector.
 83:   */
 84:   VecSet(u,1.0);
 85:   MatMult(A,u,b);

 87:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 88:                 Create the linear solver and set various options
 89:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 90:   KSPCreate(PETSC_COMM_WORLD,&ksp);

 92:   /*
 93:      Set operators. Here the matrix that defines the linear system
 94:      also serves as the matrix that defines the preconditioner.
 95:   */
 96:   KSPSetOperators(ksp,A,A);

 98:   /*
 99:      Set linear solver defaults for this problem (optional).
100:      - By extracting the KSP and PC contexts from the KSP context,
101:        we can then directly call any KSP and PC routines to set
102:        various options.
103:      - The following four statements are optional; all of these
104:        parameters could alternatively be specified at runtime via
105:        KSPSetFromOptions();
106:   */
107:   KSPGetPC(ksp,&pc);
108:   PCSetType(pc,PCJACOBI);
109:   KSPSetTolerances(ksp,1.e-5,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);

111:   /*
112:     Set runtime options, e.g.,
113:         -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
114:     These options will override those specified above as long as
115:     KSPSetFromOptions() is called _after_ any other customization
116:     routines.
117:   */
118:   KSPSetFromOptions(ksp);

120:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121:                       Solve the linear system
122:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
123:   KSPSolve(ksp,b,x);

125:   /*
126:      View solver info; we could instead use the option -ksp_view to
127:      print this info to the screen at the conclusion of KSPSolve().
128:   */
129:   KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD);

131:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132:                       Check the solution and clean up
133:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
134:   VecAXPY(x,-1.0,u);
135:   VecNorm(x,NORM_2,&norm);
136:   KSPGetIterationNumber(ksp,&its);
137:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %D\n",(double)norm,its);

139:   /*
140:      Free work space.  All PETSc objects should be destroyed when they
141:      are no longer needed.
142:   */
143:   VecDestroy(&x)); PetscCall(VecDestroy(&u);
144:   VecDestroy(&b)); PetscCall(MatDestroy(&A);
145:   KSPDestroy(&ksp);

147:   /*
148:      Always call PetscFinalize() before exiting a program.  This routine
149:        - finalizes the PETSc libraries as well as MPI
150:        - provides summary and diagnostic information if certain runtime
151:          options are chosen (e.g., -log_view).
152:   */
153:   PetscFinalize();
154:   return 0;
155: }

157: /*TEST

159:    test:
160:       args: -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always

162:    test:
163:       suffix: 2
164:       args: -pc_type sor -pc_sor_symmetric -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always

166:    test:
167:       suffix: 2_aijcusparse
168:       requires: cuda
169:       args: -pc_type sor -pc_sor_symmetric -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always -mat_type aijcusparse -vec_type cuda

171:    test:
172:       suffix: 3
173:       args: -pc_type eisenstat -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always

175:    test:
176:       suffix: 3_aijcusparse
177:       requires: cuda
178:       args: -pc_type eisenstat -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always -mat_type aijcusparse -vec_type cuda

180:    test:
181:       suffix: aijcusparse
182:       requires: cuda
183:       args: -ksp_monitor_short -ksp_gmres_cgs_refinement_type refine_always -mat_type aijcusparse -vec_type cuda
184:       output_file: output/ex1_1_aijcusparse.out

186: TEST*/