27 use,
intrinsic :: iso_fortran_env
49 integer,
public,
parameter :: &
50 SINGULARITY_NONE = 0, &
58 integer :: coulomb_singularity = 0
59 real(real64),
allocatable :: Fk(:)
61 real(real64) :: energy
67 type(singularity_t),
intent(inout) :: this
68 type(namespace_t),
intent(in) :: namespace
69 class(space_t),
intent(in) :: space
70 type(states_elec_t),
intent(in) :: st
71 type(kpoints_t),
intent(in) :: kpoints
79 if (.not.
allocated(this%Fk))
then
80 safe_allocate(this%Fk(1:st%nik))
84 if (space%is_periodic())
then
108 if (space%dim == 2 .or. space%dim > 3) default = singularity_none
110 call parse_variable(namespace,
'HFSingularity', default, this%coulomb_singularity)
113 if(this%coulomb_singularity /= singularity_none)
then
118 if(space%dim == 2)
then
123 if (this%coulomb_singularity /= singularity_none)
then
133 type(singularity_t),
intent(inout) :: this
137 safe_deallocate_a(this%Fk)
138 this%coulomb_singularity = -1
148 class(
space_t),
intent(in) :: space
149 type(states_elec_t),
intent(in) :: st
150 type(kpoints_t),
intent(in) :: kpoints
152 integer :: ik, ik2, ikpoint, nk, nsteps
153 integer :: ikx, iky, ikz, istep, kpt_start, kpt_end
154 real(real64) :: length
155 real(real64) :: kpoint(space%dim), qpoint(space%dim)
156 real(real64) :: kvol_element
158 integer :: default_nk, default_step
159 real(real64),
parameter :: singul_cnst = 7.7955541794415_real64
160 real(real64),
parameter :: log2_minus_gamma = 0.11593151565841244881_real64
169 kpt_start = st%d%kpt%start
170 kpt_end = st%d%kpt%end
172 if (.not. st%d%kpt%parallel)
then
174 kpt_start = dist_kpt%start
175 kpt_end = dist_kpt%end
178 do ik = kpt_start, kpt_end
179 ikpoint = st%d%get_kpoint_index(ik)
180 kpoint = kpoints%get_point(ikpoint, absolute_coordinates = .false.)
184 do ik2 = 1, kpoints%full%npoints
185 qpoint = kpoint - kpoints%full%red_point(:, ik2)
188 qpoint = qpoint - anint(qpoint + 1e-5_real64)
190 if (all(abs(qpoint) < 1e-6_real64)) cycle
192 this%Fk(ik) = this%Fk(ik) +
aux_funct(qpoint) * kpoints%full%weight(ik2)
194 select case(space%dim)
196 this%Fk(ik) = this%Fk(ik)/kpoints%latt%rcell_volume
198 this%Fk(ik) = this%Fk(ik)*
m_four*
m_pi/kpoints%latt%rcell_volume
202 if (st%d%kpt%parallel)
then
221 if(space%dim == 1) default_nk = 1200
224 message(1) =
'HFSingularity_Nk must be a multiple of 3.'
238 if(space%dim == 1) default_step = 15
239 call parse_variable(namespace,
'HFSingularityNsteps', default_step, nsteps)
241 select case(space%dim)
250 qpoint(1) = ikx/(
m_two*nk)*length
252 if(abs(ikx)<=nk/3) cycle
253 this%FF = this%FF +
aux_funct(qpoint)*kvol_element
255 if(istep<nsteps)
then
257 kvol_element = kvol_element/
m_three
267 this%FF = this%FF +
m_two * (
m_pi)/kpoints%latt%rcell_volume * length &
268 * (
m_one-
log(
m_pi / kpoints%latt%rcell_volume * length) + log2_minus_gamma)
279 qpoint(1) = ikx/(
m_two*nk)*length
282 qpoint(2) = iky/(
m_two*nk)*length
285 qpoint(3) = ikz/(
m_two*nk)*length
287 if (abs(ikx) <= nk/3 .and. abs(iky) <= nk/3 .and. abs(ikz) <= nk/3) cycle
289 this%FF = this%FF +
aux_funct(qpoint)*kvol_element
293 if (istep < nsteps)
then
295 kvol_element = kvol_element / 27.0_real64
303 this%FF = this%FF + singul_cnst*(kpoints%latt%rcell_volume)**(
m_twothird)/
m_pi/kpoints%latt%rcell_volume*length
310 this%FF = 4.423758_real64*(kpoints%latt%rcell_volume*
m_four)**(
m_twothird)/
m_pi/kpoints%latt%rcell_volume
315 this%FF = singul_cnst*(kpoints%latt%rcell_volume)**(
m_twothird)/
m_pi/kpoints%latt%rcell_volume
320 do ik = st%d%kpt%start, st%d%kpt%end
321 this%energy = this%energy + this%Fk(ik)*st%kweights(ik)
324 if (st%d%kpt%parallel)
then
331 this%energy = (this%energy-this%FF)*st%qtot
333 write(
message(1),
'(a,f12.6,a,a,a)')
'Debug: Singularity energy ', &
343 real(real64) function aux_funct(qq) result(ff)
344 real(real64),
intent(in) :: qq(3)
346 real(real64) :: half_a, qq_abs(space%dim)
351 select case(space%dim)
354 ff = -
log(abs(
sin(qq(1)*
m_pi))*kpoints%latt%klattice(1,1)/(
m_two*
m_pi)) + 0.11593151565841244881_real64
358 (
m_two*
sin(qq(1)*
m_pi)*
sin(qq(1)*
m_pi)*dot_product(kpoints%latt%klattice(:,1), kpoints%latt%klattice(:,1)) &
359 +
sin(qq(1)*
m_two*
m_pi)*
sin(qq(2)*
m_two*
m_pi)*dot_product(kpoints%latt%klattice(:,1), kpoints%latt%klattice(:,2))) &
360 + (
m_two*
sin(qq(2)*
m_pi)*
sin(qq(2)*
m_pi)*dot_product(kpoints%latt%klattice(:,2), kpoints%latt%klattice(:,2)) &
361 +
sin(qq(2)*
m_two*
m_pi)*
sin(qq(3)*
m_two*
m_pi)*dot_product(kpoints%latt%klattice(:,2), kpoints%latt%klattice(:,3))) &
362 + (
m_two*
sin(qq(3)*
m_pi)*
sin(qq(3)*
m_pi)*dot_product(kpoints%latt%klattice(:,3), kpoints%latt%klattice(:,3)) &
363 +
sin(qq(3)*
m_two*
m_pi)*
sin(qq(1)*
m_two*
m_pi)*dot_product(kpoints%latt%klattice(:,3), kpoints%latt%klattice(:,1)))))
369 ff = (half_a)**2/(
m_three-
cos(qq_abs(1)*half_a)*
cos(qq_abs(2)*half_a) &
370 -
cos(qq_abs(1)*half_a)*
cos(qq_abs(3)*half_a) &
371 -
cos(qq_abs(3)*half_a)*
cos(qq_abs(2)*half_a))
double log(double __x) __attribute__((__nothrow__
double sin(double __x) __attribute__((__nothrow__
double cos(double __x) __attribute__((__nothrow__
subroutine, public distributed_end(this)
subroutine, public distributed_init(this, total, comm, tag, scalapack_compat)
Distribute N instances across M processes of communicator comm
integer, parameter, public spin_polarized
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
real(real64), parameter, public m_four
real(real64), parameter, public m_third
real(real64), parameter, public m_pi
some mathematical constants
real(real64), parameter, public m_twothird
real(real64), parameter, public m_epsilon
real(real64), parameter, public m_half
real(real64), parameter, public m_one
real(real64), parameter, public m_three
subroutine, public kpoints_to_absolute(latt, kin, kout)
subroutine, public messages_not_implemented(feature, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
type(mpi_grp_t), public mpi_world
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
subroutine singularity_correction(this, namespace, space, st, kpoints)
integer, parameter, public singularity_general
subroutine, public singularity_end(this)
integer, parameter, public singularity_sphere
subroutine, public singularity_init(this, namespace, space, st, kpoints)
integer, parameter, public singularity_gygi
This module handles spin dimensions of the states and the k-point distribution.
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
character(len=20) pure function, public units_abbrev(this)
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
real(real64) function aux_funct(qq)
Distribution of N instances over mpi_grpsize processes, for the local rank mpi_grprank....