Class fermion_thermo_tl (o2scl)

O2scl : Class List

template<class fermion_t = fermion, class fd_inte_t = fermi_dirac_integ_gsl, class be_inte_t = bessel_K_exp_integ_gsl, class root_t = root_cern<>, class func_t = funct, class fp_t = double>
class fermion_thermo_tl : public o2scl::fermion_zerot_tl<double>

Fermion with finite-temperature thermodynamics [abstract base].

This is an abstract base for the computation of finite-temperature fermionic statistics. Different children (e.g. fermion_rel_tl) use different techniques to computing the momentum integrations.

Because massless fermions at finite temperature are much simpler, there are separate member functions included in this class to handle them. The functions massless_calc_density() and massless_calc_mu() compute the thermodynamics of massless fermions at finite temperature given the density or the chemical potentials. The functions massless_pair_density() and massless_pair_mu() perform the same task, but automatically include antiparticles.

The function massless_calc_density() uses a root object to solve for the chemical potential as a function of the density. The default is an object of type root_cern. The function massless_pair_density() does not need to use the root object because of the simplification afforded by the inclusion of antiparticles.

Todo

In class fermion_thermo_tl:

  • Future: Create a Chebyshev approximation for inverting the the Fermi functions for massless_calc_density() functions?

Template types:

  • fermion_t: the type of particle object which will be passed to the methods defined by this class

  • fd_inte_t: the integration object for massless fermions

  • be_inte_t: the integration object for the nondegenerate limit

  • root_t: the solver for massless fermions

  • func_t: the function object type

  • fp_t: the floating point type

Subclassed by o2scl::fermion_rel_tl< fermion_deriv_tl< double > >, o2scl::fermion_rel_tl< fermion_tl< cpp_dec_float_25 >, fermi_dirac_integ_direct< cpp_dec_float_25, funct_cdf35, cpp_dec_float_35 >, bessel_K_exp_integ_boost< cpp_dec_float_25, cpp_dec_float_35 >, inte_double_exp_boost<>, inte_double_exp_boost<>, root_cern< funct_cdf25, cpp_dec_float_25 >, funct_cdf25, cpp_dec_float_25 >, o2scl::fermion_rel_tl< fermion_tl< long double >, fermi_dirac_integ_direct< long double, funct_cdf25, cpp_dec_float_25 >, bessel_K_exp_integ_boost< long double, cpp_dec_float_25 >, inte_double_exp_boost<>, inte_double_exp_boost<>, root_cern< funct_ld, long double >, funct_ld, long double >, o2scl::fermion_rel_tl< fermion_deriv >, o2scl::fermion_nonrel_tl< fermion_t, fd_inte_t, be_inte_t, root_t, func_t, fp_t >, o2scl::fermion_rel_tl< fermion_t, fd_inte_t, be_inte_t, nit_t, dit_t, root_t, func_t, fp_t >

Massless fermions

root_t *massless_root

A pointer to the solver for massless fermions.

root_t def_massless_root

The default solver for massless_calc_density()

We default to a solver of type root_cern here since we don’t have a bracket or a derivative.

inline void massless_calc_mu(fermion_t &f, fp_t temper)

Finite temperature massless fermions.

inline void massless_calc_density(fermion_t &f, fp_t temper)

Finite temperature massless fermions.

inline void massless_pair_mu(fermion_t &f, fp_t temper)

Finite temperature massless fermions and antifermions.

When \( \mu=0 \) (e.g. for photons), the density is zero, the energy density is

\[ \varepsilon = \frac{7g T^4}{120 \pi^2} \]
and the pressure is \( P = \varepsilon/3 \). When the chemical potential is finite, we have
\[\begin{eqnarray*} n &=& \frac{g \mu^3}{6 \pi^2} \left[ 1 + \left( \frac{\pi T}{\mu} \right)^2 \right] \varepsilon &=& \frac{g \mu^4}{8 \pi^2} \left[ 1 + 2 \left( \frac{\pi T}{\mu} \right)^2 + \frac{7}{15} \left( \frac{\pi T}{\mu} \right)^4 \right] \end{eqnarray*}\]
and the pressure is \( P = \varepsilon/3 \).

inline void massless_pair_density(fermion_t &f, fp_t temper)

Finite temperature massless fermions and antifermions.

In the cases \( n^3 \gg T \) and \( T \gg n^3 \) , expansions are used instead of the exact formulas to avoid loss of precision.

In particular, using the parameter

\[ \alpha = \frac{g^2 \pi^2 T^6}{243 n^2} \]
and defining the expression
\[ \mathrm{cbt} = \alpha^{-1/6} \left( -1 + \sqrt{1+\alpha}\right)^{1/3} \]
we can write the chemical potential as
\[ \mu = \frac{\pi T}{\sqrt{3}} \left(\frac{1}{\mathrm{cbt}} - \mathrm{cbt} \right) \]

These expressions, however, do not work well when \( \alpha \) is very large or very small, so series expansions are used whenever \( \alpha > 10^{4} \) or \( \alpha < 3 \times 10^{-4} \). For small \( \alpha \),

\[ \left(\frac{1}{\mathrm{cbt}} - \mathrm{cbt} \right) \approx \frac{2^{1/3}}{\alpha^{1/6}} - \frac{\alpha^{1/6}}{2^{1/3}} + \frac{\alpha^{5/6}}{6{\cdot}2^{2/3}} + \frac{\alpha^{7/6}}{12{\cdot}2^{1/3}} - \frac{\alpha^{11/6}}{18{\cdot}2^{2/3}} - \frac{5 \alpha^{13/6}}{144{\cdot}2^{1/3}} + \frac{77 \alpha^{17/6}}{2592{\cdot}2^{2/3}} \]
and for large \( \alpha \),
\[ \left(\frac{1}{\mathrm{cbt}} - \mathrm{cbt} \right) \approx \frac{2}{3} \sqrt{\frac{1}{\alpha}} - \frac{8}{81} \left(\frac{1}{\alpha}\right)^{3/2} + \frac{32}{729} \left(\frac{1}{\alpha}\right)^{5/2} \]

