00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 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 O2SCL_CFL_NJL_EOS_H 00024 #define O2SCL_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 00044 #include <o2scl/gsl_deriv.h> 00045 #include <o2scl/cern_deriv.h> 00046 #include <o2scl/gsl_mmin_conf.h> 00047 #include <o2scl/poly.h> 00048 00049 #include <o2scl/part.h> 00050 #include <o2scl/nambujl_eos.h> 00051 00052 #ifndef DOXYGENP 00053 namespace o2scl { 00054 #endif 00055 00056 /** \brief Nambu Jona-Lasinio model with a schematic CFL 00057 di-quark interaction at finite temperature 00058 00059 The variable B0 must be set before use. 00060 00061 The original Lagrangian is 00062 00063 \f[ 00064 {\cal L} = 00065 {\cal L}_{\mathrm{Dirac}} + 00066 {\cal L}_{\mathrm{4-fermion}} + 00067 {\cal L}_{\mathrm{6-fermion}} + 00068 {\cal L}_{CSC1} + 00069 {\cal L}_{CSC2} 00070 \f] 00071 00072 \f[ 00073 {\cal L}_{\mathrm{Dirac}} = \bar{q}_{i \alpha} 00074 \left( i {\partial} \delta_{i j} \delta_{\alpha \beta} - 00075 m_{i j} \delta_{\alpha \beta} - \mu_{i j,~\alpha \beta} \gamma^0 00076 \right) q_{j \beta} 00077 \f] 00078 00079 \f[ 00080 {\cal L}_{\mathrm{4-fermion}} = 00081 G_S \sum_{a=0}^8 \left[ 00082 \left( \bar{q} \lambda^a_f q \right)^2 + 00083 \left( \bar{q} i \gamma_5 \lambda^a_f q \right)^2 00084 \right] 00085 \f] 00086 00087 \f[ 00088 {\cal L}_{\mathrm{6-fermion}} = 00089 - G_{D} \left[ 00090 {\mathrm det}_{i j} \, \bar{q}_{i \alpha} \left( 1 + i \gamma_5 \right) 00091 q_{j \beta} + 00092 {\mathrm det}_{i j} \, \bar{q}_{i \alpha} \left( 1 - i \gamma_5 \right) 00093 q_{j \beta} \right] \delta_{\alpha \beta} 00094 \f] 00095 00096 \f[ 00097 {\cal L}_{CSC1} = 00098 G_{DIQ} \sum_k \sum_{\gamma} \left[ 00099 \left(\bar{q}_{i \alpha} \epsilon_{i j k} 00100 \epsilon_{\alpha \beta \gamma} q^C_{j \beta}\right) 00101 \left(\bar{q}_{i^{\prime} \alpha^{\prime}}^C 00102 \epsilon_{i^{\prime} j^{\prime} k} \epsilon_{\alpha^{\prime} 00103 \beta^{\prime} \gamma} q_{j^{\prime} \beta^{\prime}}\right)\right] 00104 \f] 00105 00106 \f[ 00107 {\cal L}_{CSC2} = 00108 G_{DIQ} \sum_k \sum_{\gamma} \left[ 00109 \left(\bar{q}_{i \alpha} i \gamma_5 \epsilon_{i j k} 00110 \epsilon_{\alpha \beta \gamma} q^C_{j \beta}\right) 00111 \left(\bar{q}_{i^{\prime} \alpha^{\prime}}^C i \gamma_5 00112 \epsilon_{i^{\prime} j^{\prime} k} \epsilon_{\alpha^{\prime} 00113 \beta^{\prime} \gamma} q_{j^{\prime} \beta^{\prime}}\right) \right] \,, 00114 \f] 00115 00116 where \f$ \mu \f$ is the quark number chemical potential. 00117 couplings \f$ G_S \f$, \f$ G_D \f$, and \f$ G_{DIQ} \f$ 00118 ultra-violet three-momentum cutoff, \f$ \Lambda \f$ 00119 00120 The thermodynamic potential is 00121 \f[ 00122 \Omega(\mu_i,\left<\bar{q} q\right>_i,\left< q q\right>_i,T) 00123 = \Omega_{\mathrm{vac}}+\Omega_{\mathrm{stat}} + \Omega_0 00124 \f] 00125 where \f$ i \f$ runs over all nine (three colors times three 00126 flavors) quarks. We assume that the condensates are independent 00127 of color and 00128 that the quark chemical potentials 00129 are of the form 00130 \f$ \mu_Q=\mu_{\mathrm{Flavor(Q)}}+\mu_{\mathrm{Color(Q)}} \f$ 00131 with 00132 \f[ 00133 \mu_{\mathrm{red}} = \mu_3 + \mu_8/\sqrt{3} 00134 \quad 00135 \mu_{\mathrm{green}} = -\mu_3 + \mu_8/\sqrt{3} 00136 \quad 00137 \mu_{\mathrm{blue}} = -2 \mu_8 /\sqrt{3} 00138 \f] 00139 With these assumptions, the thermodynamic potential as given 00140 by the function thd_potential(), is a function of 12 variables 00141 \f[ 00142 \Omega(\mu_u, \mu_d, \mu_s, \mu_3, \mu_8, \left<\bar{u} u\right>, 00143 \left<\bar{d} d\right>, \left<\bar{s} s\right>, 00144 \left< u d\right>, \left< u s\right>, \left< d s\right>, 00145 T) 00146 \f] 00147 00148 The individual terms are 00149 \f[ 00150 \Omega_{\mathrm{stat}} = - \frac{1}{2} \int \frac{d^3 00151 p}{\left(2 \pi\right)^3} \, \sum_{i=1}^{72} \left[ \frac{\lambda_i}{2} + 00152 T \ln{\left(1 + e^{-\lambda_i/T} \right)} \right] 00153 \f] 00154 00155 \f[ 00156 \Omega_{\mathrm{vac}} = 00157 - 2 G_S \sum_{i=u,d,s} \langle {\bar q_i} q_i \rangle^2 00158 +4 G_D \left<{\bar u} u\right> \left<{\bar d} d \right> 00159 \left<{\bar s} s\right> 00160 + \sum_k \sum_{\gamma} \frac{\left|\Delta^{k \gamma}\right|^2}{4 G_{DIQ}} 00161 \f] 00162 00163 where \f$ \lambda_i \f$ are the eigenvalues of the (72 by 72) matrix 00164 (calculated by the function eigenvalues()) 00165 \f[ 00166 D = \left[ 00167 \begin{array}{cc} 00168 - \gamma^0 \vec{\gamma} \cdot \vec{p} - M_{i} \gamma^0 + \mu_{i \alpha} 00169 & \Delta i \gamma^0 \gamma_5 C \\ 00170 i \Delta \gamma^0 C \gamma_5 00171 & - \gamma^0 \vec{\gamma}^T \cdot \vec{p} + M_{i} \gamma^0 00172 - \mu_{i \alpha}\end{array} 00173 \right] 00174 \f] 00175 and \f$ C \f$ is the charge conjugation matrix (in the Dirac 00176 representation). 00177 00178 The values of the various condensates are usually determined by 00179 the condition 00180 \f[ 00181 \frac{\partial \Omega}{\left<\bar{q} q\right>_i} = 0 00182 \quad 00183 \frac{\partial \Omega}{\left<q q\right>_i} = 0 00184 \f] 00185 00186 Note that setting fixed_mass to \c true and setting all of the 00187 gaps to zero when \c gap_limit is less than zero will reproduce 00188 an analog of the bag model with a momentum cutoff. 00189 00190 The variable nambujl_eos::fromqq is automatically set to true in 00191 the constructor, as computations with \c fromqq=false are not 00192 implemented. 00193 00194 \future This class internally mixes ovector, omatrix, gsl_vector 00195 and gsl_matrix objects in a confusing and non-optimal way. Fix this. 00196 00197 \future Allow user to change derivative object? This isn't possible 00198 right now because the stepsize parameter of the derivative 00199 object is used. 00200 00201 */ 00202 class cfl_njl_eos : public nambujl_eos { 00203 public: 00204 00205 cfl_njl_eos(); 00206 00207 virtual ~cfl_njl_eos(); 00208 00209 /** 00210 \brief Set the parameters and the bag constant 'B0' 00211 00212 This function allows the user to specify the momentum cutoff, 00213 \c lambda, the four-fermion coupling \c fourferm, the 00214 six-fermion coupling from the 't Hooft interaction \c sixferm, 00215 and the color-superconducting coupling, \c fourgap. If 0.0 is 00216 given for any of the values, then the default is used (\f$ 00217 \Lambda=602.3/(\hbar c), G=1.835/\Lambda^2, K=12.36/\Lambda^5 00218 \f$). 00219 00220 If the four-fermion coupling that produces a gap is not 00221 specified, it is automatically set to 3/4 G, which is the 00222 value obtained from the Fierz transformation. 00223 00224 The value of the shift in the bag constant nambujl_eos::B0 is 00225 automatically calculated to ensure that the vacuum has zero 00226 energy density and zero pressure. The functions set_quarks() 00227 and set_thermo() must be used before hand to specify the \ref 00228 quark and \ref thermo objects. 00229 00230 */ 00231 virtual int set_parameters(double lambda=0.0, double fourferm=0.0, 00232 double sixferm=0.0, double fourgap=0.0); 00233 00234 /** \brief Calculate equation of state as a function of chemical 00235 potentials 00236 00237 If any of the gaps given in quark::del are greater than 00238 gap_limit, then the respective gap equations are 00239 solved. If fixed_mass is false, then the quark condensates are 00240 solved for. If color_neut is true, then color neutrality is 00241 enforced. 00242 */ 00243 virtual int calc_temp_p(quark &u, quark &d, quark &s, 00244 const double temper, thermo &th); 00245 00246 /** 00247 \brief Calculate the EOS 00248 00249 Calculate the EOS from the quark condensates in \c u.qq, \c 00250 d.qq and \c s.qq. Return the mass gap equations in \c qq1, \c 00251 qq2, \c qq3, and the normal gap equations in \c gap1, \c gap2, 00252 and \c gap3. 00253 00254 Using \c fromqq=false as in nambujl_eos and nambujl_eos does not 00255 work here and will return an error. Also, the quarks must be 00256 set through quark_eos::quark_set() before use. 00257 00258 If all of the gaps are less than gap_limit, then the 00259 nambujl_eos::calc_temp_p() is used, and \c gap1, \c gap2, and 00260 \c gap3 are set to equal \c u.del, \c d.del, and \c s.del, 00261 respectively. 00262 00263 */ 00264 virtual int calc_eq_temp_p(quark &u, quark &d, quark &s, double &qq1, 00265 double &qq2, double &qq3, double &gap1, 00266 double &gap2, double &gap3, double mu3, 00267 double mu8, 00268 double &n3, double &n8, thermo &qb, 00269 const double temper); 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]. No space is allocated for the array by the 00278 function. 00279 */ 00280 virtual int eigenvalues(double mom, ovector_view &egv, double mu3, 00281 double mu8); 00282 00283 /** \brief Direct calculation of the thermodynamic potential 00284 */ 00285 virtual double thd_potential(quark &u, quark &d, quark &s, 00286 double mu3, double mu8, 00287 const double ltemper); 00288 00289 /// \name Numerical methods 00290 //@{ 00291 /// Set the minimizers for the thermodynamic potential 00292 int set_minimizers(multi_min<void *, multi_funct<void *> > &mi1, 00293 multi_min<void *, multi_funct<void *> > &mi2) { 00294 min1=&mi1; 00295 min2=&mi2; 00296 return 0; 00297 } 00298 00299 /// Set the routine for solving quartics 00300 int set_quartic(quartic_real_coeff &q) { 00301 quartic=&q; 00302 quartic_set=true; 00303 return 0; 00304 } 00305 00306 /// The first default minimizer 00307 gsl_mmin_conf<void *, multi_funct<void *> > def_min1; 00308 /// The second default minimizers 00309 gsl_mmin_conf<void *, multi_funct<void *> > def_min2; 00310 00311 /** 00312 \brief The default quartic routine 00313 00314 Slightly better accuracy (with slower execution times) can be 00315 achieved using \ref gsl_poly_real_coeff which polishes the 00316 roots of the quartics. For example 00317 00318 \code 00319 cfl_njl_eos cfl; 00320 gsl_poly_real_coeff gp; 00321 cfl.set_quartic(gp); 00322 \endcode 00323 */ 00324 cern_quartic_real_coeff def_quartic; 00325 //@} 00326 00327 /** 00328 \brief Smallest allowable gap (default 0.0) 00329 00330 If any of the gaps are below this value, then it is assumed 00331 that they are zero and the equation of state is simplified 00332 accordingly. If all of the gaps are less than gap_limit, then 00333 the results from nambujl_eos are used in 00334 calc_eq_temp_p(), calc_temp_p() and thd_potential(). 00335 */ 00336 double gap_limit; 00337 00338 /** \brief If this is true, then finite temperature 00339 corrections are ignored (default false) 00340 00341 This implements some simplifications in the 00342 momentum integration that are not possible at finite temperature. 00343 */ 00344 bool zerot; 00345 00346 /** \brief Use a fixed quark mass and ignore the quark 00347 condensates 00348 */ 00349 bool fixed_mass; 00350 00351 /// If true, then ensure color neutrality 00352 bool color_neut; 00353 00354 /** 00355 \brief Diquark coupling constant (default 3 G/4) 00356 00357 The default value is the one derived from a Fierz 00358 transformation. (\ref Buballa04) 00359 */ 00360 double GD; 00361 00362 /** 00363 \brief Faster eigenvalues 00364 00365 This assumes that the strange quark mass can be used 00366 as an effective chemical potential 00367 \f$ \mu_s^{\prime}=\mu_s-m_s^2/2/\mu_s\f$ 00368 and calculates 00369 the eigenvalues to first order in the four variables 00370 \f[ 00371 \alpha=\frac{\mu_d-\mu_u}{\mu_u} 00372 \f] 00373 \f[ 00374 \beta=\frac{\mu_s-\mu_u}{\mu_u} 00375 \f] 00376 \f[ 00377 \gamma=\frac{\Delta_{us}-\Delta_{ud}}{\Delta_{ud}} 00378 \f] 00379 \f[ 00380 \delta=\frac{\Delta_{ds}-\Delta_{ud}}{\Delta_{ud}} 00381 \f] 00382 00383 */ 00384 bool faster; 00385 00386 /// The quality of the most recently obtained solution 00387 double solver_qual; 00388 00389 /// If true, use the solver, otherwise, use the minimizers 00390 bool use_solver; 00391 00392 /// Return string denoting type ("cfl_njl_eos") 00393 virtual const char *type() { return "cfl_njl_eos"; }; 00394 00395 /** \brief The function used to take derivatives of the thermodynamic 00396 potential (used by calc_temp_p() ) 00397 00398 To use this function just to compute the thermodynamic potential, 00399 make sure \ref var is 0 (or \c NULL) and use 00400 \code 00401 tpot(0.0,0); 00402 \endcode 00403 since the arguments are both ignored. 00404 */ 00405 virtual double tpot(double var, void *&pa); 00406 00407 #ifndef DOXYGEN_INTERNAL 00408 00409 protected: 00410 00411 friend class io_tlate<cfl_njl_eos>; 00412 00413 /// Function specifying gap equations to solve 00414 int solve_fun_p(size_t nv, const ovector_view &x, ovector_view &y, 00415 void *&vp); 00416 00417 /// Outer minimization function for color neutrality 00418 int min_fun2(size_t nv, const ovector_view &x, double &y, 00419 void *&vp); 00420 00421 /// Inner minimization function for mass and gap equations 00422 int min_fun(size_t nv, const ovector_view &x, double &y, 00423 void *&vp); 00424 00425 /** \brief A function that takes the derivative of the thermodynamic 00426 potential with respect to the variable whose address is 'dvar'. 00427 */ 00428 double derivfn(double *dvar); 00429 00430 /// The integrand of the thermodynamic potential 00431 double tpot_integrand(double p, void *&pa); 00432 00433 /// The integrand of the entropy 00434 double entropy_integrand(double p, void *&pa); 00435 00436 /** 00437 \brief Treat the simply gapped quarks in all cases gracefully 00438 00439 This function uses the quarks \c q1 and \c q2 to construct the 00440 eigenvalues of the inverse propagator, properly handling the 00441 either zero or finite quark mass and either zero or finite 00442 quark gaps. In the case of finite quark mass and finite quark 00443 gaps, the quartic solver is used. 00444 00445 The chemical potentials are separated so we can add the 00446 color chemical potentials to the quark chemical potentials 00447 if necessary. 00448 00449 This function is used by eigenvalues(). It does not work for 00450 the "ur-dg-sb" set of quarks which are paired in a non-trivial 00451 way. 00452 00453 \todo Only the "ms" part of the quarks is referenced, so 00454 we should rewrite to use only double's as function arguments, 00455 and avoid passing pointers to quark objects. 00456 */ 00457 int gapped_eigenvalues(quark *q1, quark *q2, double mom, 00458 double mu1, double mu2, double delta, 00459 ovector_cx_view &ev); 00460 00461 /// Temperature 00462 double temper; 00463 00464 /// 3rd gluon chemical potential 00465 double smu3; 00466 00467 /// 8th gluon chemical potential 00468 double smu8; 00469 00470 /// \name Numerical methods 00471 //@{ 00472 /// The routine to solve quartics 00473 quartic_real_coeff *quartic; 00474 00475 /** 00476 \brief true if the quartic routine has been set 00477 00478 This is used for the I/O. 00479 */ 00480 bool quartic_set; 00481 00482 /** \brief The derivative routine 00483 */ 00484 gsl_deriv<void *,funct<void *> > de; 00485 00486 /// Pointer to inner minimizer 00487 multi_min<void *, multi_funct<void *> > *min1; 00488 00489 /// Pointer to outer minimizer 00490 multi_min<void *, multi_funct<void *> > *min2; 00491 //@} 00492 00493 /// Pointer to the variable to be differentiated by 00494 double *var; 00495 00496 /// \name For computing eigenvalues 00497 //@{ 00498 /// Inverse propagator matrix 00499 omatrix_cx m; 00500 00501 /// Inverse propagator matrix 00502 omatrix_cx m2; 00503 00504 /// Storage for eigenvalues 00505 ovector eval; 00506 00507 /// Storage for eigenvalues 00508 ovector eval2; 00509 00510 /// Workspace for eigenvalue computation 00511 gsl_eigen_herm_workspace *w; 00512 00513 /// Workspace for eigenvalue computation 00514 gsl_eigen_herm_workspace *w2; 00515 //@} 00516 00517 /// Switches for solving the gap equations 00518 //@{ 00519 bool ugap_solve; 00520 bool dgap_solve; 00521 bool sgap_solve; 00522 //@} 00523 00524 private: 00525 00526 cfl_njl_eos(const cfl_njl_eos &); 00527 cfl_njl_eos& operator=(const cfl_njl_eos&); 00528 00529 #endif 00530 00531 }; 00532 00533 template<> int io_tlate<cfl_njl_eos>::input 00534 (cinput *co, in_file_format *ins, cfl_njl_eos *nj); 00535 template<> int io_tlate<cfl_njl_eos>::output 00536 (coutput *co, out_file_format *outs, cfl_njl_eos *nj); 00537 template<> const char *io_tlate<cfl_njl_eos>::type(); 00538 00539 typedef io_tlate<cfl_njl_eos> cfl_njl_eos_io_type; 00540 00541 #ifndef DOXYGENP 00542 } 00543 #endif 00544 00545 #endif