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