Class eos_quark_cfl (o2scl)¶
-
class eos_quark_cfl : public o2scl::eos_quark_njl¶
Nambu Jona-Lasinio model with a schematic CFL di-quark interaction at finite temperature.
The variable B0 must be set before use.
The original Lagrangian is
\[ {\cal L} = {\cal L}_{\mathrm{Dirac}} + {\cal L}_{\mathrm{4-fermion}} + {\cal L}_{\mathrm{6-fermion}} + {\cal L}_{CSC1} + {\cal L}_{CSC2} \]\[ {\cal L}_{\mathrm{Dirac}} = \bar{q}_{i \alpha} \left( i {\partial} \delta_{i j} \delta_{\alpha \beta} - m_{i j} \delta_{\alpha \beta} - \mu_{i j,~\alpha \beta} \gamma^0 \right) q_{j \beta} \]\[ {\cal L}_{\mathrm{4-fermion}} = G_S \sum_{a=0}^8 \left[ \left( \bar{q} \lambda^a_f q \right)^2 + \left( \bar{q} i \gamma_5 \lambda^a_f q \right)^2 \right] \]\[ {\cal L}_{\mathrm{6-fermion}} = - G_{D} \left[ {\mathrm det}_{i j} \, \bar{q}_{i \alpha} \left( 1 + i \gamma_5 \right) q_{j \beta} + {\mathrm det}_{i j} \, \bar{q}_{i \alpha} \left( 1 - i \gamma_5 \right) q_{j \beta} \right] \delta_{\alpha \beta} \]\[ {\cal L}_{CSC1} = G_{DIQ} \sum_k \sum_{\gamma} \left[ \left(\bar{q}_{i \alpha} \epsilon_{i j k} \epsilon_{\alpha \beta \gamma} q^C_{j \beta}\right) \left(\bar{q}_{i^{\prime} \alpha^{\prime}}^C \epsilon_{i^{\prime} j^{\prime} k} \epsilon_{\alpha^{\prime} \beta^{\prime} \gamma} q_{j^{\prime} \beta^{\prime}}\right)\right] \]\[ {\cal L}_{CSC2} = G_{DIQ} \sum_k \sum_{\gamma} \left[ \left(\bar{q}_{i \alpha} i \gamma_5 \epsilon_{i j k} \epsilon_{\alpha \beta \gamma} q^C_{j \beta}\right) \left(\bar{q}_{i^{\prime} \alpha^{\prime}}^C i \gamma_5 \epsilon_{i^{\prime} j^{\prime} k} \epsilon_{\alpha^{\prime} \beta^{\prime} \gamma} q_{j^{\prime} \beta^{\prime}}\right) \right] \,, \]where \( \mu \)
is the quark number chemical potential.
couplings
\( G_S \), \( G_D \), and \( G_{DIQ} \) ultra-violet three-momentum cutoff, \( \Lambda \)The thermodynamic potential is
\[ \Omega(\mu_i,\left<\bar{q} q\right>_i,\left< q q\right>_i,T) = \Omega_{\mathrm{vac}}+\Omega_{\mathrm{stat}} + \Omega_0 \]where \( i \) runs over all nine (three colors times three flavors) quarks. We assume that the condensates are independent of color and that the quark chemical potentials are of the form \( \mu_Q=\mu_{\mathrm{Flavor(Q)}}+\mu_{\mathrm{Color(Q)}} \) with\[ \mu_{\mathrm{red}} = \mu_3 + \mu_8/\sqrt{3} \quad \mu_{\mathrm{green}} = -\mu_3 + \mu_8/\sqrt{3} \quad \mu_{\mathrm{blue}} = -2 \mu_8 /\sqrt{3} \]With these assumptions, the thermodynamic potential as given by the function thd_potential(), is a function of 12 variables\[ \Omega(\mu_u, \mu_d, \mu_s, \mu_3, \mu_8, \left<\bar{u} u\right>, \left<\bar{d} d\right>, \left<\bar{s} s\right>, \left< u d\right>, \left< u s\right>, \left< d s\right>, T) \]The individual terms are
\[ \Omega_{\mathrm{stat}} = - \frac{1}{2} \int \frac{d^3 p}{\left(2 \pi\right)^3} \, \sum_{i=1}^{72} \left[ \frac{\lambda_i}{2} + T \ln{\left(1 + e^{-\lambda_i/T} \right)} \right] \]\[ \Omega_{\mathrm{vac}} = - 2 G_S \sum_{i=u,d,s} \langle {\bar q_i} q_i \rangle^2 +4 G_D \left<{\bar u} u\right> \left<{\bar d} d \right> \left<{\bar s} s\right> + \sum_k \sum_{\gamma} \frac{\left|\Delta^{k \gamma}\right|^2}{4 G_{DIQ}} \]where \( \lambda_i \) are the eigenvalues of the (72 by 72) matrix (calculated by the function eigenvalues())
\[\begin{split} D = \left[ \begin{array}{cc} - \gamma^0 \vec{\gamma} \cdot \vec{p} - M_{i} \gamma^0 + \mu_{i \alpha} & \Delta i \gamma^0 \gamma_5 C \\ i \Delta \gamma^0 C \gamma_5 & - \gamma^0 \vec{\gamma}^T \cdot \vec{p} + M_{i} \gamma^0 - \mu_{i \alpha}\end{array} \right] \end{split}\]and \( C \) is the charge conjugation matrix (in the Dirac representation).The values of the various condensates are usually determined by the condition
\[ \frac{\partial \Omega}{\left<\bar{q} q\right>_i} = 0 \quad \frac{\partial \Omega}{\left<q q\right>_i} = 0 \]Note that setting fixed_mass to
true
and setting all of the gaps to zero whengap_limit
is less than zero will reproduce an analog of the bag model with a momentum cutoff.The variable eos_quark_njl::fromqq is automatically set to true in the constructor, as computations with
fromqq=false
are not implemented.- Idea for Future:
This class internally mixes ubvector, ubmatrix, gsl_vector and gsl_matrix objects in a confusing and non-optimal way. Fix this.
Allow user to change derivative object? This isn’t possible right now because the stepsize parameter of the derivative object is used.
References:
Created for [Steiner02].
Subclassed by o2scl::eos_quark_cfl6
Numerical methods
-
typedef boost::numeric::ublas::vector<double> ubvector¶
-
typedef boost::numeric::ublas::vector<std::complex<double>> ubvector_complex¶
-
typedef boost::numeric::ublas::matrix<std::complex<double>> ubmatrix_complex¶
-
quartic_real_coeff *quartic¶
The routine to solve quartics.
For computing eigenvalues
-
gsl_matrix_complex *iprop¶
Inverse propagator matrix.
-
gsl_matrix_complex *eivec¶
The eigenvectors.
-
ubmatrix_complex dipdgapu¶
The derivative of the inverse propagator wrt the ds gap.
-
ubmatrix_complex dipdgapd¶
The derivative of the inverse propagator wrt the us gap.
-
ubmatrix_complex dipdgaps¶
The derivative of the inverse propagator wrt the ud gap.
-
gsl_vector *eval¶
The eigenvalues.
-
gsl_eigen_hermv_workspace *w¶
Workspace for eigenvalue computation.
For the integration
-
double rescale_error(double err, double result_abs, double result_asc)¶
The error scaling function for integ_err.
-
int integ_err(double a, double b, const size_t nr, ubvector &res, double &err2)¶
A new version of inte_qng_gsl to integrate several functions at the same time.
-
eos_quark_cfl(const eos_quark_cfl&)¶
-
eos_quark_cfl &operator=(const eos_quark_cfl&)¶
Public Functions
-
eos_quark_cfl()¶
-
virtual ~eos_quark_cfl()¶
-
virtual int set_parameters_cfl(double lambda = 0.0, double fourferm = 0.0, double sixferm = 0.0, double fourgap = 0.0)¶
Set the parameters and the bag constant ‘B0’.
This function allows the user to specify the momentum cutoff,
lambda
, the four-fermion couplingfourferm
, the six-fermion coupling from the ‘t Hooft interactionsixferm
, and the color-superconducting coupling,fourgap
. If 0.0 is given for any of the values, then the default is used ( \( \Lambda=602.3/(\hbar c), G=1.835/\Lambda^2, K=12.36/\Lambda^5 \)).If the four-fermion coupling that produces a gap is not specified, it is automatically set to 3/4 G, which is the value obtained from the Fierz transformation.
The value of the shift in the bag constant eos_quark_njl::B0 is automatically calculated to ensure that the vacuum has zero energy density and zero pressure. The functions set_quarks() and set_thermo() must be used before hand to specify the quark and thermo objects.
-
virtual int calc_eq_temp_p(quark &u, quark &d, quark &s, double &qq1, double &qq2, double &qq3, double &gap1, double &gap2, double &gap3, double mu3, double mu8, double &n3, double &n8, thermo &qb, double temper)¶
Calculate the EOS.
Calculate the EOS from the quark condensates in
u.qq
,d.qq
ands.qq
. Return the mass gap equations inqq1
,qq2
,qq3
, and the normal gap equations ingap1
,gap2
, andgap3
.Using
fromqq=false
as in eos_quark_njl and eos_quark_njl does not work here and will return an error. Also, the quarks must be set through eos_quark::quark_set() before use.If all of the gaps are less than gap_limit, then the eos_quark_njl::calc_temp_p() is used, and
gap1
,gap2
, andgap3
are set to equalu.del
,d.del
, ands.del
, respectively.Todo
In eos_quark_cfl::calc_eq_temp_p(): It surprises me that n3 is not -res[11]. Is there a sign error in the color densities?
-
virtual int test_derivatives(double lmom, double mu3, double mu8, test_mgr &t)¶
Check the derivatives specified by eigenvalues()
-
virtual int eigenvalues(double lmom, double mu3, double mu8, double egv[36], double dedmuu[36], double dedmud[36], double dedmus[36], double dedmu[36], double dedmd[36], double dedms[36], double dedu[36], double dedd[36], double deds[36], double dedmu3[36], double dedmu8[36])¶
Calculate the energy eigenvalues as a function of the momentum.
Given the momentum
mom
, and the chemical potentials associated with the third and eighth gluons (mu3
andmu8
), the energy eigenvalues are computed in egv[0] … egv[35].
-
inline int set_quartic(quartic_real_coeff<> &q)¶
Set the routine for solving quartics.
-
int test_normal_eigenvalues(test_mgr &t)¶
Test the routine to compute the eigenvalues of non-superfluid fermions.
-
int test_gapped_eigenvalues(test_mgr &t)¶
Test the routine to compute the eigenvalues of superfluid fermions.
-
inline virtual const char *type()¶
Return string denoting type (“eos_quark_cfl”)
Public Members
-
double eq_limit¶
The equal mass threshold.
-
bool integ_test¶
Set to true to test the integration (default false)
-
quartic_real_coeff_cern def_quartic¶
The default quartic routine.
Slightly better accuracy (with slower execution times) can be achieved using poly_real_coeff_gsl which polishes the roots of the quartics. For example
eos_quark_cfl cfl; poly_real_coeff_gsl gp; cfl.set_quartic(gp);
-
double gap_limit¶
Smallest allowable gap (default 0.0)
If any of the gaps are below this value, then it is assumed that they are zero and the equation of state is simplified accordingly. If all of the gaps are less than gap_limit, then the results from eos_quark_njl are used in calc_eq_temp_p(), calc_temp_p() and thd_potential().
-
bool zerot¶
If this is true, then finite temperature corrections are ignored (default false)
This implements some simplifications in the momentum integration that are not possible at finite temperature.
-
bool fixed_mass¶
Use a fixed quark mass and ignore the quark condensates.
-
bool color_neut¶
If true, then ensure color neutrality.
-
double GD¶
Diquark coupling constant (default 3 G/4)
The default value is the one derived from a Fierz transformation. ([Buballa04])
-
double inte_epsabs¶
The absolute precision for the integration (default \( 10^{-4} \) )
This is analogous to gsl_inte::epsabs
-
double inte_epsrel¶
The relative precision for the integration (default \( 10^{-4} \) )
This is analogous to gsl_inte::epsrel
-
size_t inte_npoints¶
The number of points used in the last integration (default 0)
This returns 21, 43, or 87 depending on the number of function evaluations needed to obtain the desired precision. If it the routine failes to obtain the desired precision, then this variable is set to 88.
Protected Functions
-
virtual int integrands(double p, double res[])¶
The integrands.
res[0] is the thermodynamic potential, \( \Omega \)
res[1] is \( d -\Omega / d T \)
res[2] is \( d \Omega / d \mu_u \)
res[3] is \( d \Omega / d \mu_d \)
res[4] is \( d \Omega / d \mu_s \)
res[5] is \( d \Omega / d m_u \)
res[6] is \( d \Omega / d m_d \)
res[7] is \( d \Omega / d m_s \)
res[8] is \( d \Omega / d \Delta_{ds} \)
res[9] is \( d \Omega / d \Delta_{us} \)
res[10] is \( d \Omega / d \Delta_{ud} \)
res[11] is \( d \Omega / d \mu_3 \)
res[12] is \( d \Omega / d \mu_8 \)
-
int normal_eigenvalues(double m, double lmom, double mu, double lam[2], double dldmu[2], double dldm[2])¶
Compute ungapped eigenvalues and the appropriate derivatives.
-
int gapped_eigenvalues(double m1, double m2, double lmom, double mu1, double mu2, double tdelta, double lam[4], double dldmu1[4], double dldmu2[4], double dldm1[4], double dldm2[4], double dldg[4])¶
Treat the simply gapped quarks in all cases gracefully.
This function uses the quarks
q1
andq2
to construct the eigenvalues of the inverse propagator, properly handling the either zero or finite quark mass and either zero or finite quark gaps. In the case of finite quark mass and finite quark gaps, the quartic solver is used.The chemical potentials are separated so we can add the color chemical potentials to the quark chemical potentials if necessary.
This function is used by eigenvalues(). It does not work for the “ur-dg-sb” set of quarks which are paired in a non-trivial way.
Todo
In eos_quark_cfl::gapped_eigenvalues(): In the code, the equal mass case seems to be commented out. Why?