Actual source code: ex16.c

petsc-3.10.5 2019-03-28
Report Typos and Errors

  2: /* Usage:  mpiexec ex16 [-help] [all PETSc options] */

  4: static char help[] = "Solves a sequence of linear systems with different right-hand-side vectors.\n\
  5: Input parameters include:\n\
  6:   -ntimes <ntimes>  : number of linear systems to solve\n\
  7:   -view_exact_sol   : write exact solution vector to stdout\n\
  8:   -m <mesh_x>       : number of mesh points in x-direction\n\
  9:   -n <mesh_n>       : number of mesh points in y-direction\n\n";

 11: /*T
 12:    Concepts: KSP^repeatedly solving linear systems;
 13:    Concepts: KSP^Laplacian, 2d
 14:    Concepts: Laplacian, 2d
 15:    Processors: n
 16: T*/

 18: /*
 19:   Include "petscksp.h" so that we can use KSP solvers.  Note that this file
 20:   automatically includes:
 21:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 22:      petscmat.h - matrices
 23:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 24:      petscviewer.h - viewers               petscpc.h  - preconditioners
 25: */
 26:  #include <petscksp.h>

 28: int main(int argc,char **args)
 29: {
 30:   Vec            x,b,u;  /* approx solution, RHS, exact solution */
 31:   Mat            A;        /* linear system matrix */
 32:   KSP            ksp;     /* linear solver context */
 33:   PetscReal      norm;     /* norm of solution error */
 35:   PetscInt       ntimes,i,j,k,Ii,J,Istart,Iend;
 36:   PetscInt       m   = 8,n = 7,its;
 37:   PetscBool      flg = PETSC_FALSE;
 38:   PetscScalar    v,one = 1.0,rhs;

 40:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 41:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 42:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);

 44:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 45:          Compute the matrix for use in solving a series of
 46:          linear systems of the form, A x_i = b_i, for i=1,2,...
 47:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 48:   /*
 49:      Create parallel matrix, specifying only its global dimensions.
 50:      When using MatCreate(), the matrix format can be specified at
 51:      runtime. Also, the parallel partitioning of the matrix is
 52:      determined by PETSc at runtime.
 53:   */
 54:   MatCreate(PETSC_COMM_WORLD,&A);
 55:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
 56:   MatSetFromOptions(A);
 57:   MatSetUp(A);

 59:   /*
 60:      Currently, all PETSc parallel matrix formats are partitioned by
 61:      contiguous chunks of rows across the processors.  Determine which
 62:      rows of the matrix are locally owned.
 63:   */
 64:   MatGetOwnershipRange(A,&Istart,&Iend);

 66:   /*
 67:      Set matrix elements for the 2-D, five-point stencil in parallel.
 68:       - Each processor needs to insert only elements that it owns
 69:         locally (but any non-local elements will be sent to the
 70:         appropriate processor during matrix assembly).
 71:       - Always specify global rows and columns of matrix entries.
 72:    */
 73:   for (Ii=Istart; Ii<Iend; Ii++) {
 74:     v = -1.0; i = Ii/n; j = Ii - i*n;
 75:     if (i>0)   {J = Ii - n; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
 76:     if (i<m-1) {J = Ii + n; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
 77:     if (j>0)   {J = Ii - 1; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
 78:     if (j<n-1) {J = Ii + 1; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
 79:     v = 4.0; MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES);
 80:   }

 82:   /*
 83:      Assemble matrix, using the 2-step process:
 84:        MatAssemblyBegin(), MatAssemblyEnd()
 85:      Computations can be done while messages are in transition
 86:      by placing code between these two statements.
 87:   */
 88:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 89:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 91:   /*
 92:      Create parallel vectors.
 93:       - When using VecCreate(), VecSetSizes() and VecSetFromOptions(),
 94:         we specify only the vector's global
 95:         dimension; the parallel partitioning is determined at runtime.
 96:       - When solving a linear system, the vectors and matrices MUST
 97:         be partitioned accordingly.  PETSc automatically generates
 98:         appropriately partitioned matrices and vectors when MatCreate()
 99:         and VecCreate() are used with the same communicator.
100:       - Note: We form 1 vector from scratch and then duplicate as needed.
101:   */
102:   VecCreate(PETSC_COMM_WORLD,&u);
103:   VecSetSizes(u,PETSC_DECIDE,m*n);
104:   VecSetFromOptions(u);
105:   VecDuplicate(u,&b);
106:   VecDuplicate(b,&x);

108:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109:                 Create the linear solver and set various options
110:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

112:   /*
113:      Create linear solver context
114:   */
115:   KSPCreate(PETSC_COMM_WORLD,&ksp);

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

123:   /*
124:     Set runtime options, e.g.,
125:         -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
126:     These options will override those specified above as long as
127:     KSPSetFromOptions() is called _after_ any other customization
128:     routines.
129:   */
130:   KSPSetFromOptions(ksp);

132:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133:        Solve several linear systems of the form  A x_i = b_i
134:        I.e., we retain the same matrix (A) for all systems, but
135:        change the right-hand-side vector (b_i) at each step.

137:        In this case, we simply call KSPSolve() multiple times.  The
138:        preconditioner setup operations (e.g., factorization for ILU)
139:        be done during the first call to KSPSolve() only; such operations
140:        will NOT be repeated for successive solves.
141:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

143:   ntimes = 2;
144:   PetscOptionsGetInt(NULL,NULL,"-ntimes",&ntimes,NULL);
145:   for (k=1; k<ntimes+1; k++) {

147:     /*
148:        Set exact solution; then compute right-hand-side vector.  We use
149:        an exact solution of a vector with all elements equal to 1.0*k.
150:     */
151:     rhs  = one * (PetscReal)k;
152:     VecSet(u,rhs);
153:     MatMult(A,u,b);

155:     /*
156:        View the exact solution vector if desired
157:     */
158:     PetscOptionsGetBool(NULL,NULL,"-view_exact_sol",&flg,NULL);
159:     if (flg) {VecView(u,PETSC_VIEWER_STDOUT_WORLD);}

161:     KSPSolve(ksp,b,x);

163:     /*
164:        Check the error
165:     */
166:     VecAXPY(x,-1.0,u);
167:     VecNorm(x,NORM_2,&norm);
168:     KSPGetIterationNumber(ksp,&its);
169:     /*
170:        Print convergence information.  PetscPrintf() produces a single
171:        print statement from all processes that share a communicator.
172:     */
173:     PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g System %D: iterations %D\n",(double)norm,k,its);
174:   }

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

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_view).
192:   */
193:   PetscFinalize();
194:   return ierr;
195: }


198: /*TEST

200:    test:
201:       nsize: 2
202:       args: -ntimes 4 -ksp_gmres_cgs_refinement_type refine_always

204: TEST*/