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
33 use sort_oct_m
34 use unit_oct_m
37
38 implicit none
39
40 private
41 public :: &
42 smear_t, &
43 smear_init, &
44 smear_copy, &
53
54 type smear_t
55 private
56 integer, public :: method
57 real(real64), public :: dsmear
58 real(real64) :: dsmear_cond
59 real(real64), public :: e_fermi
60 real(real64) :: e_fermi_cond
61
62 integer, public :: el_per_state
63 real(real64), public :: ef_occ
64 real(real64) :: ef_occ_cond
65 logical, public :: integral_occs
66 integer, public :: MP_n
67 integer :: fermi_count
68 integer :: nik_factor
69 integer :: nspins
70
71 ! PhotonDoping
72 logical, public:: photodop
73 real(real64) :: nephotodop
74 integer :: photodop_bandmin
75
76 end type smear_t
77
78 integer, parameter, public :: &
79 SMEAR_SEMICONDUCTOR = 1, &
81 smear_cold = 3, &
83 smear_spline = 5, &
84 smear_fixed_occ = 6, &
85 smear_gaussian = 7, &
87
88 real(real64), parameter :: TOL_SMEAR = 1e-6_real64
89
90contains
91
92 !--------------------------------------------------
93 subroutine smear_init(this, namespace, ispin, fixed_occ, integral_occs, kpoints)
94 type(smear_t), intent(out) :: this
95 type(namespace_t), intent(in) :: namespace
96 integer, intent(in) :: ispin
97 logical, intent(in) :: fixed_occ
98 logical, intent(in) :: integral_occs
99 type(kpoints_t), intent(in) :: kpoints
100
101 push_sub(smear_init)
102
103 this%integral_occs = integral_occs
104
105 !%Variable SmearingFunction
106 !%Type integer
107 !%Default semiconducting
108 !%Section States
109 !%Description
110 !% This is the function used to smear the electronic occupations.
111 !% It is ignored if the <tt>Occupations</tt> block is set.
112 !%Option semiconducting 1
113 !% Semiconducting occupations, <i>i.e.</i> the lowest lying states are occupied
114 !% until no more electrons are left.
115 !%Option fermi_dirac 2
116 !% Simple Fermi-Dirac distribution. In this case, <tt>Smearing</tt> has
117 !% the meaning of an electronic temperature. DN Mermin, <i>Phys. Rev.</i> <b>137</b>, A1441 (1965).
118 !%Option cold_smearing 3
119 !% N Marzari, D Vanderbilt, A De Vita, and MC Payne, <i>Phys. Rev. Lett.</i> <b>82</b>, 3296 (1999).
120 !%Option methfessel_paxton 4
121 !% M Methfessel and AT Paxton, <i>Phys. Rev. B</i> <b>40</b>, 3616 (1989).
122 !% In this case, the variable <tt>SmearingMPOrder</tt> sets the order of the smearing.
123 !% Occupations may be negative.
124 !%Option spline_smearing 5
125 !% Nearly identical to Gaussian smearing.
126 !% JM Holender, MJ Gillan, MC Payne, and AD Simpson, <i>Phys. Rev. B</i> <b>52</b>, 967 (1995).
127 !%Option gaussian_smearing 7
128 !% Gaussian smearing.
129 !%Option lorentzian_smearing 8
130 !% Lorentzian smearing.
131 !%End
132 if (fixed_occ) then
133 this%method = smear_fixed_occ
134 else
135 call parse_variable(namespace, 'SmearingFunction', smear_semiconductor, this%method)
136 if (.not. varinfo_valid_option('SmearingFunction', this%method)) call messages_input_error(namespace, 'SmearingFunction')
137 call messages_print_var_option('SmearingFunction', this%method, namespace=namespace)
138 end if
139
140 !%Variable Smearing
141 !%Type float
142 !%Default 0.1 eV
143 !%Section States
144 !%Description
145 !% If <tt>Occupations</tt> is not set, <tt>Smearing</tt> is the
146 !% smearing width used in the <tt>SmearingFunction</tt> to distribute the electrons
147 !% among the existing states.
148 !%End
149 this%dsmear = 1e-14_real64
150 if (this%method /= smear_semiconductor .and. this%method /= smear_fixed_occ) then
151 call parse_variable(namespace, 'Smearing', 0.1_real64 / (m_two * p_ry), this%dsmear, units_inp%energy)
152 end if
154 !%Variable PhotoDopingSmearing
155 !%Type float
156 !%Default 0.1 eV
157 !%Section States
158 !%Description
159 !% If <tt>Occupations</tt> is not set, <tt>PhotoDopingSmearing</tt> is the
160 !% smearing width used in the <tt>SmearingFunction</tt> to distribute the electrons in conduction bands
161 !% among the existing states.
162 !%End
163 this%dsmear_cond = this%dsmear
164 if (this%method /= smear_semiconductor .and. this%method /= smear_fixed_occ) then
165 call parse_variable(namespace, 'PhotoDopingSmearing', 0.1_real64 / (m_two * p_ry), this%dsmear_cond, units_inp%energy)
166 end if
168 !%Variable PhotoDopingNumElectrons
169 !%Type float
170 !%Default 0
171 !%Section States
172 !%Description
173 !% This variable allows to set the number of valence electrons taken out to put
174 !% into the conduction bands. The method follows Marini and Calandra, PRB 104, 144103 (2021).
175 !%End
176 call parse_variable(namespace, 'PhotoDopingNumElectrons', m_zero, this%nephotodop)
177
178 if (this%nephotodop > m_min_occ) then
179 this%photodop=.true.
180 else
181 this%photodop=.false.
182 endif
184
185 !%Variable PhotoDopingBand
186 !%Type integer
187 !%Default 1
188 !%Section States
189 !%Description
190 !% This variable specifies where the minium (starting) conduction band index for the photodoping electrons.
191 !%End
192 call parse_variable(namespace, 'PhotoDopingBand', 1, this%photodop_bandmin)
193
194
195 call messages_obsolete_variable(namespace, 'ElectronicTemperature', 'Smearing')
196
197 if (ispin == 1) then ! unpolarized
198 this%el_per_state = 2
199 else
200 this%el_per_state = 1
201 end if
202
203 if (ispin == 2) then
204 this%nspins = 2
205 else
206 this%nspins = 1
207 end if
208
209 if (this%method == smear_semiconductor) then
210 this%nik_factor = kpoints_kweight_denominator(kpoints)
211
212 if (this%nik_factor == 0) then
213 message(1) = "k-point weights in KPoints or KPointsReduced blocks must be rational numbers for semiconducting smearing."
214 call messages_fatal(1, namespace=namespace)
215 end if
216 end if
217
218 this%MP_n = 0
219 if (this%method == smear_methfessel_paxton) then
220 !%Variable SmearingMPOrder
221 !%Type integer
222 !%Default 1
223 !%Section States
224 !%Description
225 !% Sets the order of the Methfessel-Paxton smearing function.
226 !%End
227 call parse_variable(namespace, 'SmearingMPOrder', 1, this%MP_n)
228 end if
229
230 pop_sub(smear_init)
231 end subroutine smear_init
232
233
234 !--------------------------------------------------
235 subroutine smear_copy(to, from)
236 type(smear_t), intent(out) :: to
237 type(smear_t), intent(in) :: from
238
239 push_sub(smear_copy)
240
241 to%method = from%method
242 to%dsmear = from%dsmear
243 to%e_fermi = from%e_fermi
244 to%el_per_state = from%el_per_state
245 to%ef_occ = from%ef_occ
246 to%MP_n = from%MP_n
247 to%fermi_count = from%fermi_count
248 to%nik_factor = from%nik_factor
249
250 pop_sub(smear_copy)
251 end subroutine smear_copy
252
253
254 !--------------------------------------------------
255 subroutine smear_find_fermi_energy(this, namespace, eigenvalues, occupations, &
256 qtot, nik, nst, kweights)
257 type(smear_t), intent(inout) :: this
258 type(namespace_t), intent(in) :: namespace
259 real(real64), intent(in) :: eigenvalues(:,:), occupations(:,:)
260 real(real64), intent(in) :: qtot, kweights(:)
261 integer, intent(in) :: nik, nst
262
263 real(real64), parameter :: tol = 1.0e-10_real64
264 integer :: ist, ik, iter, maxq, weight, sumq_int, sum_weight
265 real(real64) :: sumq_frac
266 logical :: conv
267 real(real64), allocatable :: eigenval_list(:)
268 integer, allocatable :: k_list(:), reorder(:)
269 integer :: fermi_count_up, fermi_count_down
270
272
273 maxq = this%el_per_state * nst * this%nspins
274 if (maxq - qtot <= -tol) then ! not enough states
275 message(1) = 'Not enough states'
276 write(message(2),'(6x,a,f12.6,a,i10)')'(total charge = ', qtot, &
277 ' max charge = ', maxq
278 call messages_fatal(2, namespace=namespace)
279 end if
280
281 conv = .true.
282 if (this%method == smear_fixed_occ) then ! Fermi energy: last of occupied states
283 ist_cycle: do ist = nst, 1, -1
284 if (any(occupations(ist, :) > m_min_occ)) then
285 this%e_fermi = eigenvalues(ist, 1)
286 this%ef_occ = occupations(ist, 1) / this%el_per_state
287 do ik = 2, nik
288 if (eigenvalues(ist, ik) > this%e_fermi .and. occupations(ist, ik) > m_min_occ) then
289 this%e_fermi = eigenvalues(ist, ik)
290 this%ef_occ = occupations(ist, ik) / this%el_per_state
291 end if
292 end do
293 exit ist_cycle
294 end if
295 end do ist_cycle
296
297 else if (this%method == smear_semiconductor) then
298 ! first we sort the eigenvalues
299 safe_allocate(eigenval_list(1:nst * nik))
300 safe_allocate( k_list(1:nst * nik))
301 safe_allocate( reorder(1:nst * nik))
302
303 iter = 1
304 do ist = 1, nst
305 do ik = 1, nik
306 eigenval_list(iter) = eigenvalues(ist, ik)
307 k_list(iter) = ik
308 reorder(iter) = iter
309 iter = iter + 1
310 end do
311 end do
312
313 call sort(eigenval_list, reorder)
314
315 sumq_int = floor(qtot) * this%nik_factor
316 sumq_frac = qtot * this%nik_factor - sumq_int
317 if ( sumq_frac + tol > this%nik_factor ) then
318 sumq_int = sumq_int + this%nik_factor
319 sumq_frac = sumq_frac - this%nik_factor
320 end if
321 if (sumq_frac < tol) sumq_frac = m_zero
322
323 do iter = 1, nst * nik
324 weight = int(kweights(k_list(reorder(iter))) * this%nik_factor + m_half)
325 if (.not. weight > 0) cycle
326 this%e_fermi = eigenval_list(iter)
327 this%ef_occ = (sumq_int + sumq_frac) / (weight * this%el_per_state)
328
329 if (sumq_int - weight * this%el_per_state <= 0) then
330 ! count how many occupied states are at the fermi level,
331 ! this is required later to fill the states
332 this%fermi_count = 1
333 fermi_count_up = 1
334 fermi_count_down = 1
335 sum_weight = weight
336 do
337 if (iter - fermi_count_down < 1) exit
338 if (abs(this%e_fermi - eigenval_list(iter - fermi_count_down)) > tol_smear) exit
339 weight = int(kweights(k_list(reorder(iter - fermi_count_down))) * this%nik_factor + m_half)
340 if (weight > m_epsilon) then
341 sumq_int = sumq_int + weight * this%el_per_state
342 sum_weight = sum_weight + weight
343 end if
344 fermi_count_down = fermi_count_down + 1
345 end do
346 do
347 if (iter + fermi_count_up > nst*nik) exit
348 if (abs(this%e_fermi - eigenval_list(iter + fermi_count_up)) > tol_smear) exit
349 weight = int(kweights(k_list(reorder(iter + fermi_count_up))) * this%nik_factor + m_half)
350 if (weight > m_epsilon) then
351 sum_weight = sum_weight + weight
352 end if
353 fermi_count_up = fermi_count_up + 1
354 end do
355 this%fermi_count = fermi_count_up + fermi_count_down - 1
356 this%ef_occ = (sumq_int + sumq_frac) / (sum_weight * this%el_per_state)
357 exit
358 end if
359
360 sumq_int = sumq_int - weight * this%el_per_state
361 end do
362 assert(this%ef_occ < m_one + 10.0_real64*m_epsilon)
363
364 safe_deallocate_a(eigenval_list)
365 safe_deallocate_a(k_list)
366 safe_deallocate_a(reorder)
367
368 else ! bisection
369 if (this%photodop) then
370 ! find fermi energy for valence electrons
371 call bisection_find_fermi_energy(this, namespace, this%dsmear, tol, &
372 eigenvalues, kweights, nik, qtot-this%nephotodop, 1, this%photodop_bandmin-1, &
373 this%e_fermi, this%ef_occ)
374
375 ! find fermi energy for conduction electrons
376 call bisection_find_fermi_energy(this, namespace, this%dsmear_cond, tol, &
377 eigenvalues, kweights, nik, this%nephotodop, this%photodop_bandmin, nst, &
378 this%e_fermi_cond, this%ef_occ_cond)
379 else
380 call bisection_find_fermi_energy(this, namespace, this%dsmear, tol, &
381 eigenvalues, kweights, nik, qtot, 1, nst, this%e_fermi, this%ef_occ)
382 end if
383 end if
384
386 end subroutine smear_find_fermi_energy
387
388 ! ---------------------------------------------------------
389 subroutine bisection_find_fermi_energy(this, namespace, dsmear_in, tol, &
390 eigenvalues, kweights, nik, q_in, start_band, end_band, e_fermi, ef_occ)
391 type(smear_t), intent(inout) :: this
392 type(namespace_t), intent(in) :: namespace
393 real(real64), intent(in) :: dsmear_in
394 real(real64), intent(in) :: tol
395 integer, intent(in) :: nik
396 real(real64), intent(in) :: eigenvalues(:,:)
397 real(real64), intent(in) :: kweights(:), q_in
398 integer, intent(in) :: start_band, end_band
399 real(real64), intent(out) :: e_fermi, ef_occ
400
401 integer, parameter :: nitmax = 200
402 integer :: ist, ik, iter
403 real(real64) :: drange, dsmear, emin, emax, xx, sumq
404 logical :: conv
405
407
408 dsmear = max(1e-14_real64, dsmear_in)
409 drange = dsmear * sqrt(-log(tol * 0.01_real64))
410
411 emin = minval(eigenvalues) - drange
412 emax = maxval(eigenvalues) + drange
413
414 do iter = 1, nitmax
415 e_fermi = m_half * (emin + emax)
416 sumq = m_zero
417
418 do ik = 1, nik
419 do ist = start_band, end_band
420 xx = (e_fermi - eigenvalues(ist, ik))/dsmear
421 sumq = sumq + kweights(ik) * this%el_per_state * &
422 smear_step_function(this, xx)
423 end do
424 end do
425
426 conv = (abs(sumq - q_in) <= tol)
427 if (conv) exit
428
429 if (sumq <= q_in) emin = e_fermi
430 if (sumq >= q_in) emax = e_fermi
431
432 ef_occ = smear_step_function(this, m_zero)
433 end do
434
435 if (.not. conv) then
436 message(1) = 'Fermi: did not converge.'
437 call messages_fatal(1, namespace=namespace)
438 end if
439
441 end subroutine bisection_find_fermi_energy
442
443 ! ---------------------------------------------------------
444 subroutine smear_fill_occupations(this, eigenvalues, occupations, nik, nst)
445 type(smear_t), intent(in) :: this
446 real(real64), intent(in) :: eigenvalues(:,:)
447 real(real64), intent(inout) :: occupations(:,:)
448 integer, intent(in) :: nik, nst
449
450 integer :: ik, ist, ifermi
451 real(real64) :: dsmear, xx, dsmear_cond
452
453 push_sub(smear_fill_occupations)
454
455 if (this%method == smear_fixed_occ) then
456 ! do nothing
457 else if (this%method == smear_semiconductor) then
458 assert(this%fermi_count > 0 .and. this%fermi_count <= nik*nst)
459
460 ifermi = 0
461 do ik = 1, nik
462 do ist = 1, nst
463 xx = eigenvalues(ist, ik) - this%e_fermi
464 if (xx < -tol_smear) then
465 occupations(ist, ik) = this%el_per_state
466 else if (abs(xx) <= tol_smear .and. ifermi < this%fermi_count) then
467 occupations(ist, ik) = this%ef_occ * this%el_per_state
468 ifermi = ifermi + 1
469 else
470 occupations(ist, ik) = m_zero
471 end if
472
473 end do
474 end do
475
476 else
477 dsmear = max(1e-14_real64, this%dsmear)
478 if (this%photodop) then
479 dsmear_cond = max(1e-14_real64, this%dsmear_cond)
480 do ik = 1, nik
481 do ist = 1, nst
482 if (ist < this%photodop_bandmin) then
483 ! valence electrons
484 xx = (this%e_fermi - eigenvalues(ist, ik))/dsmear
485 else
486 ! conduction electrons
487 xx = (this%e_fermi_cond - eigenvalues(ist, ik))/dsmear_cond
488 end if
489 occupations(ist, ik) = smear_step_function(this, xx) * this%el_per_state
490 end do
491 end do
492 else
493 do ik = 1, nik
494 do ist = 1, nst
495 xx = (this%e_fermi - eigenvalues(ist, ik))/dsmear
496 occupations(ist, ik) = smear_step_function(this, xx) * this%el_per_state
497 end do
498 end do
499 end if
500 end if
501
503 end subroutine smear_fill_occupations
504
505
506 !--------------------------------------------------
507 real(real64) function smear_calc_entropy(this, eigenvalues, &
508 nik, nst, kweights, occ) result(entropy)
509 type(smear_t), intent(inout) :: this
510 real(real64), intent(in) :: eigenvalues(:,:)
511 real(real64), intent(in) :: kweights(:)
512 integer, intent(in) :: nik, nst
513 real(real64), intent(in) :: occ(:, :)
514
515 integer :: ist, ik
516 real(real64) :: dsmear, xx, term, ff
517
518 push_sub(smear_calc_entropy)
519
520 entropy = m_zero
521
522 if (this%method == smear_fixed_occ .or. this%method == smear_semiconductor) then
523 ! Fermi-Dirac entropy, not quite the same as will be obtained with true smearing
524 ! RM Wentzcovitch, JL Martins, and PB Allen, Phys. Rev. B 45, 11372 (1992) eqn (5)
525 ! also N Marzari PhD thesis p 117, http://quasiamore.mit.edu/phd/Marzari_PhD.pdf
526 do ik = 1, nik
527 do ist = 1, nst
528 ff = occ(ist, ik) / this%el_per_state
529 if (ff > m_zero .and. ff < m_one) then
530 term = ff * log(ff) + (1 - ff) * log(1 - ff)
531 else ! we have semiconducting smearing, or perverse occupations as in Methfessel-Paxton
532 term = m_zero
533 end if
534 entropy = entropy - kweights(ik) * this%el_per_state * term
535 end do
536 end do
537 else
538 dsmear = max(1e-14_real64, this%dsmear)
540 do ik = 1, nik
541 do ist = 1, nst
542 if (eigenvalues(ist, ik) < huge(m_one)) then
543 xx = (this%e_fermi - eigenvalues(ist, ik)) / dsmear
544 entropy = entropy - kweights(ik) * this%el_per_state * &
545 smear_entropy_function(this, xx)
546 end if
547 end do
548 end do
549 end if
550
551 pop_sub(smear_calc_entropy)
552 end function smear_calc_entropy
553
554
555 ! ---------------------------------------------------------
556 real(real64) function smear_delta_function(this, xx) result(deltaf)
557 type(smear_t), intent(in) :: this
558 real(real64), intent(in) :: xx
559
560 real(real64), parameter :: maxarg = 200.0_real64
561 real(real64) :: xp, arg, hd, hp, aa
562 integer :: ii, ni
563
564 ! no PUSH_SUB, called too often
565
566 ! smear_delta_function is not defined for SMEAR_FIXED_OCC
567 assert(this%method /= smear_fixed_occ)
568
569 deltaf = m_zero
570 select case (this%method)
572 if (abs(xx) <= m_epsilon) then
573 deltaf = this%ef_occ
574 end if
575
576 case (smear_fermi_dirac)
577 if (abs(xx) <= 36.0_real64) then
578 deltaf = m_one / (m_two + exp(-xx) + exp(xx))
579 end if
580
581 case (smear_cold)
582 xp = xx - m_one / sqrt(m_two)
583 arg = min(maxarg, xp**2)
584
585 deltaf = exp(-arg) / sqrt(m_pi) * (m_two - sqrt(m_two) * xx)
586
588 arg = min(maxarg, xx**2)
589 deltaf = exp(-arg) / sqrt(m_pi)
590
591 if (this%MP_n > 0) then ! recursion
592 hd = m_zero
593 hp = exp(-arg)
594 ni = 0
595 aa = m_one / sqrt(m_pi)
596 do ii = 1, this%MP_n
597 hd = m_two * xx * hp - m_two * ni * hd
598 ni = ni + 1
599 aa = -aa / (m_four * ii)
600 hp = m_two * xx * hd - m_two * ni * hp
601 ni = ni + 1
602 deltaf = deltaf + aa * hp
603 end do
604 end if
605
606 case (smear_spline)
607 xp = abs(xx) + m_one / sqrt(m_two)
608 deltaf = sqrt(m_e) * xp * exp(-xp * xp)
609
610 case (smear_gaussian)
611 deltaf = exp(-xx * xx) / sqrt(m_pi)
612
613 case (smear_lorentzian)
614 deltaf = m_one / (m_pi * (xx * xx + m_one))
615
616 end select
617
618 end function smear_delta_function
619
620
621 ! ---------------------------------------------------------
622 real(real64) function smear_step_function(this, xx) result(stepf)
623 type(smear_t), intent(in) :: this
624 real(real64), intent(in) :: xx
625
626 real(real64), parameter :: maxarg = 200.0_real64
627 real(real64) :: xp, arg, hd, hp, aa
628 integer :: ii, ni
629
630 push_sub(smear_step_function)
631
632 ! smear_step_function is not defined for SMEAR_FIXED_OCC
633 assert(this%method /= smear_fixed_occ)
634
635 stepf = m_zero
636 select case (this%method)
638 if (xx > m_zero) then
639 stepf = m_one
640 else if (abs(xx) <= m_epsilon) then
641 stepf = this%ef_occ
642 end if
643
644 case (smear_fermi_dirac)
645 if (xx > maxarg) then
646 stepf = m_one
647 else if (xx > -maxarg) then
648 stepf = m_one / (m_one + exp(-xx))
649 end if
650
652 xp = xx - m_one / sqrt(m_two)
653 arg = min(maxarg, xp**2)
654
655 stepf = m_half * erf(xp) + &
656 m_one / sqrt(m_two * m_pi) * exp(-arg) + m_half
657
659 stepf = m_half * erfc(-xx)
660
661 if (this%MP_n > 0) then ! recursion
662 hd = m_zero
663 arg = min(maxarg, xx**2)
664 hp = exp(-arg)
665 ni = 0
666 aa = m_one / sqrt(m_pi)
667 do ii = 1, this%MP_n
668 hd = m_two * xx * hp - m_two * ni * hd
669 ni = ni + 1
670 aa = -aa / (m_four * ii)
671 stepf = stepf - aa * hd
672 hp = m_two * xx * hd - m_two * ni * hp
673 ni = ni + 1
674 end do
675 end if
676
677 case (smear_spline)
678 if (xx <= m_zero) then
679 xp = xx - m_one / sqrt(m_two)
680 stepf = m_half * sqrt(m_e) * exp(-xp * xp)
681 else
682 xp = xx + m_one / sqrt(m_two)
683 stepf = m_one - m_half * sqrt(m_e) * exp(-xp * xp)
684 end if
685
686 case (smear_gaussian)
687 stepf = m_half * erfc(-xx)
688
689 case (smear_lorentzian)
690 stepf = atan(xx) / m_pi + m_half
691
692 end select
693
694 pop_sub(smear_step_function)
695 end function smear_step_function
696
697
698 ! ---------------------------------------------------------
700 real(real64) function smear_entropy_function(this, xx) result(entropyf)
701 type(smear_t), intent(in) :: this
702 real(real64), intent(in) :: xx
703
704 real(real64), parameter :: maxarg = 200.0_real64
705 real(real64) :: xp, arg, hd, hp, hpm1, aa
706 integer :: ii, ni
707
708 push_sub(smear_entropy_function)
709
710 ! smear_entropy_function is not defined for SMEAR_FIXED_OCC
711 assert(this%method /= smear_fixed_occ)
712
713 entropyf = m_zero
714 select case (this%method)
716
718 if (abs(xx) <= 36.0_real64) then
719 xp = m_one / (m_one + exp(-xx))
720 entropyf = xp * log(xp) + (m_one - xp) * log(m_one - xp)
721 end if
722
723 case (smear_cold)
724 xp = xx - m_one / sqrt(m_two)
725 arg = min(maxarg, xp**2)
726
727 entropyf = m_one / sqrt(m_two * m_pi) * xp * exp(-arg)
728
730 arg = min(maxarg, xx**2)
731 entropyf = -m_half * exp(-arg) / sqrt(m_pi)
732
733 if (this%MP_n > 0) then ! recursion
734 hd = m_zero
735 hp = exp(-arg)
736 ni = 0
737 aa = m_one / sqrt(m_pi)
738 do ii = 1, this%MP_n
739 hd = m_two * xx * hp - m_two * ni * hd
740 ni = ni + 1
741 hpm1 = hp
742 hp = m_two * xx * hd - m_two * ni * hp
743 ni = ni + 1
744 aa = -aa / (m_four * ii)
745 entropyf = entropyf - aa * (m_half * hp + hpm1 * ni)
746 end do
747 end if
748
749 case (smear_spline)
750 xp = abs(xx) + m_one / sqrt(m_two)
751 entropyf = -sqrt(m_e) * (abs(xx) * exp(-xp * xp) / m_two + sqrt(m_pi) / m_four * erfc(xp))
752
753 case (smear_gaussian)
754 entropyf = -exp(-xx * xx) / m_two / sqrt(m_pi)
755
756 case (smear_lorentzian)
757 entropyf = m_ninf
758
759 end select
760
762 end function smear_entropy_function
763
764
765 ! ---------------------------------------------------------
766 logical pure function smear_is_semiconducting(this) result(answer)
767 type(smear_t), intent(in) :: this
768
769 answer = this%method == smear_semiconductor
770
771 end function smear_is_semiconducting
772
773 subroutine smear_write_info(this, namespace, iunit)
774 type(smear_t), intent(in) :: this
775 type(namespace_t), intent(in) :: namespace
776 integer, optional, intent(in) :: iunit
777
778
779 if (this%method /= smear_semiconductor .and. this%method /= smear_fixed_occ) then
780 if (this%photodop) then
781 write(message(1), '(a,f12.6,1x,a)') "Fermi energy (valence ) = ", &
782 units_from_atomic(units_out%energy, this%e_fermi), units_abbrev(units_out%energy)
783 write(message(2), '(a,f12.6,1x,a)') "Fermi energy (conduction) = ", &
784 units_from_atomic(units_out%energy, this%e_fermi_cond), units_abbrev(units_out%energy)
785 call messages_info(2, iunit=iunit, namespace=namespace)
786 else
787 write(message(1), '(a,f12.6,1x,a)') "Fermi energy = ", &
788 units_from_atomic(units_out%energy, this%e_fermi), units_abbrev(units_out%energy)
789 call messages_info(1, iunit=iunit, namespace=namespace)
790 end if
791 end if
792
793 end subroutine smear_write_info
794
796end module smear_oct_m
797
798!! Local Variables:
799!! mode: f90
800!! coding: utf-8
801!! 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:1650
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1023
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:691
real(real64) function, public smear_calc_entropy(this, eigenvalues, nik, nst, kweights, occ)
Definition: smear.F90:604
subroutine, public smear_fill_occupations(this, eigenvalues, occupations, nik, nst)
Definition: smear.F90:540
subroutine, public smear_find_fermi_energy(this, namespace, eigenvalues, occupations, qtot, nik, nst, kweights)
Definition: smear.F90:352
real(real64) function, public smear_entropy_function(this, xx)
This function is defined as .
Definition: smear.F90:796
real(real64) function, public smear_delta_function(this, xx)
Definition: smear.F90:652
integer, parameter, public smear_fermi_dirac
Definition: smear.F90:173
integer, parameter, public smear_methfessel_paxton
Definition: smear.F90:173
subroutine, public smear_write_info(this, namespace, iunit)
Definition: smear.F90:869
real(real64) function, public smear_step_function(this, xx)
Definition: smear.F90:718
integer, parameter, public smear_semiconductor
Definition: smear.F90:173
integer, parameter, public smear_lorentzian
Definition: smear.F90:173
subroutine, public smear_copy(to, from)
Definition: smear.F90:331
integer, parameter, public smear_fixed_occ
Definition: smear.F90:173
integer, parameter, public smear_spline
Definition: smear.F90:173
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:486
subroutine, public smear_init(this, namespace, ispin, fixed_occ, integral_occs, kpoints)
Definition: smear.F90:189
integer, parameter, public smear_cold
Definition: smear.F90:173
logical pure function, public smear_is_semiconducting(this)
Definition: smear.F90:862
integer, parameter, public smear_gaussian
Definition: smear.F90:173
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)