etoile_rot.C

00001 /*
00002  * Methods for the class Etoile_rot
00003  *
00004  * (see file etoile.h for documentation)
00005  */
00006 
00007 /*
00008  *   Copyright (c) 2000-2001 Eric Gourgoulhon
00009  *
00010  *   This file is part of LORENE.
00011  *
00012  *   LORENE is free software; you can redistribute it and/or modify
00013  *   it under the terms of the GNU General Public License as published by
00014  *   the Free Software Foundation; either version 2 of the License, or
00015  *   (at your option) any later version.
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 
00029 char etoile_rot_C[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/etoile_rot.C,v 1.4 2004/03/25 10:29:07 j_novak Exp $" ;
00030 
00031 /*
00032  * $Id: etoile_rot.C,v 1.4 2004/03/25 10:29:07 j_novak Exp $
00033  * $Log: etoile_rot.C,v $
00034  * Revision 1.4  2004/03/25 10:29:07  j_novak
00035  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
00036  *
00037  * Revision 1.3  2001/12/06 15:11:43  jl_zdunik
00038  * Introduction of the new function f_eq() in the class Etoile_rot
00039  *
00040  * Revision 1.2  2001/12/04 21:27:53  e_gourgoulhon
00041  *
00042  * All writing/reading to a binary file are now performed according to
00043  * the big endian convention, whatever the system is big endian or
00044  * small endian, thanks to the functions fwrite_be and fread_be
00045  *
00046  * Revision 1.1.1.1  2001/11/20 15:19:28  e_gourgoulhon
00047  * LORENE
00048  *
00049  * Revision 2.17  2001/10/24  15:36:20  eric
00050  * Ajout de la fonction display_poly.
00051  *
00052  * Revision 2.16  2001/10/16  14:49:02  eric
00053  * Appel de get_omega_c() pour avoir la valeur centrale de Omega.
00054  * Affichage different si rotation differentielle.
00055  *
00056  * Revision 2.15  2001/09/13  08:32:01  eric
00057  * Ajout du facteur de compacite M/R dans l'affichage.
00058  *
00059  * Revision 2.14  2001/06/20  14:20:56  novak
00060  * Appel a Etoile_rot::set_der0x0 dans del_deriv (au lieu de set_der0x0
00061  * tout court).
00062  *
00063  * Revision 2.13  2001/03/26 09:30:58  jlz
00064  * New members p_espec_isco and p_lspec_isco.
00065  *
00066  * Revision 2.12  2000/11/20  21:42:02  eric
00067  * Appel de fait_nphi() dans le constructeur par lecture de fichier.
00068  *
00069  * Revision 2.11  2000/11/18  23:18:30  eric
00070  * Modifs affichage.
00071  *
00072  * Revision 2.10  2000/11/18  21:09:57  eric
00073  * Ajout des membres  p_r_isco et p_f_isco.
00074  *
00075  * Revision 2.9  2000/11/07  16:33:08  eric
00076  * Modif affichage.
00077  *
00078  * Revision 2.8  2000/10/12  15:37:01  eric
00079  * Ajout de la fonction fait_nphi().
00080  *
00081  * Revision 2.7  2000/09/18  16:15:12  eric
00082  * Ajout du membre tkij.
00083  *
00084  * Revision 2.6  2000/08/31  15:38:00  eric
00085  * Bases spectrales standards pour bbb et b_car dans le constructeur
00086  * standard (initialisation a la metrique plate).
00087  *
00088  * Revision 2.5  2000/08/31  11:25:45  eric
00089  * Ajout des membres tnphi et ak_car.
00090  *
00091  * Revision 2.4  2000/08/25  12:28:29  eric
00092  * Modif affichage.
00093  *
00094  * Revision 2.3  2000/08/18  14:01:59  eric
00095  * Ajout de partial_display
00096  *
00097  * Revision 2.2  2000/08/17  12:40:04  eric
00098  * *** empty log message ***
00099  *
00100  * Revision 2.1  2000/07/21  16:31:26  eric
00101  * *** empty log message ***
00102  *
00103  * Revision 1.1  2000/07/20  15:32:37  eric
00104  * Initial revision
00105  *
00106  *
00107  * $Header: /cvsroot/Lorene/C++/Source/Etoile/etoile_rot.C,v 1.4 2004/03/25 10:29:07 j_novak Exp $
00108  *
00109  */
00110 
00111 // Headers C
00112 #include "math.h"
00113 
00114 // Headers Lorene
00115 #include "etoile.h"
00116 #include "eos.h"
00117 #include "nbr_spx.h"
00118 #include "utilitaires.h"
00119 #include "unites.h"     
00120 
00121                 //--------------//
00122                 // Constructors //
00123                 //--------------//
00124 
00125 // Standard constructor
00126 // --------------------
00127 Etoile_rot::Etoile_rot(Map& mpi, int nzet_i, bool relat, const Eos& eos_i)
00128                : Etoile(mpi, nzet_i, relat, eos_i), 
00129              bbb(mpi), 
00130              b_car(mpi), 
00131              nphi(mpi), 
00132              tnphi(mpi), 
00133              uuu(mpi), 
00134              logn(logn_auto), 
00135              nuf(mpi), 
00136              nuq(mpi), 
00137              dzeta(beta_auto), 
00138              tggg(mpi), 
00139              w_shift(mpi, 1, CON, mp.get_bvect_cart()), 
00140              khi_shift(mpi), 
00141              tkij(mpi, 2, COV, mp.get_bvect_cart()), 
00142              ak_car(mpi), 
00143              ssjm1_nuf(mpi), 
00144              ssjm1_nuq(mpi), 
00145              ssjm1_dzeta(mpi), 
00146              ssjm1_tggg(mpi), 
00147              ssjm1_khi(mpi), 
00148              ssjm1_wshift(mpi, 1, CON, mp.get_bvect_cart())          
00149 {
00150 
00151     // Initialization to a static state : 
00152     omega = 0 ; 
00153     uuu = 0 ; 
00154 
00155     // Initialization to a flat metric : 
00156     bbb = 1 ;
00157     bbb.set_std_base() ; 
00158     b_car = 1 ;
00159     b_car.set_std_base() ;
00160     nphi = 0 ;   
00161     tnphi = 0 ;   
00162     nuf = 0 ; 
00163     nuq = 0 ; 
00164     tggg = 0 ;   
00165 
00166     w_shift.set_etat_qcq() ; 
00167     for (int i=0; i<3; i++) {
00168     w_shift.set(i) = 0 ; 
00169     }
00170 
00171     khi_shift.set_etat_qcq() ; 
00172     khi_shift.set() = 0 ; 
00173 
00174     tkij.set_etat_zero() ; 
00175 
00176     ak_car = 0 ; 
00177 
00178     ssjm1_nuf = 0 ; 
00179     ssjm1_nuq = 0 ; 
00180     ssjm1_dzeta = 0 ; 
00181     ssjm1_tggg = 0 ; 
00182     ssjm1_khi = 0 ; 
00183 
00184     ssjm1_wshift.set_etat_qcq() ; 
00185     for (int i=0; i<3; i++) {
00186     ssjm1_wshift.set(i) = 0 ; 
00187     }
00188     
00189     // Pointers of derived quantities initialized to zero : 
00190     set_der_0x0() ;
00191     
00192 }
00193 
00194 // Copy constructor
00195 // ----------------
00196 
00197 Etoile_rot::Etoile_rot(const Etoile_rot& et)
00198                : Etoile(et), 
00199              bbb(et.bbb), 
00200              b_car(et.b_car), 
00201              nphi(et.nphi), 
00202              tnphi(et.tnphi), 
00203              uuu(et.uuu), 
00204              logn(logn_auto), 
00205              nuf(et.nuf), 
00206              nuq(et.nuq),
00207              dzeta(beta_auto), 
00208              tggg(et.tggg), 
00209              w_shift(et.w_shift), 
00210              khi_shift(et.khi_shift), 
00211              tkij(et.tkij),
00212              ak_car(et.ak_car), 
00213              ssjm1_nuf(et.ssjm1_nuf), 
00214              ssjm1_nuq(et.ssjm1_nuq), 
00215              ssjm1_dzeta(et.ssjm1_dzeta), 
00216              ssjm1_tggg(et.ssjm1_tggg), 
00217              ssjm1_khi(et.ssjm1_khi), 
00218              ssjm1_wshift(et.ssjm1_wshift)           
00219 {
00220     omega = et.omega ; 
00221 
00222     // Pointers of derived quantities initialized to zero : 
00223     set_der_0x0() ;
00224 }    
00225 
00226 
00227 // Constructor from a file
00228 // -----------------------
00229 Etoile_rot::Etoile_rot(Map& mpi, const Eos& eos_i, FILE* fich)
00230                : Etoile(mpi, eos_i, fich), 
00231              bbb(mpi), 
00232              b_car(mpi), 
00233              nphi(mpi), 
00234              tnphi(mpi), 
00235              uuu(mpi), 
00236              logn(logn_auto), 
00237              nuf(mpi), 
00238              nuq(mpi), 
00239              dzeta(beta_auto), 
00240              tggg(mpi), 
00241              w_shift(mpi, 1, CON, mp.get_bvect_cart()), 
00242              khi_shift(mpi), 
00243              tkij(mpi, 2, COV, mp.get_bvect_cart()), 
00244              ak_car(mpi), 
00245              ssjm1_nuf(mpi), 
00246              ssjm1_nuq(mpi), 
00247              ssjm1_dzeta(mpi), 
00248              ssjm1_tggg(mpi), 
00249              ssjm1_khi(mpi), 
00250              ssjm1_wshift(mpi, 1, CON, mp.get_bvect_cart())          
00251 {
00252 
00253     // Etoile parameters
00254     // -----------------
00255 
00256     // omega is read in the file:     
00257     fread_be(&omega, sizeof(double), 1, fich) ;     
00258           
00259    
00260     // Read of the saved fields:
00261     // ------------------------
00262 
00263     Tenseur nuf_file(mp, fich) ; 
00264     nuf = nuf_file ; 
00265         
00266     Tenseur nuq_file(mp, fich) ; 
00267     nuq = nuq_file ; 
00268         
00269     Tenseur tggg_file(mp, fich) ; 
00270     tggg = tggg_file ; 
00271         
00272     Tenseur w_shift_file(mp, mp.get_bvect_cart(), fich) ; 
00273     w_shift = w_shift_file ;
00274     
00275     Tenseur khi_shift_file(mp, fich) ; 
00276     khi_shift = khi_shift_file ;
00277     
00278     fait_shift() ;      // constructs shift from w_shift and khi_shift
00279     fait_nphi() ;       // constructs N^phi from (N^x,N^y,N^z)
00280     
00281     Cmp ssjm1_nuf_file(mp, *(mp.get_mg()), fich) ; 
00282     ssjm1_nuf = ssjm1_nuf_file ; 
00283 
00284     Cmp ssjm1_nuq_file(mp, *(mp.get_mg()), fich) ; 
00285     ssjm1_nuq = ssjm1_nuq_file ; 
00286 
00287     Cmp ssjm1_dzeta_file(mp, *(mp.get_mg()), fich) ; 
00288     ssjm1_dzeta = ssjm1_dzeta_file ; 
00289 
00290     Cmp ssjm1_tggg_file(mp, *(mp.get_mg()), fich) ; 
00291     ssjm1_tggg = ssjm1_tggg_file ; 
00292 
00293     Cmp ssjm1_khi_file(mp, *(mp.get_mg()), fich) ; 
00294     ssjm1_khi = ssjm1_khi_file ; 
00295 
00296     Tenseur ssjm1_wshift_file(mp, mp.get_bvect_cart(), fich) ; 
00297     ssjm1_wshift = ssjm1_wshift_file ; 
00298 
00299     // All other fields are initialized to zero : 
00300     // ----------------------------------------
00301     bbb = 0 ; 
00302     b_car = 0 ; 
00303     uuu = 0 ; 
00304 
00305     // Pointers of derived quantities initialized to zero 
00306     // --------------------------------------------------
00307     set_der_0x0() ;
00308     
00309 }
00310 
00311                 //------------//
00312                 // Destructor //
00313                 //------------//
00314 
00315 Etoile_rot::~Etoile_rot(){
00316 
00317     del_deriv() ; 
00318 
00319 }
00320 
00321         //----------------------------------//
00322         // Management of derived quantities //
00323         //----------------------------------//
00324 
00325 void Etoile_rot::del_deriv() const {
00326 
00327     Etoile::del_deriv() ; 
00328 
00329     if (p_angu_mom != 0x0) delete p_angu_mom ; 
00330     if (p_tsw != 0x0) delete p_tsw ; 
00331     if (p_grv2 != 0x0) delete p_grv2 ; 
00332     if (p_grv3 != 0x0) delete p_grv3 ; 
00333     if (p_r_circ != 0x0) delete p_r_circ ; 
00334     if (p_aplat != 0x0) delete p_aplat ; 
00335     if (p_z_eqf != 0x0) delete p_z_eqf ; 
00336     if (p_z_eqb != 0x0) delete p_z_eqb ; 
00337     if (p_z_pole != 0x0) delete p_z_pole ; 
00338     if (p_mom_quad != 0x0) delete p_mom_quad ;
00339     if (p_r_isco != 0x0) delete p_r_isco ;
00340     if (p_f_isco != 0x0) delete p_f_isco ;
00341     if (p_lspec_isco != 0x0) delete p_lspec_isco ;
00342     if (p_espec_isco != 0x0) delete p_espec_isco ;
00343     if (p_f_eq != 0x0) delete p_f_eq ;
00344     
00345     Etoile_rot::set_der_0x0() ; 
00346 }               
00347 
00348 
00349 
00350 
00351 void Etoile_rot::set_der_0x0() const {
00352 
00353     Etoile::set_der_0x0() ;
00354 
00355     p_angu_mom = 0x0 ; 
00356     p_tsw = 0x0 ;
00357     p_grv2 = 0x0 ;
00358     p_grv3 = 0x0 ;
00359     p_r_circ = 0x0 ;
00360     p_aplat = 0x0 ;
00361     p_z_eqf = 0x0 ;
00362     p_z_eqb = 0x0 ;
00363     p_z_pole = 0x0 ;
00364     p_mom_quad = 0x0 ;
00365     p_r_isco = 0x0 ;
00366     p_f_isco = 0x0 ;
00367     p_lspec_isco = 0x0 ;
00368     p_espec_isco = 0x0 ;
00369     p_f_eq = 0x0 ;
00370     
00371 }               
00372 
00373 void Etoile_rot::del_hydro_euler() {
00374 
00375     Etoile::del_hydro_euler() ; 
00376 
00377     del_deriv() ; 
00378 
00379 }               
00380 
00381 
00382                 //--------------//
00383                 //  Assignment  //
00384                 //--------------//
00385 
00386 // Assignment to another Etoile_rot
00387 // --------------------------------
00388 void Etoile_rot::operator=(const Etoile_rot& et) {
00389 
00390     // Assignment of quantities common to all the derived classes of Etoile
00391     Etoile::operator=(et) ;     
00392 
00393     // Assignement of proper quantities of class Etoile_rot
00394     omega = et.omega ; 
00395 
00396     bbb = et.bbb ; 
00397     b_car = et.b_car ; 
00398     nphi = et.nphi ; 
00399     tnphi = et.tnphi ; 
00400     uuu = et.uuu ; 
00401     nuf = et.nuf ; 
00402     nuq = et.nuq ; 
00403     tggg = et.tggg ; 
00404     w_shift = et.w_shift ;
00405     khi_shift = et.khi_shift ;
00406     tkij = et.tkij ; 
00407     ak_car = et.ak_car ;
00408     ssjm1_nuf = et.ssjm1_nuf ;
00409     ssjm1_nuq = et.ssjm1_nuq ;
00410     ssjm1_dzeta = et.ssjm1_dzeta ; 
00411     ssjm1_tggg = et.ssjm1_tggg ;
00412     ssjm1_khi = et.ssjm1_khi ;
00413     ssjm1_wshift = et.ssjm1_wshift ; 
00414     
00415     del_deriv() ;  // Deletes all derived quantities
00416 
00417 }   
00418 
00419                 //--------------//
00420                 //    Outputs   //
00421                 //--------------//
00422 
00423 // Save in a file
00424 // --------------
00425 void Etoile_rot::sauve(FILE* fich) const {
00426     
00427     Etoile::sauve(fich) ; 
00428     
00429     fwrite_be(&omega, sizeof(double), 1, fich) ;        
00430     
00431     nuf.sauve(fich) ; 
00432     nuq.sauve(fich) ; 
00433     tggg.sauve(fich) ; 
00434     w_shift.sauve(fich) ; 
00435     khi_shift.sauve(fich) ; 
00436     
00437     ssjm1_nuf.sauve(fich) ; 
00438     ssjm1_nuq.sauve(fich) ; 
00439     ssjm1_dzeta.sauve(fich) ; 
00440     ssjm1_tggg.sauve(fich) ; 
00441     ssjm1_khi.sauve(fich) ; 
00442     ssjm1_wshift.sauve(fich) ; 
00443     
00444     
00445 }
00446 
00447 // Printing
00448 // --------
00449 
00450 ostream& Etoile_rot::operator>>(ostream& ost) const {
00451     
00452   using namespace Unites ;
00453 
00454     Etoile::operator>>(ost) ; 
00455 
00456     double omega_c = get_omega_c() ; 
00457     
00458     ost << endl ; 
00459     if (omega != __infinity) {
00460     ost << "Uniformly rotating star" << endl ; 
00461     ost << "-----------------------" << endl ; 
00462     
00463     double freq = omega / (2.*M_PI) ;  
00464     ost << "Omega : " << omega * f_unit 
00465         << " rad/s     f : " << freq * f_unit << " Hz" << endl ; 
00466     ost << "Rotation period : " << 1000. / (freq * f_unit) << " ms"
00467         << endl ;
00468     
00469     }
00470     else {
00471     ost << "Differentially rotating star" << endl ; 
00472     ost << "----------------------------" << endl ; 
00473     
00474     double freq = omega_c / (2.*M_PI) ;  
00475     ost << "Central value of Omega : " << omega_c * f_unit 
00476         << " rad/s     f : " << freq * f_unit << " Hz" << endl ; 
00477     ost << "Central rotation period : " << 1000. / (freq * f_unit) << " ms"
00478         << endl ;
00479     
00480     }
00481     
00482        
00483     double nphi_c = nphi()(0, 0, 0, 0) ;
00484     if ( (omega_c==0) && (nphi_c==0) ) {
00485         ost << "Central N^phi :               " << nphi_c << endl ;
00486     }
00487     else{
00488         ost << "Central N^phi/Omega :         " << nphi_c / omega_c << endl ;
00489     }
00490         
00491     ost << "Error on the virial identity GRV2 : " << endl ; 
00492     ost << "GRV2 = " << grv2() << endl ; 
00493     ost << "Error on the virial identity GRV3 : " << endl ; 
00494     double xgrv3 = grv3(&ost) ; 
00495     ost << "GRV3 = " << xgrv3 << endl ; 
00496 
00497     double mom_quad_38si = mom_quad() * rho_unit * (pow(r_unit, double(5.)) 
00498                             / double(1.e38) ) ;
00499     ost << "Quadrupole moment Q : " << mom_quad_38si << " 10^38 kg m^2"
00500          << endl ; 
00501     ost << "Q / (M R_circ^2) :       " 
00502      << mom_quad() / ( mass_g() * pow( r_circ(), 2. ) ) << endl ; 
00503     ost << "c^4 Q / (G^2 M^3) :      " 
00504      << mom_quad() / ( pow(qpig/(4*M_PI), 2.) * pow(mass_g(), 3.) ) 
00505      << endl ; 
00506 
00507     ost << "Angular momentum J :      " 
00508      << angu_mom()/( qpig / (4* M_PI) * msol*msol) << " G M_sol^2 / c"
00509      << endl ; 
00510     ost << "c J / (G M^2) :           " 
00511      << angu_mom()/( qpig / (4* M_PI) * pow(mass_g(), 2.) ) << endl ;           
00512 
00513     if (omega != __infinity) {
00514     double mom_iner = angu_mom() / omega ; 
00515     double mom_iner_38si = mom_iner * rho_unit * (pow(r_unit, double(5.)) 
00516         / double(1.e38) ) ; 
00517     ost << "Moment of inertia:       " << mom_iner_38si << " 10^38 kg m^2"
00518         << endl ; 
00519     }
00520 
00521     ost << "Ratio T/W :            " << tsw() << endl ; 
00522     ost << "Circumferential equatorial radius R_circ :     " 
00523     << r_circ()/km << " km" << endl ;  
00524     ost << "Coordinate equatorial radius r_eq : " << ray_eq()/km << " km" 
00525      << endl ;  
00526     ost << "Flattening r_pole/r_eq :        " << aplat() << endl ; 
00527 
00528     double compact = qpig/(4.*M_PI) * mass_g() / r_circ() ; 
00529     ost << "Compaction parameter M_g / R_circ : " << compact << endl ; 
00530 
00531     int lsurf = nzet - 1; 
00532     int nt = mp.get_mg()->get_nt(lsurf) ; 
00533     int nr = mp.get_mg()->get_nr(lsurf) ; 
00534     ost << "Equatorial value of the velocity U:         " 
00535      << uuu()(nzet-1, 0, nt-1, nr-1) << " c" << endl ; 
00536     ost << "Redshift at the equator (forward) : " << z_eqf() << endl ; 
00537     ost << "Redshift at the equator (backward): " << z_eqb() << endl ; 
00538     ost << "Redshift at the pole              : " << z_pole() << endl ; 
00539 
00540 
00541     ost << "Central value of log(N)        : " 
00542     << logn()(0, 0, 0, 0) << endl ; 
00543 
00544     ost << "Central value of dzeta=log(AN) : " 
00545     << dzeta()(0, 0, 0, 0) << endl ; 
00546 
00547     if ( (omega_c==0) && (nphi_c==0) ) {
00548         ost << "Central N^phi :               " << nphi_c << endl ;
00549     }
00550     else{
00551         ost << "Central N^phi/Omega :         " << nphi_c / omega_c << endl ;
00552     }
00553 
00554     ost << "  ... w_shift (NB: components in the star Cartesian frame) [c] :  "
00555     << endl  
00556     << w_shift(0)(0, 0, 0, 0) << "  " 
00557     << w_shift(1)(0, 0, 0, 0) << "  " 
00558     << w_shift(2)(0, 0, 0, 0) << endl ; 
00559 
00560     ost << "Central value of khi_shift [km c] : " 
00561         << khi_shift()(0, 0, 0, 0) / km << endl ; 
00562 
00563     ost << "Central value of B^2 : " << b_car()(0,0,0,0) <<  endl ; 
00564 
00565     Tbl diff_a_b = diffrel( a_car(), b_car() ) ;
00566     ost << 
00567     "Relative discrepancy in each domain between the metric coef. A^2 and B^2 : "
00568        << endl ;
00569     for (int l=0; l<diff_a_b.get_taille(); l++) {
00570         ost << diff_a_b(l) << "  " ;
00571     }
00572     ost << endl ;
00573 
00574     // Approximate formula for R_isco = 6 R_g (1-(2/3)^1.5 j )
00575     // up to the first order in j
00576     double jdimless = angu_mom() / ( ggrav * pow(mass_g(), 2.) ) ;
00577     double r_grav = ggrav * mass_g() ;
00578     double r_isco_appr = 6. * r_grav * ( 1. - pow(2./3.,1.5) * jdimless ) ;
00579 
00580     // Approximate formula for the ISCO frequency
00581     //   freq_ms = 6^{-1.5}/2pi/R_g (1+11*6^(-1.5) j )
00582     // up to the first order in j
00583     double f_isco_appr = ( 1. + 11. /6. /sqrt(6.) * jdimless ) / r_grav /
00584                             (12. * M_PI ) / sqrt(6.) ;
00585 
00586     ost << endl << "Innermost stable circular orbit (ISCO) : " << endl ;
00587     double xr_isco = r_isco(&ost) ;
00588     ost <<"    circumferential radius r_isco = "
00589         << xr_isco / km << " km" << endl ;
00590     ost <<"     (approx. 6M + 1st order in j : "
00591         << r_isco_appr / km << " km)" << endl ;
00592     ost <<"                      (approx. 6M : "
00593         << 6. * r_grav / km << " km)" << endl ;
00594     ost <<"    orbital frequency f_isco = "
00595         << f_isco() * f_unit << " Hz" << endl ;
00596     ost <<"     (approx. 1st order in j : "
00597         << f_isco_appr * f_unit << " Hz)" << endl ;
00598         
00599 
00600     return ost ;
00601     
00602 }
00603 
00604 
00605 void Etoile_rot::partial_display(ostream& ost) const {
00606     
00607   using namespace Unites ;
00608 
00609     double omega_c = get_omega_c() ; 
00610     double freq = omega_c / (2.*M_PI) ;  
00611     ost << "Central Omega : " << omega_c * f_unit 
00612         << " rad/s     f : " << freq * f_unit << " Hz" << endl ; 
00613     ost << "Rotation period : " << 1000. / (freq * f_unit) << " ms"
00614      << endl ;
00615     ost << endl << "Central enthalpy : " << ent()(0,0,0,0) << " c^2" << endl ; 
00616     ost << "Central proper baryon density : " << nbar()(0,0,0,0) 
00617     << " x 0.1 fm^-3" << endl ; 
00618     ost << "Central proper energy density : " << ener()(0,0,0,0) 
00619     << " rho_nuc c^2" << endl ; 
00620     ost << "Central pressure : " << press()(0,0,0,0) 
00621     << " rho_nuc c^2" << endl ; 
00622 
00623     ost << "Central value of log(N)        : " 
00624     << logn()(0, 0, 0, 0) << endl ; 
00625     ost << "Central lapse N :      " << nnn()(0,0,0,0) <<  endl ; 
00626     ost << "Central value of dzeta=log(AN) : " 
00627     << dzeta()(0, 0, 0, 0) << endl ; 
00628     ost << "Central value of A^2 : " << a_car()(0,0,0,0) <<  endl ; 
00629     ost << "Central value of B^2 : " << b_car()(0,0,0,0) <<  endl ; 
00630 
00631     double nphi_c = nphi()(0, 0, 0, 0) ;
00632     if ( (omega_c==0) && (nphi_c==0) ) {
00633         ost << "Central N^phi :               " << nphi_c << endl ;
00634     }
00635     else{
00636         ost << "Central N^phi/Omega :         " << nphi_c / omega_c 
00637                             << endl ;
00638     }
00639         
00640 
00641     int lsurf = nzet - 1; 
00642     int nt = mp.get_mg()->get_nt(lsurf) ; 
00643     int nr = mp.get_mg()->get_nr(lsurf) ; 
00644     ost << "Equatorial value of the velocity U:         " 
00645      << uuu()(nzet-1, 0, nt-1, nr-1) << " c" << endl ; 
00646 
00647     ost << endl 
00648     << "Coordinate equatorial radius r_eq =    " 
00649     << ray_eq()/km << " km" << endl ;  
00650     ost << "Flattening r_pole/r_eq :        " << aplat() << endl ; 
00651 
00652 }
00653 
00654 
00655 double Etoile_rot::get_omega_c() const {
00656     
00657     return omega ; 
00658     
00659 }
00660 
00661 
00662 // display_poly
00663 // ------------
00664 
00665 void Etoile_rot::display_poly(ostream& ost) const {
00666 
00667   using namespace Unites ;
00668 
00669   const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( &eos ) ;    
00670 
00671     if (p_eos_poly != 0x0) {
00672 
00673     double kappa = p_eos_poly->get_kap() ; 
00674     double gamma = p_eos_poly->get_gam() ;  ; 
00675 
00676     // kappa^{n/2}
00677     double kap_ns2 = pow( kappa,  0.5 /(gamma-1) ) ; 
00678     
00679     // Polytropic unit of length in terms of r_unit : 
00680     double r_poly = kap_ns2 / sqrt(ggrav) ; 
00681     
00682     // Polytropic unit of time in terms of t_unit :
00683     double t_poly = r_poly ; 
00684 
00685     // Polytropic unit of mass in terms of m_unit :
00686     double m_poly = r_poly / ggrav ; 
00687     
00688     // Polytropic unit of angular momentum in terms of j_unit :
00689     double j_poly = r_poly * r_poly / ggrav ; 
00690     
00691     // Polytropic unit of density in terms of rho_unit :
00692     double rho_poly = 1. / (ggrav * r_poly * r_poly) ;  
00693 
00694     ost.precision(10) ; 
00695     ost << endl << "Quantities in polytropic units : " << endl ; 
00696     ost  << "==============================" << endl ; 
00697     ost << " ( r_poly = " << r_poly / km << " km )" << endl ; 
00698     ost << "  n_c     : " << nbar()(0, 0, 0, 0) / rho_poly << endl ; 
00699     ost << "  e_c     : " << ener()(0, 0, 0, 0) / rho_poly << endl ; 
00700     ost << "  Omega_c : " << get_omega_c() * t_poly << endl ; 
00701     ost << "  P_c     : " << 2.*M_PI / get_omega_c() / t_poly << endl ; 
00702     ost << "  M_bar   : " << mass_b() / m_poly << endl ; 
00703     ost << "  M       : " << mass_g() / m_poly << endl ; 
00704     ost << "  J   : " << angu_mom() / j_poly << endl ; 
00705     ost << "  r_eq    : " << ray_eq() / r_poly << endl ; 
00706     ost << "  R_circ  : " << r_circ() / r_poly << endl ; 
00707     
00708     
00709     }
00710     
00711 
00712 } 
00713 
00714 
00715 
00716 
00717 
00718 
00719                 //-------------------------//
00720                 //  Computational routines //
00721                 //-------------------------//
00722                 
00723 void Etoile_rot::fait_shift() {
00724 
00725     Tenseur d_khi = khi_shift.gradient() ; 
00726     
00727     if (d_khi.get_etat() == ETATQCQ) { 
00728     d_khi.dec2_dzpuis() ;   // divide by r^2 in the external compactified
00729                 // domain
00730     }
00731  
00732     // x_k dW^k/dx_i
00733     
00734     Tenseur x_d_w = skxk( w_shift.gradient() ) ;
00735     x_d_w.dec_dzpuis() ;
00736     
00737     double lambda = double(1) / double(3) ; 
00738 
00739     // The final computation is done component by component because
00740     // d_khi and x_d_w are covariant comp. whereas w_shift is
00741     // contravariant
00742     
00743     shift.set_etat_qcq() ; 
00744     
00745     for (int i=0; i<3; i++) {
00746     shift.set(i) = (lambda+2)/2./(lambda+1) * w_shift(i)
00747         - (lambda/2./(lambda+1))  * (d_khi(i) + x_d_w(i)) ;      
00748     }
00749     
00750     shift.set_triad( *(w_shift.get_triad()) ) ; 
00751     
00752 } 
00753 
00754 
00755 
00756 void Etoile_rot::fait_nphi() {
00757 
00758     if ( shift.get_etat() == ETATZERO ) {
00759     tnphi = 0 ; 
00760     nphi = 0 ;
00761     return ;  
00762     }
00763 
00764     assert( shift.get_etat() == ETATQCQ ) ; 
00765     
00766     // Computation of tnphi
00767     // --------------------
00768     tnphi.set_etat_qcq() ; 
00769     
00770     mp.comp_p_from_cartesian( shift(0), shift(1), tnphi.set() ) ; 
00771 
00772     // Computation of nphi
00773     // -------------------
00774     
00775     nphi = tnphi ; 
00776     (nphi.set()).div_rsint() ; 
00777     
00778 }

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