Octopus
singularity.F90
Go to the documentation of this file.
1!! Copyright (C) 2019-2022 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
19#include "global.h"
20
22 use comm_oct_m
23 use debug_oct_m
26 use global_oct_m
27 use, intrinsic :: iso_fortran_env
31 use mpi_oct_m
32 use parser_oct_m
35 use space_oct_m
38 use unit_oct_m
40
41 implicit none
42
43 private
44 public :: &
48
49 integer, public, parameter :: &
50 SINGULARITY_NONE = 0, &
52 singularity_gygi = 2, &
54
55
56 type singularity_t
57 !For the treatment of the singularity in solids
58 integer :: coulomb_singularity = 0
59 real(real64), allocatable :: Fk(:)
60 real(real64) :: FF
61 real(real64) :: energy
62 end type singularity_t
63
64contains
65
66 subroutine singularity_init(this, namespace, space, st, kpoints)
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
72
73 integer :: default
74
75 push_sub(singularity_init)
76
77 this%energy = m_zero
78
79 if (.not. allocated(this%Fk)) then
80 safe_allocate(this%Fk(1:st%nik))
81 this%Fk(1:st%nik) = m_zero
82 this%FF = m_zero
83
84 if (space%is_periodic()) then
85
86 !%Variable HFSingularity
87 !%Type integer
88 !%Default general
89 !%Section Hamiltonian::XC
90 !%Description
91 !% (Experimental) This variable selects the method used for the treatment of the
92 !% singularity of the Coulomb potential in Hatree-Fock and hybrid-functional DFT calculations.
93 !% This shoulbe be only applied for periodic systems and is only
94 !% used for FFT kernels of the Poisson solvers.
95 !%Option none 0
96 !% The singularity is replaced by zero.
97 !%Option general 1
98 !% The general treatment of the singularity, as described in Carrier et al, PRB 75 205126 (2007).
99 !% This is the default option
100 !%Option fcc 2
101 !% The treatment of the singulariy as described in Gygi and Baldereschi, PRB 34, 4405 (1986).
102 !% This is formally valid for cubic systems only.
103 !%Option spherical_bz 3
104 !% The divergence in q=0 is treated analytically assuming a spherical Brillouin zone
105 !%End
106
107 default = singularity_general
108 if (space%dim == 2 .or. space%dim > 3) default = singularity_none
109
110 call parse_variable(namespace, 'HFSingularity', default, this%coulomb_singularity)
111 call messages_print_var_option('HFSingularity', this%coulomb_singularity, namespace=namespace)
112
113 if(this%coulomb_singularity /= singularity_none) then
114 if(this%coulomb_singularity /= singularity_general .and. space%dim == 1) then
115 call messages_not_implemented("HFSingularity /= general for 1D")
116 end if
117
118 if(space%dim == 2) then
119 call messages_not_implemented("HFSingularity /= none for 2D")
120 end if
121 end if
122
123 if (this%coulomb_singularity /= singularity_none) then
124 call singularity_correction(this, namespace, space, st, kpoints)
125 end if
126 end if
127 end if
128
129 pop_sub(singularity_init)
130 end subroutine singularity_init
131
132 subroutine singularity_end(this)
133 type(singularity_t), intent(inout) :: this
134
135 push_sub(singularity_end)
136
137 safe_deallocate_a(this%Fk)
138 this%coulomb_singularity = -1
139
140 pop_sub(singularity_end)
141 end subroutine singularity_end
143 !This routine implements the general tratment of the singularity for periodic solids,
144 !as described in Carrier et al. PRB 75, 205126 (2007)
145 subroutine singularity_correction(this, namespace, space, st, kpoints)
146 type(singularity_t), intent(inout) :: this
147 type(namespace_t), intent(in) :: namespace
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
157 type(distributed_t) :: dist_kpt
158 integer :: default_nk, default_step
159 real(real64), parameter :: singul_cnst = 7.7955541794415_real64 !The constant is 4*pi*(3/(4*pi))^1/3
160 real(real64), parameter :: log2_minus_gamma = 0.11593151565841244881_real64 !The constant is log(2) - gamma, where gamma is the Euler-Mascheroni constant
161 push_sub(singularity_correction)
162
163 call profiling_in("COULOMB_SINGULARITY")
164
165 !At the moment this is only implemented in 3D and in 1D.
166 assert(space%dim == 3 .or. this%coulomb_singularity == singularity_general)
167
168 kpt_start = st%d%kpt%start
169 kpt_end = st%d%kpt%end
170
171 if (.not. st%d%kpt%parallel) then
172 call distributed_init(dist_kpt, st%nik, mpi_world%comm, "singularity")
173 kpt_start = dist_kpt%start
174 kpt_end = dist_kpt%end
175 end if
176
177 do ik = kpt_start, kpt_end
178 ikpoint = st%d%get_kpoint_index(ik)
179 kpoint = kpoints%get_point(ikpoint, absolute_coordinates = .false.)
180
181 this%Fk(ik) = m_zero
182
183 do ik2 = 1, kpoints%full%npoints
184 qpoint = kpoint - kpoints%full%red_point(:, ik2)
185
186 !We remove potential umklapp
187 qpoint = qpoint - anint(qpoint + 1e-5_real64)
188
189 if (all(abs(qpoint) < 1e-6_real64)) cycle
190
191 this%Fk(ik) = this%Fk(ik) + aux_funct(qpoint) * kpoints%full%weight(ik2)
192 end do
193 select case(space%dim)
194 case(1)
195 this%Fk(ik) = this%Fk(ik)/kpoints%latt%rcell_volume
196 case(3)
197 this%Fk(ik) = this%Fk(ik)*m_four*m_pi/kpoints%latt%rcell_volume
198 end select
199 end do
200
201 if (st%d%kpt%parallel) then
202 call comm_allreduce(st%d%kpt%mpi_grp, this%Fk)
203 else
204 call comm_allreduce(dist_kpt%mpi_grp, this%Fk)
205 call distributed_end(dist_kpt)
206 end if
207
208 if (this%coulomb_singularity == singularity_general) then
209 !%Variable HFSingularityNk
210 !%Type integer
211 !%Default 60 in 3D, 1200 in 1D
212 !%Section Hamiltonian::XC
213 !%Description
214 !% Number of k-point used (total number of k-points) is (2*Nk+1)^3) in the numerical integration
215 !% of the auxiliary function f(q). See PRB 75, 205126 (2007) for more details.
216 !% Only for HFSingularity=general.
217 !% Also used in 1D.
218 !%End
219 default_nk = 60
220 if(space%dim == 1) default_nk = 1200
221 call parse_variable(namespace, 'HFSingularityNk', default_nk, nk)
222 if (nk <= 0) then
223 message(1) = 'HFSingularityNk must be > 0.'
224 call messages_fatal(1, namespace=namespace)
225 end if
226 if (mod(nk, 3) /= 0) then
227 message(1) = 'HFSingularityNk must be a multiple of 3.'
228 call messages_fatal(1, namespace=namespace)
229 end if
230
231 !%Variable HFSingularityNsteps
232 !%Type integer
233 !%Default 7 in 3D, 15 in 1D
234 !%Section Hamiltonian::XC
235 !%Description
236 !% Number of grid refinement steps in the numerical integration of the auxiliary function f(q).
237 !% See PRB 75, 205126 (2007) for more details. Only for HFSingularity=general.
238 !% Also used in 1D.
239 !%End
240 default_step = 7
241 if(space%dim == 1) default_step = 15
242 call parse_variable(namespace, 'HFSingularityNsteps', default_step, nsteps)
243
244 select case(space%dim)
245 case(1)
246 this%FF = m_zero
247 length = m_one
248 kvol_element = (m_one/(m_two*nk+m_one))*((m_two*m_pi))/kpoints%latt%rcell_volume
249 qpoint = m_zero
250 do istep = 1, nsteps
251
252 do ikx = 0, nk
253 qpoint(1) = ikx/(m_two*nk)*length
254
255 if(abs(ikx)<=nk/3) cycle
256 this%FF = this%FF + aux_funct(qpoint)*kvol_element
257 end do
258 if(istep<nsteps) then
259 length = length/m_three
260 kvol_element = kvol_element/m_three
261 end if
262 end do
263
264 !We have a factor two because we used the fact that f(q)=f(-q)
265 !We multiply by 1/((2*pi)^1)
266 this%FF = this%FF*m_two/((m_two*m_pi))
267
268 !We add the remaining part
269 !The constant is log(2) - gamma, where gamma is the Euler-Mascheroni constant
270 this%FF = this%FF + m_two * (m_pi)/kpoints%latt%rcell_volume * length &
271 * (m_one-log(m_pi / kpoints%latt%rcell_volume * length) + log2_minus_gamma)
272 case(2)
273 call messages_not_implemented("General Coulomb singularity in the 2D case")
274
275 case(3)
276 this%FF = m_zero
277 length = m_one
278 kvol_element = (m_one/(m_two*nk+m_one))**3*((m_two*m_pi)**3)/kpoints%latt%rcell_volume
279 do istep = 1, nsteps
280
281 do ikx = 0, nk
282 qpoint(1) = ikx/(m_two*nk)*length
283
284 do iky = -nk, nk
285 qpoint(2) = iky/(m_two*nk)*length
286
287 do ikz = -nk, nk
288 qpoint(3) = ikz/(m_two*nk)*length
289
290 if (abs(ikx) <= nk/3 .and. abs(iky) <= nk/3 .and. abs(ikz) <= nk/3) cycle
291
292 this%FF = this%FF + aux_funct(qpoint)*kvol_element
293 end do
294 end do
295 end do
296 if (istep < nsteps) then
297 length = length / m_three
298 kvol_element = kvol_element / 27.0_real64
299 end if
300 end do
301
302 !We have a factor two because we used the fact that f(q)=f(-q)
303 !We multiply by 4*pi/((2*pi)^3)
304 this%FF = this%FF*8.0_real64*m_pi/((m_two*m_pi)**3)
305 !The remaining part is treated as a spherical BZ
306 this%FF = this%FF + singul_cnst*(kpoints%latt%rcell_volume)**(m_twothird)/m_pi/kpoints%latt%rcell_volume*length
307
308 end select
309
310 else if (this%coulomb_singularity == singularity_gygi) then
311 !See Eq. (7) of PRB 34, 4405 (1986)
312 !Here we use the fact that the fcc volume is a^3/4
313 this%FF = 4.423758_real64*(kpoints%latt%rcell_volume*m_four)**(m_twothird)/m_pi/kpoints%latt%rcell_volume
314
315 else
316 !The constant is 4*pi*(3/(4*pi))^1/3
317 !We multiply by 4*pi/(2*pi^3)
318 this%FF = singul_cnst*(kpoints%latt%rcell_volume)**(m_twothird)/m_pi/kpoints%latt%rcell_volume
319 end if
320
321
322 this%energy = m_zero
323 do ik = st%d%kpt%start, st%d%kpt%end
324 this%energy = this%energy + this%Fk(ik)*st%kweights(ik)
325 end do
326
327 if (st%d%kpt%parallel) then
328 call comm_allreduce(st%d%kpt%mpi_grp, this%energy)
329 end if
330
331 ! In the case, the number of k-point is doubled, so the energy is twice larger
332 if (st%d%ispin == spin_polarized) this%energy = m_half * this%energy
333
334 this%energy = (this%energy-this%FF)*st%qtot
335
336 write(message(1), '(a,f12.6,a,a,a)') 'Debug: Singularity energy ', &
337 units_from_atomic(units_out%energy, this%energy), &
338 ' [',trim(units_abbrev(units_out%energy)),']'
339 call messages_info(1, namespace=namespace, debug_only=.true.)
340
341 call profiling_out("COULOMB_SINGULARITY")
343
344 contains
345
346 real(real64) function aux_funct(qq) result(ff)
347 real(real64), intent(in) :: qq(space%dim)
348
349 real(real64) :: half_a, qq_abs(space%dim)
350
351 ! no PUSH/POP as called too often
352
353 if (this%coulomb_singularity == singularity_general) then
354 select case(space%dim)
355 case(1)
356 !The constant is -log(2) + gamma, where gamma is the Euler-Mascheroni constant
357 ff = -log(abs(sin(qq(1)*m_pi))*kpoints%latt%klattice(1,1)/(m_two*m_pi)) + 0.11593151565841244881_real64
358 case(3)
359 !See Eq. (16) of PRB 75, 205126 (2007)
360 ff = (m_two*m_pi)**2/(m_two*( &
361 (m_two*sin(qq(1)*m_pi)*sin(qq(1)*m_pi)*dot_product(kpoints%latt%klattice(:,1), kpoints%latt%klattice(:,1)) &
362 + sin(qq(1)*m_two*m_pi)*sin(qq(2)*m_two*m_pi)*dot_product(kpoints%latt%klattice(:,1), kpoints%latt%klattice(:,2))) &
363 + (m_two*sin(qq(2)*m_pi)*sin(qq(2)*m_pi)*dot_product(kpoints%latt%klattice(:,2), kpoints%latt%klattice(:,2)) &
364 + sin(qq(2)*m_two*m_pi)*sin(qq(3)*m_two*m_pi)*dot_product(kpoints%latt%klattice(:,2), kpoints%latt%klattice(:,3))) &
365 + (m_two*sin(qq(3)*m_pi)*sin(qq(3)*m_pi)*dot_product(kpoints%latt%klattice(:,3), kpoints%latt%klattice(:,3)) &
366 + sin(qq(3)*m_two*m_pi)*sin(qq(1)*m_two*m_pi)*dot_product(kpoints%latt%klattice(:,3), kpoints%latt%klattice(:,1)))))
367 end select
368 else
369 half_a = m_half*(kpoints%latt%rcell_volume*m_four)**(m_third)
370 call kpoints_to_absolute(kpoints%latt, qq, qq_abs)
371 !See Eq. (6) of PRB 34, 4405 (1986)
372 ff = (half_a)**2/(m_three-cos(qq_abs(1)*half_a)*cos(qq_abs(2)*half_a) &
373 -cos(qq_abs(1)*half_a)*cos(qq_abs(3)*half_a) &
374 -cos(qq_abs(3)*half_a)*cos(qq_abs(2)*half_a))
375 end if
376
377 end function aux_funct
378
379 end subroutine singularity_correction
380
381end module singularity_oct_m
382
383!! Local Variables:
384!! mode: f90
385!! coding: utf-8
386!! End:
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
Definition: global.F90:190
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_four
Definition: global.F90:192
real(real64), parameter, public m_third
Definition: global.F90:195
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:186
real(real64), parameter, public m_twothird
Definition: global.F90:196
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_one
Definition: global.F90:189
real(real64), parameter, public m_three
Definition: global.F90:191
subroutine, public kpoints_to_absolute(latt, kin, kout)
Definition: kpoints.F90:1031
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1114
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
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 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.
Definition: unit.F90:132
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:223
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....
int true(void)