Octopus
smear.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
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! Some pieces glued from PW/w(0,1)gauss.f90 from PWSCF
20
21#include "global.h"
22
23module smear_oct_m
24 use debug_oct_m
25 use global_oct_m
26 use, intrinsic :: iso_fortran_env
31 use parser_oct_m
34 use sort_oct_m
35 use unit_oct_m
38
39 implicit none
40
41 private
42 public :: &
43 smear_t, &
44 smear_init, &
45 smear_copy, &
54
55 type smear_t
56 private
57 integer, public :: method
58 real(real64), public :: dsmear
59 real(real64) :: dsmear_cond
60 real(real64), public :: e_fermi
61 real(real64) :: e_fermi_cond
62
63 integer, public :: el_per_state
64 real(real64), public :: ef_occ
65 real(real64) :: ef_occ_cond
66 logical, public :: integral_occs
67 integer, public :: MP_n
68 integer :: fermi_count
69 integer :: nik_factor
70 integer :: nspins
71
72 type(simplex_t), pointer :: tetra_simplex => null()
73
74 ! PhotonDoping
75 logical, public:: photodop
76 real(real64) :: nephotodop
77 integer :: photodop_bandmin
78
79 end type smear_t
80
81 integer, parameter, public :: &
82 SMEAR_SEMICONDUCTOR = 1, &
84 smear_cold = 3, &
86 smear_spline = 5, &
87 smear_fixed_occ = 6, &
88 smear_gaussian = option__smearingfunction__gaussian_smearing, &
89 smear_lorentzian = option__smearingfunction__lorentzian_smearing, &
90 smear_tetrahedra = option__smearingfunction__tetrahedra, &
91 smear_tetrahedra_opt = option__smearingfunction__tetrahedra_opt
92
93 real(real64), parameter :: TOL_SMEAR = 1e-6_real64
94
95contains
96
97 !--------------------------------------------------
98 subroutine smear_init(this, namespace, ispin, fixed_occ, integral_occs, kpoints)
99 type(smear_t), intent(out) :: this
100 type(namespace_t), intent(in) :: namespace
101 integer, intent(in) :: ispin
102 logical, intent(in) :: fixed_occ
103 logical, intent(in) :: integral_occs
104 type(kpoints_t), intent(in) :: kpoints
105
106 push_sub(smear_init)
107
108 this%integral_occs = integral_occs
109
110 !%Variable SmearingFunction
111 !%Type integer
112 !%Default semiconducting
113 !%Section States
114 !%Description
115 !% This is the function used to smear the electronic occupations.
116 !% It is ignored if the <tt>Occupations</tt> block is set.
117 !%Option semiconducting 1
118 !% Semiconducting occupations, <i>i.e.</i> the lowest lying states are occupied
119 !% until no more electrons are left.
120 !%Option fermi_dirac 2
121 !% Simple Fermi-Dirac distribution. In this case, <tt>Smearing</tt> has
122 !% the meaning of an electronic temperature. DN Mermin, <i>Phys. Rev.</i> <b>137</b>, A1441 (1965).
123 !%Option cold_smearing 3
124 !% N Marzari, D Vanderbilt, A De Vita, and MC Payne, <i>Phys. Rev. Lett.</i> <b>82</b>, 3296 (1999).
125 !%Option methfessel_paxton 4
126 !% M Methfessel and AT Paxton, <i>Phys. Rev. B</i> <b>40</b>, 3616 (1989).
127 !% In this case, the variable <tt>SmearingMPOrder</tt> sets the order of the smearing.
128 !% Occupations may be negative.
129 !%Option spline_smearing 5
130 !% Nearly identical to Gaussian smearing.
131 !% JM Holender, MJ Gillan, MC Payne, and AD Simpson, <i>Phys. Rev. B</i> <b>52</b>, 967 (1995).
132 !%Option gaussian_smearing
133 !% Gaussian smearing.
134 !%Option lorentzian_smearing
135 !% Lorentzian smearing.
136 !%Option tetrahedra
137 !% Linear tetrahedron method as in P. E. Bloechl, et al., <i>Phys. Rev. B</i> <b>49</b>, 16223 (1994).
138 !% Requires a regular Monkhorst-Pack generated by Octopus.
139 !%Option tetrahedra_opt
140 !% Improved tetrahedron method as in M. Kawamura, et al., <i>Phys. Rev. B</i> <b>89</b>, 094515 (2014).
141 !% Requires a regular Monkhorst-Pack generated by Octopus.
142 !%End
143 if (fixed_occ) then
144 this%method = smear_fixed_occ
145 else
146 call parse_variable(namespace, 'SmearingFunction', smear_semiconductor, this%method)
147 if (.not. varinfo_valid_option('SmearingFunction', this%method)) call messages_input_error(namespace, 'SmearingFunction')
148 call messages_print_var_option('SmearingFunction', this%method, namespace=namespace)
149 end if
151 !%Variable Smearing
152 !%Type float
153 !%Default 0.1 eV
154 !%Section States
155 !%Description
156 !% If <tt>Occupations</tt> is not set, <tt>Smearing</tt> is the
157 !% smearing width used in the <tt>SmearingFunction</tt> to distribute the electrons
158 !% among the existing states.
159 !%End
160 this%dsmear = 1e-14_real64
161 if (this%method /= smear_semiconductor .and. this%method /= smear_fixed_occ) then
162 call parse_variable(namespace, 'Smearing', 0.1_real64 / (m_two * p_ry), this%dsmear, units_inp%energy)
163 end if
165 !%Variable PhotoDopingSmearing
166 !%Type float
167 !%Default 0.1 eV
168 !%Section States
169 !%Description
170 !% If <tt>Occupations</tt> is not set, <tt>PhotoDopingSmearing</tt> is the
171 !% smearing width used in the <tt>SmearingFunction</tt> to distribute the electrons in conduction bands
172 !% among the existing states.
173 !%End
174 this%dsmear_cond = this%dsmear
175 if (this%method /= smear_semiconductor .and. this%method /= smear_fixed_occ) then
176 call parse_variable(namespace, 'PhotoDopingSmearing', 0.1_real64 / (m_two * p_ry), this%dsmear_cond, units_inp%energy)
177 end if
178
179 !%Variable PhotoDopingNumElectrons
180 !%Type float
181 !%Default 0
182 !%Section States
183 !%Description
184 !% This variable allows to set the number of valence electrons taken out to put
185 !% into the conduction bands. The method follows Marini and Calandra, PRB 104, 144103 (2021).
186 !%End
187 call parse_variable(namespace, 'PhotoDopingNumElectrons', m_zero, this%nephotodop)
189 if (this%nephotodop > m_min_occ) then
190 this%photodop=.true.
191 else
192 this%photodop=.false.
193 endif
194
195
196 !%Variable PhotoDopingBand
197 !%Type integer
198 !%Default 1
199 !%Section States
200 !%Description
201 !% This variable specifies where the minium (starting) conduction band index for the photodoping electrons.
202 !%End
203 call parse_variable(namespace, 'PhotoDopingBand', 1, this%photodop_bandmin)
204
205
206 call messages_obsolete_variable(namespace, 'ElectronicTemperature', 'Smearing')
207
208 if (ispin == 1) then ! unpolarized
209 this%el_per_state = 2
210 else
211 this%el_per_state = 1
212 end if
213
214 if (ispin == 2) then
215 this%nspins = 2
216 else
217 this%nspins = 1
218 end if
219
220 if (this%method == smear_semiconductor) then
221 this%nik_factor = kpoints_kweight_denominator(kpoints)
222
223 if (this%nik_factor == 0) then
224 message(1) = "k-point weights in KPoints or KPointsReduced blocks must be rational numbers for semiconducting smearing."
225 call messages_fatal(1, namespace=namespace)
226 end if
227 else if (this%method == smear_tetrahedra .or. this%method == smear_tetrahedra_opt) then
228 if (.not. associated(kpoints%reduced%simplex)) then
229 message(1) = "Tetrahedron smearing requires a simplex mesh for the k-point grid."
230 call messages_fatal(1, namespace=namespace)
231 end if
232
233 this%tetra_simplex => kpoints%reduced%simplex
234 end if
235
236 this%MP_n = 0
237 if (this%method == smear_methfessel_paxton) then
238 !%Variable SmearingMPOrder
239 !%Type integer
240 !%Default 1
241 !%Section States
242 !%Description
243 !% Sets the order of the Methfessel-Paxton smearing function.
244 !%End
245 call parse_variable(namespace, 'SmearingMPOrder', 1, this%MP_n)
246 end if
247
248 pop_sub(smear_init)
249 end subroutine smear_init
250
251
252 !--------------------------------------------------
253 subroutine smear_copy(to, from)
254 type(smear_t), intent(out) :: to
255 type(smear_t), intent(in) :: from
256
257 push_sub(smear_copy)
258
259 to%method = from%method
260 to%dsmear = from%dsmear
261 to%e_fermi = from%e_fermi
262 to%el_per_state = from%el_per_state
263 to%ef_occ = from%ef_occ
264 to%MP_n = from%MP_n
265 to%fermi_count = from%fermi_count
266 to%nik_factor = from%nik_factor
267 to%tetra_simplex => from%tetra_simplex
268
269 pop_sub(smear_copy)
270 end subroutine smear_copy
271
272
273 !--------------------------------------------------
274 subroutine smear_find_fermi_energy(this, namespace, eigenvalues, occupations, &
275 qtot, nik, nst, kweights)
276 type(smear_t), intent(inout) :: this
277 type(namespace_t), intent(in) :: namespace
278 real(real64), intent(in) :: eigenvalues(:,:), occupations(:,:)
279 real(real64), intent(in) :: qtot, kweights(:)
280 integer, intent(in) :: nik, nst
281
282 real(real64), parameter :: tol = 1.0e-10_real64
283 integer :: ist, ik, iter, maxq, weight, sumq_int, sum_weight
284 real(real64) :: sumq_frac
285 logical :: conv
286 real(real64), allocatable :: eigenval_list(:)
287 integer, allocatable :: k_list(:), reorder(:)
288 integer :: fermi_count_up, fermi_count_down
289
291
292 maxq = this%el_per_state * nst * this%nspins
293 if (maxq - qtot <= -tol) then ! not enough states
294 message(1) = 'Not enough states'
295 write(message(2),'(6x,a,f12.6,a,i10)')'(total charge = ', qtot, &
296 ' max charge = ', maxq
297 call messages_fatal(2, namespace=namespace)
298 end if
299
300 conv = .true.
301 if (this%method == smear_fixed_occ) then ! Fermi energy: last of occupied states
302 ist_cycle: do ist = nst, 1, -1
303 if (any(occupations(ist, :) > m_min_occ)) then
304 this%e_fermi = eigenvalues(ist, 1)
305 this%ef_occ = occupations(ist, 1) / this%el_per_state
306 do ik = 2, nik
307 if (eigenvalues(ist, ik) > this%e_fermi .and. occupations(ist, ik) > m_min_occ) then
308 this%e_fermi = eigenvalues(ist, ik)
309 this%ef_occ = occupations(ist, ik) / this%el_per_state
310 end if
311 end do
312 exit ist_cycle
313 end if
314 end do ist_cycle
315
316 else if (this%method == smear_semiconductor) then
317 ! first we sort the eigenvalues
318 safe_allocate(eigenval_list(1:nst * nik))
319 safe_allocate( k_list(1:nst * nik))
320 safe_allocate( reorder(1:nst * nik))
321
322 iter = 1
323 do ist = 1, nst
324 do ik = 1, nik
325 eigenval_list(iter) = eigenvalues(ist, ik)
326 k_list(iter) = ik
327 reorder(iter) = iter
328 iter = iter + 1
329 end do
330 end do
331
332 call sort(eigenval_list, reorder)
333
334 sumq_int = floor(qtot) * this%nik_factor
335 sumq_frac = qtot * this%nik_factor - sumq_int
336 if ( sumq_frac + tol > this%nik_factor ) then
337 sumq_int = sumq_int + this%nik_factor
338 sumq_frac = sumq_frac - this%nik_factor
339 end if
340 if (sumq_frac < tol) sumq_frac = m_zero
341
342 do iter = 1, nst * nik
343 weight = int(kweights(k_list(reorder(iter))) * this%nik_factor + m_half)
344 if (.not. weight > 0) cycle
345 this%e_fermi = eigenval_list(iter)
346 this%ef_occ = (sumq_int + sumq_frac) / (weight * this%el_per_state)
347
348 if (sumq_int - weight * this%el_per_state <= 0) then
349 ! count how many occupied states are at the fermi level,
350 ! this is required later to fill the states
351 this%fermi_count = 1
352 fermi_count_up = 1
353 fermi_count_down = 1
354 sum_weight = weight
355 do
356 if (iter - fermi_count_down < 1) exit
357 if (abs(this%e_fermi - eigenval_list(iter - fermi_count_down)) > tol_smear) exit
358 weight = int(kweights(k_list(reorder(iter - fermi_count_down))) * this%nik_factor + m_half)
359 if (weight > m_epsilon) then
360 sumq_int = sumq_int + weight * this%el_per_state
361 sum_weight = sum_weight + weight
362 end if
363 fermi_count_down = fermi_count_down + 1
364 end do
365 do
366 if (iter + fermi_count_up > nst*nik) exit
367 if (abs(this%e_fermi - eigenval_list(iter + fermi_count_up)) > tol_smear) exit
368 weight = int(kweights(k_list(reorder(iter + fermi_count_up))) * this%nik_factor + m_half)
369 if (weight > m_epsilon) then
370 sum_weight = sum_weight + weight
371 end if
372 fermi_count_up = fermi_count_up + 1
373 end do
374 this%fermi_count = fermi_count_up + fermi_count_down - 1
375 this%ef_occ = (sumq_int + sumq_frac) / (sum_weight * this%el_per_state)
376 exit
377 end if
378
379 sumq_int = sumq_int - weight * this%el_per_state
380 end do
381 assert(this%ef_occ < m_one + 10.0_real64*m_epsilon)
382
383 safe_deallocate_a(eigenval_list)
384 safe_deallocate_a(k_list)
385 safe_deallocate_a(reorder)
386
387 else if (this%method == smear_tetrahedra .or. this%method == smear_tetrahedra_opt) then
388 if (this%photodop) then
389 call messages_not_implemented('Photodoping for SmearingFunction = tetrahedra/tetrahedra_opt')
390 end if
391 call bisection_find_fermi_energy_tetra(this, namespace, tol, eigenvalues, nik, nst, qtot, this%e_fermi)
392 else ! bisection
393 if (this%photodop) then
394 ! find fermi energy for valence electrons
395 call bisection_find_fermi_energy(this, namespace, this%dsmear, tol, &
396 eigenvalues, kweights, nik, qtot-this%nephotodop, 1, this%photodop_bandmin-1, &
397 this%e_fermi, this%ef_occ)
398
399 ! find fermi energy for conduction electrons
400 call bisection_find_fermi_energy(this, namespace, this%dsmear_cond, tol, &
401 eigenvalues, kweights, nik, this%nephotodop, this%photodop_bandmin, nst, &
402 this%e_fermi_cond, this%ef_occ_cond)
403 else
404 call bisection_find_fermi_energy(this, namespace, this%dsmear, tol, &
405 eigenvalues, kweights, nik, qtot, 1, nst, this%e_fermi, this%ef_occ)
406 end if
407 end if
408
410 end subroutine smear_find_fermi_energy
411
412 ! ---------------------------------------------------------
413 subroutine bisection_find_fermi_energy(this, namespace, dsmear_in, tol, &
414 eigenvalues, kweights, nik, q_in, start_band, end_band, e_fermi, ef_occ)
415 type(smear_t), intent(inout) :: this
416 type(namespace_t), intent(in) :: namespace
417 real(real64), intent(in) :: dsmear_in
418 real(real64), intent(in) :: tol
419 integer, intent(in) :: nik
420 real(real64), intent(in) :: eigenvalues(:,:)
421 real(real64), intent(in) :: kweights(:), q_in
422 integer, intent(in) :: start_band, end_band
423 real(real64), intent(out) :: e_fermi, ef_occ
424
425 integer, parameter :: nitmax = 200
426 integer :: ist, ik, iter
427 real(real64) :: drange, dsmear, emin, emax, xx, sumq
428 logical :: conv
429
431
432 dsmear = max(1e-14_real64, dsmear_in)
433 drange = dsmear * sqrt(-log(tol * 0.01_real64))
434
435 emin = minval(eigenvalues) - drange
436 emax = maxval(eigenvalues) + drange
437
438 do iter = 1, nitmax
439 e_fermi = m_half * (emin + emax)
440 sumq = m_zero
441
442 do ik = 1, nik
443 do ist = start_band, end_band
444 xx = (e_fermi - eigenvalues(ist, ik))/dsmear
445 sumq = sumq + kweights(ik) * this%el_per_state * &
446 smear_step_function(this, xx)
447 end do
448 end do
449
450 conv = (abs(sumq - q_in) <= tol)
451 if (conv) exit
452
453 if (sumq <= q_in) emin = e_fermi
454 if (sumq >= q_in) emax = e_fermi
455
456 ef_occ = smear_step_function(this, m_zero)
457 end do
458
459 if (.not. conv) then
460 message(1) = 'Fermi: did not converge.'
461 call messages_fatal(1, namespace=namespace)
462 end if
463
465 end subroutine bisection_find_fermi_energy
466
467 ! ---------------------------------------------------------
487 subroutine bisection_find_fermi_energy_tetra(this, namespace, tol, eigenvalues, nik, nst, q_in, e_fermi)
488 type(smear_t), intent(inout) :: this
489 type(namespace_t), intent(in) :: namespace
490 real(real64), intent(in) :: tol
491 integer, intent(in) :: nik, nst
492 real(real64), intent(in) :: eigenvalues(:,:)
493 real(real64), intent(in) :: q_in
494 real(real64), intent(out) :: e_fermi
495
496 integer, parameter :: nitmax = 200
497 integer :: iter
498 real(real64) :: emin, emax, sumq
499 logical :: conv
500
501 integer :: ist, ii, ll, is, ik, ns, nik_red, nvertices
502 real(real64) :: e_simplex(20), w_simplex(4), dos_simplex(4)
503
504 push_sub_with_profile(bisection_find_fermi_energy_tetra)
505
506 assert(associated(this%tetra_simplex))
507 assert(this%tetra_simplex%n_simplices > 0)
508 assert(this%tetra_simplex%sdim > 0)
509
510 ns = this%nspins
511 assert(ns > 0)
512 assert(mod(nik, ns) == 0)
513
514 nik_red = nik / ns
515 assert(maxval(this%tetra_simplex%simplices) <= nik_red)
516
517 nvertices = this%tetra_simplex%rdim + 1
518 assert(nvertices > 0)
519 assert(nvertices <= this%tetra_simplex%sdim)
520
521 emin = minval(eigenvalues)
522 emax = maxval(eigenvalues)
523 if (.not. (emin < emax)) then
524 write(message(1), '(a)') 'Tetra bisection: invalid energy bracket from eigenvalues.'
525 call messages_fatal(1, namespace=namespace)
526 end if
527
528 do iter = 1, nitmax
529 e_fermi = m_half * (emin + emax)
530 sumq = m_zero
531 !$omp parallel do collapse(2) private(ist, ii, is, ll, ik, e_simplex, w_simplex, dos_simplex) reduction(+:sumq)
532 do ist = 1, nst
533 do ii = 1, this%tetra_simplex%n_simplices
534 do is = 0, ns - 1
535 do ll = 1, this%tetra_simplex%sdim
536 ik = ns * (this%tetra_simplex%simplices(ii, ll) - 1) + 1
537 e_simplex(ll) = eigenvalues(ist, ik + is)
538 end do
539 w_simplex(:) = m_zero
540
541 call simplex_weights(this%tetra_simplex%rdim, e_simplex(1:this%tetra_simplex%sdim), e_fermi, &
542 w_simplex(1:nvertices), dos_simplex(1:nvertices))
543
544 do ll = 1, nvertices
545 sumq = sumq + this%el_per_state * w_simplex(ll) / this%tetra_simplex%n_points
546 end do
547 end do
548 end do
549 end do
550 !$omp end parallel do
551
552 conv = (abs(sumq - q_in) <= tol)
553 if (conv) exit
554
555 if (sumq <= q_in) emin = e_fermi
556 if (sumq >= q_in) emax = e_fermi
557 end do
558
559 if (.not. conv) then
560 message(1) = 'Fermi (tetrahedra): did not converge.'
561 call messages_fatal(1, namespace=namespace)
562 end if
563
564 pop_sub_with_profile(bisection_find_fermi_energy_tetra)
566
567 ! ---------------------------------------------------------
568 subroutine smear_fill_occupations(this, eigenvalues, occupations, kweights, nik, nst)
569 type(smear_t), intent(in) :: this
570 real(real64), intent(in) :: eigenvalues(:,:)
571 real(real64), intent(inout) :: occupations(:,:)
572 real(real64), intent(in) :: kweights(:)
573 integer, intent(in) :: nik, nst
574
575 integer :: ik, ist, ifermi
576 real(real64) :: dsmear, xx, dsmear_cond
577
578 push_sub(smear_fill_occupations)
579
580 if (this%method == smear_fixed_occ) then
581 ! do nothing
582 else if (this%method == smear_semiconductor) then
583 assert(this%fermi_count > 0 .and. this%fermi_count <= nik*nst)
584
585 ifermi = 0
586 do ik = 1, nik
587 do ist = 1, nst
588 xx = eigenvalues(ist, ik) - this%e_fermi
589 if (xx < -tol_smear) then
590 occupations(ist, ik) = this%el_per_state
591 else if (abs(xx) <= tol_smear .and. ifermi < this%fermi_count) then
592 occupations(ist, ik) = this%ef_occ * this%el_per_state
593 ifermi = ifermi + 1
594 else
595 occupations(ist, ik) = m_zero
596 end if
597
598 end do
599 end do
600 else if (this%method == smear_tetrahedra .or. this%method == smear_tetrahedra_opt) then
601 if (this%photodop) then
602 call messages_not_implemented('Photodoping for SmearingFunction = tetrahedra/tetrahedra_opt')
603 end if
604 block
605 integer :: ii, ll, is, ns, nik_red, nvertices
606 real(real64) :: e_simplex(20), w_simplex(4), dos_simplex(4)
607
608 assert(associated(this%tetra_simplex))
609 assert(this%tetra_simplex%n_simplices > 0)
610 assert(this%tetra_simplex%sdim > 0)
611
612 ns = this%nspins
613 assert(ns > 0)
614 assert(mod(nik, ns) == 0)
615
616 nik_red = nik / ns
617 assert(maxval(this%tetra_simplex%simplices) <= nik_red)
618
619 nvertices = this%tetra_simplex%rdim + 1
620 assert(nvertices > 0)
621 assert(nvertices <= this%tetra_simplex%sdim)
622
623 occupations(:, :) = m_zero
624 !$omp parallel do collapse(2) private(ist, ii, is, ll, ik, e_simplex, w_simplex, dos_simplex) reduction(+:occupations)
625 do ist = 1, nst
626 do ii = 1, this%tetra_simplex%n_simplices
627 do is = 0, ns - 1
628 do ll = 1, this%tetra_simplex%sdim
629 ik = ns * (this%tetra_simplex%simplices(ii, ll) - 1) + 1
630 e_simplex(ll) = eigenvalues(ist, ik + is)
631 end do
632 w_simplex(:) = m_zero
633
634 call simplex_weights(this%tetra_simplex%rdim, e_simplex(1:this%tetra_simplex%sdim), this%e_fermi, &
635 w_simplex(1:nvertices), dos_simplex(1:nvertices))
636
637 do ll = 1, nvertices
638 ik = ns * (this%tetra_simplex%simplices(ii, ll) - 1) + 1
639 assert(kweights(ik + is) > m_epsilon)
640 occupations(ist, ik + is) = occupations(ist, ik + is) + this%el_per_state * &
641 w_simplex(ll) / (this%tetra_simplex%n_points * kweights(ik + is))
642 end do
643 end do
644 end do
645 end do
646 !$omp end parallel do
647 end block
648 else
649 dsmear = max(1e-14_real64, this%dsmear)
650 if (this%photodop) then
651 dsmear_cond = max(1e-14_real64, this%dsmear_cond)
652 do ik = 1, nik
653 do ist = 1, nst
654 if (ist < this%photodop_bandmin) then
655 ! valence electrons
656 xx = (this%e_fermi - eigenvalues(ist, ik))/dsmear
657 else
658 ! conduction electrons
659 xx = (this%e_fermi_cond - eigenvalues(ist, ik))/dsmear_cond
660 end if
661 occupations(ist, ik) = smear_step_function(this, xx) * this%el_per_state
662 end do
663 end do
664 else
665 do ik = 1, nik
666 do ist = 1, nst
667 xx = (this%e_fermi - eigenvalues(ist, ik))/dsmear
668 occupations(ist, ik) = smear_step_function(this, xx) * this%el_per_state
669 end do
670 end do
671 end if
672 end if
673
675 end subroutine smear_fill_occupations
676
677 !--------------------------------------------------
678 real(real64) function smear_calc_entropy(this, eigenvalues, &
679 nik, nst, kweights, occ) result(entropy)
680 type(smear_t), intent(inout) :: this
681 real(real64), intent(in) :: eigenvalues(:,:)
682 real(real64), intent(in) :: kweights(:)
683 integer, intent(in) :: nik, nst
684 real(real64), intent(in) :: occ(:, :)
685
686 integer :: ist, ik
687 real(real64) :: dsmear, xx, term, ff
688
689 push_sub(smear_calc_entropy)
690
691 entropy = m_zero
692
693 if (this%method == smear_fixed_occ .or. this%method == smear_semiconductor) then
694 ! Fermi-Dirac entropy, not quite the same as will be obtained with true smearing
695 ! RM Wentzcovitch, JL Martins, and PB Allen, Phys. Rev. B 45, 11372 (1992) eqn (5)
696 ! also N Marzari PhD thesis p 117, http://quasiamore.mit.edu/phd/Marzari_PhD.pdf
697 do ik = 1, nik
698 do ist = 1, nst
699 ff = occ(ist, ik) / this%el_per_state
700 if (ff > m_zero .and. ff < m_one) then
701 term = ff * log(ff) + (1 - ff) * log(1 - ff)
702 else ! we have semiconducting smearing, or perverse occupations as in Methfessel-Paxton
703 term = m_zero
704 end if
705 entropy = entropy - kweights(ik) * this%el_per_state * term
706 end do
707 end do
708 else
709 dsmear = max(1e-14_real64, this%dsmear)
710
711 do ik = 1, nik
712 do ist = 1, nst
713 if (eigenvalues(ist, ik) < huge(m_one)) then
714 xx = (this%e_fermi - eigenvalues(ist, ik)) / dsmear
715 entropy = entropy - kweights(ik) * this%el_per_state * &
716 smear_entropy_function(this, xx)
717 end if
718 end do
719 end do
720 end if
721
722 pop_sub(smear_calc_entropy)
723 end function smear_calc_entropy
724
725
726 ! ---------------------------------------------------------
727 real(real64) function smear_delta_function(this, xx) result(deltaf)
728 type(smear_t), intent(in) :: this
729 real(real64), intent(in) :: xx
730
731 real(real64), parameter :: maxarg = 200.0_real64
732 real(real64) :: xp, arg, hd, hp, aa
733 integer :: ii, ni
734
735 ! no PUSH_SUB, called too often
736
737 ! smear_delta_function is not defined for SMEAR_FIXED_OCC
738 assert(this%method /= smear_fixed_occ)
739
740 deltaf = m_zero
741 select case (this%method)
743 if (abs(xx) <= m_epsilon) then
744 deltaf = this%ef_occ
745 end if
746
747 case (smear_fermi_dirac)
748 if (abs(xx) <= 36.0_real64) then
749 deltaf = m_one / (m_two + exp(-xx) + exp(xx))
750 end if
751
752 case (smear_cold)
753 xp = xx - m_one / sqrt(m_two)
754 arg = min(maxarg, xp**2)
755
756 deltaf = exp(-arg) / sqrt(m_pi) * (m_two - sqrt(m_two) * xx)
757
759 arg = min(maxarg, xx**2)
760 deltaf = exp(-arg) / sqrt(m_pi)
761
762 if (this%MP_n > 0) then ! recursion
763 hd = m_zero
764 hp = exp(-arg)
765 ni = 0
766 aa = m_one / sqrt(m_pi)
767 do ii = 1, this%MP_n
768 hd = m_two * xx * hp - m_two * ni * hd
769 ni = ni + 1
770 aa = -aa / (m_four * ii)
771 hp = m_two * xx * hd - m_two * ni * hp
772 ni = ni + 1
773 deltaf = deltaf + aa * hp
774 end do
775 end if
776
777 case (smear_spline)
778 xp = abs(xx) + m_one / sqrt(m_two)
779 deltaf = sqrt(m_e) * xp * exp(-xp * xp)
780
781 case (smear_gaussian)
782 deltaf = exp(-xx * xx) / sqrt(m_pi)
783
784 case (smear_lorentzian)
785 deltaf = m_one / (m_pi * (xx * xx + m_one))
786
787 end select
788
789 end function smear_delta_function
790
791
792 ! ---------------------------------------------------------
793 real(real64) function smear_step_function(this, xx) result(stepf)
794 type(smear_t), intent(in) :: this
795 real(real64), intent(in) :: xx
796
797 real(real64), parameter :: maxarg = 200.0_real64
798 real(real64) :: xp, arg, hd, hp, aa
799 integer :: ii, ni
800
801 push_sub(smear_step_function)
802
803 ! smear_step_function is not defined for SMEAR_FIXED_OCC
804 assert(this%method /= smear_fixed_occ)
805
806 stepf = m_zero
807 select case (this%method)
809 if (xx > m_zero) then
810 stepf = m_one
811 else if (abs(xx) <= m_epsilon) then
812 stepf = this%ef_occ
813 end if
814
815 case (smear_fermi_dirac)
816 if (xx > maxarg) then
817 stepf = m_one
818 else if (xx > -maxarg) then
819 stepf = m_one / (m_one + exp(-xx))
820 end if
821
823 xp = xx - m_one / sqrt(m_two)
824 arg = min(maxarg, xp**2)
825
826 stepf = m_half * erf(xp) + &
827 m_one / sqrt(m_two * m_pi) * exp(-arg) + m_half
828
830 stepf = m_half * erfc(-xx)
831
832 if (this%MP_n > 0) then ! recursion
833 hd = m_zero
834 arg = min(maxarg, xx**2)
835 hp = exp(-arg)
836 ni = 0
837 aa = m_one / sqrt(m_pi)
838 do ii = 1, this%MP_n
839 hd = m_two * xx * hp - m_two * ni * hd
840 ni = ni + 1
841 aa = -aa / (m_four * ii)
842 stepf = stepf - aa * hd
843 hp = m_two * xx * hd - m_two * ni * hp
844 ni = ni + 1
845 end do
846 end if
847
848 case (smear_spline)
849 if (xx <= m_zero) then
850 xp = xx - m_one / sqrt(m_two)
851 stepf = m_half * sqrt(m_e) * exp(-xp * xp)
852 else
853 xp = xx + m_one / sqrt(m_two)
854 stepf = m_one - m_half * sqrt(m_e) * exp(-xp * xp)
855 end if
856
857 case (smear_gaussian)
858 stepf = m_half * erfc(-xx)
859
860 case (smear_lorentzian)
861 stepf = atan(xx) / m_pi + m_half
862
863 end select
864
865 pop_sub(smear_step_function)
866 end function smear_step_function
867
868
869 ! ---------------------------------------------------------
871 real(real64) function smear_entropy_function(this, xx) result(entropyf)
872 type(smear_t), intent(in) :: this
873 real(real64), intent(in) :: xx
874
875 real(real64), parameter :: maxarg = 200.0_real64
876 real(real64) :: xp, arg, hd, hp, hpm1, aa
877 integer :: ii, ni
878
879 push_sub(smear_entropy_function)
880
881 ! smear_entropy_function is not defined for SMEAR_FIXED_OCC
882 assert(this%method /= smear_fixed_occ)
883
884 entropyf = m_zero
885 select case (this%method)
887
889 if (abs(xx) <= 36.0_real64) then
890 xp = m_one / (m_one + exp(-xx))
891 entropyf = xp * log(xp) + (m_one - xp) * log(m_one - xp)
892 end if
893
894 case (smear_cold)
895 xp = xx - m_one / sqrt(m_two)
896 arg = min(maxarg, xp**2)
897
898 entropyf = m_one / sqrt(m_two * m_pi) * xp * exp(-arg)
899
901 arg = min(maxarg, xx**2)
902 entropyf = -m_half * exp(-arg) / sqrt(m_pi)
903
904 if (this%MP_n > 0) then ! recursion
905 hd = m_zero
906 hp = exp(-arg)
907 ni = 0
908 aa = m_one / sqrt(m_pi)
909 do ii = 1, this%MP_n
910 hd = m_two * xx * hp - m_two * ni * hd
911 ni = ni + 1
912 hpm1 = hp
913 hp = m_two * xx * hd - m_two * ni * hp
914 ni = ni + 1
915 aa = -aa / (m_four * ii)
916 entropyf = entropyf - aa * (m_half * hp + hpm1 * ni)
917 end do
918 end if
919
920 case (smear_spline)
921 xp = abs(xx) + m_one / sqrt(m_two)
922 entropyf = -sqrt(m_e) * (abs(xx) * exp(-xp * xp) / m_two + sqrt(m_pi) / m_four * erfc(xp))
923
924 case (smear_gaussian)
925 entropyf = -exp(-xx * xx) / m_two / sqrt(m_pi)
926
927 case (smear_lorentzian)
928 entropyf = m_ninf
929
930 end select
931
933 end function smear_entropy_function
934
935
936 ! ---------------------------------------------------------
937 logical pure function smear_is_semiconducting(this) result(answer)
938 type(smear_t), intent(in) :: this
939
940 answer = this%method == smear_semiconductor
941
942 end function smear_is_semiconducting
943
944 subroutine smear_write_info(this, namespace, iunit)
945 type(smear_t), intent(in) :: this
946 type(namespace_t), intent(in) :: namespace
947 integer, optional, intent(in) :: iunit
948
949
950 if (this%method /= smear_semiconductor .and. this%method /= smear_fixed_occ) then
951 if (this%photodop) then
952 write(message(1), '(a,f12.6,1x,a)') "Fermi energy (valence ) = ", &
953 units_from_atomic(units_out%energy, this%e_fermi), units_abbrev(units_out%energy)
954 write(message(2), '(a,f12.6,1x,a)') "Fermi energy (conduction) = ", &
955 units_from_atomic(units_out%energy, this%e_fermi_cond), units_abbrev(units_out%energy)
956 call messages_info(2, iunit=iunit, namespace=namespace)
957 else
958 write(message(1), '(a,f12.6,1x,a)') "Fermi energy = ", &
959 units_from_atomic(units_out%energy, this%e_fermi), units_abbrev(units_out%energy)
960 call messages_info(1, iunit=iunit, namespace=namespace)
961 end if
962 end if
963
964 end subroutine smear_write_info
965
967end module smear_oct_m
968
969!! Local Variables:
970!! mode: f90
971!! coding: utf-8
972!! End:
This is the common interface to a sorting routine. It performs the shell algorithm,...
Definition: sort.F90:156
double log(double __x) __attribute__((__nothrow__
double exp(double __x) __attribute__((__nothrow__
double atan(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
double floor(double __x) __attribute__((__nothrow__
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public p_ry
Definition: global.F90:230
real(real64), parameter, public m_epsilon
Definition: global.F90:207
real(real64), parameter, public m_half
Definition: global.F90:197
real(real64), parameter, public m_one
Definition: global.F90:192
real(real64), parameter, public m_min_occ
Minimal occupation that is considered to be non-zero.
Definition: global.F90:218
integer function, public kpoints_kweight_denominator(this)
Definition: kpoints.F90:1802
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1067
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:999
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:409
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:690
subroutine, public simplex_weights(rdim, esimplex, eF, weights, dos)
Get the weights and DOS contribution of a single simplex.
Definition: simplex.F90:426
integer, parameter, public smear_tetrahedra_opt
Definition: smear.F90:176
real(real64) function, public smear_calc_entropy(this, eigenvalues, nik, nst, kweights, occ)
Definition: smear.F90:775
integer, parameter, public smear_tetrahedra
Definition: smear.F90:176
subroutine, public smear_find_fermi_energy(this, namespace, eigenvalues, occupations, qtot, nik, nst, kweights)
Definition: smear.F90:371
real(real64) function, public smear_entropy_function(this, xx)
This function is defined as .
Definition: smear.F90:967
real(real64) function, public smear_delta_function(this, xx)
Definition: smear.F90:823
integer, parameter, public smear_fermi_dirac
Definition: smear.F90:176
subroutine bisection_find_fermi_energy_tetra(this, namespace, tol, eigenvalues, nik, nst, q_in, e_fermi)
Finds the Fermi energy using the tetrahedron method via bisection.
Definition: smear.F90:583
integer, parameter, public smear_methfessel_paxton
Definition: smear.F90:176
subroutine, public smear_write_info(this, namespace, iunit)
Definition: smear.F90:1040
real(real64) function, public smear_step_function(this, xx)
Definition: smear.F90:889
integer, parameter, public smear_semiconductor
Definition: smear.F90:176
subroutine, public smear_fill_occupations(this, eigenvalues, occupations, kweights, nik, nst)
Definition: smear.F90:664
integer, parameter, public smear_lorentzian
Definition: smear.F90:176
subroutine, public smear_copy(to, from)
Definition: smear.F90:349
integer, parameter, public smear_fixed_occ
Definition: smear.F90:176
integer, parameter, public smear_spline
Definition: smear.F90:176
subroutine bisection_find_fermi_energy(this, namespace, dsmear_in, tol, eigenvalues, kweights, nik, q_in, start_band, end_band, e_fermi, ef_occ)
Definition: smear.F90:510
subroutine, public smear_init(this, namespace, ispin, fixed_occ, integral_occs, kpoints)
Definition: smear.F90:194
integer, parameter, public smear_cold
Definition: smear.F90:176
logical pure function, public smear_is_semiconducting(this)
Definition: smear.F90:1033
integer, parameter, public smear_gaussian
Definition: smear.F90:176
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:119
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
This module defines the unit system, used for input and output.
type(unit_system_t), public units_inp
the units systems for reading and writing
int true(void)