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