72 integer,
parameter :: MAX_NPROJECTIONS = 4
73 integer,
parameter :: MAX_L = 5
89 integer :: nprojections
90 integer,
public :: lmax
91 integer,
public :: lloc
95 type(submesh_t),
public :: sphere
100 type(hgh_projector_t),
allocatable,
public :: hgh_p(:, :)
101 type(kb_projector_t),
allocatable,
public :: kb_p(:, :)
102 type(rkb_projector_t),
allocatable,
public :: rkb_p(:, :)
103 complex(real64),
allocatable,
public :: phase(:, :, :)
109 logical elemental function projector_is_null(p)
110 type(projector_t),
intent(in) :: p
112 projector_is_null = (p%type ==
proj_none)
116 logical elemental function projector_is(p, type)
117 type(projector_t),
intent(in) :: p
118 integer,
intent(in) :: type
124 type(projector_t),
intent(inout) :: p
125 type(pseudopotential_t),
target,
intent(in) :: pseudo
126 type(namespace_t),
intent(in) :: namespace
127 integer,
intent(in) :: dim
128 integer,
intent(in) :: reltype
130 type(ps_t),
pointer :: ps
147 p%type = ps%projector_type
149 if (p%type == proj_kb .and. reltype == 1)
then
150 if (ps%relativistic_treatment == proj_j_dependent)
then
153 call messages_write(
"Spin-orbit coupling for species '"//trim(pseudo%get_label())//
" is not available.")
154 call messages_warning(namespace=namespace)
159 case (proj_kb, proj_rkb)
160 p%nprojections = ps%kbc
172 integer,
intent(in) :: dim
173 type(states_elec_dim_t),
intent(in) :: std
174 type(boundaries_t),
intent(in) :: bnd
175 type(kpoints_t),
intent(in) :: kpoints
176 real(real64),
optional,
allocatable,
intent(in) :: vec_pot(:)
177 real(real64),
optional,
allocatable,
intent(in) :: vec_pot_var(:, :)
179 integer :: ns, iq, is, ikpoint
180 real(real64) :: kr, kpoint(dim)
181 integer :: nphase, iphase
182 real(real64),
allocatable :: diff(:,:)
188 if (bnd%spiralBC) nphase = 3
190 if (.not.
allocated(this%phase) .and. ns > 0)
then
191 safe_allocate(this%phase(1:ns, 1:nphase, std%kpt%start:std%kpt%end))
198 safe_allocate(diff(1:dim, 1:ns))
202 diff(:, is) = this%sphere%rel_x(:,is) + this%sphere%center - this%sphere%mesh%x(this%sphere%map(is), :)
205 do iq = std%kpt%start, std%kpt%end
206 ikpoint = std%get_kpoint_index(iq)
209 assert(ikpoint <= kpoints_number(kpoints))
212 kpoint(1:dim) = kpoints%get_point(ikpoint)
214 do iphase = 1, nphase
220 kr = sum(kpoint(1:dim)*diff(1:dim, is))
222 if (
present(vec_pot))
then
223 if (
allocated(vec_pot)) kr = kr + sum(vec_pot(1:dim)*diff(1:dim, is))
226 if (
present(vec_pot_var))
then
227 if (
allocated(vec_pot_var)) kr = kr + sum(vec_pot_var(1:dim, this%sphere%map(is)) &
228 *(this%sphere%rel_x(:, is)+this%sphere%center))
231 if (bnd%spiralBC .and. iphase > 1)
then
232 kr = kr + (2*(iphase-1)-3)*sum(bnd%spiral_q(1:dim)*diff(1:dim, is))
235 this%phase(is, iphase, iq) =
exp(-m_zi*kr)
242 safe_deallocate_a(diff)
251 class(pseudopotential_t),
intent(in) :: ps
252 real(real64),
intent(in) :: so_strength
261 safe_allocate(p%hgh_p(0:p%lmax, -p%lmax:p%lmax))
263 if (ll == p%lloc) cycle
265 call hgh_projector_init(p%hgh_p(ll, mm), p%sphere, p%reltype, ps, ll, mm, so_strength)
270 safe_allocate(p%kb_p(0:p%lmax, -p%lmax:p%lmax))
272 if (ll == p%lloc) cycle
274 call kb_projector_init(p%kb_p(ll, mm), p%sphere, ps, ll, mm)
279 safe_allocate(p%rkb_p(1:p%lmax, -p%lmax:p%lmax))
281 if (ll == p%lloc) cycle
283 call rkb_projector_init(p%rkb_p(ll, mm), p%sphere, ps, ll, mm, so_strength)
287 if (p%lloc /= 0)
then
288 safe_allocate(p%kb_p(1, 1))
289 call kb_projector_init(p%kb_p(1, 1), p%sphere, ps, 0, 0)
305 call submesh_end(p%sphere)
310 if (ll == p%lloc) cycle
312 call hgh_projector_end(p%hgh_p(ll, mm))
315 safe_deallocate_a(p%hgh_p)
319 if (ll == p%lloc) cycle
321 call kb_projector_end(p%kb_p(ll, mm))
324 safe_deallocate_a(p%kb_p)
328 if (ll == p%lloc) cycle
330 call rkb_projector_end(p%rkb_p(ll, mm))
333 safe_deallocate_a(p%rkb_p)
334 if (p%lloc /= 0)
then
335 call kb_projector_end(p%kb_p(1, 1))
336 safe_deallocate_a(p%kb_p)
343 safe_deallocate_a(p%phase)
350#include "projector_inc.F90"
353#include "complex.F90"
354#include "projector_inc.F90"
double exp(double __x) __attribute__((__nothrow__
This module implements batches of mesh functions.
This module implements common operations on batches of mesh functions.
Module implementing boundary conditions in Octopus.
This module implements the underlying real-space grid.
This module defines the meshes, which are used in Octopus.
real(real64) function, public dprojector_matrix_element(pj, bnd, dim, ik, psia, psib)
dprojector_matrix_element calculates <psia|projector|psib>
subroutine, public projector_build(p, ps, so_strength)
subroutine, public dprojector_commute_r(pj, mesh, bnd, dim, idir, ik, psi, cpsi)
This function calculates |cpsi> += [x, V_nl] |psi>
logical elemental function, public projector_is(p, type)
subroutine, public dproject_psi(mesh, bnd, pj, npj, dim, psi, ppsi, ik)
dproject_psi calculates the action of a projector on the psi wavefunction. The result is summed up to...
subroutine, public projector_init(p, pseudo, namespace, dim, reltype)
subroutine, public projector_init_phases(this, dim, std, bnd, kpoints, vec_pot, vec_pot_var)
subroutine, public dproject_psi_batch(mesh, bnd, pj, npj, dim, psib, ppsib)
To optimize the application of the non-local operator in parallel, the projectors are applied in step...
complex(real64) function, public zprojector_matrix_element(pj, bnd, dim, ik, psia, psib)
zprojector_matrix_element calculates <psia|projector|psib>
subroutine, public zproject_sphere(mesh, pj, dim, psi, ppsi)
subroutine, public zprojector_commute_r_allatoms_alldir(pj, ions, mesh, dim, bnd, ik, psi, cpsi)
This function calculates |cpsi> += [x, V_nl] |psi>
subroutine, public projector_end(p)
subroutine, public zproject_psi(mesh, bnd, pj, npj, dim, psi, ppsi, ik)
zproject_psi calculates the action of a projector on the psi wavefunction. The result is summed up to...
subroutine, public zproject_psi_batch(mesh, bnd, pj, npj, dim, psib, ppsib)
To optimize the application of the non-local operator in parallel, the projectors are applied in step...
subroutine, public dproject_sphere(mesh, pj, dim, psi, ppsi)
subroutine, public dprojector_commute_r_allatoms_alldir(pj, ions, mesh, dim, bnd, ik, psi, cpsi)
This function calculates |cpsi> += [x, V_nl] |psi>
subroutine, public zprojector_commute_r(pj, mesh, bnd, dim, idir, ik, psi, cpsi)
This function calculates |cpsi> += [x, V_nl] |psi>
logical elemental function, public projector_is_null(p)
integer, parameter, public proj_none
This module handles spin dimensions of the states and the k-point distribution.
The projector data type is intended to hold the local and non-local parts of the pseudopotentials....