30 #include <o2scl/test_mgr.h> 31 #include <o2scl/eos_quark_cfl.h> 33 #ifndef DOXYGEN_NO_O2NS 214 double &qq1,
double &qq2,
double &qq3,
215 double &gap1,
double &gap2,
double &gap3,
216 double mu3,
double mu8,
217 double &n3,
double &n8,
thermo &qb,
221 virtual int integrands(
double p,
double res[]);
238 virtual int eigenvalues6(
double lmom,
double mu3,
double mu8,
239 double egv[36],
double dedmuu[36],
240 double dedmud[36],
double dedmus[36],
241 double dedmu[36],
double dedmd[36],
242 double dedms[36],
double dedu[36],
243 double dedd[36],
double deds[36],
244 double dedmu3[36],
double dedmu8[36]);
253 virtual int make_matrices(
double lmom,
double mu3,
double mu8,
254 double egv[36],
double dedmuu[36],
255 double dedmud[36],
double dedmus[36],
256 double dedmu[36],
double dedmd[36],
257 double dedms[36],
double dedu[36],
258 double dedd[36],
double deds[36],
259 double dedmu3[36],
double dedmu8[36]);
265 virtual const char *
type() {
return "eos_quark_cfl6"; };
279 #ifndef DOXYGEN_INTERNAL 315 gsl_eigen_hermv_workspace *
w6;
326 #ifndef DOXYGEN_NO_O2NS Nambu Jona-Lasinio model with a schematic CFL di-quark interaction at finite temperature.
ubmatrix_complex dipdqqd
The derivative wrt the down quark condensate.
virtual int integrands(double p, double res[])
The momentum integrands.
ubmatrix_complex dipdgaps
The derivative wrt the ud gap.
ubmatrix_complex dipdgapd
The derivative wrt the us gap.
static const int mat_size
The size of the matrix to be diagonalized.
gsl_matrix_complex * eivec6
The eigenvectors.
gsl_eigen_hermv_workspace * w6
GSL workspace for the eigenvalue computation.
virtual int test_derivatives(double lmom, double mu3, double mu8, test_mgr &t)
Check the derivatives specified by eigenvalues()
ubmatrix_complex dipdqqs
The derivative wrt the strange quark condensate.
An EOS like eos_quark_cfl but with a color-superconducting 't Hooft interaction.
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.
gsl_vector * eval6
Storage for the eigenvalues.
gsl_matrix_complex * iprop6
Storage for the inverse propagator.
double temper
Temperature.
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.
double kdlimit
The absolute value below which the CSC 't Hooft coupling is ignored(default )
int set_masses()
Set the quark effective masses from the gaps and the condensates.
ubmatrix_complex dipdqqu
The derivative wrt the up quark condensate.
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.
virtual const char * type()
Return string denoting type ("eos_quark_cfl6")
double KD
The color superconducting 't Hooft coupling (default 0)
ubmatrix_complex dipdgapu
The derivative wrt the ds gap.