et_equil_spher_regu.C

00001 /*
00002  * Method of class Etoile to compute a static spherical configuration
00003  *  by regularizing source.
00004  *
00005  * (see file etoile.h for documentation).
00006  *
00007  */
00008 
00009 /*
00010  *   Copyright (c) 2000-2001 Eric Gourgoulhon
00011  *   Copyright (c) 2000-2001 Keisuke Taniguchi
00012  *
00013  *   This file is part of LORENE.
00014  *
00015  *   LORENE is free software; you can redistribute it and/or modify
00016  *   it under the terms of the GNU General Public License as published by
00017  *   the Free Software Foundation; either version 2 of the License, or
00018  *   (at your option) any later version.
00019  *
00020  *   LORENE is distributed in the hope that it will be useful,
00021  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00022  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00023  *   GNU General Public License for more details.
00024  *
00025  *   You should have received a copy of the GNU General Public License
00026  *   along with LORENE; if not, write to the Free Software
00027  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00028  *
00029  */
00030 
00031 
00032 char et_equil_spher_regu_C[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/et_equil_spher_regu.C,v 1.2 2004/03/25 10:29:04 j_novak Exp $" ;
00033 
00034 /*
00035  * $Id: et_equil_spher_regu.C,v 1.2 2004/03/25 10:29:04 j_novak Exp $
00036  * $Log: et_equil_spher_regu.C,v $
00037  * Revision 1.2  2004/03/25 10:29:04  j_novak
00038  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
00039  *
00040  * Revision 1.1.1.1  2001/11/20 15:19:28  e_gourgoulhon
00041  * LORENE
00042  *
00043  * Revision 2.14  2001/03/06  16:34:04  keisuke
00044  * Change the regularization degree k_div=1 to the arbitrary one.
00045  *
00046  * Revision 2.13  2001/02/07  09:46:11  eric
00047  * unsgam1 est desormais donne par Eos::der_nbar_ent (cas newtonien)
00048  *   ou Eos::der_ener_ent (cas relativiste).
00049  *
00050  * Revision 2.12  2000/09/26  15:42:35  keisuke
00051  * Correction erreur: la triade de duu_div doit etre celle de mp et non celle
00052  *  de l'objet temporaire mpaff.
00053  *
00054  * Revision 2.11  2000/09/26  06:56:04  keisuke
00055  * Suppress "int" from the declaration of k_div.
00056  *
00057  * Revision 2.10  2000/09/22  15:53:36  keisuke
00058  * d_logn_auto_div est desormais un membre de la classe Etoile.
00059  *
00060  * Revision 2.9  2000/09/08  12:23:17  keisuke
00061  * Correct a typological error in the equation.
00062  *
00063  * Revision 2.8  2000/09/07  15:37:23  keisuke
00064  * Add a new argument in poisson_regular
00065  *  and suppress logn_auto_total.
00066  *
00067  * Revision 2.7  2000/09/06  12:47:35  keisuke
00068  * Suppress #include "graphique.h".
00069  *
00070  * Revision 2.6  2000/09/06  12:39:59  keisuke
00071  * Save the map in the every iterative step after operating "homothetie".
00072  *
00073  * Revision 2.5  2000/09/04  16:15:05  keisuke
00074  * Change the argument of Map_af::poisson_regular.
00075  *
00076  * Revision 2.4  2000/08/31  16:05:26  keisuke
00077  * Modify the arguments of Map_af::poisson_regular.
00078  *
00079  * Revision 2.3  2000/08/29  11:38:25  eric
00080  * Ajout des membres k_div et logn_auto_div a la classe Etoile.
00081  *
00082  * Revision 2.2  2000/08/28  16:10:39  keisuke
00083  * Add "nzet" in the argumant of poisson_regular.
00084  *
00085  * Revision 2.1  2000/08/25  14:58:15  keisuke
00086  * Modif (Virial theorem).
00087  *
00088  * Revision 2.0  2000/08/25  09:01:33  keisuke
00089  * *** empty log message ***
00090  *
00091  *
00092  * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_equil_spher_regu.C,v 1.2 2004/03/25 10:29:04 j_novak Exp $
00093  *
00094  */
00095 
00096 // Headers C
00097 #include <math.h>
00098 
00099 // Headers Lorene
00100 #include "etoile.h"
00101 #include "eos.h"
00102 #include "param.h"
00103 #include "unites.h"
00104 
00105 //********************************************************************
00106 
00107 void Etoile::equil_spher_regular(double ent_c, double precis){
00108 
00109     // Fundamental constants and units
00110     // -------------------------------
00111   using namespace Unites ;
00112 
00113     // Initializations
00114     // ---------------
00115 
00116     //    k_div = 1 ;       // Regularization index
00117 
00118     cout << "Input the regularization degree (k_div) : " ;
00119     cin >> k_div ;          // Regularization index
00120 
00121     const Mg3d* mg = mp.get_mg() ;
00122     int nz = mg->get_nzone() ;      // total number of domains
00123 
00124     // Index of the point at phi=0, theta=pi/2 at the surface of the star:
00125     int l_b = nzet - 1 ;
00126     int i_b = mg->get_nr(l_b) - 1 ;
00127     int j_b = mg->get_nt(l_b) - 1 ;
00128     int k_b = 0 ;
00129 
00130     // Value of the enthalpy defining the surface of the star
00131     double ent_b = 0 ;
00132 
00133     // Initialization of the enthalpy field to the constant value ent_c :
00134 
00135     ent = ent_c ;
00136     ent.annule(nzet, nz-1) ;
00137 
00138     // Corresponding profiles of baryon density, energy density and pressure
00139 
00140     equation_of_state() ;
00141 
00142     // Initial metric 
00143     a_car = 1 ;      // this value will remain unchanged in the Newtonian case
00144     beta_auto = 0 ;  // this value will remain unchanged in the Newtonian case
00145 
00146     // Auxiliary quantities
00147     // --------------------
00148 
00149     // Affine mapping for solving the Poisson equations
00150     Map_af mpaff(mp) ;  
00151 
00152     Param par_nul   ;    // Param (null) for Map_af::poisson.
00153 
00154     Tenseur ent_jm1(mp) ;   // Enthalpy at previous step
00155     ent_jm1 = 0 ;
00156 
00157     Tenseur source(mp) ;
00158     Tenseur logn(mp) ;
00159     Tenseur logn_quad(mp) ;
00160     logn = 0 ;
00161     logn_quad = 0 ;
00162 
00163     Tenseur dlogn_auto_regu(mp, 1, COV, mp.get_bvect_spher()) ;
00164 
00165     Cmp source_regu(mp) ;
00166     Cmp source_div(mp) ;
00167 
00168     Cmp dlogn(mp) ;
00169     dlogn = 0 ;
00170     Cmp dbeta(mp) ;
00171 
00172     Tenseur dlogn_auto(mp, 1, COV, mp.get_bvect_spher()) ;
00173     Tenseur dlogn_quad(mp) ;
00174     dlogn_quad = 0 ;
00175 
00176     double diff_ent = 1 ;     
00177     int mermax = 200 ;      // Max number of iterations
00178 
00179     double alpha_r = 1 ;
00180 
00181     //=========================================================================
00182     //          Start of iteration
00183     //=========================================================================
00184 
00185     for(int mer=0 ; (diff_ent > precis) && (mer<mermax) ; mer++ ) {
00186 
00187       cout << "-----------------------------------------------" << endl ;
00188       cout << "step: " << mer << endl ;
00189       cout << "alpha_r: " << alpha_r << endl ;
00190       cout << "diff_ent = " << diff_ent << endl ;    
00191 
00192       //-----------------------------------------------------
00193       // Resolution of Poisson equation for ln(N)
00194       //-----------------------------------------------------
00195 
00196       // Matter part of ln(N)
00197       // --------------------
00198 
00199       double unsgam1 ;    // effective power of H in the source 
00200               // close to the surface
00201 
00202       if (relativistic) {
00203     
00204     source = a_car * (ener + 3*press) ;
00205 
00206     // 1/(gam-1) = dln(e)/dln(H) close to the surface : 
00207     unsgam1 = eos.der_ener_ent_p(ent_b + 1e-10*(ent_c-ent_b)) ; 
00208       }
00209       else {
00210     
00211     source = nbar ;
00212     
00213     // 1/(gam-1) = dln(n)/dln(H) close to the surface : 
00214     unsgam1 = eos.der_nbar_ent_p(ent_b + 1e-10*(ent_c-ent_b)) ; 
00215       }
00216 
00217       (source.set()).set_dzpuis(4) ; 
00218 
00219       source.set_std_base() ;       // Sets the standard spectral bases. 
00220 
00221       logn_auto.set_etat_qcq() ;
00222       logn_auto_regu.set_etat_qcq() ;
00223       logn_auto_div.set_etat_qcq() ;
00224 
00225       source_regu.std_base_scal() ;
00226       source_div.std_base_scal() ;
00227 
00228       mpaff.poisson_regular(source(), k_div, nzet, unsgam1, par_nul,
00229                 logn_auto.set(), logn_auto_regu.set(),
00230                 logn_auto_div.set(),
00231                 d_logn_auto_div, source_regu, source_div) ;
00232 
00233       dlogn_auto_regu = logn_auto_regu.gradient_spher() ;
00234 
00235       dlogn_auto = dlogn_auto_regu + d_logn_auto_div ;
00236 
00237       // NB: at this stage logn_auto is in machine units, not in c^2
00238 
00239       // Quadratic part of ln(N)
00240       // -----------------------
00241 
00242       if (relativistic) {
00243 
00244     mpaff.dsdr(beta_auto(), dbeta) ;
00245 
00246     source = - dlogn * dbeta ;
00247 
00248     logn_quad.set_etat_qcq() ;
00249 
00250     mpaff.poisson(source(), par_nul, logn_quad.set()) ; 
00251 
00252     dlogn_quad.set_etat_qcq() ;
00253 
00254     mpaff.dsdr(logn_quad(), dlogn_quad.set()) ;
00255 
00256       }
00257 
00258       //-----------------------------------------------------
00259       // Computation of the new radial scale
00260       //-----------------------------------------------------
00261 
00262       // alpha_r (r = alpha_r r') is determined so that the enthalpy
00263       // takes the requested value ent_b at the stellar surface
00264 
00265       double nu_mat0_b  = logn_auto()(l_b, k_b, j_b, i_b) ; 
00266       double nu_mat0_c  = logn_auto()(0, 0, 0, 0) ; 
00267 
00268       double nu_quad0_b  = logn_quad()(l_b, k_b, j_b, i_b) ; 
00269       double nu_quad0_c  = logn_quad()(0, 0, 0, 0) ; 
00270 
00271       double alpha_r2 = ( ent_c - ent_b - nu_quad0_b + nu_quad0_c )
00272     / ( qpig*(nu_mat0_b - nu_mat0_c) ) ;
00273 
00274       alpha_r = sqrt(alpha_r2) ;
00275 
00276       // New radial scale
00277       // -----------------
00278       mpaff.homothetie( alpha_r ) ;
00279 
00280       // The mapping is transfered to that of the star:
00281       // ----------------------------------------------
00282       mp = mpaff ;
00283 
00284       d_logn_auto_div.set_triad( mp.get_bvect_spher() ) ;  // Absolutely necessary !!!
00285 
00286       //--------------------
00287       // First integral
00288       //--------------------
00289 
00290       // Gravitation potential in units c^2 :
00291       logn_auto = alpha_r2*qpig * logn_auto ;
00292       logn = logn_auto + logn_quad ;
00293 
00294       // Enthalpy in all space
00295       double logn_c = logn()(0, 0, 0, 0) ;
00296       ent = ent_c - logn() + logn_c ;
00297 
00298       //---------------------
00299       // Equation of state
00300       //---------------------
00301 
00302       equation_of_state() ;
00303 
00304       // derivative of gravitation potential in units c^2 :
00305       dlogn_auto = alpha_r*qpig * dlogn_auto ;
00306       dlogn = dlogn_auto(0) + dlogn_quad() ;
00307 
00308       if (relativistic) {
00309 
00310     //----------------------------
00311     // Equation for beta = ln(AN)
00312     //----------------------------
00313 
00314     mpaff.dsdr(beta_auto(), dbeta) ;
00315 
00316     source = 3 * qpig * a_car * press ;
00317 
00318     source = source() 
00319       - 0.5 * (  dlogn * dlogn + dbeta * dbeta ) ;
00320 
00321     source.set_std_base() ;    // Sets the standard spectral bases.
00322 
00323     beta_auto.set_etat_qcq() ; 
00324 
00325     mpaff.poisson(source(), par_nul, beta_auto.set()) ;
00326 
00327     // Metric coefficient A^2 update
00328 
00329     a_car = exp(2*(beta_auto - logn)) ;
00330 
00331       }
00332 
00333       // Relative difference with enthalpy at the previous step
00334       // ------------------------------------------------------
00335 
00336       diff_ent = norme( diffrel(ent(), ent_jm1()) ) / nzet ; 
00337 
00338       // Next step
00339       // ---------
00340 
00341       ent_jm1 = ent ;
00342 
00343 
00344     }  // End of iteration loop 
00345 
00346     //=========================================================================
00347     //          End of iteration
00348     //=========================================================================
00349 
00350 
00351     // Sets value to all the Tenseur's of the star
00352     // -------------------------------------------
00353 
00354     // ... hydro 
00355     ent.annule(nzet, nz-1) ;    // enthalpy set to zero at the exterior of 
00356                 // the star
00357     ener_euler = ener ;
00358     s_euler = 3 * press ;
00359     gam_euler = 1 ;
00360     u_euler = 0 ;
00361 
00362     // ... metric
00363     nnn = exp( unsurc2 * logn ) ;
00364     shift = 0 ;
00365 
00366     // Info printing
00367     // -------------
00368 
00369     cout << endl
00370      << "Characteristics of the star obtained by Etoile::equil_spher_regular : "
00371      << endl
00372      << "-----------------------------------------------------------------" 
00373      << endl ;
00374 
00375     double ray = mp.val_r(l_b, 1., M_PI/2., 0) ;
00376     cout << "Coordinate radius  :       " << ray / km << " km" << endl ;
00377 
00378     double rcirc = ray * sqrt( a_car()(l_b, k_b, j_b, i_b) ) ;
00379 
00380     double compact = qpig/(4.*M_PI) * mass_g() / rcirc ;
00381 
00382     cout << "Circumferential radius R : " << rcirc/km  << " km" << endl ;
00383     cout << "Baryon mass :       " << mass_b()/msol << " Mo" << endl ;
00384     cout << "Gravitational mass M :  " << mass_g()/msol << " Mo" << endl ;
00385     cout << "Compacity parameter GM/(c^2 R) : " << compact << endl ;
00386 
00387 
00388     //-----------------
00389     // Virial theorem
00390     //-----------------
00391 
00392     //... Pressure term
00393 
00394     source = qpig * a_car * sqrt(a_car) * s_euler ;
00395     source.set_std_base() ;
00396     double vir_mat = source().integrale() ;
00397 
00398     //... Gravitational term
00399 
00400     Cmp tmp = beta_auto().dsdr() - dlogn ;
00401 
00402     source =  - ( dlogn * dlogn - 0.5 * tmp * tmp ) * sqrt(a_car())  ;
00403 
00404     source.set_std_base() ;
00405     double vir_grav = source().integrale() ;
00406 
00407     //... Relative error on the virial identity GRV3
00408 
00409     double grv3 = ( vir_mat + vir_grav ) / vir_mat ;
00410 
00411     cout << "Virial theorem GRV3 : " << endl ;
00412     cout << "     3P term    : " << vir_mat << endl ;
00413     cout << "     grav. term : " << vir_grav << endl ;
00414     cout << "     relative error : " << grv3 << endl ;
00415 
00416 
00417 }

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