Octopus
ion_electron_local_potential.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2021 Nicolas Tancogne-Dejean
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
23 use atom_oct_m
24 use comm_oct_m
25 use debug_oct_m
27 use epot_oct_m
28 use global_oct_m
31 use ions_oct_m
32 use, intrinsic :: iso_fortran_env
35 use mesh_oct_m
41 use ps_oct_m
44 use space_oct_m
50
51 implicit none
52
53 private
54 public :: &
56
58 private
59
60 type(poisson_t), pointer :: psolver
61 type(mesh_t), pointer :: mesh
62 type(space_t), pointer :: space
63
64 ! This information should be copied, and not obtained thanks to pointer
65 ! This is a temporary change here
66 type(distributed_t), pointer, public :: atoms_dist => null()
67 type(atom_t), pointer, public :: atom(:) => null()
68 real(real64), pointer, public :: pos(:,:) => null()
69 type(lattice_vectors_t), pointer :: latt => null()
70
71 ! Temporary pointer to namespace
72 type(namespace_t), pointer :: namespace
73
74 logical :: have_density
75
76 contains
77 procedure :: init => ion_electron_local_potential_init
78 procedure :: bind => ion_electron_local_potential_bind
79 procedure :: calculate => ion_electron_local_potential_calculate
80 procedure :: calculate_energy => ion_electron_local_potential_calculate_energy
81 procedure :: end => ion_electron_local_potential_end
84
88
89contains
90
91 ! ---------------------------------------------------------
92 function ion_electron_local_potential_constructor(partner) result(this)
93 class(interaction_partner_t), target, intent(inout) :: partner
94 class(ion_electron_local_potential_t), pointer :: this
95
97
98 allocate(this)
99
100 this%label = "ion-electron local"
101
102 this%partner => partner
103
106
107 ! ---------------------------------------------------------
108 subroutine ion_electron_local_potential_init(this, mesh, psolver, ions, namespace)
109 class(ion_electron_local_potential_t), intent(inout) :: this
110 class(mesh_t), target, intent(in) :: mesh
111 type(poisson_t), target, intent(in) :: psolver
112 type(ions_t), target, intent(in) :: ions
113 type(namespace_t), target, intent(in) :: namespace
114
115 integer :: ia
116
118
119 this%mesh => mesh
120 this%psolver => psolver
121
122 safe_allocate(this%potential(1:mesh%np, 1))
123
124 this%have_density = .false.
125 do ia = 1, ions%nspecies
126 if (local_potential_has_density(ions%space, ions%species(ia)%s)) then
127 this%have_density = .true.
128 exit
129 end if
130 end do
131
132 this%atoms_dist => ions%atoms_dist
133 this%atom => ions%atom
134 this%space => ions%space
135 this%pos => ions%pos
136 this%latt => ions%latt
137
138 this%namespace => namespace
139
142
143 ! ---------------------------------------------------------
146 subroutine ion_electron_local_potential_bind(this, psolver, ions)
147 class(ion_electron_local_potential_t), intent(inout) :: this
148 type(poisson_t), target, intent(in) :: psolver
149 type(ions_t), target, intent(in) :: ions
150
153 this%psolver => psolver
154 this%atoms_dist => ions%atoms_dist
155 this%atom => ions%atom
156 this%space => ions%space
157 this%pos => ions%pos
158 this%latt => ions%latt
159
163 ! ---------------------------------------------------------
164 ! Note: this code is adapted from epot_local_potential. There are code duplication at the moment
166 class(ion_electron_local_potential_t), intent(inout) :: this
168 real(real64), allocatable :: density(:), rho(:), vl(:)
169 type(submesh_t) :: sphere
170 type(ps_t), pointer :: ps
171 integer :: ia, ip
172 real(real64) :: radius
176 call profiling_in("ION_ELEC_LOC_INTER")
178 this%potential(:,1) = m_zero
179
180 if (this%have_density) then
181 safe_allocate(density(1:this%mesh%np))
182 density = m_zero
183 end if
184
185 do ia = this%atoms_dist%start, this%atoms_dist%end
186 ! Local potential, we can get it by solving the Poisson equation
187 ! (for all-electron species or pseudopotentials in periodic
188 ! systems) or by applying it directly to the grid
189
190 if (local_potential_has_density(this%space, this%atom(ia)%species)) then
191
192 safe_allocate(rho(1:this%mesh%np))
193 call species_get_long_range_density(this%atom(ia)%species, this%namespace, this%space, this%latt, &
194 this%pos(:, ia), this%mesh, rho, sphere)
195 call lalg_axpy(this%mesh%np, m_one, rho, density)
196 safe_deallocate_a(rho)
197
198 else
199
200 safe_allocate(vl(1:this%mesh%np))
201 call species_get_local(this%atom(ia)%species, this%namespace, this%space, this%latt, &
202 this%pos(:, ia), this%mesh, vl)
203 call lalg_axpy(this%mesh%np, m_one, vl, this%potential(:,1))
204 safe_deallocate_a(vl)
205
206 end if
207
208 !then we add the localized part
209 select type(spec=>this%atom(ia)%species)
210 type is(pseudopotential_t)
211
212 ps => spec%ps
213
214 radius = min(ps%vl%x_threshold*1.05_real64, spline_range_max(ps%vl))
215 if ( .not. submesh_compatible(sphere,radius,this%pos(:,ia), minval(this%mesh%spacing(1:this%space%dim))) ) then
216 call submesh_init(sphere, this%space, this%mesh, this%latt, this%pos(:, ia), radius)
217 end if
218 safe_allocate(vl(1:sphere%np))
219
220 do ip = 1, sphere%np
221 if(sphere%r(ip) <= radius) then
222 vl(ip) = spline_eval(ps%vl, sphere%r(ip))
223 else
224 vl(ip) = 0
225 end if
226 end do
227
228 call submesh_add_to_mesh(sphere, vl, this%potential(:,1))
229
230 safe_deallocate_a(vl)
231 nullify(ps)
232 end select
233 call submesh_end(sphere)
234 end do
235
236 ! reduce over atoms if required
237 if (this%atoms_dist%parallel) then
238 call comm_allreduce(this%atoms_dist%mpi_grp, this%potential(:,1), dim = this%mesh%np)
239 if (this%have_density) then
240 call comm_allreduce(this%atoms_dist%mpi_grp, density, dim = this%mesh%np)
241 end if
242 end if
243
244 ! now we solve the poisson equation with the density of all nodes
245 if (this%have_density) then
246
247 safe_allocate(vl(1:this%mesh%np_part))
248 call dpoisson_solve(this%psolver, this%namespace, vl, density)
249 call lalg_axpy(this%mesh%np, m_one, vl, this%potential(:,1))
250 safe_deallocate_a(vl)
251
252 end if
253 safe_deallocate_a(density)
254
255 call profiling_out("ION_ELEC_LOC_INTER")
256
259
260 ! ---------------------------------------------------------
261 subroutine ion_electron_local_potential_end(this)
262 class(ion_electron_local_potential_t), intent(inout) :: this
263
265
266 safe_deallocate_a(this%potential)
267 nullify(this%mesh)
268 nullify(this%psolver)
269
270 call interaction_end(this)
271
274
275
276 ! ---------------------------------------------------------
278 type(ion_electron_local_potential_t), intent(inout) :: this
279
281
283
286
287 ! ---------------------------------------------------------
289 class(ion_electron_local_potential_t), intent(inout) :: this
290
292
294 assert(.false.)
295
297
299
300
302
303!! Local Variables:
304!! mode: f90
305!! coding: utf-8
306!! End:
constant times a vector plus a vector
Definition: lalg_basic.F90:173
logical pure function, public local_potential_has_density(space, species)
Definition: epot.F90:503
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_one
Definition: global.F90:192
This module defines the abstract interaction_t class, and some auxiliary classes for interactions.
subroutine, public interaction_end(this)
This module defines classes and functions for interaction partners.
class(ion_electron_local_potential_t) function, pointer ion_electron_local_potential_constructor(partner)
subroutine ion_electron_local_potential_bind(this, psolver, ions)
Rebind ion-electron local potential helper pointers after a Hamiltonian copy. Ensure the local potent...
subroutine ion_electron_local_potential_init(this, mesh, psolver, ions, namespace)
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public dpoisson_solve(this, namespace, pot, rho, all_nodes, kernel, reset)
Calculates the Poisson equation. Given the density returns the corresponding potential.
Definition: poisson.F90:875
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
Definition: ps.F90:116
This module defines the quantity_t class and the IDs for quantities, which can be exposed by a system...
Definition: quantity.F90:140
subroutine, public species_get_local(species, namespace, space, latt, pos, mesh, vl)
used when the density is not available, or otherwise the Poisson eqn would be used instead
subroutine, public species_get_long_range_density(species, namespace, space, latt, pos, mesh, rho, sphere_inout, nlr_x)
real(real64) function, public spline_eval(spl, x)
Definition: splines.F90:443
real(real64) pure function, public spline_range_max(this)
Definition: splines.F90:1165
subroutine, public submesh_end(this)
Definition: submesh.F90:733
logical function, public submesh_compatible(this, radius, center, dx)
Definition: submesh.F90:715
subroutine, public submesh_init(this, space, mesh, latt, center, rc)
Definition: submesh.F90:282
A type storing the information and data about a pseudopotential.
Definition: ps.F90:188
A submesh is a type of mesh, used for the projectors in the pseudopotentials It contains points on a ...
Definition: submesh.F90:177
int true(void)