et_bin_vel_pot.C

00001 /*
00002  * Method of class Etoile_bin to compute the velocity scalar potential $\psi$
00003  * by solving the continuity equation.
00004  *
00005  * (see file etoile.h for documentation).
00006  *
00007  */
00008 
00009 /*
00010  *   Copyright (c) 2000-2001 Eric Gourgoulhon
00011  *
00012  *   This file is part of LORENE.
00013  *
00014  *   LORENE is free software; you can redistribute it and/or modify
00015  *   it under the terms of the GNU General Public License as published by
00016  *   the Free Software Foundation; either version 2 of the License, or
00017  *   (at your option) any later version.
00018  *
00019  *   LORENE is distributed in the hope that it will be useful,
00020  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00021  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00022  *   GNU General Public License for more details.
00023  *
00024  *   You should have received a copy of the GNU General Public License
00025  *   along with LORENE; if not, write to the Free Software
00026  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00027  *
00028  */
00029 
00030 
00031 char et_bin_vel_pot_C[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_vel_pot.C,v 1.14 2007/10/18 14:26:43 e_gourgoulhon Exp $" ;
00032 
00033 /*
00034  * $Id: et_bin_vel_pot.C,v 1.14 2007/10/18 14:26:43 e_gourgoulhon Exp $
00035  * $Log: et_bin_vel_pot.C,v $
00036  * Revision 1.14  2007/10/18 14:26:43  e_gourgoulhon
00037  * Changed the call to Eos::der_nbar_ent in order to allow for MEos type
00038  * of equation of state.
00039  *
00040  * Revision 1.13  2007/10/16 21:56:26  e_gourgoulhon
00041  * Can deal with more than one domain into the star,
00042  * thanks to the new function Map_radial::poisson_compact.
00043  *
00044  * Revision 1.12  2005/10/18 13:12:33  p_grandclement
00045  * update of the mixted binary codes
00046  *
00047  * Revision 1.11  2004/05/25 15:38:38  f_limousin
00048  * Minor modifs.
00049  *
00050  * Revision 1.10  2004/05/10 10:17:27  f_limousin
00051  * Add a new member ssjm1_psi of class Etoile for the resolution of the
00052  * oisson_interne equation
00053  *
00054  * Revision 1.9  2004/04/19 11:26:17  f_limousin
00055  * Add a new function Etoile_bin::velocity_potential( , , , ) for the
00056  * case of strange stars
00057  *
00058  * Revision 1.8  2004/04/08 17:02:00  f_limousin
00059  * Modif to avoid an error in the compilation
00060  *
00061  * Revision 1.7  2004/04/08 16:52:58  f_limousin
00062  * Minor change
00063  *
00064  * Revision 1.6  2004/04/08 16:36:36  f_limousin
00065  * Implement the resolution of the continuity equation for strange
00066  * stars.
00067  *
00068  * Revision 1.5  2003/10/24 11:43:57  e_gourgoulhon
00069  * beta is now computed as ln(AN) in the case beta_auto
00070  * is undefined (for instance, if the companion is black hole).
00071  *
00072  * Revision 1.4  2003/01/17 13:38:56  f_limousin
00073  * Add comments
00074  *
00075  * Revision 1.3  2003/01/13 15:31:50  e_gourgoulhon
00076  * Suppressed the desaliasing
00077  *  (did not worked due to missing basis in ylm).
00078  *
00079  * Revision 1.2  2002/12/10 14:44:21  k_taniguchi
00080  * Change the multiplication "*" to "%"
00081  *   and flat_scalar_prod to flat_scalar_prod_desal.
00082  *
00083  * Revision 1.1.1.1  2001/11/20 15:19:28  e_gourgoulhon
00084  * LORENE
00085  *
00086  * Revision 2.9  2001/02/23  15:18:59  eric
00087  * Modification du calcul de zeta_h pour eviter division par zero
00088  *   dans les domaines externes a l'etoile.
00089  *
00090  * Revision 2.8  2001/02/07  09:47:42  eric
00091  * zeta_h est desormais donne par Eos::der_nbar_ent.
00092  *
00093  * Revision 2.7  2000/12/22  13:10:03  eric
00094  * Prolongement C^1 de dpsi en dehors de l'etoile.
00095  *
00096  * Revision 2.6  2000/03/22  12:56:44  eric
00097  * Nouveau prototype d'Etoile_bin::velocity_potential : l'erreur est
00098  * retournee en double.
00099  *
00100  * Revision 2.5  2000/02/25  17:35:29  eric
00101  * Annulation de la source dans les zones externes avant l'appel a
00102  * poisson_compact.
00103  *
00104  * Revision 2.4  2000/02/22  11:42:55  eric
00105  * Test resolution de l'equation.
00106  *
00107  * Revision 2.3  2000/02/22  10:42:25  eric
00108  * Correction erreur dans les termes sources: multiplication par unsurc2 de
00109  *  termes relativistes.
00110  *
00111  * Revision 2.2  2000/02/21  15:05:50  eric
00112  * Traitement du cas psi0 = 0 .
00113  *
00114  * Revision 2.1  2000/02/21  13:59:39  eric
00115  * Remplacement du membre psi par psi0.
00116  * Modif calcul de d_psi a la fin.
00117  *
00118  * Revision 2.0  2000/02/17  18:50:44  eric
00119  * *** empty log message ***
00120  *
00121  *
00122  * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_vel_pot.C,v 1.14 2007/10/18 14:26:43 e_gourgoulhon Exp $
00123  *
00124  */
00125 
00126 // Headers Lorene
00127 #include "scalar.h"
00128 #include "metrique.h" 
00129 #include "etoile.h"
00130 #include "eos.h"
00131 #include "param.h"
00132 #include "et_bin_nsbh.h"
00133 #include "utilitaires.h"
00134 
00135 // Local prototype
00136 Cmp raccord_c1(const Cmp& uu, int l1) ; 
00137 
00138 double Etoile_bin::velocity_potential(int mermax, double precis, double relax) {
00139   
00140   // Which star is that ?
00141   const Et_bin_nsbh* pnsbh = dynamic_cast<const Et_bin_nsbh*>(this) ;
00142   
00143   if (eos.identify() == 5 || eos.identify() == 4 || 
00144       eos.identify() == 3) {
00145     
00146     // Routine used for binary strange stars.
00147 
00148     int nzm1 = mp.get_mg()->get_nzone() - 1 ;    
00149 
00150     //----------------------------------
00151     // Specific relativistic enthalpy           ---> hhh
00152     //----------------------------------
00153    
00154     Tenseur hhh = exp(unsurc2 * ent) ;  // = 1 at the Newtonian limit
00155     hhh.set_std_base() ;
00156 
00157     //----------------------------------------------
00158     //  Computation of W^i = - A^2 h Gamma_n B^i/N
00159     // See Eq (62) from Gourgoulhon et al. (2001)
00160     //----------------------------------------------
00161 
00162     Tenseur www = - a_car * hhh * gam_euler * bsn ; 
00163     
00164     www.change_triad( mp.get_bvect_cart() ) ;   // components on the mapping
00165                         // Cartesian basis
00166     
00167     //-------------------------------------------------
00168     // Constant value of W^i at the center of the star
00169     //-------------------------------------------------
00170     
00171     Tenseur v_orb(mp, 1, CON, mp.get_bvect_cart()) ; 
00172     
00173     v_orb.set_etat_qcq() ; 
00174     for (int i=0; i<3; i++) {
00175     v_orb.set(i) = www(i)(0, 0, 0, 0) ; 
00176     }
00177 
00178     v_orb.annule(nzm1, nzm1) ;  // set to zero in the ZEC
00179 
00180 
00181     v_orb.set_triad( *(www.get_triad()) ) ;  
00182     v_orb.set_std_base() ;
00183 
00184               
00185     //-------------------------------------------------
00186     // Source and coefficients a,b for poisson_compact (idenpendent from psi0)
00187     //-------------------------------------------------
00188     
00189     Cmp dndh_log = eos.der_nbar_ent(ent(), nzet) ; 
00190 
00191     // In order to avoid any division by zero in the computation of zeta_h
00192     //  the value of dndh_log is set to 1 in the external domains:
00193     for (int l=nzet; l <= nzm1; l++) {
00194     dndh_log.set(l) = 1 ; 
00195     }
00196     
00197     double erreur ;             
00198    
00199     Tenseur zeta_h( ent() / dndh_log ) ;
00200     zeta_h.set_std_base() ;
00201     
00202     Scalar zeta_h_scalar (zeta_h()) ;
00203     zeta_h_scalar.set_outer_boundary(0, (ent() / dndh_log)(0,0,0,0)) ;
00204     for (int l=1; l<=nzm1; l++)
00205     zeta_h_scalar.set_domain(l) = 1 ;
00206     
00207     Cmp zeta_h_cmp (zeta_h_scalar) ;
00208     zeta_h.set() = zeta_h_cmp ;
00209     zeta_h.set_std_base() ;
00210 
00211     
00212 
00213     Tenseur beta(mp) ;
00214     
00215     if (pnsbh!=0x0) {
00216     beta = log( sqrt(a_car) * nnn ) ; 
00217     beta.set_std_base() ;
00218     }
00219     else {
00220     beta = beta_auto + beta_comp ; 
00221     }
00222 
00223     Tenseur tmp_zeta = 1 - unsurc2 * zeta_h ;
00224     tmp_zeta.set_std_base() ;
00225     
00226     Tenseur bb = tmp_zeta * ent.gradient_spher()
00227     + unsurc2 * zeta_h * beta.gradient_spher() ;
00228     
00229     Tenseur entmb = ent - beta ; 
00230     
00231     Tenseur grad_ent (ent.gradient()) ;
00232     grad_ent.change_triad(mp.get_bvect_spher()) ;
00233 
00234     // Source for the poisson equation 
00235     // See Eq (63) from Gourgoulhon et al. (2001)
00236     Tenseur source = flat_scalar_prod( www - v_orb, ent.gradient() )
00237     + unsurc2 * zeta_h * (
00238         flat_scalar_prod( v_orb, entmb.gradient() )
00239         + flat_scalar_prod( www, gam_euler.gradient() )
00240         / gam_euler ) ; 
00241 
00242     for (int l=1; l<=nzm1; l++)
00243     source.set().annule(l) ;
00244    
00245     source = (source - flat_scalar_prod(bb, psi0.gradient_spher())) 
00246     / zeta_h ;
00247     source.annule(nzet, nzm1) ; 
00248 
00249     Param par ; 
00250     int niter ; 
00251     par.add_int(mermax) ; 
00252     par.add_double(precis, 0) ; 
00253     par.add_double(relax, 1) ; 
00254     par.add_int_mod(niter) ; 
00255 
00256     par.add_cmp_mod(ssjm1_psi, 0) ;
00257     
00258     if (psi0.get_etat() == ETATZERO) {
00259     psi0.set_etat_qcq() ; 
00260     psi0.set() = 0 ; 
00261     }
00262 
00263     int nr = mp.get_mg()->get_nr(0);
00264     int nt = mp.get_mg()->get_nt(0);
00265     int np = mp.get_mg()->get_np(0);
00266 
00267     cout << "nr = " << nr << "   nt = " << nt << "   np = " << np << endl ;
00268 
00269     cout << "psi0" << endl << norme(psi0()/(nr*nt*np)) << endl ;
00270     cout << "d(psi)/dr" << endl << norme(psi0.set().dsdr()/(nr*nt*np)) << endl ;
00271 
00272     Valeur lim(mp.get_mg()->get_angu()) ;
00273     lim.annule_hard() ;
00274 
00275     Tenseur normal (mp, 1, CON, mp.get_bvect_cart()) ;
00276     Tenseur normal2 (mp, 1, COV, mp.get_bvect_cart()) ;
00277     normal.set_etat_qcq() ;
00278     normal2.set_etat_qcq() ;
00279 
00280     const Coord& rr0 = mp.r ;
00281     Tenseur rr(mp) ;
00282     rr.set_etat_qcq() ;
00283     rr.set() = rr0 ;
00284     rr.set_std_base() ;
00285 
00286     Tenseur_sym plat(mp, 2, COV, mp.get_bvect_cart() ) ;
00287     plat.set_etat_qcq() ;
00288     for (int i=0; i<3; i++) {
00289     for (int j=0; j<i; j++) {
00290         plat.set(i,j) = 0 ;
00291     }
00292     plat.set(i,i) = 1 ;
00293     }
00294     plat.set_std_base() ;
00295     
00296     Metrique flat(plat, true) ; 
00297     Tenseur dcov_r = rr.derive_cov(flat) ;
00298 
00299 
00300     for (int i=0; i<3; i++) {
00301     normal.set(i) = dcov_r(i) ;
00302     normal2.set(i) = dcov_r(i) ;
00303     }
00304 
00305     normal.change_triad(mp.get_bvect_spher()) ;
00306     normal2.change_triad(mp.get_bvect_spher()) ;
00307  
00308 
00309 
00310     Tenseur bsn0 (bsn) ;
00311     bsn0.change_triad(mp.get_bvect_cart()) ;
00312     Tenseur aa (mp, 1, CON, mp.get_bvect_cart()) ;
00313     aa = - v_orb - a_car * gam_euler * hhh * bsn0 ;
00314     aa.change_triad(mp.get_bvect_spher()) ;
00315     
00316 
00317     Tenseur dcov_psi = psi0.derive_cov(flat) ;
00318     dcov_psi.change_triad(mp.get_bvect_spher()) ;
00319 
00320     Cmp limite (mp) ;
00321     limite =  ( - dcov_psi(1) * normal(1) - dcov_psi(2) * normal(2)
00322         + contract(aa, 0, normal2, 0)()) 
00323     /normal(0) ;
00324     
00325     for (int j=0; j<nt; j++)
00326     for (int k=0; k<np; k++)
00327         lim.set(0, k, j, 0) = limite(0, k, j, nr-1) ;
00328     
00329 //  cout << "lim" << endl << lim << endl ;
00330 
00331     lim.std_base_scal() ;
00332     Cmp resu (psi0()) ;
00333     source().poisson_neumann_interne(lim, par, resu) ;
00334     psi0 = resu ;
00335 
00336 /*
00337     resu.va.ylm() ;
00338     Scalar psi00(resu) ;
00339     psi00.spectral_display("psi00") ;
00340 
00341     cout << "value of d(psi)/dr at the surface after poisson" << endl ;
00342     for (int j=0; j<nt; j++)
00343     for (int k=0; k<np; k++)
00344     cout << "j = " << j << " ; k = " << k << " : " << 
00345     psi0.set().dsdr()(0, k, j, nr-1) << endl ;
00346 */
00347     for (int l=1; l<=nzm1; l++)
00348     psi0.set().annule(l) ;
00349     
00350 
00351     //---------------------------------------------------
00352     // Check of the solution  
00353     //---------------------------------------------------
00354     
00355     Cmp laplacien_psi0 =  psi0().laplacien() ; 
00356     
00357     erreur = diffrel(laplacien_psi0, source())(0) ; 
00358 
00359     cout << "Check of the resolution of the continuity equation for strange stars: " 
00360      << endl ; 
00361     cout << "norme(source) : " << norme(source())(0) << endl 
00362      << "Error in the solution : " << erreur << endl ; 
00363     
00364     //--------------------------------
00365     // Computation of grad(psi)
00366     //--------------------------------
00367     
00368     // The computation is done component by component because psi0.gradient()
00369     // is a covariant vector, whereas v_orb is a contravariant one. 
00370     
00371     d_psi.set_etat_qcq() ; 
00372     
00373     for (int i=0; i<3; i++) {
00374     d_psi.set(i) = (psi0.gradient())(i) + v_orb(i) ; 
00375     }
00376 
00377     d_psi.set_triad(  *(v_orb.get_triad()) ) ; 
00378    
00379     // C^1 continuation of d_psi outside the star
00380     //  (to ensure a smooth enthalpy field accross the stellar surface)
00381     // ----------------------------------------------------------------
00382     
00383     d_psi.annule(nzet, nzm1) ;    
00384     for (int i=0; i<3; i++) {
00385     d_psi.set(i) = raccord_c1(d_psi(i), nzet) ; 
00386     }
00387     
00388     assert( d_psi.get_triad() == &(mp.get_bvect_cart()) ) ; 
00389 
00390     d_psi.change_triad(ref_triad) ; 
00391     
00392     return erreur ; 
00393  
00394 
00395   }  // End of strange stars case
00396 
00397 //=============================================================================
00398 
00399   else {
00400     
00401     int nzm1 = mp.get_mg()->get_nzone() - 1 ;    
00402 
00403     //----------------------------------
00404     // Specific relativistic enthalpy           ---> hhh
00405     //----------------------------------
00406     
00407     Tenseur hhh = exp(unsurc2 * ent) ;  // = 1 at the Newtonian limit
00408     hhh.set_std_base() ;
00409 
00410     //----------------------------------------------
00411     //  Computation of W^i = - A^2 h Gamma_n B^i/N
00412     // See Eq (62) from Gourgoulhon et al. (2001)
00413     //----------------------------------------------
00414 
00415     Tenseur www = - a_car * hhh * gam_euler * bsn ; 
00416     
00417     www.change_triad( mp.get_bvect_cart() ) ;   // components on the mapping
00418     // Cartesian basis
00419     
00420     //-------------------------------------------------
00421     // Constant value of W^i at the center of the star
00422     //-------------------------------------------------
00423     
00424     Tenseur v_orb(mp, 1, CON, mp.get_bvect_cart()) ; 
00425     
00426     v_orb.set_etat_qcq() ; 
00427     for (int i=0; i<3; i++) {
00428     v_orb.set(i) = www(i)(0, 0, 0, 0) ; 
00429     }
00430 
00431     v_orb.set_triad( *(www.get_triad()) ) ;  
00432     v_orb.set_std_base() ;
00433 
00434               
00435     //-------------------------------------------------
00436     // Source and coefficients a,b for poisson_compact (idenpendent from psi0)
00437     //-------------------------------------------------
00438     
00439     Cmp dndh_log(mp) ; 
00440     dndh_log = 0 ; 
00441 
00442     for (int l=0; l<nzet; l++) {
00443 
00444         Param par ;       // Paramater for multi-domain equation of state
00445         par.add_int(l) ;
00446 
00447         dndh_log = dndh_log +  eos.der_nbar_ent(ent(), 1, l, &par) ;
00448 
00449     }
00450     
00451     // Cmp dndh_log = eos.der_nbar_ent(ent(), nzet) ; 
00452 
00453     // In order to avoid any division by zero in the computation of zeta_h
00454     //  the value of dndh_log is set to 1 in the external domains:
00455     for (int l=nzet; l <= nzm1; l++) {
00456     dndh_log.set(l) = 1 ; 
00457     }
00458     
00459     Tenseur zeta_h( ent() / dndh_log ) ;
00460     zeta_h.set_std_base() ;
00461     
00462     Tenseur beta(mp) ; 
00463     
00464     if (pnsbh!=0x0) {
00465     beta = log( sqrt(a_car) * nnn ) ;
00466     beta.set_std_base() ;
00467     }
00468     else {
00469     beta = beta_auto + beta_comp ; 
00470     }
00471     
00472     Tenseur tmp_zeta = 1 - unsurc2 * zeta_h ;
00473     tmp_zeta.set_std_base() ;
00474     
00475     Tenseur bb = tmp_zeta * ent.gradient_spher()
00476     + unsurc2 * zeta_h * beta.gradient_spher() ;
00477     
00478     Tenseur entmb = ent - beta ; 
00479     
00480     // See Eq (63) from Gourgoulhon et al. (2001)
00481     Tenseur source = flat_scalar_prod( www - v_orb, ent.gradient() )
00482     + unsurc2 * zeta_h * (
00483         flat_scalar_prod( v_orb, entmb.gradient() )
00484         + flat_scalar_prod( www, gam_euler.gradient() )
00485         / gam_euler ) ; 
00486     
00487     
00488     source.annule(nzet, nzm1) ; 
00489     
00490     //---------------------------------------------------
00491     // Resolution by means of Map_radial::poisson_compact 
00492     //---------------------------------------------------
00493     
00494     Param par ; 
00495     int niter ; 
00496     par.add_int(mermax) ; 
00497     par.add_double(precis, 0) ; 
00498     par.add_double(relax, 1) ; 
00499     par.add_int_mod(niter) ; 
00500     
00501     
00502     if (psi0.get_etat() == ETATZERO) {
00503     psi0.set_etat_qcq() ; 
00504     psi0.set() = 0 ; 
00505     }
00506     
00507     source.set().va.ylm() ;
00508     
00509     mp.poisson_compact(nzet, source(), zeta_h(), bb, par, psi0.set() ) ;
00510 
00511     //---------------------------------------------------
00512     // Check of the solution  
00513     //---------------------------------------------------
00514     
00515     Tenseur bb_dpsi0 = flat_scalar_prod( bb, psi0.gradient_spher() ) ;
00516     
00517     Cmp oper = zeta_h() * psi0().laplacien() + bb_dpsi0() ; 
00518     
00519     source.set().va.ylm_i() ;
00520     
00521     cout << "Check of the resolution of the continuity equation : "  << endl ; 
00522     Tbl terr = diffrel(oper, source()) ;
00523     double erreur = 0 ; 
00524     for (int l=0; l<nzet; l++) { 
00525         double err = terr(l) ; 
00526         cout << " domain " << l << " : norme(source) : " << norme(source())(l) 
00527             << "    relative error : " << err << endl ;
00528         if (err > erreur) erreur = err ; 
00529     } 
00530     // arrete() ; 
00531     
00532    //--------------------------------
00533    // Computation of grad(psi)
00534    //--------------------------------
00535     
00536    // The computation is done component by component because psi0.gradient()
00537    // is a covariant vector, whereas v_orb is a contravariant one. 
00538     
00539     d_psi.set_etat_qcq() ; 
00540     
00541     for (int i=0; i<3; i++) {
00542     d_psi.set(i) = (psi0.gradient())(i) + v_orb(i) ; 
00543     }
00544 
00545     d_psi.set_triad(  *(v_orb.get_triad()) ) ; 
00546    
00547    // C^1 continuation of d_psi outside the star
00548    //  (to ensure a smooth enthalpy field accross the stellar surface)
00549    // ----------------------------------------------------------------
00550     
00551     d_psi.annule(nzet, nzm1) ;    
00552     for (int i=0; i<3; i++) {
00553     d_psi.set(i) = raccord_c1(d_psi(i), nzet) ; 
00554     }
00555     
00556     assert( d_psi.get_triad() == &(mp.get_bvect_cart()) ) ; 
00557 
00558     d_psi.change_triad(ref_triad) ; 
00559     
00560     return erreur ; 
00561  
00562 
00563   }
00564 }

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