Octopus
poisson_multigrid.F90
Go to the documentation of this file.
1!! Copyright (C) 2005-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch, X. Andrade
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
23 use debug_oct_m
25 use global_oct_m
26 use, intrinsic :: iso_fortran_env
28 use mesh_oct_m
35 use parser_oct_m
38 use space_oct_m
40
41 implicit none
42
43 private
44 public :: &
49
51 private
52
53 type(poisson_corr_t) :: corrector
54 type(multigrid_t) :: mgrid ! multigrid object
55 type(mg_solver_t) :: mg_solver
56
57 integer :: multigrid_cycle
58
59 end type poisson_mg_solver_t
60
61contains
62
63 ! ---------------------------------------------------------
64 subroutine poisson_multigrid_init(this, namespace, space, mesh, der, stencil, mc, ml, thr)
65 type(poisson_mg_solver_t), intent(out) :: this
66 type(namespace_t), intent(in) :: namespace
67 class(space_t), intent(in) :: space
68 type(mesh_t), intent(inout) :: mesh
69 type(derivatives_t), intent(in) :: der
70 type(stencil_t), intent(in) :: stencil
71 type(multicomm_t), intent(in) :: mc
72 integer, intent(in) :: ml
73 real(real64), intent(in) :: thr
74
75 integer :: i
76
78
79 call poisson_corrections_init(this%corrector, namespace, space, ml, mesh)
80
81 call multigrid_init(this%mgrid, namespace, space, mesh, der, stencil, mc)
82
83 ! For the multigrid solver to work, we need to set a pointer from one operator
84 ! to the corresponding one on the coarser grid.
85 do i = 1, this%mgrid%n_levels
86 this%mgrid%level(i - 1)%der%lapl%coarser => this%mgrid%level(i)%der%lapl
87 end do
88
89 call multigrid_solver_init(this%mg_solver, namespace, space, mesh, thr)
90
91 !%Variable PoissonMultigridCycle
92 !%Type integer
93 !%Section Hamiltonian::Poisson::Multigrid
94 !%Description
95 !% The flavor of multigrid cycle
96 !%Option v_shape 1
97 !% V-shape cycle
98 !%Option w_shape 2
99 !% W-shape cycle
100 !%Option fmg 3
101 !% Full multigrid solver
102 !%End
103 call parse_variable(namespace, 'PoissonMultigridCycle', mg_v_shape, this%multigrid_cycle)
104
106 end subroutine poisson_multigrid_init
107
108
109 ! ---------------------------------------------------------
110 subroutine poisson_multigrid_end(this)
111 type(poisson_mg_solver_t), intent(inout) :: this
112
113 push_sub(poisson_multigrid_end)
115 call poisson_corrections_end(this%corrector)
116
117 call multigrid_end(this%mgrid)
118
119 pop_sub(poisson_multigrid_end)
120 end subroutine poisson_multigrid_end
121
122
123 ! ---------------------------------------------------------
125 subroutine poisson_multigrid_solver(this, namespace, der, pot, rho)
126 type(poisson_mg_solver_t), intent(in) :: this
127 type(namespace_t), intent(in) :: namespace
128 type(derivatives_t), intent(in) :: der
129 real(real64), intent(out) :: pot(:)
130 real(real64), contiguous, intent(in) :: rho(:)
131
132 integer :: ip
133 real(real64), allocatable :: vh_correction(:), res(:), cor(:)
134
136
137 ! correction for treating boundaries
138 safe_allocate(vh_correction(1:der%mesh%np_part))
139 safe_allocate(res(1:der%mesh%np_part))
140 safe_allocate(cor(1:der%mesh%np_part))
141
142 call correct_rho(this%corrector, der, rho, res, vh_correction)
143 call lalg_scal(der%mesh%np, -m_four*m_pi, res)
144
145 select case (this%multigrid_cycle)
147 cor = m_zero
148 call multigrid_iterative_solver(this%mg_solver, namespace, this%mgrid%level(0)%der, this%mgrid%level(0)%der%lapl, &
149 cor, res, this%multigrid_cycle)
150 case(mg_fmg)
151 call multigrid_fmg_solver(this%mg_solver, namespace, this%mgrid%level(0)%der, this%mgrid%level(0)%der%lapl, cor, res)
152 case default
153 assert(.false.)
154 end select
155
156 do ip = 1, der%mesh%np
157 pot(ip) = cor(ip) + vh_correction(ip)
158 end do
159
160 safe_deallocate_a(vh_correction)
161 safe_deallocate_a(res)
162 safe_deallocate_a(cor)
163
165 end subroutine poisson_multigrid_solver
166
168
169!! Local Variables:
170!! mode: f90
171!! coding: utf-8
172!! End:
scales a vector by a constant
Definition: lalg_basic.F90:156
Module implementing boundary conditions in Octopus.
Definition: boundaries.F90:122
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_four
Definition: global.F90:191
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:185
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
subroutine, public multigrid_end(mgrid)
Definition: multigrid.F90:502
subroutine, public multigrid_init(mgrid, namespace, space, mesh, der, stencil, mc, nlevels)
Definition: multigrid.F90:185
This modules provides the routines for solving Ax=b using the V-shaped multigrid method.
recursive subroutine, public multigrid_fmg_solver(this, namespace, der, op, sol, rhs)
Full multigrid (FMG) solver.
subroutine, public multigrid_iterative_solver(this, namespace, der, op, sol, rhs, multigrid_shape)
An iterative multigrid solver.
subroutine, public multigrid_solver_init(this, namespace, space, mesh, thr)
integer, parameter, public mg_fmg
integer, parameter, public mg_v_shape
integer, parameter, public mg_w_shape
subroutine, public poisson_corrections_end(this)
subroutine, public poisson_corrections_init(this, namespace, space, ml, mesh)
subroutine, public correct_rho(this, der, rho, rho_corrected, vh_correction)
subroutine, public poisson_multigrid_solver(this, namespace, der, pot, rho)
A multigrid Poisson solver with corrections at the boundaries.
subroutine, public poisson_multigrid_end(this)
subroutine, public poisson_multigrid_init(this, namespace, space, mesh, der, stencil, mc, ml, thr)
This module defines stencils used in Octopus.
Definition: stencil.F90:135