Equation of State Sub-Library: Version 0.910
cfl_njl_eos.h
00001 /*
00002   -------------------------------------------------------------------
00003   
00004   Copyright (C) 2006-2012, 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   /** \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       \future This class internally mixes ovector, omatrix, gsl_vector
00197       and gsl_matrix objects in a confusing and non-optimal way. Fix this. 
00198       \future 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 
00202       \hline
00203       <b>References:</b>
00204 
00205       Created for \ref Steiner02.
00206   */
00207   class cfl_njl_eos : public nambujl_eos {
00208   public:
00209     
00210     cfl_njl_eos();
00211     
00212     virtual ~cfl_njl_eos();
00213     
00214     /** \brief Set the parameters and the bag constant 'B0'
00215         
00216         This function allows the user to specify the momentum cutoff,
00217         \c lambda, the four-fermion coupling \c fourferm, the
00218         six-fermion coupling from the 't Hooft interaction \c sixferm,
00219         and the color-superconducting coupling, \c fourgap. If 0.0 is
00220         given for any of the values, then the default is used (\f$
00221         \Lambda=602.3/(\hbar c), G=1.835/\Lambda^2, K=12.36/\Lambda^5
00222         \f$).
00223         
00224         If the four-fermion coupling that produces a gap is not
00225         specified, it is automatically set to 3/4 G, which is the
00226         value obtained from the Fierz transformation.
00227         
00228         The value of the shift in the bag constant nambujl_eos::B0 is
00229         automatically calculated to ensure that the vacuum has zero
00230         energy density and zero pressure. The functions set_quarks()
00231         and set_thermo() must be used before hand to specify the \ref
00232         quark and \ref thermo objects.
00233         
00234     */
00235     virtual int set_parameters(double lambda=0.0, double fourferm=0.0, 
00236                                double sixferm=0.0, double fourgap=0.0);
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, 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     /** \brief Calculate the energy eigenvalues as a function of the
00268         momentum
00269         
00270         Given the momentum \c mom, and the chemical potentials
00271         associated with the third and eighth gluons (\c mu3 and \c
00272         mu8), the energy eigenvalues are computed in egv[0]
00273         ... egv[35].
00274     */
00275     virtual int eigenvalues(double lmom, double mu3, double mu8, 
00276                             double egv[36], double dedmuu[36], 
00277                             double dedmud[36], double dedmus[36], 
00278                             double dedmu[36], double dedmd[36], 
00279                             double dedms[36], double dedu[36], 
00280                             double dedd[36], double deds[36], 
00281                             double dedmu3[36], double dedmu8[36]);
00282     
00283     /// Set the routine for solving quartics
00284     int set_quartic(quartic_real_coeff &q) { 
00285       quartic=&q; 
00286       return 0;
00287     }
00288     
00289     /// The equal mass threshold
00290     double eq_limit;
00291     
00292     /// Set to true to test the integration (default false)
00293     bool integ_test;
00294     
00295     /// Test the integration routines
00296     int test_integration(test_mgr &t);
00297 
00298     //
00299     //gsl_quadratic_real_coeff quad;
00300 
00301     /** \brief The default quartic routine
00302         
00303         Slightly better accuracy (with slower execution times) can be
00304         achieved using \ref gsl_poly_real_coeff which polishes the
00305         roots of the quartics. For example
00306 
00307         \code
00308         cfl_njl_eos cfl;
00309         gsl_poly_real_coeff gp;
00310         cfl.set_quartic(gp);
00311         \endcode
00312      */
00313     cern_quartic_real_coeff def_quartic;
00314     //@}
00315     
00316     /** \brief Test the routine to compute the eigenvalues of 
00317         non-superfluid fermions
00318     */
00319     int test_normal_eigenvalues(test_mgr &t);
00320 
00321     /** \brief Test the routine to compute the eigenvalues of 
00322         superfluid fermions
00323     */
00324     int test_gapped_eigenvalues(test_mgr &t);
00325     
00326     /** \brief Smallest allowable gap (default 0.0)
00327         
00328         If any of the gaps are below this value, then it is assumed
00329         that they are zero and the equation of state is simplified
00330         accordingly. If all of the gaps are less than gap_limit, then
00331         the results from nambujl_eos are used in
00332         calc_eq_temp_p(), calc_temp_p() and thd_potential().
00333     */
00334     double gap_limit;
00335     
00336     /** \brief If this is true, then finite temperature
00337         corrections are ignored (default false)
00338         
00339         This implements some simplifications in the 
00340         momentum integration that are not possible at finite temperature.
00341     */
00342     bool zerot;
00343     
00344     /** \brief Use a fixed quark mass and ignore the quark
00345         condensates
00346     */
00347     bool fixed_mass;
00348 
00349     /// If true, then ensure color neutrality
00350     bool color_neut;
00351     
00352     /** \brief Diquark coupling constant (default 3 G/4)
00353         
00354         The default value is the one derived from a Fierz
00355         transformation. (\ref Buballa04)
00356     */
00357     double GD;
00358 
00359     /** \brief The absolute precision for the integration
00360         (default \f$ 10^{-4} \f$ )
00361         
00362         This is analogous to gsl_inte::epsabs
00363     */
00364     double inte_epsabs;
00365 
00366     /** \brief The relative precision for the integration 
00367         (default \f$ 10^{-4} \f$ )
00368 
00369         This is analogous to gsl_inte::epsrel
00370     */
00371     double inte_epsrel;
00372     
00373     /** \brief The number of points used in the last integration
00374         (default 0)
00375 
00376         This returns 21, 43, or 87 depending on the number of function
00377         evaluations needed to obtain the desired precision. If it
00378         the routine failes to obtain the desired precision, then
00379         this variable is set to 88.
00380     */
00381     size_t inte_npoints;
00382 
00383     /// Return string denoting type ("cfl_njl_eos")
00384     virtual const char *type() { return "cfl_njl_eos"; };
00385     
00386 #ifndef DOXYGEN_INTERNAL
00387     
00388   protected:
00389     
00390     /** \brief The integrands
00391         
00392         - res[0] is the thermodynamic potential, \f$ \Omega \f$
00393         - res[1] is \f$ d -\Omega / d T \f$
00394         - res[2] is \f$ d \Omega / d \mu_u \f$
00395         - res[3] is \f$ d \Omega / d \mu_d \f$
00396         - res[4] is \f$ d \Omega / d \mu_s \f$
00397         - res[5] is \f$ d \Omega / d m_u \f$
00398         - res[6] is \f$ d \Omega / d m_d \f$
00399         - res[7] is \f$ d \Omega / d m_s \f$
00400         - res[8] is \f$ d \Omega / d \Delta_{ds} \f$
00401         - res[9] is \f$ d \Omega / d \Delta_{us} \f$
00402         - res[10] is \f$ d \Omega / d \Delta_{ud} \f$
00403         - res[11] is \f$ d \Omega / d \mu_3 \f$
00404         - res[12] is \f$ d \Omega / d \mu_8 \f$
00405     */
00406     virtual int integrands(double p, double res[]);
00407     
00408     /// Compute ungapped eigenvalues and the appropriate derivatives
00409     int normal_eigenvalues(double m, double lmom, double mu, 
00410                            double lam[2], double dldmu[2], 
00411                            double dldm[2]);
00412 
00413     /** \brief Treat the simply gapped quarks in all cases gracefully
00414 
00415         This function uses the quarks \c q1 and \c q2 to construct the
00416         eigenvalues of the inverse propagator, properly handling the
00417         either zero or finite quark mass and either zero or finite
00418         quark gaps. In the case of finite quark mass and finite quark
00419         gaps, the quartic solver is used.
00420         
00421         The chemical potentials are separated so we can add the 
00422         color chemical potentials to the quark chemical potentials
00423         if necessary.
00424 
00425         This function is used by eigenvalues(). It does not work for
00426         the "ur-dg-sb" set of quarks which are paired in a non-trivial
00427         way.
00428 
00429         \todo In the code, the equal mass case seems to be commented
00430         out. Why?
00431     */
00432     int gapped_eigenvalues(double m1, double m2, double lmom,
00433                            double mu1, double mu2, double tdelta,
00434                            double lam[4], double dldmu1[4], 
00435                            double dldmu2[4], double dldm1[4],
00436                            double dldm2[4], double dldg[4]);    
00437     
00438     /// Temperature
00439     double temper;
00440     
00441     /// 3rd gluon chemical potential
00442     double smu3;
00443     
00444     /// 8th gluon chemical potential
00445     double smu8;
00446 
00447     /// \name Numerical methods
00448     //@{
00449     /// The routine to solve quartics
00450     quartic_real_coeff *quartic;
00451     
00452     /// \name For computing eigenvalues
00453     //@{
00454     /// Inverse propagator matrix
00455     omatrix_cx iprop;
00456     
00457     /// The eigenvectors
00458     omatrix_cx eivec;
00459 
00460     /// The derivative of the inverse propagator wrt the ds gap
00461     omatrix_cx dipdgapu;
00462     /// The derivative of the inverse propagator wrt the us gap
00463     omatrix_cx dipdgapd;
00464     /// The derivative of the inverse propagator wrt the ud gap
00465     omatrix_cx dipdgaps;
00466 
00467     /// The eigenvalues
00468     ovector eval;
00469 
00470     /// Workspace for eigenvalue computation
00471     gsl_eigen_hermv_workspace *w;
00472     //@}
00473 
00474     /// \name For the integration
00475     //@{
00476     /// The error scaling function for integ_err
00477     double rescale_error(double err, double result_abs, 
00478                          double result_asc);
00479     
00480     /** \brief A new version of gsl_inte_qng to integrate several
00481         functions at the same time
00482     */
00483     int integ_err(double a, double b, const size_t nr,
00484                   ovector &res, double &err2);
00485     //@}
00486 
00487   private:
00488     
00489     cfl_njl_eos(const cfl_njl_eos &);
00490     cfl_njl_eos& operator=(const cfl_njl_eos&);
00491     
00492 #endif
00493     
00494   };
00495   
00496 #ifndef DOXYGENP
00497 }
00498 #endif
00499 
00500 #endif
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Friends

Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).

Get Object-oriented Scientific Computing
Lib at SourceForge.net. Fast, secure and Free Open Source software
downloads.