00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 2008, Andrew W. Steiner 00005 00006 This file is part of O2scl. 00007 00008 O2scl is free software; you can redistribute it and/or modify 00009 it under the terms of the GNU General Public License as published by 00010 the Free Software Foundation; either version 3 of the License, or 00011 (at your option) any later version. 00012 00013 O2scl is distributed in the hope that it will be useful, 00014 but WITHOUT ANY WARRANTY; without even the implied warranty of 00015 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00016 GNU General Public License for more details. 00017 00018 You should have received a copy of the GNU General Public License 00019 along with O2scl. If not, see <http://www.gnu.org/licenses/>. 00020 00021 ------------------------------------------------------------------- 00022 */ 00023 #ifndef CFL_NJL_EOS_H 00024 #define CFL_NJL_EOS_H 00025 00026 #include <iostream> 00027 #include <complex> 00028 00029 #include <gsl/gsl_math.h> 00030 #include <gsl/gsl_eigen.h> 00031 #include <gsl/gsl_complex_math.h> 00032 #include <gsl/gsl_complex.h> 00033 #include <gsl/gsl_sort_vector.h> 00034 #include <gsl/gsl_poly.h> 00035 00036 #include <o2scl/constants.h> 00037 #include <o2scl/err_hnd.h> 00038 #include <o2scl/ovector_cx_tlate.h> 00039 #include <o2scl/omatrix_tlate.h> 00040 #include <o2scl/omatrix_cx_tlate.h> 00041 #include <o2scl/multi_funct.h> 00042 #include <o2scl/mm_funct.h> 00043 #include <o2scl/cx_arith.h> 00044 #include <o2scl/vec_arith.h> 00045 #include <o2scl/gsl_inte_qng.h> 00046 #include <o2scl/poly.h> 00047 #include <o2scl/test_mgr.h> 00048 #include <o2scl/columnify.h> 00049 00050 #include <o2scl/part.h> 00051 #include <o2scl/nambujl_eos.h> 00052 00053 #ifndef DOXYGENP 00054 namespace o2scl { 00055 #endif 00056 00057 /** 00058 \brief Nambu Jona-Lasinio model with a schematic CFL 00059 di-quark interaction at finite temperature 00060 00061 The variable B0 must be set before use. 00062 00063 The original Lagrangian is 00064 00065 \f[ 00066 {\cal L} = 00067 {\cal L}_{\mathrm{Dirac}} + 00068 {\cal L}_{\mathrm{4-fermion}} + 00069 {\cal L}_{\mathrm{6-fermion}} + 00070 {\cal L}_{CSC1} + 00071 {\cal L}_{CSC2} 00072 \f] 00073 00074 \f[ 00075 {\cal L}_{\mathrm{Dirac}} = \bar{q}_{i \alpha} 00076 \left( i {\partial} \delta_{i j} \delta_{\alpha \beta} - 00077 m_{i j} \delta_{\alpha \beta} - \mu_{i j,~\alpha \beta} \gamma^0 00078 \right) q_{j \beta} 00079 \f] 00080 00081 \f[ 00082 {\cal L}_{\mathrm{4-fermion}} = 00083 G_S \sum_{a=0}^8 \left[ 00084 \left( \bar{q} \lambda^a_f q \right)^2 + 00085 \left( \bar{q} i \gamma_5 \lambda^a_f q \right)^2 00086 \right] 00087 \f] 00088 00089 \f[ 00090 {\cal L}_{\mathrm{6-fermion}} = 00091 - G_{D} \left[ 00092 {\mathrm det}_{i j} \, \bar{q}_{i \alpha} \left( 1 + i \gamma_5 \right) 00093 q_{j \beta} + 00094 {\mathrm det}_{i j} \, \bar{q}_{i \alpha} \left( 1 - i \gamma_5 \right) 00095 q_{j \beta} \right] \delta_{\alpha \beta} 00096 \f] 00097 00098 \f[ 00099 {\cal L}_{CSC1} = 00100 G_{DIQ} \sum_k \sum_{\gamma} \left[ 00101 \left(\bar{q}_{i \alpha} \epsilon_{i j k} 00102 \epsilon_{\alpha \beta \gamma} q^C_{j \beta}\right) 00103 \left(\bar{q}_{i^{\prime} \alpha^{\prime}}^C 00104 \epsilon_{i^{\prime} j^{\prime} k} \epsilon_{\alpha^{\prime} 00105 \beta^{\prime} \gamma} q_{j^{\prime} \beta^{\prime}}\right)\right] 00106 \f] 00107 00108 \f[ 00109 {\cal L}_{CSC2} = 00110 G_{DIQ} \sum_k \sum_{\gamma} \left[ 00111 \left(\bar{q}_{i \alpha} i \gamma_5 \epsilon_{i j k} 00112 \epsilon_{\alpha \beta \gamma} q^C_{j \beta}\right) 00113 \left(\bar{q}_{i^{\prime} \alpha^{\prime}}^C i \gamma_5 00114 \epsilon_{i^{\prime} j^{\prime} k} \epsilon_{\alpha^{\prime} 00115 \beta^{\prime} \gamma} q_{j^{\prime} \beta^{\prime}}\right) \right] \,, 00116 \f] 00117 00118 where \f$ \mu \f$ is the quark number chemical potential. 00119 couplings \f$ G_S \f$, \f$ G_D \f$, and \f$ G_{DIQ} \f$ 00120 ultra-violet three-momentum cutoff, \f$ \Lambda \f$ 00121 00122 The thermodynamic potential is 00123 \f[ 00124 \Omega(\mu_i,\left<\bar{q} q\right>_i,\left< q q\right>_i,T) 00125 = \Omega_{\mathrm{vac}}+\Omega_{\mathrm{stat}} + \Omega_0 00126 \f] 00127 where \f$ i \f$ runs over all nine (three colors times three 00128 flavors) quarks. We assume that the condensates are independent 00129 of color and 00130 that the quark chemical potentials 00131 are of the form 00132 \f$ \mu_Q=\mu_{\mathrm{Flavor(Q)}}+\mu_{\mathrm{Color(Q)}} \f$ 00133 with 00134 \f[ 00135 \mu_{\mathrm{red}} = \mu_3 + \mu_8/\sqrt{3} 00136 \quad 00137 \mu_{\mathrm{green}} = -\mu_3 + \mu_8/\sqrt{3} 00138 \quad 00139 \mu_{\mathrm{blue}} = -2 \mu_8 /\sqrt{3} 00140 \f] 00141 With these assumptions, the thermodynamic potential as given 00142 by the function thd_potential(), is a function of 12 variables 00143 \f[ 00144 \Omega(\mu_u, \mu_d, \mu_s, \mu_3, \mu_8, \left<\bar{u} u\right>, 00145 \left<\bar{d} d\right>, \left<\bar{s} s\right>, 00146 \left< u d\right>, \left< u s\right>, \left< d s\right>, 00147 T) 00148 \f] 00149 00150 The individual terms are 00151 \f[ 00152 \Omega_{\mathrm{stat}} = - \frac{1}{2} \int \frac{d^3 00153 p}{\left(2 \pi\right)^3} \, \sum_{i=1}^{72} \left[ \frac{\lambda_i}{2} + 00154 T \ln{\left(1 + e^{-\lambda_i/T} \right)} \right] 00155 \f] 00156 00157 \f[ 00158 \Omega_{\mathrm{vac}} = 00159 - 2 G_S \sum_{i=u,d,s} \langle {\bar q_i} q_i \rangle^2 00160 +4 G_D \left<{\bar u} u\right> \left<{\bar d} d \right> 00161 \left<{\bar s} s\right> 00162 + \sum_k \sum_{\gamma} \frac{\left|\Delta^{k \gamma}\right|^2}{4 G_{DIQ}} 00163 \f] 00164 00165 where \f$ \lambda_i \f$ are the eigenvalues of the (72 by 72) matrix 00166 (calculated by the function eigenvalues()) 00167 \f[ 00168 D = \left[ 00169 \begin{array}{cc} 00170 - \gamma^0 \vec{\gamma} \cdot \vec{p} - M_{i} \gamma^0 + \mu_{i \alpha} 00171 & \Delta i \gamma^0 \gamma_5 C \\ 00172 i \Delta \gamma^0 C \gamma_5 00173 & - \gamma^0 \vec{\gamma}^T \cdot \vec{p} + M_{i} \gamma^0 00174 - \mu_{i \alpha}\end{array} 00175 \right] 00176 \f] 00177 and \f$ C \f$ is the charge conjugation matrix (in the Dirac 00178 representation). 00179 00180 The values of the various condensates are usually determined by 00181 the condition 00182 \f[ 00183 \frac{\partial \Omega}{\left<\bar{q} q\right>_i} = 0 00184 \quad 00185 \frac{\partial \Omega}{\left<q q\right>_i} = 0 00186 \f] 00187 00188 Note that setting fixed_mass to \c true and setting all of the 00189 gaps to zero when \c gap_limit is less than zero will reproduce 00190 an analog of the bag model with a momentum cutoff. 00191 00192 The variable nambujl_eos::fromqq is automatically set to true in 00193 the constructor, as computations with \c fromqq=false are not 00194 implemented. 00195 00196 \todo Improve documentation (Note that it appears to report 00197 the member functions twice?) 00198 \todo Allow user to change derivative object? This isn't possible 00199 right now because the stepsize parameter of the derivative 00200 object is used. 00201 \todo This class internally mixes ovector, omatrix, gsl_vector 00202 and gsl_matrix objects in a confusing and non-optimal way. Fix this. 00203 00204 */ 00205 class cfl_njl_eos : public nambujl_eos { 00206 public: 00207 00208 cfl_njl_eos(); 00209 00210 virtual ~cfl_njl_eos(); 00211 00212 /** 00213 \brief Set the parameters and the bag constant 'B0' 00214 00215 This function allows the user to specify the momentum cutoff, 00216 \c lambda, the four-fermion coupling \c fourferm, the 00217 six-fermion coupling from the 't Hooft interaction \c sixferm, 00218 and the color-superconducting coupling, \c fourgap. If 0.0 is 00219 given for any of the values, then the default is used (\f$ 00220 \Lambda=602.3/(\hbar c), G=1.835/\Lambda^2, K=12.36/\Lambda^5 00221 \f$). 00222 00223 If the four-fermion coupling that produces a gap is not 00224 specified, it is automatically set to 3/4 G, which is the 00225 value obtained from the Fierz transformation. 00226 00227 The value of the shift in the bag constant nambujl_eos::B0 is 00228 automatically calculated to ensure that the vacuum has zero 00229 energy density and zero pressure. The functions set_quarks() 00230 and set_thermo() must be used before hand to specify the \ref 00231 quark and \ref thermo objects. 00232 00233 */ 00234 virtual int set_parameters(double lambda=0.0, double fourferm=0.0, 00235 double sixferm=0.0, double fourgap=0.0); 00236 00237 /** 00238 \brief Calculate the EOS 00239 00240 Calculate the EOS from the quark condensates in \c u.qq, \c 00241 d.qq and \c s.qq. Return the mass gap equations in \c qq1, \c 00242 qq2, \c qq3, and the normal gap equations in \c gap1, \c gap2, 00243 and \c gap3. 00244 00245 Using \c fromqq=false as in nambujl_eos and nambujl_eos does not 00246 work here and will return an error. Also, the quarks must be 00247 set through quark_eos::quark_set() before use. 00248 00249 If all of the gaps are less than gap_limit, then the 00250 nambujl_eos::calc_temp_p() is used, and \c gap1, \c gap2, and 00251 \c gap3 are set to equal \c u.del, \c d.del, and \c s.del, 00252 respectively. 00253 00254 \todo It surprises me that n3 is not -res[11]. Is there 00255 a sign error in the color densities? 00256 */ 00257 virtual int calc_eq_temp_p(quark &u, quark &d, quark &s, double &qq1, 00258 double &qq2, double &qq3, double &gap1, 00259 double &gap2, double &gap3, double mu3, 00260 double mu8, double &n3, double &n8, 00261 thermo &qb, const double temper); 00262 00263 /// Check the derivatives specified by eigenvalues() 00264 virtual int test_derivatives(double lmom, double mu3, double mu8, 00265 test_mgr &t); 00266 00267 /** 00268 \brief Calculate the energy eigenvalues as a function of the 00269 momentum 00270 00271 Given the momentum \c mom, and the chemical potentials 00272 associated with the third and eighth gluons (\c mu3 and \c 00273 mu8), the energy eigenvalues are computed in egv[0] 00274 ... egv[35]. 00275 */ 00276 virtual int eigenvalues(double lmom, double mu3, double mu8, 00277 double egv[36], double dedmuu[36], 00278 double dedmud[36], double dedmus[36], 00279 double dedmu[36], double dedmd[36], 00280 double dedms[36], double dedu[36], 00281 double dedd[36], double deds[36], 00282 double dedmu3[36], double dedmu8[36]); 00283 00284 /// Set the routine for solving quartics 00285 int set_quartic(quartic_real_coeff &q) { 00286 quartic=&q; 00287 return 0; 00288 } 00289 00290 /// The equal mass threshold 00291 double eq_limit; 00292 00293 /// Set to true to test the integration (default false) 00294 bool integ_test; 00295 00296 /// Test the integration routines 00297 int test_integration(test_mgr &t); 00298 00299 // 00300 //gsl_quadratic_real_coeff quad; 00301 00302 /** 00303 \brief The default quartic routine 00304 00305 Slightly better accuracy (with slower execution times) can be 00306 achieved using \ref gsl_poly_real_coeff which polishes the 00307 roots of the quartics. For example 00308 00309 \code 00310 cfl_njl_eos cfl; 00311 gsl_poly_real_coeff gp; 00312 cfl.set_quartic(gp); 00313 \endcode 00314 */ 00315 cern_quartic_real_coeff def_quartic; 00316 //@} 00317 00318 /** 00319 \brief Test the routine to compute the eigenvalues of 00320 non-superfluid fermions 00321 */ 00322 int test_normal_eigenvalues(test_mgr &t); 00323 00324 /** 00325 \brief Test the routine to compute the eigenvalues of 00326 superfluid fermions 00327 */ 00328 int test_gapped_eigenvalues(test_mgr &t); 00329 00330 /** 00331 \brief Smallest allowable gap (default 0.0) 00332 00333 If any of the gaps are below this value, then it is assumed 00334 that they are zero and the equation of state is simplified 00335 accordingly. If all of the gaps are less than gap_limit, then 00336 the results from nambujl_eos are used in 00337 calc_eq_temp_p(), calc_temp_p() and thd_potential(). 00338 */ 00339 double gap_limit; 00340 00341 /** 00342 \brief If this is true, then finite temperature 00343 corrections are ignored (default false) 00344 00345 This implements some simplifications in the 00346 momentum integration that are not possible at finite temperature. 00347 */ 00348 bool zerot; 00349 00350 /** 00351 \brief Use a fixed quark mass and ignore the quark 00352 condensates 00353 */ 00354 bool fixed_mass; 00355 00356 /// If true, then ensure color neutrality 00357 bool color_neut; 00358 00359 /** 00360 \brief Diquark coupling constant (default 3 G/4) 00361 00362 The default value is the one derived from a Fierz 00363 transformation. (\ref Buballa04) 00364 */ 00365 double GD; 00366 00367 /** 00368 \brief The absolute precision for the integration 00369 (default \f$ 10^{-4} \f$ ) 00370 00371 This is analogous to gsl_inte::epsabs 00372 */ 00373 double inte_epsabs; 00374 00375 /** 00376 \brief The relative precision for the integration 00377 (default \f$ 10^{-4} \f$ ) 00378 00379 This is analogous to gsl_inte::epsrel 00380 */ 00381 double inte_epsrel; 00382 00383 /** 00384 \brief The number of points used in the last integration 00385 (default 0) 00386 00387 This returns 21, 43, or 87 depending on the number of function 00388 evaluations needed to obtain the desired precision. If it 00389 the routine failes to obtain the desired precision, then 00390 this variable is set to 88. 00391 */ 00392 size_t inte_npoints; 00393 00394 /// Return string denoting type ("cfl_njl_eos") 00395 virtual const char *type() { return "cfl_njl_eos"; }; 00396 00397 #ifndef DOXYGEN_INTERNAL 00398 00399 protected: 00400 00401 friend class io_tlate<cfl_njl_eos>; 00402 00403 /** 00404 \brief The integrands 00405 00406 - res[0] is the thermodynamic potential, \f$ \Omega \f$ 00407 - res[1] is \f$ d -\Omega / d T \f$ 00408 - res[2] is \f$ d \Omega / d \mu_u \f$ 00409 - res[3] is \f$ d \Omega / d \mu_d \f$ 00410 - res[4] is \f$ d \Omega / d \mu_s \f$ 00411 - res[5] is \f$ d \Omega / d m_u \f$ 00412 - res[6] is \f$ d \Omega / d m_d \f$ 00413 - res[7] is \f$ d \Omega / d m_s \f$ 00414 - res[8] is \f$ d \Omega / d \Delta_{ds} \f$ 00415 - res[9] is \f$ d \Omega / d \Delta_{us} \f$ 00416 - res[10] is \f$ d \Omega / d \Delta_{ud} \f$ 00417 - res[11] is \f$ d \Omega / d \mu_3 \f$ 00418 - res[12] is \f$ d \Omega / d \mu_8 \f$ 00419 */ 00420 virtual int integrands(double p, double res[]); 00421 00422 /// Compute ungapped eigenvalues and the appropriate derivatives 00423 int normal_eigenvalues(double m, double lmom, double mu, 00424 double lam[2], double dldmu[2], 00425 double dldm[2]); 00426 00427 /** 00428 \brief Treat the simply gapped quarks in all cases gracefully 00429 00430 This function uses the quarks \c q1 and \c q2 to construct the 00431 eigenvalues of the inverse propagator, properly handling the 00432 either zero or finite quark mass and either zero or finite 00433 quark gaps. In the case of finite quark mass and finite quark 00434 gaps, the quartic solver is used. 00435 00436 The chemical potentials are separated so we can add the 00437 color chemical potentials to the quark chemical potentials 00438 if necessary. 00439 00440 This function is used by eigenvalues(). It does not work for 00441 the "ur-dg-sb" set of quarks which are paired in a non-trivial 00442 way. 00443 00444 \todo Only the "ms" part of the quarks is referenced, so 00445 we should rewrite to use only double's as function arguments, 00446 and avoid passing pointers to quark objects. 00447 */ 00448 int gapped_eigenvalues(double m1, double m2, double lmom, 00449 double mu1, double mu2, double tdelta, 00450 double lam[4], double dldmu1[4], 00451 double dldmu2[4], double dldm1[4], 00452 double dldm2[4], double dldg[4]); 00453 00454 /// Temperature 00455 double temper; 00456 00457 /// 3rd gluon chemical potential 00458 double smu3; 00459 00460 /// 8th gluon chemical potential 00461 double smu8; 00462 00463 /// \name Numerical methods 00464 //@{ 00465 /// The routine to solve quartics 00466 quartic_real_coeff *quartic; 00467 00468 /// \name For computing eigenvalues 00469 //@{ 00470 /// Inverse propagator matrix 00471 omatrix_cx iprop; 00472 00473 /// The eigenvectors 00474 omatrix_cx eivec; 00475 00476 /// The derivative of the inverse propagator wrt the ds gap 00477 omatrix_cx dipdgapu; 00478 /// The derivative of the inverse propagator wrt the us gap 00479 omatrix_cx dipdgapd; 00480 /// The derivative of the inverse propagator wrt the ud gap 00481 omatrix_cx dipdgaps; 00482 00483 /// The eigenvalues 00484 ovector eval; 00485 00486 /// Workspace for eigenvalue computation 00487 gsl_eigen_hermv_workspace *w; 00488 //@} 00489 00490 /// \name For the integration 00491 //@{ 00492 /// The error scaling function for integ_err 00493 double rescale_error(double err, const double result_abs, 00494 const double result_asc); 00495 00496 /** 00497 \brief A new version of gsl_inte_qng to integrate several 00498 functions at the same time 00499 */ 00500 int integ_err(double a, double b, const size_t nr, 00501 ovector &res, double &err2); 00502 //@} 00503 00504 private: 00505 00506 cfl_njl_eos(const cfl_njl_eos &); 00507 cfl_njl_eos& operator=(const cfl_njl_eos&); 00508 00509 #endif 00510 00511 }; 00512 00513 #ifndef DOXYGENP 00514 } 00515 #endif 00516 00517 #endif
Documentation generated with Doxygen and provided under the GNU Free Documentation License. See License Information for details.
Project hosting provided by
,
O2scl Sourceforge Project Page