Octopus
poisson_cg.F90
Go to the documentation of this file.
1!! Cropyright (C) 2004-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
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
22 use debug_oct_m
24 use global_oct_m
25 use, intrinsic :: iso_fortran_env
27 use mesh_oct_m
33
34 implicit none
35
36 private
37 public :: &
42
43 real(real64), public :: threshold
44 integer, public :: maxiter
45
46contains
47
48
49 ! ---------------------------------------------------------
50 subroutine poisson_cg_init(thr, itr)
51 integer, intent(in) :: itr
52 real(real64), intent(in) :: thr
53
54 push_sub(poisson_cg_init)
55 threshold = thr
56 maxiter = itr
57 pop_sub(poisson_cg_init)
58 end subroutine poisson_cg_init
59
60
61 ! ---------------------------------------------------------
62 subroutine poisson_cg_end
63
64 end subroutine poisson_cg_end
65
66
67 ! ---------------------------------------------------------
68 subroutine poisson_cg1(namespace, der, corrector, pot, rho)
69 type(namespace_t), intent(in) :: namespace
70 type(derivatives_t), target, intent(in) :: der
71 type(poisson_corr_t), intent(in) :: corrector
72 real(real64), intent(inout) :: pot(:)
73 real(real64), intent(in) :: rho(:)
74
75 integer :: iter
76 real(real64) :: res
77 real(real64), allocatable :: wk(:), lwk(:), zk(:), pk(:)
78
79 push_sub(poisson_cg1)
80
81 safe_allocate( wk(1:der%mesh%np_part))
82 safe_allocate(lwk(1:der%mesh%np_part))
83 safe_allocate( zk(1:der%mesh%np_part))
84 safe_allocate( pk(1:der%mesh%np_part))
85
86 ! build initial guess for the potential
87 wk(1:der%mesh%np) = pot(1:der%mesh%np)
88 call poisson_boundary_conditions(corrector, der%mesh, rho, wk)
89 call dderivatives_lapl(der, wk, lwk, ghost_update = .true., set_bc = .false.)
90
91 zk(1:der%mesh%np) = -m_four*m_pi*rho(1:der%mesh%np) - lwk(1:der%mesh%np)
92 safe_deallocate_a(wk)
93 safe_deallocate_a(lwk) ! they are no longer needed
94
95 der_pointer => der
96 mesh_pointer => der%mesh
97 call lalg_copy(der%mesh%np, zk, pk)
98 pk(der%mesh%np + 1:) = m_zero
99 iter = maxiter
100 call dconjugate_gradients(der%mesh%np, pk, zk, internal_laplacian_op, internal_dotp, iter, res, threshold)
101 if (res >= threshold) then
102 message(1) = 'Conjugate-gradients Poisson solver did not converge.'
103 write(message(2), '(a,i8)') ' Iter = ',iter
104 write(message(3), '(a,e14.6)') ' Res = ', res
105 call messages_warning(3, namespace=namespace)
106 end if
107 nullify(der_pointer, mesh_pointer)
108 pot(1:der%mesh%np) = pot(1:der%mesh%np) + pk(1:der%mesh%np)
109
110 safe_deallocate_a(zk)
111 safe_deallocate_a(pk)
112 pop_sub(poisson_cg1)
113 end subroutine poisson_cg1
115
116 ! ---------------------------------------------------------
117 subroutine poisson_cg2(namespace, der, pot, rho)
118 type(namespace_t), intent(in) :: namespace
119 type(derivatives_t), target, intent(in) :: der
120 real(real64), contiguous, intent(inout) :: pot(:)
121 real(real64), intent(in) :: rho(:)
122
123 integer :: iter, ip
124 real(real64), allocatable :: potc(:), rhs(:)
125 real(real64) :: res
126
127 push_sub(poisson_cg2)
128
129 iter = maxiter
130 der_pointer => der
131 mesh_pointer => der%mesh
132
133 safe_allocate(rhs(1:der%mesh%np))
134 safe_allocate(potc(1:der%mesh%np_part))
135
136 do ip = 1, der%mesh%np
137 rhs(ip) = -4.0_real64*m_pi*rho(ip)
138 end do
139 call lalg_copy(der%mesh%np, pot, potc)
140
141 call dconjugate_gradients(der%mesh%np, potc, rhs, internal_laplacian_op, internal_dotp, iter, res, threshold)
142
143 if (res >= threshold) then
144 message(1) = 'Conjugate-gradients Poisson solver did not converge.'
145 write(message(2), '(a,i8)') ' Iter = ', iter
146 write(message(3), '(a,e14.6)') ' Res = ', res
147 call messages_warning(3, namespace=namespace)
148 end if
149
150 call lalg_copy(der%mesh%np, potc, pot)
151
152 nullify(der_pointer, mesh_pointer)
153
154 safe_deallocate_a(rhs)
155 safe_deallocate_a(potc)
156
157 pop_sub(poisson_cg2)
158 end subroutine poisson_cg2
159
160end module poisson_cg_oct_m
162!! Local Variables:
163!! mode: f90
164!! coding: utf-8
165!! End:
Copies a vector x, to a vector y.
Definition: lalg_basic.F90:185
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public dderivatives_lapl(der, ff, op_ff, ghost_update, set_bc, factor)
apply the Laplacian to a mesh 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 the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:543
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public poisson_cg2(namespace, der, pot, rho)
Definition: poisson_cg.F90:211
subroutine, public poisson_cg1(namespace, der, corrector, pot, rho)
Definition: poisson_cg.F90:162
subroutine, public poisson_cg_init(thr, itr)
Definition: poisson_cg.F90:144
subroutine, public poisson_cg_end
Definition: poisson_cg.F90:156
subroutine, public internal_laplacian_op(xx, yy)
real(real64) function, public internal_dotp(xx, yy)
subroutine, public poisson_boundary_conditions(this, mesh, rho, pot)
type(mesh_t), pointer, public mesh_pointer
type(derivatives_t), pointer, public der_pointer
This module is intended to contain "only mathematical" functions and procedures.
Definition: solvers.F90:115
int true(void)