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