Class eos_quark_cfl6 (o2scl)

O2scl : Class List

class eos_quark_cfl6 : public o2scl::eos_quark_cfl

An EOS like eos_quark_cfl but with a color-superconducting ‘t Hooft interaction.

Beginning with the Lagrangian:

\[ {\cal L} = {\cal L}_{Dirac} + {\cal L}_{NJL} + {\cal L}_{'t Hooft} + {\cal L}_{SC} + {\cal L}_{SC6} \]
\[ {\cal L}_{Dirac} = {\bar q} \left( i \partial -m - \mu \gamma^0 \right) q \]
\[ {\cal L}_{NJL} = G_S \sum_{a=0}^8 \left[ \left( {\bar q} \lambda^a q \right)^2 - \left( {\bar q} \lambda^a \gamma^5 q \right)^2 \right] \]
\[ {\cal L}_{'t Hooft} = G_D \left[ \mathrm{det}_f {\bar q} \left(1-\gamma^5 \right) q +\mathrm{det}_f {\bar q} \left(1+\gamma^5 \right) q \right] \]
\[ {\cal L}_{SC} = G_{DIQ} \left( {\bar q}_{i \alpha} i \gamma^5 \varepsilon^{i j k} \varepsilon^{\alpha \beta \gamma} q^C_{j \beta} \right) \left( {\bar q}_{\ell \delta} i \gamma^5 \epsilon^{\ell m k} \epsilon^{\delta \varepsilon \gamma} q^C_{m \varepsilon} \right) \]
\[ {\cal L}_{SC6} = K_D \left( {\bar q}_{i \alpha} i \gamma^5 \varepsilon^{i j k} \varepsilon^{\alpha \beta \gamma} q^C_{j \beta} \right) \left( {\bar q}_{\ell \delta} i \gamma^5 \epsilon^{\ell m n} \epsilon^{\delta \varepsilon \eta} q^C_{m \varepsilon} \right) \left( {\bar q}_{k \gamma} q_{n \eta} \right) \]

We can simplify the relevant terms in \({\cal L}_{NJL}\):

\[ {\cal L}_{NJL} = G_S \left[ \left({\bar u} u\right)^2+ \left({\bar d} d\right)^2+ \left({\bar s} s\right)^2 \right] \]
and in \({\cal L}_{'t Hooft}\):
\[ {\cal L}_{NJL} = G_D \left( {\bar u} u {\bar d} d {\bar s} s \right) \]

Using the definition:

\[ \Delta^{k \gamma} = \left< {\bar q} i \gamma^5 \epsilon \epsilon q^C_{} \right> \]
and the ansatzes:
\[ ({\bar q}_1 q_2) ({\bar q}_3 q_4) \rightarrow {\bar q}_1 q_2 \left< {\bar q}_3 q_4 \right> +{\bar q}_3 q_4 \left< {\bar q}_1 q_2 \right> -\left< {\bar q}_1 q_2 \right> \left< {\bar q}_3 q_4 \right> \]
\[ ({\bar q}_1 q_2) ({\bar q}_3 q_4) ({\bar q}_5 q_6) \rightarrow {\bar q}_1 q_2 \left< {\bar q}_3 q_4 \right> \left< {\bar q}_5 q_6 \right> +{\bar q}_3 q_4 \left< {\bar q}_1 q_2 \right> \left< {\bar q}_5 q_6 \right> +{\bar q}_5 q_6 \left< {\bar q}_1 q_2 \right> \left< {\bar q}_3 q_4 \right> -2\left< {\bar q}_1 q_2 \right> \left< {\bar q}_3 q_4 \right> \left< {\bar q}_5 q_6 \right> \]
for the mean field approximation, we can rewrite the Lagrangian
\[ {\cal L}_{NJL} = 2 G_S \left[ \left( {\bar u} u \right) \left< {\bar u} u \right> +\left( {\bar d} d \right) \left< {\bar d} d \right> +\left( {\bar s} s \right) \left< {\bar s} s \right> - \left< {\bar u} u \right>^2 - \left< {\bar d} d \right>^2 - \left< {\bar s} s \right>^2 \right] \]
\[ {\cal L}_{'t Hooft} = - 2 G_D \left[ \left( {\bar u} u \right) \left< {\bar u} u \right> \left< {\bar s} s \right> + \left( {\bar d} d \right) \left< {\bar u} u \right> \left< {\bar s} s \right> + \left( {\bar s} s \right) \left< {\bar u} u \right> \left< {\bar d} d \right> - 2 \left< {\bar u} u \right>\left< {\bar d} d \right> \left< {\bar s} s \right> \right] \]
\[ {\cal L}_{SC} = G_{DIQ} \left[ \Delta^{k \gamma} \left( {\bar q}_{\ell \delta} i \gamma^5 \epsilon^{\ell m k} \epsilon^{\delta \varepsilon \gamma} q^C_{m \varepsilon} \right) + \left( {\bar q}_{i \alpha} i \gamma^5 \varepsilon^{i j k} \varepsilon^{\alpha \beta \gamma} q^C_{j \beta} \right) \Delta^{k \gamma \dagger} - \Delta^{k \gamma} \Delta^{k \gamma \dagger} \right] \]
\[ {\cal L}_{SC6} = K_D \left[ \left( {\bar q}_{m \varepsilon} q_{n \eta} \right) \Delta^{k \gamma} \Delta^{m \varepsilon \dagger} + \left( {\bar q}_{i \alpha} i \gamma^5 \varepsilon^{i j k} \varepsilon^{\alpha \beta \gamma} q^C_{j \beta} \right) \Delta^{m \varepsilon \dagger} \left< {\bar q}_{m \varepsilon} q_{n \eta} \right> \right] \]
\[ + K_D \left[\Delta^{k \gamma} \left( {\bar q}_{\ell \delta} i \gamma^5 \epsilon^{\ell m n} \epsilon^{\delta \varepsilon \eta} q^C_{m \varepsilon} \right) \left< {\bar q}_{m \varepsilon} q_{n \eta} \right> -2 \Delta^{k \gamma} \Delta^{m \varepsilon \dagger} \left< {\bar q}_{m \varepsilon} q_{n \eta} \right> \right] \]

If we make the definition \( {\tilde \Delta} = 2 G_{DIQ} \Delta \)

References:

Created for [Steiner05].

Public Types

typedef boost::numeric::ublas::vector<std::complex<double>> ubvector_complex
typedef boost::numeric::ublas::matrix<std::complex<double>> ubmatrix_complex

Public Functions

eos_quark_cfl6()
virtual ~eos_quark_cfl6()
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. Return the mass gap equations in qq1, qq2, qq3, and the normal gap equations in gap1, gap2, and gap3.

Using fromqq=true as in eos_quark_njl and nambujl_temp_eos does not work here and will return an error.

If all of the gaps are less than gap_limit, then the nambujl_temp_eos::calc_temp_p() is used, and gap1, gap2, and gap3 are set to equal u.del, d.del, and s.del, respectively.

virtual int integrands(double p, double res[])

The momentum integrands.

virtual int test_derivatives(double lmom, double mu3, double mu8, test_mgr &t)

Check the derivatives specified by eigenvalues()

virtual int eigenvalues6(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 and their derivatives.

Given the momentum mom, and the chemical potentials associated with the third and eighth gluons (mu3 and mu8), this computes the eigenvalues of the inverse propagator and the assocated derivatives.

Note that this is not the same as eos_quark_cfl::eigenvalues() which returns dedmu rather dedqqu.

virtual int make_matrices(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])

Construct the matrices, but don’t solve the eigenvalue problem.

This is used by check_derivatives() to make sure that the derivative entries are right.

inline virtual const char *type()

Return string denoting type (“eos_quark_cfl6”)

Public Members

double KD

The color superconducting ‘t Hooft coupling (default 0)

double kdlimit

The absolute value below which the CSC ‘t Hooft coupling is ignored(default \( 10^{-6} \))

Protected Functions

int set_masses()

Set the quark effective masses from the gaps and the condensates.

Protected Attributes

gsl_matrix_complex *iprop6

Storage for the inverse propagator.

gsl_matrix_complex *eivec6

The eigenvectors.

ubmatrix_complex dipdgapu

The derivative wrt the ds gap.

ubmatrix_complex dipdgapd

The derivative wrt the us gap.

ubmatrix_complex dipdgaps

The derivative wrt the ud gap.

ubmatrix_complex dipdqqu

The derivative wrt the up quark condensate.

ubmatrix_complex dipdqqd

The derivative wrt the down quark condensate.

ubmatrix_complex dipdqqs

The derivative wrt the strange quark condensate.

gsl_vector *eval6

Storage for the eigenvalues.

gsl_eigen_hermv_workspace *w6

GSL workspace for the eigenvalue computation.

Protected Static Attributes

static const int mat_size = 36

The size of the matrix to be diagonalized.

Private Functions

eos_quark_cfl6(const eos_quark_cfl6&)
eos_quark_cfl6 &operator=(const eos_quark_cfl6&)