Octopus
atomic_orbital.F90
Go to the documentation of this file.
1!! Copyright (C) 2016 N. Tancogne-Dejean
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#include "global.h"
19
25 use debug_oct_m
26 use global_oct_m
27 use, intrinsic :: iso_fortran_env
32 use math_oct_m
33 use mesh_oct_m
38 use ps_oct_m
40 use space_oct_m
45
46 implicit none
47
48 private
49
50 public :: &
59
60 character(len=1), parameter :: &
61 l_notation(0:3) = (/ 's', 'p', 'd', 'f' /)
62
63 integer, public, parameter :: MAX_L = 4
64
65contains
66
67 ! ---------------------------------------------------------
68 real(real64) function atomic_orbital_get_radius(species, mesh, iorb, ispin, truncation, threshold) result(radius)
69 class(species_t), intent(in) :: species
70 type(mesh_t), intent(in) :: mesh
71 integer, intent(in) :: iorb, ispin
72 integer(int64), intent(in) :: truncation
73 real(real64), intent(in) :: threshold
74
75 integer :: ii, ll, mm
76
77 push_sub(atomic_orbital_get_radius)
78
79 call species%get_iwf_ilm(iorb, ispin, ii, ll, mm)
80
81 if (truncation == option__aotruncation__ao_full) then
82 radius = species%get_iwf_radius(ii, ispin, threshold)
83 else
84 radius = species%get_iwf_radius(ii, ispin)
85
86 if (truncation == option__aotruncation__ao_box) then
87 ! if the orbital is larger than the size of the box, we restrict it to this size,
88 ! otherwise the orbital will overlap more than one time with the simulation box.
89 ! This would induces phase problem if the complete mesh is used instead of the sphere
90 radius = min(radius, minval(mesh%box%bounding_box_l(1:mesh%box%dim)-mesh%spacing(1:mesh%box%dim)*1.01_real64))
91 else
92 !If asked, we truncate the orbital to the radius on the projector spheres
93 !of the NL part of the pseudopotential.
94 !This is a way to garanty no overlap between orbitals of different atoms.
95 select type(species)
96 type is(pseudopotential_t)
97 radius = min(radius, species%get_radius())
98 end select
99 end if
100
101 end if
102 ! make sure that if the spacing is too large, the orbitals fit in a few points at least
103 radius = max(radius, m_two*maxval(mesh%spacing))
104
105 pop_sub(atomic_orbital_get_radius)
106 end function atomic_orbital_get_radius
107
108
109#include "undef.F90"
110#include "real.F90"
111#include "atomic_orbital_inc.F90"
112
113#include "undef.F90"
114#include "complex.F90"
115#include "atomic_orbital_inc.F90"
116
117end module atomic_orbital_oct_m
subroutine, public zget_atomic_orbital(namespace, space, latt, pos, species, mesh, sm, ii, ll, jj, os, orbind, radius, d_dim, use_mesh, normalize)
This routine returns the atomic orbital basis – provided by the pseudopotential structure in geo.
subroutine, public dget_atomic_orbital(namespace, space, latt, pos, species, mesh, sm, ii, ll, jj, os, orbind, radius, d_dim, use_mesh, normalize)
This routine returns the atomic orbital basis – provided by the pseudopotential structure in geo.
subroutine, public zatomic_orbital_get_submesh_safe(species, submesh, ii, ll, mm, ispin, phi)
character(len=1), dimension(0:3), parameter, public l_notation
real(real64) function, public atomic_orbital_get_radius(species, mesh, iorb, ispin, truncation, threshold)
subroutine, public datomic_orbital_get_submesh(species, submesh, ii, ll, mm, ispin, phi, derivative)
subroutine, public datomic_orbital_get_submesh_safe(species, submesh, ii, ll, mm, ispin, phi)
subroutine, public zatomic_orbital_get_submesh(species, submesh, ii, ll, mm, ispin, phi, derivative)
real(real64), parameter, public m_two
Definition: global.F90:189
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
Definition: ps.F90:114