This approach works to within about 1 part in \( 10^{12} \), and is tested in fermion_ts.cpp.

Todo

In function massless_pair_density()

  • Future: This could be improved by including more terms in the expansions.

inline void set_massless_root(root_t &rp)

Set the solver for use in massless_calc_density()

inline virtual const char *type()

Return string denoting type (“fermion_thermo”)

inline void ndeg_terms(size_t j, fp_t tt, fp_t xx, fp_t m, bool inc_rest_mass, bool inc_antip, fp_t &pterm, fp_t &nterm, fp_t &enterm, fp_t &edterm)

Compute a term in the nondegenerate expansion.

This function uses be_integ to perform the Bessel integrals.

inline fp_t massless_solve_fun(fp_t x, fermion_t &f, fp_t temper)

Solve for the chemical potential for massless fermions.

Public Functions

inline fermion_thermo_tl()
inline virtual ~fermion_thermo_tl()
inline bool calc_mu_ndeg(fermion_t &f, fp_t temper, fp_t prec = 1.0e-18, bool inc_antip = false, int verbose = 0)

Non-degenerate expansion for fermions.

Attempts to evaluate thermodynamics of a non-degenerate fermion. If the result is accurate to within the requested precision, this function returns true, and otherwise this function returns false and the values in stored in the pr, n, en, and ed field are meaningless.

If \( \mu \) is negative and sufficiently far from zero, then the thermodynamic quantities are smaller than the smallest representable double-precision number. In this case, this function will return true and report all quantities as zero.

The following uses the notation of [Johns96].

Defining \( \psi \equiv (\mu-m)/T \), \( t \equiv T/m \), and \( d \equiv g~m^4/(2 \pi^2) \) the pressure in the non-degenerate limit ( \( \psi \rightarrow - \infty \)) is

\[ P = d \sum_{n=1}^{\infty} P_n \]
where
\[ P_n \equiv \left(-1\right)^{n+1} \left(\frac{t^2}{n^2}\right) e^{n \left(\psi+1/t\right)} K_2 \left( \frac{n}{t} \right) \]
The density is then
\[ n = d \sum_{n=1}^{\infty} \frac{n P_n}{T} \]
and the entropy density is
\[ s = \frac{d}{m} \sum_{n=1}^{\infty} \left\{ \frac{2 P_n}{t} -\frac{n P_n}{t^2}+ \frac{\left(-1\right)^{n+1}}{2 n} e^{n \left(\psi+1/t\right)} \left[ K_1 \left( \frac{n}{t} \right)+K_3 \left( \frac{n}{t} \right) \right] \right\} \]

This function is accurate over a wide range of conditions when \( \psi < -4 \).

The ratio of the nth term to the first term in the pressure series is

\[ R_n \equiv \frac{P_{n}}{P_{1}} = \frac{(-1)^{n+1} e^{(n-1)(\psi+1/t)} K_2(n/t) }{n^2 K_2(1/t)} \]
This function currently uses 20 terms in the series and immediately returns false if \( |R_{20}| \) is greater than prec

In the nondegenerate and nonrelativistic ( \( t \rightarrow 0 \)) limit, the argument to the Bessel functions and the exponential becomes too large. In this case, it’s better to use the expansions, e.g. for \( x \equiv n/t \rightarrow \infty \),

\[ \sqrt{\frac{2 x}{\pi}} e^{x} K_2(x) \approx 1 + \frac{3}{8 x} - \frac{15}{128 x^2} + ... \]

inline bool pair_den_ndeg(fermion_t &f, fp_t temper, fp_t prec = 1.0e-18)

Compute the net density including antiparticles in the nondegenerate approximation.

This function is experimental.

inline bool calc_mu_deg(fermion_t &f, fp_t temper, fp_t prec = 1.0e-18)

Degenerate expansion for fermions.

Attempts to evaulate thermodynamics of a degenerate fermion. If the result is accurate to within the requested precision, this function returns true, and otherwise this function returns false and the values in stored in the pr, n, en, and ed field are meaningless.

The pressure, density, and energy density, should be accurate to the requested precision, but the first term in the series expansion for the entropy is zero, so the entropy is one order lower in accuracy.

Todo

In function calc_mu_deg()

  • Future: Make a function like this for dndm, dsdT, etc. for fermion_deriv .

virtual int calc_mu(fermion_t &f, fp_t temper) = 0

Calculate properties as function of chemical potential.

virtual int calc_density(fermion_t &f, fp_t temper) = 0

Calculate properties as function of density.

Note

This function returns an integer value, in contrast to calc_mu(), because of the potential for non-convergence.

virtual int pair_mu(fermion_t &f, fp_t temper) = 0

Calculate properties with antiparticles as function of chemical potential.

virtual int pair_density(fermion_t &f, fp_t temper) = 0

Calculate properties with antiparticles as function of density.

Note

This function returns an integer value, in contrast to pair_mu(), because of the potential for non-convergence.

Public Members

fd_inte_t fd_integ

Object for Fermi-Dirac integrals (for massless fermions)

be_inte_t be_integ

Object for Bessel-exp integrals (for nondegenerate expansion)

Protected Attributes

o2scl::cubic_real_coeff_gsl2<fp_t> crcg2

Cubic solver for massless fermions.