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