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