et_bin_hydro.C

00001 /*
00002  * Methods of the class Etoile_bin for computing hydro quantities
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 et_bin_hydro_C[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_hydro.C,v 1.5 2005/08/29 15:10:16 p_grandclement Exp $" ;
00030 
00031 /*
00032  * $Id: et_bin_hydro.C,v 1.5 2005/08/29 15:10:16 p_grandclement Exp $
00033  * $Log: et_bin_hydro.C,v $
00034  * Revision 1.5  2005/08/29 15:10:16  p_grandclement
00035  * Addition of things needed :
00036  *   1) For BBH with different masses
00037  *   2) Provisory files for the mixted binaries (Bh and NS) : THIS IS NOT
00038  *   WORKING YET !!!
00039  *
00040  * Revision 1.4  2003/03/03 19:46:09  f_limousin
00041  * Set standard bases for s_euler.
00042  *
00043  * Revision 1.3  2003/01/17 13:34:56  f_limousin
00044  * Replace A^2*flat_scalar_prod by sprod
00045  *
00046  * Revision 1.2  2002/12/10 15:29:29  k_taniguchi
00047  * Change the multiplication "*" to "%"
00048  *   and flat_scalar_prod to flat_scalar_prod_desal.
00049  *
00050  * Revision 1.1.1.1  2001/11/20 15:19:28  e_gourgoulhon
00051  * LORENE
00052  *
00053  * Revision 2.11  2000/12/22  13:10:32  eric
00054  * Modif des annulations en dehors de l'etoile.
00055  *
00056  * Revision 2.10  2000/03/29  11:47:53  eric
00057  * Suppression affichage.
00058  *
00059  * Revision 2.9  2000/03/01  16:12:09  eric
00060  * Appel de set_std_base sur u_euler dans le cas irrotationnel.
00061  *
00062  * Revision 2.8  2000/02/24  09:34:10  keisuke
00063  * Ajout de l'appel a set_std_base() sur wit_w.
00064  *
00065  * Revision 2.7  2000/02/21  13:59:09  eric
00066  * Appel de set_dzpuis(0) sur loggam.
00067  *
00068  * Revision 2.6  2000/02/21  08:43:17  keisuke
00069  * Ajout de l'appel a set_std_base() sur loggam.
00070  *
00071  * Revision 2.5  2000/02/17  14:42:45  eric
00072  * Ajout de l'appel a set_std_base() sur gam_euler.
00073  *
00074  * Revision 2.4  2000/02/14  11:06:15  eric
00075  * Suppression de l'appel a annule(nzet.nzm1) sur gam_euler dans le cas en
00076  *   corotation.
00077  * Ajout de l'appel a annule(nzet,nzm1) sur wit_w.
00078  *
00079  * Revision 2.3  2000/02/12  18:36:23  eric
00080  * Appel de set_std_base sur loggam.
00081  *
00082  * Revision 2.2  2000/02/08  19:29:03  eric
00083  * Calcul sur les tenseurs.
00084  * wit_w est calcule.
00085  *
00086  * Revision 2.1  2000/02/04  16:37:28  eric
00087  * Premiere version implementee (non testee !)/
00088  *
00089  * Revision 2.0  2000/01/31  15:57:30  eric
00090  * *** empty log message ***
00091  *
00092  *
00093  * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_hydro.C,v 1.5 2005/08/29 15:10:16 p_grandclement Exp $
00094  *
00095  */
00096 
00097 // Headers C
00098 
00099 // Headers Lorene
00100 #include "etoile.h"
00101 
00102 void Etoile_bin::hydro_euler(){
00103 
00104     int nz = mp.get_mg()->get_nzone() ; 
00105     int nzm1 = nz - 1 ; 
00106 
00107     //----------------------------------
00108     // Specific relativistic enthalpy           ---> hhh
00109     //----------------------------------
00110     
00111     Tenseur hhh = exp(unsurc2 * ent) ;  // = 1 at the Newtonian limit
00112     hhh.set_std_base() ;
00113 
00114     //---------------------------------------------------
00115     // Lorentz factor between the co-orbiting             ---> gam0
00116     // observer and the Eulerian one
00117     // See Eq (23) and (24) from Gourgoulhon et al. (2001)
00118     //---------------------------------------------------
00119 
00120     Tenseur gam0 = 1/sqrt(1-unsurc2*sprod(bsn,bsn)) ;
00121     gam0.set_std_base() ;
00122 
00123     //------------------------------------------
00124     // Lorentz factor and 3-velocity of the fluid 
00125     //  with respect to the Eulerian observer
00126     //------------------------------------------
00127     
00128     if (irrotational) { 
00129 
00130 //##    cout << "Etoile_bin::hydro_euler : ### warning : " << endl ; 
00131 //  cout << "   d_psi.d_psi would be better computed in spher. coord. !"
00132 //##         << endl ; 
00133 
00134         d_psi.set_std_base() ;
00135 
00136     // See Eq (32) from Gourgoulhon et al. (2001)
00137     gam_euler = sqrt( 1 + unsurc2 * sprod(d_psi, d_psi)
00138                     / (hhh%hhh) ) ; 
00139 
00140     gam_euler.set_std_base() ;  // sets the standard spectral bases for
00141                     // a scalar field
00142 
00143 
00144     // See Eq (31) from Gourgoulhon et al. (2001)
00145     Tenseur vtmp = d_psi / ( hhh % gam_euler % a_car ) ; 
00146 
00147     // The assignment of u_euler is performed component by component 
00148     //  because u_euler is contravariant and d_psi is covariant
00149     if (vtmp.get_etat() == ETATZERO) {
00150         u_euler.set_etat_zero() ; 
00151     }
00152     else{
00153         assert(vtmp.get_etat() == ETATQCQ) ; 
00154         u_euler.set_etat_qcq() ; 
00155         for (int i=0; i<3; i++) {
00156         u_euler.set(i) = vtmp(i) ;
00157         }
00158         u_euler.set_triad( *(vtmp.get_triad()) ) ; 
00159     }
00160     
00161     u_euler.set_std_base() ; 
00162 
00163     }
00164     else {      // Rigid rotation
00165             // --------------
00166 
00167     gam_euler = gam0 ; 
00168     gam_euler.set_std_base() ;  // sets the standard spectral bases for
00169                     // a scalar field
00170 
00171     u_euler = -bsn ; 
00172 
00173     }
00174     
00175     //------------------------------------
00176     //  Energy density E with respect to the Eulerian observer
00177     // See Eq (53) from Gourgoulhon et al. (2001)  
00178     //------------------------------------
00179     
00180     ener_euler = gam_euler % gam_euler % ( ener + press ) - press ; 
00181     
00182     //------------------------------------
00183     // Trace of the stress tensor with respect to the Eulerian observer
00184     // See Eq (54) from Gourgoulhon et al. (2001)  
00185     //------------------------------------
00186 
00187     //Tenseur tmp00 = sprod(u_euler, u_euler) ; 
00188     //cout << hex << u_euler(0).va.base.b[0] << endl ; 
00189     //cout << hex << u_euler(1).va.base.b[0] << endl ; 
00190     //cout << hex << u_euler(2).va.base.b[0] << endl ; 
00191     //cout << hex << tmp00().va.base.b[0] << endl ; 
00192 
00193     s_euler = 3 * press  +  ( ener_euler + press ) *
00194       sprod(u_euler, u_euler) ;
00195     s_euler.set_std_base() ; 
00196 
00197 
00198     //-------------------------------------------
00199     //  Lorentz factor between the fluid and        ---> gam
00200     //  co-orbiting observers
00201     // See Eq (58) from Gourgoulhon et al. (2001)  
00202     //--------------------------------------------
00203     
00204     Tenseur tmp = ( 1+unsurc2*sprod(bsn,u_euler) ) ;
00205     tmp.set_std_base() ;
00206     Tenseur gam = gam0 % gam_euler % tmp ;
00207 
00208     //-------------------------------------------
00209     //  Spatial projection of the fluid 3-velocity
00210     //  with respect to the co-orbiting observer
00211     //--------------------------------------------
00212 
00213     wit_w = gam_euler / gam % u_euler + gam0 % bsn ; 
00214 
00215     wit_w.set_std_base() ;  // set the bases for spectral expansions
00216 
00217     wit_w.annule(nzm1) ;    // zero in the ZEC
00218 
00219 
00220     //-------------------------------------------
00221     //  Logarithm of the Lorentz factor between 
00222     //  the fluid and co-orbiting observers
00223     //--------------------------------------------
00224 
00225     if (relativistic) {
00226 
00227     loggam = log( gam ) ;
00228     
00229         loggam.set_std_base() ;   // set the bases for spectral expansions
00230     }
00231     else {
00232   
00233     loggam = 0.5 * flat_scalar_prod_desal(wit_w, wit_w) ;
00234 
00235         loggam.set_std_base() ;   // set the bases for spectral expansions
00236 
00237 //### Forcage a zero des termes en sin(m*phi) : 
00238 //  loggam.coef() ; 
00239 //  int np = mgrille->np[0] ; 
00240 //  int nt = mgrille->nt[0] ;
00241 //  int nr = mgrille->nr[0] ;
00242 //  int ntnr = nt * nr ; 
00243 //  for (int k = 1; k < np+2; k+=2) {
00244 //      for (int j = 0; j < nt; j++) {
00245 //      for (int i = 0; i<nr; i++) {
00246 //          loggam.c_cf[0]->t[0]->t[ntnr*k + nr*j + i] = 0 ;
00247 //      }
00248 //      }
00249 //  }
00250 //  loggam.c_ajx[0] = 0 ; 
00251 //  loggam.c_aj = 0 ; 
00252 //  loggam.coef_i() ; 
00253 //###
00254     }
00255 
00256 
00257     //-------------------------------------------------
00258     // Velocity fields set to zero in external domains
00259     //-------------------------------------------------
00260 
00261     loggam.annule(nzm1) ;       // zero in the ZEC only
00262 
00263     wit_w.annule(nzm1) ;        // zero outside the star     
00264 
00265     u_euler.annule(nzm1) ;  // zero outside the star     
00266 
00267 
00268     if (loggam.get_etat() != ETATZERO) {
00269     (loggam.set()).set_dzpuis(0) ; 
00270     }
00271     
00272     //################
00273     // verification: test on gam_euler
00274 
00275     // if (irrotational) {
00276     
00277     // Tenseur gam_test = 1. / sqrt( 1 - unsurc2 * sprod(u_euler, u_euler) ) ;
00278     
00279     // cout << "Etoile_bin::hydro_euler: test of gam_euler : " << endl ; 
00280     // cout << "  maximum error : " << endl ;
00281     // cout << max(gam_test() - gam_euler()) << endl ; 
00282     //cout << "  relative error : " << endl ; 
00283     // cout << diffrel(gam_test(), gam_euler()) << endl ; 
00284     
00285     // }
00286         
00287     //##################
00288 
00289 
00290     //### Test
00291 
00292     if (!irrotational) {
00293 
00294     //  Tenseur diff = gam - 1 ; 
00295     //  cout << "Etoile_bin::hydro_euler: rigid motion: difference gam <-> 1 : " 
00296     //       << endl ; 
00297     //  cout << norme(diff()) / norme(gam()) << endl ; 
00298     //
00299     //  cout << "Etoile_bin::hydro_euler: rigid motion: difference wit_w <-> 0 : " 
00300     //       << endl ; 
00301     //  cout << "Component x : " << endl << norme(wit_w(0)) << endl ; 
00302     //  cout << "Component y : " << endl << norme(wit_w(1)) << endl ; 
00303     //  cout << "Component z : " << endl << norme(wit_w(2)) << endl ; 
00304 
00305 //####
00306     gam = 1 ; 
00307     loggam = 0 ; 
00308     wit_w = 0 ; 
00309     }
00310     
00311     // The derived quantities are obsolete
00312     // -----------------------------------
00313     
00314     del_deriv() ;                
00315     
00316     
00317 }

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