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