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 integer, parameter :: nitmax = 200
258 real(real64), parameter :: tol = 1.0e-10_real64
259 integer :: ist, ik, iter, maxq, weight, sumq_int
260 real(real64) :: sumq_frac, sum_weight
261 logical :: conv
262 real(real64), allocatable :: eigenval_list(:)
263 integer, allocatable :: k_list(:), reorder(:)
264 integer :: fermi_count_up, fermi_count_down
265
267
268 maxq = this%el_per_state * nst * this%nspins
269 if (maxq - qtot <= -tol) then ! not enough states
270 message(1) = 'Not enough states'
271 write(message(2),'(6x,a,f12.6,a,i10)')'(total charge = ', qtot, &
272 ' max charge = ', maxq
273 call messages_fatal(2, namespace=namespace)
274 end if
275
276 conv = .true.
277 if (this%method == smear_fixed_occ) then ! Fermi energy: last of occupied states
278 ist_cycle: do ist = nst, 1, -1
279 do ik = 1, nik
280 if (occupations(ist, ik) > 1e-5_real64) then
281 this%e_fermi = eigenvalues(ist, ik)
282 this%ef_occ = occupations(ist, ik) / this%el_per_state
283 exit ist_cycle
284 end if
285
286 end do
287 end do ist_cycle
288
289 else if (this%method == smear_semiconductor) then
290 ! first we sort the eigenvalues
291 safe_allocate(eigenval_list(1:nst * nik))
292 safe_allocate( k_list(1:nst * nik))
293 safe_allocate( reorder(1:nst * nik))
294
295 iter = 1
296 do ist = 1, nst
297 do ik = 1, nik
298 eigenval_list(iter) = eigenvalues(ist, ik)
299 k_list(iter) = ik
300 reorder(iter) = iter
301 iter = iter + 1
302 end do
303 end do
304
305 call sort(eigenval_list, reorder)
306
307 sumq_int = int(qtot) * this%nik_factor
308 sumq_frac = qtot * this%nik_factor - sumq_int
309
310 do iter = 1, nst * nik
311 weight = int(kweights(k_list(reorder(iter))) * this%nik_factor + m_half)
312 if (.not. weight > 0) cycle
313 this%e_fermi = eigenval_list(iter)
314 this%ef_occ = (sumq_int + sumq_frac) / (weight * this%el_per_state)
315
316 if (sumq_int - weight * this%el_per_state <= 0) then
317 ! count how many occupied states are at the fermi level,
318 ! this is required later to fill the states
319 this%fermi_count = 1
320 fermi_count_up = 1
321 fermi_count_down = 1
322 sum_weight = weight
323 do
324 if (iter - fermi_count_down < 1) exit
325 if (abs(this%e_fermi - eigenval_list(iter - fermi_count_down)) > tol_smear) exit
326 weight = int(kweights(k_list(reorder(iter - fermi_count_down))) * this%nik_factor + m_half)
327 if (weight > m_epsilon) then
328 sumq_int = sumq_int + weight * this%el_per_state
329 sum_weight = sum_weight + weight
330 end if
331 fermi_count_down = fermi_count_down + 1
332 end do
333 do
334 if (iter + fermi_count_up > nst*nik) exit
335 if (abs(this%e_fermi - eigenval_list(iter + fermi_count_up)) > tol_smear) exit
336 weight = int(kweights(k_list(reorder(iter + fermi_count_up))) * this%nik_factor + m_half)
337 if (weight > m_epsilon) then
338 sum_weight = sum_weight + weight
339 end if
340 fermi_count_up = fermi_count_up + 1
341 end do
342 this%fermi_count = fermi_count_up + fermi_count_down - 1
343 this%ef_occ = (sumq_int + sumq_frac) / (sum_weight * this%el_per_state)
344 exit
345 end if
346
347 sumq_int = sumq_int - weight * this%el_per_state
348 end do
349 assert(this%ef_occ < m_one + m_epsilon)
350
351 safe_deallocate_a(eigenval_list)
352 safe_deallocate_a(k_list)
353 safe_deallocate_a(reorder)
354
355 else ! bisection
356 if (this%photodop) then
357 ! find fermi energy for valence electrons
358 call bisection_find_fermi_energy(this, namespace, this%dsmear, tol, &
359 eigenvalues, kweights, nik, qtot-this%nephotodop, 1, this%photodop_bandmin-1, &
360 this%e_fermi, this%ef_occ)
361
362 ! find fermi energy for conduction electrons
363 call bisection_find_fermi_energy(this, namespace, this%dsmear_cond, tol, &
364 eigenvalues, kweights, nik, this%nephotodop, this%photodop_bandmin, nst, &
365 this%e_fermi_cond, this%ef_occ_cond)
366 else
367 call bisection_find_fermi_energy(this, namespace, this%dsmear, tol, &
368 eigenvalues, kweights, nik, qtot, 1, nst, this%e_fermi, this%ef_occ)
369 end if
370 end if
371
373 end subroutine smear_find_fermi_energy
374
375 ! ---------------------------------------------------------
376 subroutine bisection_find_fermi_energy(this, namespace, dsmear_in, tol, &
377 eigenvalues, kweights, nik, q_in, start_band, end_band, e_fermi, ef_occ)
378 type(smear_t), intent(inout) :: this
379 type(namespace_t), intent(in) :: namespace
380 real(real64), intent(in) :: dsmear_in
381 real(real64), intent(in) :: tol
382 integer, intent(in) :: nik
383 real(real64), intent(in) :: eigenvalues(:,:)
384 real(real64), intent(in) :: kweights(:), q_in
385 integer, intent(in) :: start_band, end_band
386 real(real64), intent(out) :: e_fermi, ef_occ
387
388 integer, parameter :: nitmax = 200
389 integer :: ist, ik, iter
390 real(real64) :: drange, dsmear, emin, emax, xx, sumq
391 logical :: conv
392
394
395 dsmear = max(1e-14_real64, dsmear_in)
396 drange = dsmear * sqrt(-log(tol * 0.01_real64))
397
398 emin = minval(eigenvalues) - drange
399 emax = maxval(eigenvalues) + drange
400
401 do iter = 1, nitmax
402 e_fermi = m_half * (emin + emax)
403 sumq = m_zero
404
405 do ik = 1, nik
406 do ist = start_band, end_band
407 xx = (e_fermi - eigenvalues(ist, ik))/dsmear
408 sumq = sumq + kweights(ik) * this%el_per_state * &
409 smear_step_function(this, xx)
410 end do
411 end do
412
413 conv = (abs(sumq - q_in) <= tol)
414 if (conv) exit
415
416 if (sumq <= q_in) emin = e_fermi
417 if (sumq >= q_in) emax = e_fermi
418
419 ef_occ = smear_step_function(this, m_zero)
420 end do
421
422 if (.not. conv) then
423 message(1) = 'Fermi: did not converge.'
424 call messages_fatal(1, namespace=namespace)
425 end if
426
428 end subroutine bisection_find_fermi_energy
429
430 ! ---------------------------------------------------------
431 subroutine smear_fill_occupations(this, eigenvalues, occupations, nik, nst)
432 type(smear_t), intent(in) :: this
433 real(real64), intent(in) :: eigenvalues(:,:)
434 real(real64), intent(inout) :: occupations(:,:)
435 integer, intent(in) :: nik, nst
436
437 integer :: ik, ist, ifermi
438 real(real64) :: dsmear, xx, dsmear_cond
439
440 push_sub(smear_fill_occupations)
441
442 if (this%method == smear_fixed_occ) then
443 ! do nothing
444 else if (this%method == smear_semiconductor) then
445 assert(this%fermi_count > 0 .and. this%fermi_count <= nik*nst)
446
447 ifermi = 0
448 do ik = 1, nik
449 do ist = 1, nst
450 xx = eigenvalues(ist, ik) - this%e_fermi
451 if (xx < -tol_smear) then
452 occupations(ist, ik) = this%el_per_state
453 else if (abs(xx) <= tol_smear .and. ifermi < this%fermi_count) then
454 occupations(ist, ik) = this%ef_occ * this%el_per_state
455 ifermi = ifermi + 1
456 else
457 occupations(ist, ik) = m_zero
458 end if
459
460 end do
461 end do
462
463 else
464 dsmear = max(1e-14_real64, this%dsmear)
465 if (this%photodop) then
466 dsmear_cond = max(1e-14_real64, this%dsmear_cond)
467 do ik = 1, nik
468 do ist = 1, nst
469 if (ist < this%photodop_bandmin) then
470 ! valence electrons
471 xx = (this%e_fermi - eigenvalues(ist, ik))/dsmear
472 else
473 ! conduction electrons
474 xx = (this%e_fermi_cond - eigenvalues(ist, ik))/dsmear_cond
475 end if
476 occupations(ist, ik) = smear_step_function(this, xx) * this%el_per_state
477 end do
478 end do
479 else
480 do ik = 1, nik
481 do ist = 1, nst
482 xx = (this%e_fermi - eigenvalues(ist, ik))/dsmear
483 occupations(ist, ik) = smear_step_function(this, xx) * this%el_per_state
484 end do
485 end do
486 end if
487 end if
488
490 end subroutine smear_fill_occupations
491
492
493 !--------------------------------------------------
494 real(real64) function smear_calc_entropy(this, eigenvalues, &
495 nik, nst, kweights, occ) result(entropy)
496 type(smear_t), intent(inout) :: this
497 real(real64), intent(in) :: eigenvalues(:,:)
498 real(real64), intent(in) :: kweights(:)
499 integer, intent(in) :: nik, nst
500 real(real64), intent(in) :: occ(:, :)
501
502 integer :: ist, ik
503 real(real64) :: dsmear, xx, term, ff
504
505 push_sub(smear_calc_entropy)
506
507 entropy = m_zero
508
509 if (this%method == smear_fixed_occ .or. this%method == smear_semiconductor) then
510 ! Fermi-Dirac entropy, not quite the same as will be obtained with true smearing
511 ! RM Wentzcovitch, JL Martins, and PB Allen, Phys. Rev. B 45, 11372 (1992) eqn (5)
512 ! also N Marzari PhD thesis p 117, http://quasiamore.mit.edu/phd/Marzari_PhD.pdf
513 do ik = 1, nik
514 do ist = 1, nst
515 ff = occ(ist, ik) / this%el_per_state
516 if (ff > m_zero .and. ff < m_one) then
517 term = ff * log(ff) + (1 - ff) * log(1 - ff)
518 else ! we have semiconducting smearing, or perverse occupations as in Methfessel-Paxton
519 term = m_zero
520 end if
521 entropy = entropy - kweights(ik) * this%el_per_state * term
522 end do
523 end do
524 else
525 dsmear = max(1e-14_real64, this%dsmear)
526
527 do ik = 1, nik
528 do ist = 1, nst
529 if (eigenvalues(ist, ik) < huge(m_one)) then
530 xx = (this%e_fermi - eigenvalues(ist, ik)) / dsmear
531 entropy = entropy - kweights(ik) * this%el_per_state * &
532 smear_entropy_function(this, xx)
533 end if
534 end do
535 end do
536 end if
537
538 pop_sub(smear_calc_entropy)
539 end function smear_calc_entropy
540
541
542 ! ---------------------------------------------------------
543 real(real64) function smear_delta_function(this, xx) result(deltaf)
544 type(smear_t), intent(in) :: this
545 real(real64), intent(in) :: xx
546
547 real(real64), parameter :: maxarg = 200.0_real64
548 real(real64) :: xp, arg, hd, hp, aa
549 integer :: ii, ni
550
551 ! no PUSH_SUB, called too often
552
553 ! smear_delta_function is not defined for SMEAR_FIXED_OCC
554 assert(this%method /= smear_fixed_occ)
555
556 deltaf = m_zero
557 select case (this%method)
559 if (abs(xx) <= m_epsilon) then
560 deltaf = this%ef_occ
561 end if
562
563 case (smear_fermi_dirac)
564 if (abs(xx) <= 36.0_real64) then
565 deltaf = m_one / (m_two + exp(-xx) + exp(xx))
566 end if
567
568 case (smear_cold)
569 xp = xx - m_one / sqrt(m_two)
570 arg = min(maxarg, xp**2)
571
572 deltaf = exp(-arg) / sqrt(m_pi) * (m_two - sqrt(m_two) * xx)
573
575 arg = min(maxarg, xx**2)
576 deltaf = exp(-arg) / sqrt(m_pi)
577
578 if (this%MP_n > 0) then ! recursion
579 hd = m_zero
580 hp = exp(-arg)
581 ni = 0
582 aa = m_one / sqrt(m_pi)
583 do ii = 1, this%MP_n
584 hd = m_two * xx * hp - m_two * ni * hd
585 ni = ni + 1
586 aa = -aa / (m_four * ii)
587 hp = m_two * xx * hd - m_two * ni * hp
588 ni = ni + 1
589 deltaf = deltaf + aa * hp
590 end do
591 end if
592
593 case (smear_spline)
594 xp = abs(xx) + m_one / sqrt(m_two)
595 deltaf = sqrt(m_e) * xp * exp(-xp * xp)
596
597 end select
598
599 end function smear_delta_function
600
601
602 ! ---------------------------------------------------------
603 real(real64) function smear_step_function(this, xx) result(stepf)
604 type(smear_t), intent(in) :: this
605 real(real64), intent(in) :: xx
606
607 real(real64), parameter :: maxarg = 200.0_real64
608 real(real64) :: xp, arg, hd, hp, aa
609 integer :: ii, ni
610
611 push_sub(smear_step_function)
612
613 ! smear_step_function is not defined for SMEAR_FIXED_OCC
614 assert(this%method /= smear_fixed_occ)
615
616 stepf = m_zero
617 select case (this%method)
619 if (xx > m_zero) then
620 stepf = m_one
621 else if (abs(xx) <= m_epsilon) then
622 stepf = this%ef_occ
623 end if
624
625 case (smear_fermi_dirac)
626 if (xx > maxarg) then
627 stepf = m_one
628 else if (xx > -maxarg) then
629 stepf = m_one / (m_one + exp(-xx))
630 end if
631
632 case (smear_cold)
633 xp = xx - m_one / sqrt(m_two)
634 arg = min(maxarg, xp**2)
635
636 stepf = m_half * loct_erf(xp) + &
637 m_one / sqrt(m_two * m_pi) * exp(-arg) + m_half
638
640 stepf = m_half * loct_erfc(-xx)
641
642 if (this%MP_n > 0) then ! recursion
643 hd = m_zero
644 arg = min(maxarg, xx**2)
645 hp = exp(-arg)
646 ni = 0
647 aa = m_one / sqrt(m_pi)
648 do ii = 1, this%MP_n
649 hd = m_two * xx * hp - m_two * ni * hd
650 ni = ni + 1
651 aa = -aa / (m_four * ii)
652 stepf = stepf - aa * hd
653 hp = m_two * xx * hd - m_two * ni * hp
654 ni = ni + 1
655 end do
656 end if
657
658 case (smear_spline)
659 if (xx <= m_zero) then
660 xp = xx - m_one / sqrt(m_two)
661 stepf = m_half * sqrt(m_e) * exp(-xp * xp)
662 else
663 xp = xx + m_one / sqrt(m_two)
664 stepf = m_one - m_half * sqrt(m_e) * exp(-xp * xp)
665 end if
666
667 end select
668
669 pop_sub(smear_step_function)
670 end function smear_step_function
671
672
673 ! ---------------------------------------------------------
675 real(real64) function smear_entropy_function(this, xx) result(entropyf)
676 type(smear_t), intent(in) :: this
677 real(real64), intent(in) :: xx
678
679 real(real64), parameter :: maxarg = 200.0_real64
680 real(real64) :: xp, arg, hd, hp, hpm1, aa
681 integer :: ii, ni
682
683 push_sub(smear_entropy_function)
684
685 ! smear_entropy_function is not defined for SMEAR_FIXED_OCC
686 assert(this%method /= smear_fixed_occ)
687
688 entropyf = m_zero
689 select case (this%method)
691
692 case (smear_fermi_dirac)
693 if (abs(xx) <= 36.0_real64) then
694 xp = m_one / (m_one + exp(-xx))
695 entropyf = xp * log(xp) + (m_one - xp) * log(m_one - xp)
696 end if
697
698 case (smear_cold)
699 xp = xx - m_one / sqrt(m_two)
700 arg = min(maxarg, xp**2)
701
702 entropyf = m_one / sqrt(m_two * m_pi) * xp * exp(-arg)
703
705 arg = min(maxarg, xx**2)
706 entropyf = -m_half * exp(-arg) / sqrt(m_pi)
707
708 if (this%MP_n > 0) then ! recursion
709 hd = m_zero
710 hp = exp(-arg)
711 ni = 0
712 aa = m_one / sqrt(m_pi)
713 do ii = 1, this%MP_n
714 hd = m_two * xx * hp - m_two * ni * hd
715 ni = ni + 1
716 hpm1 = hp
717 hp = m_two * xx * hd - m_two * ni * hp
718 ni = ni + 1
719 aa = -aa / (m_four * ii)
720 entropyf = entropyf - aa * (m_half * hp + hpm1 * ni)
721 end do
722 end if
723
724 case (smear_spline)
725 xp = abs(xx) + m_one / sqrt(m_two)
726 entropyf = -sqrt(m_e) * (abs(xx) * exp(-xp * xp) / m_two + sqrt(m_pi) / m_four * loct_erfc(xp))
727
728 end select
729
731 end function smear_entropy_function
732
733
734 ! ---------------------------------------------------------
735 logical pure function smear_is_semiconducting(this) result(answer)
736 type(smear_t), intent(in) :: this
737
738 answer = this%method == smear_semiconductor
739
740 end function smear_is_semiconducting
741
742 subroutine smear_write_info(this, namespace, iunit)
743 type(smear_t), intent(in) :: this
744 type(namespace_t), intent(in) :: namespace
745 integer, optional, intent(in) :: iunit
746
747
748 if (this%method /= smear_semiconductor .and. this%method /= smear_fixed_occ) then
749 if (this%photodop) then
750 write(message(1), '(a,f12.6,1x,a)') "Fermi energy (valence ) = ", &
751 units_from_atomic(units_out%energy, this%e_fermi), units_abbrev(units_out%energy)
752 write(message(2), '(a,f12.6,1x,a)') "Fermi energy (conduction) = ", &
753 units_from_atomic(units_out%energy, this%e_fermi_cond), units_abbrev(units_out%energy)
754 call messages_info(2, iunit=iunit, namespace=namespace)
755 else
756 write(message(1), '(a,f12.6,1x,a)') "Fermi energy = ", &
757 units_from_atomic(units_out%energy, this%e_fermi), units_abbrev(units_out%energy)
758 call messages_info(1, iunit=iunit, namespace=namespace)
759 end if
760 end if
761
762 end subroutine smear_write_info
763
764
765end module smear_oct_m
766
767!! Local Variables:
768!! mode: f90
769!! coding: utf-8
770!! 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__
real(real64), parameter, public m_two
Definition: global.F90:189
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public p_ry
Definition: global.F90:220
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_min_occ
Minimal occupation that is considered to be non-zero.
Definition: global.F90:210
integer function, public kpoints_kweight_denominator(this)
Definition: kpoints.F90:1609
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1057
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
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:723
real(real64) function, public smear_calc_entropy(this, eigenvalues, nik, nst, kweights, occ)
Definition: smear.F90:589
subroutine, public smear_fill_occupations(this, eigenvalues, occupations, nik, nst)
Definition: smear.F90:525
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:769
real(real64) function, public smear_delta_function(this, xx)
Definition: smear.F90:637
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:836
real(real64) function, public smear_step_function(this, xx)
Definition: smear.F90:697
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:471
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:829
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)