Actual source code: ex16.c


  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_y>       : 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 */
 34:   PetscInt       ntimes,i,j,k,Ii,J,Istart,Iend;
 35:   PetscInt       m   = 8,n = 7,its;
 36:   PetscBool      flg = PETSC_FALSE;
 37:   PetscScalar    v,one = 1.0,rhs;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

160:     KSPSolve(ksp,b,x);

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

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

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

196: /*TEST

198:    test:
199:       nsize: 2
200:       args: -ntimes 4 -ksp_gmres_cgs_refinement_type refine_always

202: TEST*/