star_rot.C

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

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