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 :: calculate => ion_electron_local_potential_calculate
79 procedure :: calculate_energy => ion_electron_local_potential_calculate_energy
80 procedure :: end => ion_electron_local_potential_end
83
87
88contains
89
90 ! ---------------------------------------------------------
91 function ion_electron_local_potential_constructor(partner) result(this)
92 class(interaction_partner_t), target, intent(inout) :: partner
93 class(ion_electron_local_potential_t), pointer :: this
94
96
97 allocate(this)
98
99 this%label = "ion-electron local"
100
101 this%partner => partner
102
105
106 ! ---------------------------------------------------------
107 subroutine ion_electron_local_potential_init(this, mesh, psolver, ions, namespace)
108 class(ion_electron_local_potential_t), intent(inout) :: this
109 class(mesh_t), target, intent(in) :: mesh
110 type(poisson_t), target, intent(in) :: psolver
111 type(ions_t), target, intent(in) :: ions
112 type(namespace_t), target, intent(in) :: namespace
113
114 integer :: ia
117
118 this%mesh => mesh
119 this%psolver => psolver
120
121 safe_allocate(this%potential(1:mesh%np, 1))
122
123 this%have_density = .false.
124 do ia = 1, ions%nspecies
125 if (local_potential_has_density(ions%space, ions%species(ia)%s)) then
126 this%have_density = .true.
127 exit
128 end if
129 end do
130
131 this%atoms_dist => ions%atoms_dist
132 this%atom => ions%atom
133 this%space => ions%space
134 this%pos => ions%pos
135 this%latt => ions%latt
136
137 this%namespace => namespace
138
141
142 ! ---------------------------------------------------------
143 ! Note: this code is adapted from epot_local_potential. There are code duplication at the moment
145 class(ion_electron_local_potential_t), intent(inout) :: this
146
147 real(real64), allocatable :: density(:), rho(:), vl(:)
148 type(submesh_t) :: sphere
149 type(ps_t), pointer :: ps
150 integer :: ia, ip
151 real(real64) :: radius
152
155 call profiling_in("ION_ELEC_LOC_INTER")
156
157 this%potential(:,1) = m_zero
158
159 if (this%have_density) then
160 safe_allocate(density(1:this%mesh%np))
161 density = m_zero
162 end if
163
164 do ia = this%atoms_dist%start, this%atoms_dist%end
165 ! Local potential, we can get it by solving the Poisson equation
166 ! (for all-electron species or pseudopotentials in periodic
167 ! systems) or by applying it directly to the grid
168
169 if (local_potential_has_density(this%space, this%atom(ia)%species)) then
171 safe_allocate(rho(1:this%mesh%np))
172 call species_get_long_range_density(this%atom(ia)%species, this%namespace, this%space, this%latt, &
173 this%pos(:, ia), this%mesh, rho, sphere)
174 call lalg_axpy(this%mesh%np, m_one, rho, density)
175 safe_deallocate_a(rho)
176
177 else
178
179 safe_allocate(vl(1:this%mesh%np))
180 call species_get_local(this%atom(ia)%species, this%namespace, this%space, this%latt, &
181 this%pos(:, ia), this%mesh, vl)
182 call lalg_axpy(this%mesh%np, m_one, vl, this%potential(:,1))
183 safe_deallocate_a(vl)
185 end if
186
187 !then we add the localized part
188 select type(spec=>this%atom(ia)%species)
189 type is(pseudopotential_t)
190
191 ps => spec%ps
192
193 radius = ps%vl%x_threshold*1.05_real64
194 if ( .not. submesh_compatible(sphere,radius,this%pos(:,ia), minval(this%mesh%spacing(1:this%space%dim))) ) then
195 call submesh_init(sphere, this%space, this%mesh, this%latt, this%pos(:, ia), radius)
196 end if
197 safe_allocate(vl(1:sphere%np))
198
199 do ip = 1, sphere%np
200 if(sphere%r(ip) <= radius) then
201 vl(ip) = spline_eval(ps%vl, sphere%r(ip))
202 else
203 vl(ip) = 0
204 end if
205 end do
206
207 call submesh_add_to_mesh(sphere, vl, this%potential(:,1))
208
209 safe_deallocate_a(vl)
210 nullify(ps)
211 end select
212 call submesh_end(sphere)
213 end do
214
215 ! reduce over atoms if required
216 if (this%atoms_dist%parallel) then
217 call comm_allreduce(this%atoms_dist%mpi_grp, this%potential(:,1), dim = this%mesh%np)
218 if (this%have_density) then
219 call comm_allreduce(this%atoms_dist%mpi_grp, density, dim = this%mesh%np)
220 end if
221 end if
222
223 ! now we solve the poisson equation with the density of all nodes
224 if (this%have_density) then
225
226 safe_allocate(vl(1:this%mesh%np_part))
227 if (poisson_solver_is_iterative(this%psolver)) then
228 vl(1:this%mesh%np) = m_zero
229 end if
230 call dpoisson_solve(this%psolver, this%namespace, vl, density)
231 call lalg_axpy(this%mesh%np, m_one, vl, this%potential(:,1))
232 safe_deallocate_a(vl)
233
234 end if
235 safe_deallocate_a(density)
236
237 call profiling_out("ION_ELEC_LOC_INTER")
238
241
242 ! ---------------------------------------------------------
243 subroutine ion_electron_local_potential_end(this)
244 class(ion_electron_local_potential_t), intent(inout) :: this
245
247
248 safe_deallocate_a(this%potential)
249 nullify(this%mesh)
250 nullify(this%psolver)
251
252 call interaction_end(this)
253
256
257
258 ! ---------------------------------------------------------
260 type(ion_electron_local_potential_t), intent(inout) :: this
261
263
265
268
269 ! ---------------------------------------------------------
271 class(ion_electron_local_potential_t), intent(inout) :: this
272
274
276 assert(.false.)
277
279
281
282
284
285!! Local Variables:
286!! mode: f90
287!! coding: utf-8
288!! End:
constant times a vector plus a vector
Definition: lalg_basic.F90:170
logical pure function, public local_potential_has_density(space, species)
Definition: epot.F90:481
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_one
Definition: global.F90:188
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_init(this, mesh, psolver, ions, namespace)
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
logical pure function, public poisson_solver_is_iterative(this)
Definition: poisson.F90:1300
subroutine, public dpoisson_solve(this, namespace, pot, rho, all_nodes, kernel)
Calculates the Poisson equation. Given the density returns the corresponding potential.
Definition: poisson.F90:892
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
Definition: ps.F90:114
This module defines the quantity_t class and the IDs for quantities, which can be exposed by a system...
Definition: quantity.F90:137
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:441
subroutine, public submesh_end(this)
Definition: submesh.F90:735
logical function, public submesh_compatible(this, radius, center, dx)
Definition: submesh.F90:717
subroutine, public submesh_init(this, space, mesh, latt, center, rc)
Definition: submesh.F90:280
int true(void)