star_bhns_spher.C

00001 /*
00002  *  Method of class Star_bhns to compute a spherical star configuration
00003  *
00004  *    (see file star_bhns.h for documentation).
00005  *
00006  */
00007 
00008 /*
00009  *   Copyright (c) 2005,2007 Keisuke Taniguchi
00010  *
00011  *   This file is part of LORENE.
00012  *
00013  *   LORENE is free software; you can redistribute it and/or modify
00014  *   it under the terms of the GNU General Public License version 2
00015  *   as published by the Free Software Foundation.
00016  *
00017  *   LORENE is distributed in the hope that it will be useful,
00018  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00019  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00020  *   GNU General Public License for more details.
00021  *
00022  *   You should have received a copy of the GNU General Public License
00023  *   along with LORENE; if not, write to the Free Software
00024  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00025  *
00026  */
00027 
00028 char star_bhns_spher_C[] = "$Header: /cvsroot/Lorene/C++/Source/Star_bhns/star_bhns_spher.C,v 1.2 2008/05/15 19:16:54 k_taniguchi Exp $" ;
00029 
00030 /*
00031  * $Id: star_bhns_spher.C,v 1.2 2008/05/15 19:16:54 k_taniguchi Exp $
00032  * $Log: star_bhns_spher.C,v $
00033  * Revision 1.2  2008/05/15 19:16:54  k_taniguchi
00034  * Change of some parameters.
00035  *
00036  * Revision 1.1  2007/06/22 01:32:19  k_taniguchi
00037  * *** empty log message ***
00038  *
00039  *
00040  * $Header: /cvsroot/Lorene/C++/Source/Star_bhns/star_bhns_spher.C,v 1.2 2008/05/15 19:16:54 k_taniguchi Exp $
00041  *
00042  */
00043 
00044 // C++ headers
00045 //#include <>
00046 
00047 // C headers
00048 #include <math.h>
00049 
00050 // Lorene headers
00051 #include "star_bhns.h"
00052 #include "param.h"
00053 #include "cmp.h"
00054 #include "tenseur.h"
00055 #include "unites.h"
00056 
00057 void Star_bhns::equil_spher_bhns(double ent_c, double precis) {
00058 
00059     // Fundamental constants and units
00060     // -------------------------------
00061     using namespace Unites ;
00062 
00063     // Initializations
00064     // ---------------
00065 
00066     const Mg3d* mg = mp.get_mg() ;
00067     int nz = mg->get_nzone() ;
00068 
00069     // Index of the point at phi=0, theta=pi/2 at the surface of the star:
00070     int l_b = nzet - 1 ;
00071     int i_b = mg->get_nr(l_b) - 1 ;
00072     int j_b = mg->get_nt(l_b) - 1 ;
00073     int k_b = 0 ;
00074 
00075     // Value of the enthalpy defining the surface of the star
00076     double ent_b = 0. ;
00077 
00078     // Initialization of the enthalpy field to the constant value ent_c :
00079     ent = ent_c ;
00080     ent.annule(nzet, nz-1) ;
00081 
00082     // Corresponding profiles of baryon density, energy density and pressure
00083     equation_of_state() ;
00084 
00085     // Initial metric
00086     lapconf_auto = 1. ;
00087     lapconf_auto.std_spectral_base() ;
00088     confo_auto = 1. ;
00089     confo_auto.std_spectral_base() ;
00090 
00091     // Auxiliary quantities
00092     // --------------------
00093 
00094     // Affine mapping for solving the Poisson equations
00095     Map_af mpaff(mp) ;
00096 
00097     Param par_nul ;  // Param (null) for Map_af::poisson
00098 
00099     Scalar ent_jm1(mp) ;  // Enthalpy at previous step
00100     ent_jm1 = ent_c ;
00101     ent_jm1.std_spectral_base() ;
00102     ent_jm1.annule(nzet, nz-1) ;
00103 
00104     Scalar source_lapconf(mp) ;
00105     Scalar source_confo(mp) ;
00106     Scalar lapconf_auto_m1(mp) ;
00107     Scalar confo_auto_m1(mp) ;
00108     lapconf_auto_m1 = 0. ;
00109     confo_auto_m1 = 0. ;
00110 
00111     double diff_ent = 1. ;
00112     int mermax = 200 ;  // Maximum number of iterations
00113 
00114     double alpha_r = 1. ;
00115 
00116     //==========================================================//
00117     //                    Start of iteration                    //
00118     //==========================================================//
00119 
00120     for (int mer=0 ; (diff_ent > precis) && (mer<mermax) ; mer++ ) {
00121 
00122       cout << "-----------------------------------------------" << endl ;
00123       cout << "step: " << mer << endl ;
00124       cout << "alpha_r: " << alpha_r << endl ;
00125       cout << "diff_ent = " << diff_ent << endl ;
00126 
00127       //----------------------------------------------------
00128       // Resolution of Poisson equation for lapconf function
00129       //----------------------------------------------------
00130 
00131       // Matter part of lapconf
00132       // ----------------------
00133       source_lapconf = qpig * lapconf_auto * pow(confo_auto,4.)
00134     * (0.5*ener + 3.*press) ;
00135 
00136       source_lapconf.inc_dzpuis(4-source_lapconf.get_dzpuis()) ;
00137       source_lapconf.std_spectral_base() ;
00138 
00139       Cmp sou_lapconf(source_lapconf) ;
00140       Cmp lapconf_cmp(lapconf_auto_m1) ;
00141       lapconf_cmp.set_etat_qcq() ;
00142 
00143       mpaff.poisson(sou_lapconf, par_nul, lapconf_cmp) ;
00144 
00145       // Re-construction of a scalar
00146       lapconf_auto_m1 = lapconf_cmp ;
00147 
00148       //-------------------------------------
00149       // Computation of the new radial scale
00150       //-------------------------------------
00151 
00152       double exp_ent_c = exp(ent_c) ;
00153       double exp_ent_b = exp(ent_b) ;
00154 
00155       double lap_auto_c = lapconf_auto_m1.val_grid_point(0,0,0,0) ;
00156       double lap_auto_b = lapconf_auto_m1.val_grid_point(l_b,k_b,j_b,i_b) ;
00157 
00158       double confo_c = confo_auto.val_grid_point(0,0,0,0) ;
00159       double confo_b = confo_auto.val_grid_point(l_b,k_b,j_b,i_b) ;
00160 
00161       double alpha_r2 = (exp_ent_b*confo_c - exp_ent_c*confo_b)
00162     / ( exp_ent_c*confo_b*lap_auto_c
00163         - exp_ent_b*confo_c*lap_auto_b )  ;
00164 
00165       alpha_r = sqrt(alpha_r2) ;
00166 
00167       // New radial scale
00168       mpaff.homothetie( alpha_r ) ;
00169 
00170       //----------------
00171       // First integral
00172       //----------------
00173 
00174       // Lapconf function
00175       lapconf_auto_m1 = alpha_r2 * lapconf_auto_m1 ;
00176       lapconf_auto = lapconf_auto_m1 + 1. ;
00177 
00178       // Enthalpy in all space
00179       double lapconfo_c = lapconf_auto.val_grid_point(0,0,0,0) ;
00180       confo_c = confo_auto.val_grid_point(0,0,0,0) ;
00181       ent = ent_c + log(lapconfo_c) - log(confo_c)
00182     - log(lapconf_auto) + log(confo_auto) ;
00183       ent.std_spectral_base() ;
00184 
00185       //-------------------
00186       // Equation of state
00187       //-------------------
00188 
00189       equation_of_state() ;
00190 
00191       //-----------------------------------------------------
00192       // Resolution of Poisson equation for conformal factor
00193       //-----------------------------------------------------
00194 
00195       source_confo = - 0.5 * qpig * pow(confo_auto,5.) * ener ;
00196       source_confo.inc_dzpuis(4-source_confo.get_dzpuis()) ;
00197       source_confo.std_spectral_base() ;
00198 
00199       Cmp sou_confo(source_confo) ;
00200       Cmp cnf_auto_cmp(confo_auto_m1) ;
00201       cnf_auto_cmp.set_etat_qcq() ;
00202 
00203       mpaff.poisson(sou_confo, par_nul, cnf_auto_cmp) ;
00204 
00205       // Re-construction of a scalr
00206       confo_auto_m1 = cnf_auto_cmp ;
00207 
00208       confo_auto = confo_auto_m1 + 1. ;
00209 
00210       // Relative difference with enthalphy at the previous step
00211       // -------------------------------------------------------
00212 
00213       diff_ent = norme( diffrel(ent, ent_jm1) ) / nzet ;
00214 
00215       // Next step
00216       // ---------
00217 
00218       ent_jm1 = ent ;
00219 
00220     } // End of iteration loop
00221 
00222     //========================================================//
00223     //                    End of iteration                    //
00224     //========================================================//
00225 
00226     // The mapping is transfered to that of the star
00227     // ---------------------------------------------
00228     mp = mpaff ;
00229 
00230     // Sets values
00231     // -----------
00232 
00233     // ... hydro
00234     ent.annule(nzet, nz-1) ;
00235 
00236     ener_euler = ener ;
00237     s_euler = 3. * press ;
00238     gam_euler = 1. ;
00239     for(int i=1; i<=3; i++)
00240       u_euler.set(i) = 0 ;
00241 
00242     // ... metric
00243     lapconf_tot = lapconf_auto ;
00244     lapse_auto = lapconf_auto / confo_auto ;
00245     lapse_tot = lapse_auto ;
00246     confo_tot = confo_auto ;
00247     psi4 = pow(confo_auto, 4.) ;
00248     for (int i=1; i<=3; i++)
00249       shift_auto.set(i) = 0. ;
00250 
00251     // Info printing
00252     // -------------
00253 
00254     cout << endl
00255      << "Characteristics of the star obtained by Star_bhns::equil_spher_bhns : "
00256      << endl
00257      << "-------------------------------------------------------------------   "
00258      << endl ;
00259 
00260     cout.precision(16) ;
00261     double ray = mp.val_r(l_b, 1., M_PI/2., 0) ;
00262     cout << "Coordinate radius :               "
00263      << ray / km << " [km]" << endl ;
00264 
00265     double rcirc = ray * sqrt(psi4.val_grid_point(l_b, k_b, j_b, i_b) ) ;
00266     double compact = qpig/(4.*M_PI) * mass_g_bhns() / rcirc ;
00267 
00268     cout << "Circumferential radius R :        "
00269      << rcirc/km  << " [km]" << endl ;
00270     cout << "Baryon mass :                     "
00271      << mass_b_bhns(0,0.,1.)/msol << " [Mo]" << endl ;
00272     cout << "Gravitational mass M :            "
00273      << mass_g_bhns()/msol << " [Mo]" << endl ;
00274     cout << "Compaction parameter GM/(c^2 R) : " << compact << endl ;
00275 
00276     //----------------
00277     // Virial theorem
00278     //----------------
00279 
00280     //... Pressure term
00281     Scalar source(mp) ;
00282     source = qpig * pow(confo_auto,6.) * s_euler ;
00283     source.std_spectral_base() ;
00284     double vir_mat = source.integrale() ;
00285 
00286     //... Gravitational term
00287     Scalar tmp1(mp) ;
00288     tmp1 = confo_auto.dsdr() ;
00289     tmp1.std_spectral_base() ;
00290 
00291     Scalar tmp2(mp) ;
00292     tmp2 = confo_auto * lapconf_auto.dsdr() / lapconf_auto  - tmp1 ;
00293     tmp2.std_spectral_base() ;
00294 
00295     source = 2. * tmp1 * tmp1 - tmp2 * tmp2 ;
00296 
00297     /*
00298     Scalar tmp1(mp) ;
00299     tmp1 = log(lapse_auto) ;
00300     tmp1.std_spectral_base() ;
00301 
00302     Scalar tmp2(mp) ;
00303     tmp2 = log(confo_auto) ;
00304     tmp2.std_spectral_base() ;
00305 
00306     source = confo_auto * confo_auto
00307       * ( 2. * tmp2.dsdr() * tmp2.dsdr() - tmp1.dsdr() * tmp1.dsdr() ) ;
00308     */
00309     source.std_spectral_base() ;        
00310     double vir_grav = source.integrale() ;
00311 
00312     //... Relative error on the virial identity GRV3
00313     double grv3 = ( vir_mat + vir_grav ) / vir_mat ;
00314 
00315     cout << "Virial theorem GRV3 : " << endl ;
00316     cout << "     3P term :        " << vir_mat << endl ;
00317     cout << "     grav. term :     " << vir_grav << endl ;
00318     cout << "     relative error : " << grv3 << endl ;
00319 
00320 }

Generated on Tue Feb 7 01:35:19 2012 for LORENE by  doxygen 1.4.6