xref: /petsc/src/ksp/pc/tests/ex8f.F90 (revision 4820e4ea99a084ae862a8c395f732bc7c0e1a6d0)
1d8606c27SBarry Smith!
2d8606c27SBarry Smith!   Tests PCMGSetResidual
3d8606c27SBarry Smith!
4d8606c27SBarry Smith! -----------------------------------------------------------------------
5d8606c27SBarry Smith
6d8606c27SBarry Smithprogram main
7d8606c27SBarry Smith#include <petsc/finclude/petscksp.h>
8d8606c27SBarry Smith  use petscksp
9d8606c27SBarry Smith  implicit none
10d8606c27SBarry Smith
11d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
12d8606c27SBarry Smith!                   Variable declarations
13d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
14d8606c27SBarry Smith!
15d8606c27SBarry Smith!  Variables:
16d8606c27SBarry Smith!     ksp     - linear solver context
17dd8e379bSPierre Jolivet!     x, b, u  - approx solution, right-hand side, exact solution vectors
18d8606c27SBarry Smith!     A        - matrix that defines linear system
19d8606c27SBarry Smith!     its      - iterations for convergence
20d8606c27SBarry Smith!     norm     - norm of error in solution
21d8606c27SBarry Smith!     rctx     - random number context
22d8606c27SBarry Smith!
23d8606c27SBarry Smith
24d8606c27SBarry Smith  Mat A
25d8606c27SBarry Smith  Vec x, b, u
26d8606c27SBarry Smith  PC pc
27d8606c27SBarry Smith  PetscInt n, dim, istart, iend
28d8606c27SBarry Smith  PetscInt i, j, jj, ii, one, zero
29d8606c27SBarry Smith  PetscErrorCode ierr
30d8606c27SBarry Smith  PetscScalar v
31d8606c27SBarry Smith  external MyResidual
32d8606c27SBarry Smith  PetscScalar pfive
33d8606c27SBarry Smith  KSP ksp
34d8606c27SBarry Smith
35d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36d8606c27SBarry Smith!                 Beginning of program
37d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
38d8606c27SBarry Smith
39d8606c27SBarry Smith  PetscCallA(PetscInitialize(ierr))
40d8606c27SBarry Smith  pfive = .5
41d8606c27SBarry Smith  n = 6
42d8606c27SBarry Smith  dim = n*n
43d8606c27SBarry Smith  one = 1
44d8606c27SBarry Smith  zero = 0
45d8606c27SBarry Smith
46d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47d8606c27SBarry Smith!      Compute the matrix and right-hand-side vector that define
48d8606c27SBarry Smith!      the linear system, Ax = b.
49d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
50d8606c27SBarry Smith
51d8606c27SBarry Smith!  Create parallel matrix, specifying only its global dimensions.
52d8606c27SBarry Smith!  When using MatCreate(), the matrix format can be specified at
53d8606c27SBarry Smith!  runtime. Also, the parallel partitioning of the matrix is
54d8606c27SBarry Smith!  determined by PETSc at runtime.
55d8606c27SBarry Smith
56d8606c27SBarry Smith  PetscCallA(MatCreate(PETSC_COMM_WORLD, A, ierr))
57d8606c27SBarry Smith  PetscCallA(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, dim, dim, ierr))
58d8606c27SBarry Smith  PetscCallA(MatSetFromOptions(A, ierr))
59d8606c27SBarry Smith  PetscCallA(MatSetUp(A, ierr))
60d8606c27SBarry Smith
61d8606c27SBarry Smith!  Currently, all PETSc parallel matrix formats are partitioned by
62d8606c27SBarry Smith!  contiguous chunks of rows across the processors.  Determine which
63d8606c27SBarry Smith!  rows of the matrix are locally owned.
64d8606c27SBarry Smith
65d8606c27SBarry Smith  PetscCallA(MatGetOwnershipRange(A, Istart, Iend, ierr))
66d8606c27SBarry Smith
67d8606c27SBarry Smith!  Set matrix elements in parallel.
68d8606c27SBarry Smith!   - Each processor needs to insert only elements that it owns
69d8606c27SBarry Smith!     locally (but any non-local elements will be sent to the
70d8606c27SBarry Smith!     appropriate processor during matrix assembly).
71d8606c27SBarry Smith!   - Always specify global rows and columns of matrix entries.
72d8606c27SBarry Smith
73d8606c27SBarry Smith  do 10, II = Istart, Iend - 1
74d8606c27SBarry Smith    v = -1.0
75d8606c27SBarry Smith    i = II/n
76d8606c27SBarry Smith    j = II - i*n
77*4820e4eaSBarry Smith    if (i > 0) then
78d8606c27SBarry Smith      JJ = II - n
795d83a8b1SBarry Smith      PetscCallA(MatSetValues(A, one, [II], one, [JJ], [v], ADD_VALUES, ierr))
80d8606c27SBarry Smith    end if
81*4820e4eaSBarry Smith    if (i < n - 1) then
82d8606c27SBarry Smith      JJ = II + n
835d83a8b1SBarry Smith      PetscCallA(MatSetValues(A, one, [II], one, [JJ], [v], ADD_VALUES, ierr))
84d8606c27SBarry Smith    end if
85*4820e4eaSBarry Smith    if (j > 0) then
86d8606c27SBarry Smith      JJ = II - 1
875d83a8b1SBarry Smith      PetscCallA(MatSetValues(A, one, [II], one, [JJ], [v], ADD_VALUES, ierr))
88d8606c27SBarry Smith    end if
89*4820e4eaSBarry Smith    if (j < n - 1) then
90d8606c27SBarry Smith      JJ = II + 1
915d83a8b1SBarry Smith      PetscCallA(MatSetValues(A, one, [II], one, [JJ], [v], ADD_VALUES, ierr))
92d8606c27SBarry Smith    end if
93d8606c27SBarry Smith    v = 4.0
945d83a8b1SBarry Smith    PetscCallA(MatSetValues(A, one, [II], one, [II], [v], ADD_VALUES, ierr))
95d8606c27SBarry Smith10  continue
96d8606c27SBarry Smith
97d8606c27SBarry Smith!  Assemble matrix, using the 2-step process:
98d8606c27SBarry Smith!       MatAssemblyBegin(), MatAssemblyEnd()
99d8606c27SBarry Smith!  Computations can be done while messages are in transition
100d8606c27SBarry Smith!  by placing code between these two statements.
101d8606c27SBarry Smith
102d8606c27SBarry Smith    PetscCallA(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY, ierr))
103d8606c27SBarry Smith    PetscCallA(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY, ierr))
104d8606c27SBarry Smith
105d8606c27SBarry Smith!  Create parallel vectors.
106d8606c27SBarry Smith!   - Here, the parallel partitioning of the vector is determined by
107d8606c27SBarry Smith!     PETSc at runtime.  We could also specify the local dimensions
108d8606c27SBarry Smith!     if desired.
109d8606c27SBarry Smith!   - Note: We form 1 vector from scratch and then duplicate as needed.
110d8606c27SBarry Smith
111d8606c27SBarry Smith    PetscCallA(VecCreate(PETSC_COMM_WORLD, u, ierr))
112d8606c27SBarry Smith    PetscCallA(VecSetSizes(u, PETSC_DECIDE, dim, ierr))
113d8606c27SBarry Smith    PetscCallA(VecSetFromOptions(u, ierr))
114d8606c27SBarry Smith    PetscCallA(VecDuplicate(u, b, ierr))
115d8606c27SBarry Smith    PetscCallA(VecDuplicate(b, x, ierr))
116d8606c27SBarry Smith
117d8606c27SBarry Smith!  Set exact solution; then compute right-hand-side vector.
118d8606c27SBarry Smith
119d8606c27SBarry Smith    PetscCallA(VecSet(u, pfive, ierr))
120d8606c27SBarry Smith    PetscCallA(MatMult(A, u, b, ierr))
121d8606c27SBarry Smith
122d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123d8606c27SBarry Smith!         Create the linear solver and set various options
124d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125d8606c27SBarry Smith
126d8606c27SBarry Smith!  Create linear solver context
127d8606c27SBarry Smith
128d8606c27SBarry Smith    PetscCallA(KSPCreate(PETSC_COMM_WORLD, ksp, ierr))
129d8606c27SBarry Smith    PetscCallA(KSPGetPC(ksp, pc, ierr))
130d8606c27SBarry Smith    PetscCallA(PCSetType(pc, PCMG, ierr))
131d8606c27SBarry Smith    PetscCallA(PCMGSetLevels(pc, one, PETSC_NULL_MPI_COMM, ierr))
132d8606c27SBarry Smith    PetscCallA(PCMGSetResidual(pc, zero, MyResidual, A, ierr))
133d8606c27SBarry Smith
134d8606c27SBarry Smith!  Set operators. Here the matrix that defines the linear system
1357addb90fSBarry Smith!  also serves as the matrix used to construct the preconditioner.
136d8606c27SBarry Smith
137d8606c27SBarry Smith    PetscCallA(KSPSetOperators(ksp, A, A, ierr))
138d8606c27SBarry Smith
139d8606c27SBarry Smith    PetscCallA(KSPDestroy(ksp, ierr))
140d8606c27SBarry Smith    PetscCallA(VecDestroy(u, ierr))
141d8606c27SBarry Smith    PetscCallA(VecDestroy(x, ierr))
142d8606c27SBarry Smith    PetscCallA(VecDestroy(b, ierr))
143d8606c27SBarry Smith    PetscCallA(MatDestroy(A, ierr))
144d8606c27SBarry Smith
145d8606c27SBarry Smith    PetscCallA(PetscFinalize(ierr))
146d8606c27SBarry Smith  end
147d8606c27SBarry Smith
148d8606c27SBarry Smith  subroutine MyResidual(A, b, x, r, ierr)
149d8606c27SBarry Smith    use petscksp
150d8606c27SBarry Smith    implicit none
151d8606c27SBarry Smith    Mat A
152d8606c27SBarry Smith    Vec b, x, r
153d8606c27SBarry Smith    integer ierr
154d8606c27SBarry Smith  end
155d8606c27SBarry Smith
156d8606c27SBarry Smith!/*TEST
157d8606c27SBarry Smith!
158d8606c27SBarry Smith!   test:
159d8606c27SBarry Smith!      nsize: 1
1603886731fSPierre Jolivet!      output_file: output/empty.out
161d8606c27SBarry Smith!TEST*/
162