Octopus
epot.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-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
21module epot_oct_m
22 use comm_oct_m
23 use debug_oct_m
27 use global_oct_m
28 use grid_oct_m
30 use ions_oct_m
32 use, intrinsic :: iso_fortran_env
35 use mesh_oct_m
37 use mpi_oct_m
40 use parser_oct_m
44 use ps_oct_m
46 use space_oct_m
55 use unit_oct_m
58 use xc_oct_m
59
60 implicit none
61
62 private
63 public :: &
64 epot_t, &
65 epot_init, &
66 epot_end, &
72
73 integer, public, parameter :: &
74 NOREL = 0, &
75 spin_orbit = 1, &
78
79 type epot_t
80 ! Components are public by default
81
82 ! Ions
83 real(real64), allocatable :: vpsl(:)
84 ! !< plus the potential from static electric fields
85 type(projector_t), allocatable :: proj(:)
86 logical :: non_local
87 integer :: natoms
88
89 ! External e-m fields
90 real(real64), allocatable :: e_field(:)
91 real(real64), allocatable :: v_ext(:)
92 real(real64), allocatable :: b_field(:)
93 real(real64), allocatable :: a_static(:,:)
94 integer :: reltype
95
98 real(real64) :: gyromagnetic_ratio
99
101 real(real64) :: so_strength
102
104 real(real64) :: eii
105 real(real64), allocatable :: fii(:, :)
106 real(real64), allocatable :: vdw_forces(:, :)
107 real(real64), allocatable :: photon_forces(:)
108
110 real(real64) :: vdw_stress(3, 3)
111
112 real(real64), allocatable, private :: local_potential(:,:)
113 logical, private :: local_potential_precalculated
115 logical, private :: have_density
116 type(poisson_t), pointer, private :: poisson_solver
117
118 logical :: nlcc = .false.
119 end type epot_t
120
121contains
122
123 ! ---------------------------------------------------------
124 subroutine epot_init(ep, namespace, gr, ions, psolver, ispin, xc_family, kpoints)
125 type(epot_t), intent(out) :: ep
126 type(namespace_t), intent(in) :: namespace
127 type(grid_t), intent(in) :: gr
128 type(ions_t), intent(inout) :: ions
129 type(poisson_t), target, intent(in) :: psolver
130 integer, intent(in) :: ispin
131 integer, intent(in) :: xc_family
132 type(kpoints_t), intent(in) :: kpoints
133
134
135 integer :: ispec, ia
136 integer :: filter
137
138 push_sub(epot_init)
139
140 !%Variable FilterPotentials
141 !%Type integer
142 !%Default filter_ts
143 !%Section Hamiltonian
144 !%Description
145 !% <tt>Octopus</tt> can filter the pseudopotentials so that they no
146 !% longer contain Fourier components larger than the mesh itself. This is
147 !% very useful to decrease the egg-box effect, and so should be used in
148 !% all instances where atoms move (<i>e.g.</i> geometry optimization,
149 !% molecular dynamics, and vibrational modes).
150 !% No filtering is applied to HGH pseudopotentials, as these are already
151 !% smooth, analytical potentials in both real and Fourier space.
152 !%Option filter_none 0
153 !% Do not filter.
154 !%Option filter_TS 2
155 !% The filter of M. Tafipolsky and R. Schmid, <i>J. Chem. Phys.</i> <b>124</b>, 174102 (2006).
156 !%Option filter_BSB 3
157 !% The filter of E. L. Briggs, D. J. Sullivan, and J. Bernholc, <i>Phys. Rev. B</i> <b>54</b>, 14362 (1996).
158 !%End
159 call parse_variable(namespace, 'FilterPotentials', ps_filter_ts, filter)
160 if (.not. varinfo_valid_option('FilterPotentials', filter)) call messages_input_error(namespace, 'FilterPotentials')
161 call messages_print_var_option("FilterPotentials", filter, namespace=namespace)
162
163 if (family_is_mgga(xc_family) .and. filter /= ps_filter_none) then
164 call messages_not_implemented("FilterPotentials different from filter_none with MGGA", namespace=namespace)
165 end if
167 if (filter == ps_filter_ts) call spline_filter_mask_init()
168 do ispec = 1, ions%nspecies
169 call ions%species(ispec)%s%init_potential(namespace, mesh_gcutoff(gr), filter)
170 end do
171
172 safe_allocate(ep%vpsl(1:gr%np))
173
174 ep%vpsl(1:gr%np) = m_zero
175
176 ! No more "UserDefinedTDPotential" from this version on.
177 call messages_obsolete_variable(namespace, 'UserDefinedTDPotential', 'TDExternalFields')
179 call messages_obsolete_variable(namespace, 'ClassicalPotential')
181 !%Variable GyromagneticRatio
182 !%Type float
183 !%Default 2.0023193043768
184 !%Section Hamiltonian
185 !%Description
186 !% The gyromagnetic ratio of the electron. This is of course a physical
187 !% constant, and the default value is the exact one that you should not
188 !% touch, unless:
189 !% (i) You want to disconnect the anomalous Zeeman term in the Hamiltonian
190 !% (then set it to zero; this number only affects that term);
191 !% (ii) You are using an effective Hamiltonian, as is the case when
192 !% you calculate a 2D electron gas, in which case you have an effective
193 !% gyromagnetic factor that depends on the material.
194 !%End
195 call parse_variable(namespace, 'GyromagneticRatio', p_g, ep%gyromagnetic_ratio)
196
197 !%Variable RelativisticCorrection
198 !%Type integer
199 !%Default non_relativistic
200 !%Section Hamiltonian
201 !%Description
202 !% The default value means that <i>no</i> relativistic correction is used. To
203 !% include spin-orbit coupling turn <tt>RelativisticCorrection</tt> to <tt>spin_orbit</tt>
204 !% (this will only work if <tt>SpinComponents</tt> has been set to <tt>non_collinear</tt>, which ensures
205 !% the use of spinors).
206 !%Option non_relativistic 0
207 !% No relativistic corrections.
208 !%Option spin_orbit 1
209 !% Spin-orbit.
210 !%Option scalar_relativistic_zora 2
211 !% scalar relativistic ZORA Hamiltonian
212 !%Option fully_relativistic_zora 3
213 !% fully relativistic spin-orbit ZORA Hamiltonian
214 !% including SR and SO terms
215 !%End
216 call parse_variable(namespace, 'RelativisticCorrection', norel, ep%reltype)
217 if (.not. varinfo_valid_option('RelativisticCorrection', ep%reltype)) then
218 call messages_input_error(namespace, 'RelativisticCorrection')
219 end if
220 if (ispin /= spinors .and. ( ep%reltype == spin_orbit .or. ep%reltype == fully_relativistic_zora ) ) then
221 message(1) = "The spin-orbit term can only be applied when using spinors."
222 call messages_fatal(1, namespace=namespace)
223 end if
224
225 if((ep%reltype == spin_orbit .or. ep%reltype == fully_relativistic_zora) .and. kpoints%use_symmetries) then
226 call messages_not_implemented("Spin-orbit coupling and k-point symmetries", namespace=namespace)
227 end if
228
229 call messages_print_var_option("RelativisticCorrection", ep%reltype, namespace=namespace)
230
231 !%Variable SOStrength
232 !%Type float
233 !%Default 1.0
234 !%Section Hamiltonian
235 !%Description
236 !% Tuning of the spin-orbit coupling strength: setting this value to zero turns off spin-orbit terms in
237 !% the Hamiltonian, and setting it to one corresponds to full spin-orbit.
238 !%End
239 if (ep%reltype == spin_orbit .or. ep%reltype == fully_relativistic_zora) then
240 call parse_variable(namespace, 'SOStrength', m_one, ep%so_strength)
241 else
242 ep%so_strength = m_one
243 end if
244
245 safe_allocate(ep%proj(1:ions%natoms))
246
247 ep%natoms = ions%natoms
248 ep%non_local = .false.
249
250 ep%eii = m_zero
251 safe_allocate(ep%fii(1:ions%space%dim, 1:ions%natoms))
252 ep%fii = m_zero
253
254 safe_allocate(ep%vdw_forces(1:ions%space%dim, 1:ions%natoms))
255 ep%vdw_forces = m_zero
256
257 safe_allocate(ep%photon_forces(1:ions%space%dim))
258 ep%photon_forces = m_zero
259
260 ep%local_potential_precalculated = .false.
261
262
263 ep%have_density = .false.
264 do ia = 1, ions%nspecies
265 if (local_potential_has_density(ions%space, ions%species(ia)%s)) then
266 ep%have_density = .true.
267 exit
268 end if
269 end do
270
271 if (ep%have_density) then
272 ep%poisson_solver => psolver
273 else
274 nullify(ep%poisson_solver)
275 end if
276
277 ! find out if we need non-local core corrections
278 ep%nlcc = .false.
279 do ia = 1, ions%nspecies
280 ep%nlcc = (ep%nlcc .or. ions%species(ia)%s%is_ps_with_nlcc())
281 end do
282
283 pop_sub(epot_init)
284 end subroutine epot_init
285
286 ! ---------------------------------------------------------
287 subroutine epot_end(ep)
288 type(epot_t), intent(inout) :: ep
289
290 integer :: iproj
291
292 push_sub(epot_end)
293
294 if (ep%have_density) then
295 nullify(ep%poisson_solver)
296 end if
297
298 safe_deallocate_a(ep%local_potential)
299 safe_deallocate_a(ep%fii)
300 safe_deallocate_a(ep%vdw_forces)
301 safe_deallocate_a(ep%vpsl)
302 safe_deallocate_a(ep%photon_forces)
303
304 ! the macroscopic fields
305 safe_deallocate_a(ep%e_field)
306 safe_deallocate_a(ep%v_ext)
307 safe_deallocate_a(ep%b_field)
308 safe_deallocate_a(ep%a_static)
309
310 do iproj = 1, ep%natoms
311 if (projector_is_null(ep%proj(iproj))) cycle
312 call projector_end(ep%proj(iproj))
313 end do
314
315 assert(allocated(ep%proj))
316 safe_deallocate_a(ep%proj)
317
318 pop_sub(epot_end)
319
320 end subroutine epot_end
321
322 ! ---------------------------------------------------------
323 subroutine epot_generate(ep, namespace, mesh, ions, st_d)
324 type(epot_t), intent(inout) :: ep
325 type(namespace_t), intent(in) :: namespace
326 class(mesh_t), target, intent(in) :: mesh
327 type(ions_t), target, intent(inout) :: ions
328 type(states_elec_dim_t), intent(inout) :: st_d
329
330 integer :: ia
331 type(ps_t), pointer :: ps
332
333 call profiling_in("EPOT_GENERATE")
334 push_sub(epot_generate)
335
336 ! Local part
337 ep%vpsl = m_zero
338
339 ! we assume that we need to recalculate the ion-ion energy
340 call ion_interaction_calculate(ions%ion_interaction, ions%space, ions%latt, ions%atom, &
341 ions%natoms, ions%pos, mesh%box%bounding_box_l, ep%eii, ep%fii)
342
343 ! the pseudopotential part.
344 do ia = 1, ions%natoms
345 select type(spec=>ions%atom(ia)%species)
346 type is(pseudopotential_t)
347 call projector_end(ep%proj(ia))
348 call projector_init(ep%proj(ia), spec, namespace, st_d%dim, ep%reltype)
349 end select
350 end do
351
352 do ia = ions%atoms_dist%start, ions%atoms_dist%end
353 if (ep%proj(ia)%type == proj_none) cycle
354 select type(spec=>ions%atom(ia)%species)
355 type is(pseudopotential_t)
356 ps => spec%ps
357 call submesh_init(ep%proj(ia)%sphere, ions%space, mesh, ions%latt, ions%pos(:, ia), ps%rc_max)
358 end select
359 end do
360
361 if (ions%atoms_dist%parallel) then
362 do ia = 1, ions%natoms
363 if (ep%proj(ia)%type == proj_none) cycle
364 select type(spec=>ions%atom(ia)%species)
365 type is(pseudopotential_t)
366 ps => spec%ps
367 call submesh_broadcast(ep%proj(ia)%sphere, ions%space, mesh, ions%pos(:, ia), ps%rc_max, &
368 ions%atoms_dist%node(ia), ions%atoms_dist%mpi_grp)
369 end select
370 end do
371 end if
372
373 do ia = 1, ions%natoms
374 select type(spec=>ions%atom(ia)%species)
375 type is(pseudopotential_t)
376 call projector_build(ep%proj(ia), spec, ep%so_strength)
377 if (.not. projector_is(ep%proj(ia), proj_none)) ep%non_local = .true.
378 end select
379 end do
381 pop_sub(epot_generate)
382 call profiling_out("EPOT_GENERATE")
383 end subroutine epot_generate
384
385 ! ---------------------------------------------------------
386
387 logical pure function local_potential_has_density(space, species) result(has_density)
388 class(space_t), intent(in) :: space
389 class(species_t), intent(in) :: species
390
391 has_density = species%has_density .or. (species%is_ps() .and. space%is_periodic())
392
393 end function local_potential_has_density
394
395 ! ---------------------------------------------------------
396 subroutine epot_local_potential(ep, namespace, space, latt, mesh, species, pos, iatom, vpsl)
397 type(epot_t), intent(in) :: ep
398 type(namespace_t), intent(in) :: namespace
399 class(space_t), intent(in) :: space
400 type(lattice_vectors_t), intent(in) :: latt
401 class(mesh_t), intent(in) :: mesh
402 class(species_t), target, intent(in) :: species
403 real(real64), intent(in) :: pos(1:space%dim)
404 integer, intent(in) :: iatom
405 real(real64), contiguous, intent(inout) :: vpsl(:)
406
407 integer :: ip
408 real(real64) :: radius
409 real(real64), allocatable :: vl(:), rho(:)
410 type(submesh_t) :: sphere
411 type(ps_t), pointer :: ps
412
413 push_sub(epot_local_potential)
414 call profiling_in("EPOT_LOCAL")
415
416 if (ep%local_potential_precalculated) then
417
418 do ip = 1, mesh%np
419 vpsl(ip) = vpsl(ip) + ep%local_potential(ip, iatom)
420 end do
421 else
422
423 !Local potential, we can get it by solving the Poisson equation
424 !(for all-electron species or pseudopotentials in periodic
425 !systems) or by applying it directly to the grid
426 safe_allocate(vl(1:mesh%np))
427
428 if (local_potential_has_density(space, species)) then
429 safe_allocate(rho(1:mesh%np))
430
431 call species_get_long_range_density(species, namespace, space, latt, pos, mesh, rho, sphere)
432
433 if (poisson_solver_is_iterative(ep%poisson_solver)) then
434 ! vl has to be initialized before entering routine
435 ! and our best guess for the potential is zero
436 vl(1:mesh%np) = m_zero
437 end if
438
439 call dpoisson_solve(ep%poisson_solver, namespace, vl, rho, all_nodes = .false.)
440
441 safe_deallocate_a(rho)
442
443 else
444
445 call species_get_local(species, namespace, space, latt, pos, mesh, vl)
446
447 end if
448
449 call lalg_axpy(mesh%np, m_one, vl, vpsl)
450 safe_deallocate_a(vl)
451
452 !the localized part
453 select type(species)
454 class is(pseudopotential_t)
455
456 ps => species%ps
457
458 radius = ps%vl%x_threshold*1.05_real64
459 if (.not. submesh_compatible(sphere, radius, pos, minval(mesh%spacing(1:space%dim)))) then
460 call submesh_end(sphere)
461 call submesh_init(sphere, space, mesh, latt, pos, radius)
462 end if
463 safe_allocate(vl(1:sphere%np))
464 vl = m_zero
465
466 do ip = 1, sphere%np
467 if(sphere%r(ip) <= radius) then
468 vl(ip) = spline_eval(ps%vl, sphere%r(ip))
469 end if
470 end do
471
472 call submesh_add_to_mesh(sphere, vl, vpsl)
473
474 safe_deallocate_a(vl)
475 nullify(ps)
476
477 end select
478 call submesh_end(sphere)
479
480 end if
481
482 call profiling_out("EPOT_LOCAL")
483 pop_sub(epot_local_potential)
484 end subroutine epot_local_potential
485
486 ! ---------------------------------------------------------
487 subroutine epot_precalc_local_potential(ep, namespace, gr, ions)
488 type(epot_t), intent(inout) :: ep
489 type(namespace_t), intent(in) :: namespace
490 type(grid_t), intent(in) :: gr
491 type(ions_t), intent(in) :: ions
492
493 integer :: iatom
494
496
497 if (.not. allocated(ep%local_potential)) then
498 safe_allocate(ep%local_potential(1:gr%np, 1:ions%natoms))
499 end if
500
501 ep%local_potential_precalculated = .false.
502
503 do iatom = 1, ions%natoms
504 ep%local_potential(1:gr%np, iatom) = m_zero
505 call epot_local_potential(ep, namespace, ions%space, ions%latt, gr, ions%atom(iatom)%species, &
506 ions%pos(:, iatom), iatom, ep%local_potential(1:gr%np, iatom))!, time)
507 end do
508 ep%local_potential_precalculated = .true.
509
511 end subroutine epot_precalc_local_potential
512
513 ! ---------------------------------------------------------
514
515 logical function epot_have_external_potentials(ep)
516 type(epot_t), intent(in) :: ep
517
519
520 epot_have_external_potentials = allocated(ep%e_field)
521
523
525
526end module epot_oct_m
527
528!! Local Variables:
529!! mode: f90
530!! coding: utf-8
531!! End:
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
integer, parameter, public spinors
logical function, public epot_have_external_potentials(ep)
Definition: epot.F90:609
logical pure function, public local_potential_has_density(space, species)
Definition: epot.F90:481
integer, parameter, public spin_orbit
Definition: epot.F90:166
integer, parameter, public scalar_relativistic_zora
Definition: epot.F90:166
subroutine, public epot_end(ep)
Definition: epot.F90:381
integer, parameter, public fully_relativistic_zora
Definition: epot.F90:166
subroutine, public epot_precalc_local_potential(ep, namespace, gr, ions)
Definition: epot.F90:581
subroutine, public epot_local_potential(ep, namespace, space, latt, mesh, species, pos, iatom, vpsl)
Definition: epot.F90:490
subroutine, public epot_init(ep, namespace, gr, ions, psolver, ispin, xc_family, kpoints)
Definition: epot.F90:218
subroutine, public epot_generate(ep, namespace, mesh, ions, st_d)
Definition: epot.F90:417
real(real64), parameter, public p_g
Definition: global.F90:225
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_one
Definition: global.F90:188
This module implements the underlying real-space grid.
Definition: grid.F90:117
subroutine, public ion_interaction_calculate(this, space, latt, atom, natoms, pos, lsize, energy, force, energy_components, force_components)
Top level routine for computing electrostatic energies and forces between ions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
real(real64) function, public mesh_gcutoff(mesh)
mesh_gcutoff returns the "natural" band limitation of the grid mesh, in terms of the maximum G vector...
Definition: mesh.F90:444
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1125
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1057
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:723
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
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
subroutine, public projector_build(p, ps, so_strength)
Definition: projector.F90:343
logical elemental function, public projector_is(p, type)
Definition: projector.F90:210
subroutine, public projector_init(p, pseudo, namespace, dim, reltype)
Definition: projector.F90:217
subroutine, public projector_end(p)
Definition: projector.F90:392
logical elemental function, public projector_is_null(p)
Definition: projector.F90:203
Definition: ps.F90:114
integer, parameter, public ps_filter_ts
Definition: ps.F90:160
integer, parameter, public ps_filter_none
Definition: ps.F90:160
integer, parameter, public proj_none
Definition: ps.F90:165
subroutine, public spline_filter_mask_init()
This module handles spin dimensions of the states and the k-point distribution.
subroutine, public submesh_broadcast(this, space, mesh, center, radius, root, mpi_grp)
Definition: submesh.F90:658
subroutine, public submesh_init(this, space, mesh, latt, center, rc)
Definition: submesh.F90:280
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
Definition: xc.F90:114
pure logical function, public family_is_mgga(family, only_collinear)
Is the xc function part of the mGGA family.
Definition: xc.F90:579
Describes mesh distribution to nodes.
Definition: mesh.F90:186
An abstract class for species. Derived classes include jellium, all electron, and pseudopotential spe...
Definition: species.F90:143
class for organizing spins and k-points
int true(void)