binary_global.C

00001 /*
00002  * Methods of class Binary to compute global quantities
00003  *
00004  * (see file binary.h for documentation)
00005  */
00006 
00007 /*
00008  *   Copyright (c) 2004 Francois Limousin
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 binary_global_C[] = "$Header: /cvsroot/Lorene/C++/Source/Binary/binary_global.C,v 1.15 2006/08/01 14:26:50 f_limousin Exp $" ;
00030 
00031 /*
00032  * $Id: binary_global.C,v 1.15 2006/08/01 14:26:50 f_limousin Exp $
00033  * $Log: binary_global.C,v $
00034  * Revision 1.15  2006/08/01 14:26:50  f_limousin
00035  * Small changes
00036  *
00037  * Revision 1.14  2006/04/11 14:25:15  f_limousin
00038  * New version of the code : improvement of the computation of some
00039  * critical sources, estimation of the dirac gauge, helical symmetry...
00040  *
00041  * Revision 1.13  2005/09/18 13:13:41  f_limousin
00042  * Extension of vphi in the compactified domain for the computation
00043  * of J_ADM by a volume integral.
00044  *
00045  * Revision 1.12  2005/09/15 14:41:04  e_gourgoulhon
00046  * The total angular momentum is now computed via a volume integral.
00047  *
00048  * Revision 1.11  2005/09/13 19:38:31  f_limousin
00049  * Reintroduction of the resolution of the equations in cartesian coordinates.
00050  *
00051  * Revision 1.10  2005/04/08 12:36:45  f_limousin
00052  * Just to avoid warnings...
00053  *
00054  * Revision 1.9  2005/02/17 17:35:00  f_limousin
00055  * Change the name of some quantities to be consistent with other classes
00056  * (for instance nnn is changed to nn, shift to beta, beta to lnq...)
00057  *
00058  * Revision 1.8  2004/07/21 11:46:24  f_limousin
00059  * Add function mass_adm_vol() to compute the ADM mass of the system
00060  * with a volume integral instead of a surface one.
00061  *
00062  * Revision 1.7  2004/05/25 14:25:53  f_limousin
00063  * Add the virial theorem for conformally flat configurations.
00064  *
00065  * Revision 1.6  2004/03/31 12:44:54  f_limousin
00066  * Minor modifs.
00067  *
00068  * Revision 1.5  2004/03/25 10:29:01  j_novak
00069  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
00070  *
00071  * Revision 1.4  2004/02/27 10:25:30  f_limousin
00072  * Modif. to avoid an error in compilation.
00073  *
00074  * Revision 1.3  2004/02/27 10:03:04  f_limousin
00075  * The computation of mass_adm() and mass_komar() is now OK !
00076  *
00077  * Revision 1.2  2004/01/20 15:21:36  f_limousin
00078  * First version
00079  *
00080  *
00081  * $Header: /cvsroot/Lorene/C++/Source/Binary/binary_global.C,v 1.15 2006/08/01 14:26:50 f_limousin Exp $
00082  *
00083  */
00084 
00085 
00086 // Headers C
00087 #include "math.h"
00088 
00089 // Headers Lorene
00090 #include "nbr_spx.h"
00091 #include "binary.h"
00092 #include "unites.h"
00093 #include "metric.h"
00094 
00095             //---------------------------------//
00096             //      ADM mass           //
00097             //---------------------------------//
00098 
00099 double Binary::mass_adm() const {
00100     
00101   using namespace Unites ;
00102   if (p_mass_adm == 0x0) {      // a new computation is requireed
00103     
00104     p_mass_adm = new double ; 
00105         
00106     *p_mass_adm = 0 ; 
00107     
00108     const Map_af map0 (et[0]->get_mp()) ;
00109     const Metric& flat = (et[0]->get_flat()) ;
00110 
00111     Vector dpsi(0.5*(et[0]->get_lnq() - 
00112                  et[0]->get_logn()).derive_cov(flat)) ;
00113   
00114     Vector ww (0.125*(contract(et[0]->get_hij().derive_cov(flat), 1, 2) 
00115                   - (et[0]->get_hij().trace(flat)).derive_con(flat))) ;
00116     
00117     dpsi.change_triad(map0.get_bvect_spher()) ;
00118     ww.change_triad(map0.get_bvect_spher()) ;
00119 
00120     // ww = 0 in Dirac gauge (Eq 174 of BGGN)
00121     Scalar integrand (dpsi(1) + 0*ww(1)) ;
00122 
00123     *p_mass_adm = map0.integrale_surface_infini (integrand) / (-qpig/2.) ;
00124     
00125     }   // End of the case where a new computation was necessary
00126     
00127     return *p_mass_adm ; 
00128     
00129 }
00130 
00131 double Binary::mass_adm_vol() const {
00132 
00133   using namespace Unites ;
00134 
00135   double massadm ;
00136   massadm = 0. ;
00137 
00138   for (int i=0; i<=1; i++) {        // loop on the stars
00139 
00140     // Declaration of all fields
00141       const Scalar& psi4 = et[i]->get_psi4() ;
00142       Scalar psi (pow(psi4, 0.25)) ;
00143       psi.std_spectral_base() ;
00144       const Scalar& ener_euler = et[i]->get_ener_euler() ;
00145       const Scalar& kcar_auto = et[i]->get_kcar_auto() ;
00146       const Scalar& kcar_comp = et[i]->get_kcar_comp() ;
00147       const Metric& gtilde = et[i]->get_gtilde() ;
00148       const Metric& flat = et[i]->get_flat() ;
00149       const Sym_tensor& hij = et[i]->get_hij() ;
00150       const Sym_tensor& hij_auto = et[i]->get_hij_auto() ;
00151       const Vector& dcov_logn = et[i]->get_dcov_logn() ;
00152       const Vector& dcov_phi = et[i]->get_dcov_phi() ;
00153       const Vector& dcov_lnq = 2*dcov_phi + dcov_logn ;
00154       const Scalar& lnq_auto = et[i]->get_lnq_auto() ;
00155       const Scalar& logn_auto = et[i]->get_logn_auto() ;
00156       const Scalar& phi_auto = 0.5 * (lnq_auto - logn_auto) ;
00157 
00158       const Tensor& dcov_hij_auto = hij_auto.derive_cov(flat) ;
00159       const Tensor& dcov_gtilde = gtilde.cov().derive_cov(flat) ;
00160       const Tensor& dcov_phi_auto = phi_auto.derive_cov(flat) ;
00161       const Tensor& dcov_logn_auto = logn_auto.derive_cov(flat) ;
00162       const Tensor& dcov_lnq_auto = lnq_auto.derive_cov(flat) ;
00163       Tensor dcovdcov_lnq_auto = lnq_auto.derive_cov(flat).derive_cov(flat) ;
00164       dcovdcov_lnq_auto.inc_dzpuis() ;
00165       Tensor dcovdcov_logn_auto = logn_auto.derive_cov(flat).derive_cov(flat) ;
00166       dcovdcov_logn_auto.inc_dzpuis() ;
00167  
00168       // Source in IWM approximation 
00169       Scalar source = - psi4 % (qpig*ener_euler + (kcar_auto + kcar_comp)/4.) 
00170     - 0*2*contract(contract(gtilde.con(), 0, dcov_phi, 0), 
00171                0, dcov_phi_auto, 0, true) ;
00172       
00173       // Source = 0 in IWM 
00174       source += 4*contract(hij, 0, 1, dcov_logn * dcov_phi_auto, 0, 1) +
00175     2*contract(hij, 0, 1, dcov_phi * dcov_phi_auto, 0, 1) +
00176     0.0625 * contract(gtilde.con(), 0, 1, contract(
00177                dcov_hij_auto, 0, 1, dcov_gtilde, 0, 1), 0, 1) - 
00178            0.125 * contract(gtilde.con(), 0, 1, contract(dcov_hij_auto, 
00179                   0, 1, dcov_gtilde, 0, 2), 0, 1) -
00180      contract(hij,0,1,dcovdcov_lnq_auto + dcov_lnq_auto*dcov_lnq,0,1) +
00181      contract(hij,0,1,dcovdcov_logn_auto + dcov_logn_auto*dcov_logn,0,1) ;
00182 
00183       source = source * psi ;
00184 
00185       source.std_spectral_base() ;
00186 
00187       massadm += - source.integrale()/qpig ;
00188   }
00189 
00190   return massadm ;
00191 }
00192 
00193             //---------------------------------//
00194             //      Komar mass         //
00195             //---------------------------------//
00196 
00197 double Binary::mass_kom() const {
00198     
00199   using namespace Unites ;
00200 
00201   if (p_mass_kom == 0x0) {      // a new computation is requireed
00202     
00203     p_mass_kom = new double ; 
00204       
00205     *p_mass_kom = 0 ; 
00206     
00207     const Tensor& logn = et[0]->get_logn() ;
00208     const Metric& flat = (et[0]->get_flat()) ;
00209     const Sym_tensor&  hij = (et[0]->get_hij()) ;
00210     Map_af map0 (et[0]->get_mp()) ; 
00211     
00212     Vector vect = logn.derive_con(flat) + 
00213                        contract(hij, 1, logn.derive_cov(flat), 0) ;
00214     vect.change_triad(map0.get_bvect_spher()) ;
00215     Scalar integrant (vect(1)) ;
00216     
00217     *p_mass_kom = map0.integrale_surface_infini (integrant) / qpig ;
00218     
00219   } // End of the case where a new computation was necessary
00220     
00221   return *p_mass_kom ; 
00222     
00223 }
00224 
00225 double Binary::mass_kom_vol() const {
00226     
00227   using namespace Unites ;
00228 
00229   double masskom ;
00230   masskom = 0. ;
00231 
00232   for (int i=0; i<=1; i++) {        // loop on the stars
00233 
00234      // Declaration of all fields
00235       const Scalar& psi4 = et[i]->get_psi4() ;
00236       const Scalar& ener_euler = et[i]->get_ener_euler() ;
00237       const Scalar& s_euler = et[i]->get_s_euler() ;
00238       const Scalar& kcar_auto = et[i]->get_kcar_auto() ;
00239       const Scalar& kcar_comp = et[i]->get_kcar_comp() ;
00240       const Metric& gtilde = et[i]->get_gtilde() ;
00241       const Metric& flat = et[i]->get_flat() ;
00242       const Sym_tensor& hij = et[i]->get_hij() ;
00243       const Scalar& logn = et[i]->get_logn_auto() + et[i]->get_logn_comp() ;
00244       const Scalar& logn_auto = et[i]->get_logn_auto() ;
00245       Scalar nn = exp(logn) ;
00246       nn.std_spectral_base() ;
00247       
00248       const Tensor& dcov_logn_auto = logn_auto.derive_cov(flat) ;
00249       const Vector& dcov_logn = et[i]->get_dcov_logn() ;
00250       const Vector& dcon_logn = et[i]->get_dcon_logn() ;
00251       const Vector& dcov_phi = et[i]->get_dcov_phi() ;
00252       Tensor dcovdcov_logn_auto = (logn_auto.derive_cov(flat))
00253     .derive_cov(flat) ;
00254       dcovdcov_logn_auto.inc_dzpuis() ;
00255 
00256       Scalar source = qpig * psi4 % (ener_euler + s_euler) ;
00257       source += psi4 % (kcar_auto + kcar_comp) ;
00258       source += - 0*contract(dcov_logn_auto, 0, dcon_logn, 0, true) 
00259       - 2. * contract(contract(gtilde.con(), 0, dcov_phi, 0), 0, 
00260               dcov_logn_auto, 0, true) ;
00261       source += - contract(hij, 0, 1, dcovdcov_logn_auto + 
00262                dcov_logn_auto*dcov_logn, 0, 1) ;
00263 
00264       source = source / qpig * nn  ;
00265   
00266       source.std_spectral_base() ;
00267 
00268       masskom += source.integrale() ;
00269       
00270   }
00271 
00272   return masskom ;
00273 
00274 }
00275 
00276 
00277             //---------------------------------//
00278             //   Total angular momentum        //
00279             //---------------------------------//
00280 
00281 const Tbl& Binary::angu_mom() const {
00282 
00283   using namespace Unites ;
00284 
00285   /*
00286     if (p_angu_mom == 0x0) {        // a new computation is requireed
00287     
00288       p_angu_mom = new Tbl(3) ; 
00289       
00290       p_angu_mom->annule_hard() ;   // fills the double array with zeros
00291   
00292       const Sym_tensor& kij_auto = et[0]->get_tkij_auto() ;
00293       const Sym_tensor& kij_comp = et[0]->get_tkij_comp() ;
00294       const Tensor& psi4 = et[0]->get_psi4() ;
00295       const Map_af map0 (kij_auto.get_mp()) ;
00296 
00297       Sym_tensor kij = (kij_auto + kij_comp) / psi4 ;
00298       kij.change_triad(map0.get_bvect_cart()) ;
00299  
00300       // X component
00301       // -----------
00302 
00303       Vector vect_x(et[0]->get_mp(), CON, map0.get_bvect_cart()) ;      
00304        
00305       for (int i=1; i<=3; i++) {
00306 
00307       Scalar kij_1 = kij(3, i) ;
00308       Scalar kij_2 = kij(2, i) ;
00309           
00310       kij_1.mult_rsint() ;
00311       Valeur vtmp = kij_1.get_spectral_va().mult_sp() ;
00312       kij_1.set_spectral_va() = vtmp ; 
00313       
00314       kij_2.mult_r() ;
00315       vtmp = kij_2.get_spectral_va().mult_ct() ;
00316       kij_2.set_spectral_va() = vtmp ; 
00317 
00318       vect_x.set(i) = kij_1 - kij_2 ;
00319       }
00320  
00321       vect_x.change_triad(map0.get_bvect_spher()) ;
00322       Scalar integrant_x (vect_x(1)) ;
00323       
00324       p_angu_mom->set(0) = map0.integrale_surface_infini (integrant_x) 
00325                       / (8*M_PI) ;
00326       
00327       // Y component
00328       // -----------
00329       
00330       Vector vect_y(et[0]->get_mp(), CON, map0.get_bvect_cart()) ;      
00331    
00332       for (int i=1; i<=3; i++) {
00333 
00334       Scalar kij_1 = kij(1, i) ;
00335       Scalar kij_2 = kij(3, i) ;      
00336       
00337       kij_1.mult_r() ;
00338       Valeur vtmp = kij_1.get_spectral_va().mult_ct() ;
00339       kij_1.set_spectral_va() = vtmp ; 
00340       
00341       kij_2.mult_rsint() ;
00342       vtmp = kij_2.get_spectral_va().mult_cp() ;
00343       kij_2.set_spectral_va() = vtmp ; 
00344      
00345       vect_y.set(i) = kij_1 - kij_2 ;
00346       }
00347 
00348       vect_y.change_triad(map0.get_bvect_spher()) ;
00349       Scalar integrant_y (vect_y(1)) ;
00350       
00351       p_angu_mom->set(1) = map0.integrale_surface_infini (integrant_y) 
00352                       / (8*M_PI) ;
00353       
00354       // Z component
00355       // -----------
00356 
00357       Vector vect_z(et[0]->get_mp(), CON, map0.get_bvect_cart()) ;      
00358 
00359       for (int i=1; i<=3; i++) {
00360 
00361       Scalar kij_1 = kij(2, i) ;
00362       Scalar kij_2 = kij(1, i) ;          
00363 
00364       kij_1.mult_rsint() ;
00365       Valeur vtmp = kij_1.get_spectral_va().mult_cp() ;
00366       kij_1.set_spectral_va() = vtmp ; 
00367 
00368       kij_2.mult_rsint() ;
00369       vtmp =  kij_2.get_spectral_va().mult_sp() ;
00370       kij_2.set_spectral_va() = vtmp ;
00371 
00372       vect_z.set(i) = kij_1 - kij_2 ;
00373       }
00374        
00375       vect_z.change_triad(map0.get_bvect_spher()) ;
00376       Scalar integrant_z (vect_z(1)) ;
00377       
00378       p_angu_mom->set(2) = map0.integrale_surface_infini (integrant_z) 
00379                      ;// (8*M_PI) ;
00380       
00381       
00382     }   // End of the case where a new computation was necessary
00383   */
00384   
00385 
00386     /*  
00387   if (p_angu_mom == 0x0) {      // a new computation is requireed
00388     p_angu_mom = new Tbl(3) ; 
00389     p_angu_mom->annule_hard() ; // fills the double array with zeros
00390     p_angu_mom->set(0) = 0. ;
00391     p_angu_mom->set(1) = 0. ;
00392 
00393     // Alignement 
00394     double orientation_un = et[0]->get_mp().get_rot_phi() ;
00395     assert ((orientation_un==0) || (orientation_un==M_PI)) ;
00396     double orientation_deux = et[1]->get_mp().get_rot_phi() ;
00397     assert ((orientation_deux==0) || (orientation_deux==M_PI)) ;
00398     int same_orient = (orientation_un == orientation_deux) ? 1 : -1 ;
00399     
00400     // Construction of an auxiliar grid and mapping
00401     int nzones = et[0]->get_mp().get_mg()->get_nzone() ;
00402     double* bornes = new double [nzones+1] ;
00403     double courant = (et[0]->get_mp().get_ori_x()-et[0]->get_mp().get_ori_x())+1 ;
00404     for (int i=nzones-1 ; i>0 ; i--) {
00405       bornes[i] = courant ;
00406       courant /= 2. ;
00407     }
00408     bornes[0] = 0 ;
00409     bornes[nzones] = __infinity ;
00410     
00411     Map_af mapping (*(et[0]->get_mp().get_mg()), bornes) ;
00412     
00413     delete [] bornes ; 
00414     
00415     // Construction of k_total
00416     Sym_tensor k_total (mapping, CON, mapping.get_bvect_cart()) ;
00417     
00418     Vector shift_un (mapping, CON, mapping.get_bvect_cart()) ;
00419     Vector shift_deux (mapping, CON, mapping.get_bvect_cart()) ;
00420     
00421     Vector beta_un (et[0]->get_beta_auto()) ;
00422     Vector beta_deux (et[1]->get_beta_auto()) ;
00423     beta_un.change_triad(et[0]->get_mp().get_bvect_cart()) ;
00424     beta_deux.change_triad(et[1]->get_mp().get_bvect_cart()) ;
00425     beta_un.std_spectral_base() ;
00426     beta_deux.std_spectral_base() ;
00427     
00428     shift_un.set(1).import(beta_un(1)) ;
00429     shift_un.set(2).import(beta_un(2)) ;
00430     shift_un.set(3).import(beta_un(3)) ;
00431  
00432     shift_deux.set(1).import(same_orient*beta_deux(1)) ;
00433     shift_deux.set(2).import(same_orient*beta_deux(2)) ;
00434     shift_deux.set(3).import(beta_deux(3)) ;
00435     
00436     Vector shift_tot (shift_un+shift_deux) ;
00437     shift_tot.std_spectral_base() ;
00438     shift_tot.annule(0, nzones-2) ;
00439     
00440     
00441     // Substract the residuals
00442     shift_tot.inc_dzpuis(2) ;
00443     shift_tot.dec_dzpuis(2) ;
00444     
00445     
00446     Sym_tensor temp_gamt (et[0]->get_gtilde().cov()) ;
00447     temp_gamt.change_triad(mapping.get_bvect_cart()) ;
00448     Metric gamt_cart (temp_gamt) ;
00449     
00450     k_total = shift_tot.ope_killing_conf(gamt_cart) / 2. ;
00451     
00452     for (int lig=1 ; lig<=3 ; lig++)
00453     for (int col=lig ; col<=3 ; col++)
00454       k_total.set(lig, col).mult_r_ced() ;
00455     
00456     Vector vecteur_un (mapping, CON, mapping.get_bvect_cart()) ;
00457     for (int i=1 ; i<=3 ; i++)
00458       vecteur_un.set(i) = k_total(1, i) ;
00459     vecteur_un.change_triad (mapping.get_bvect_spher()) ;
00460     Scalar integrant_un (vecteur_un(1)) ;
00461     
00462     Vector vecteur_deux (mapping, CON, mapping.get_bvect_cart()) ;
00463     for (int i=1 ; i<=3 ; i++)
00464       vecteur_deux.set(i) = k_total(2, i) ;
00465     vecteur_deux.change_triad (mapping.get_bvect_spher()) ;
00466     Scalar integrant_deux (vecteur_deux(1)) ;
00467     
00468     // Multiplication by y and x :
00469     integrant_un.set_spectral_va() = integrant_un.get_spectral_va()
00470       .mult_st() ;
00471     integrant_un.set_spectral_va() = integrant_un.get_spectral_va()
00472       .mult_sp() ;
00473     
00474     integrant_deux.set_spectral_va() = integrant_deux.get_spectral_va()
00475       .mult_st() ;
00476     integrant_deux.set_spectral_va() = integrant_deux.get_spectral_va()
00477       .mult_cp() ;
00478     
00479     p_angu_mom->set(2) = mapping.integrale_surface_infini (-integrant_un
00480                      +integrant_deux) / (2*qpig) ;
00481 
00482   }
00483 
00484     */
00485     
00486     if (p_angu_mom == 0x0) {        // a new computation is requireed
00487     
00488     p_angu_mom = new Tbl(3) ; 
00489     p_angu_mom->annule_hard() ; // fills the double array with zeros
00490 
00491     // Reference Cartesian vector basis of the Absolute frame
00492     Base_vect_cart bvect_ref(0.) ;  // 0. = parallel to the Absolute frame
00493     
00494     for (int i=0; i<=1; i++) {      // loop on the stars
00495 
00496         const Map& mp = et[i]->get_mp() ; 
00497         int nzm1 = mp.get_mg()->get_nzone() - 1 ; 
00498         
00499         // Function exp(-(r-r_0)^2/sigma^2)
00500         // --------------------------------
00501         
00502         double r0 = mp.val_r(nzm1-1, 1, 0, 0) ;
00503         double sigma = 1.*r0 ;
00504         
00505         Scalar rr (mp) ;
00506         rr = mp.r ;
00507         
00508         Scalar ff (mp) ;
00509         ff = exp( -(rr - r0)*(rr - r0)/sigma/sigma ) ;
00510         for (int ii=0; ii<nzm1; ii++)
00511           ff.set_domain(ii) = 1. ;
00512         ff.set_outer_boundary(nzm1, 0) ;
00513         ff.std_spectral_base() ;
00514 
00515         // Azimuthal vector d/dphi 
00516         Vector vphi(mp, CON, bvect_ref) ;       
00517         Scalar yya (mp) ;
00518         yya = mp.ya ;
00519         Scalar xxa (mp) ;
00520         xxa = mp.xa ;
00521         vphi.set(1) = - yya * ff ;  // phi^X
00522         vphi.set(2) = xxa * ff ; 
00523         vphi.set(3) = 0 ;  
00524 
00525         vphi.set(1).set_outer_boundary(nzm1, 0) ;
00526         vphi.set(2).set_outer_boundary(nzm1, 0) ;
00527     
00528         vphi.std_spectral_base() ; 
00529         vphi.change_triad(mp.get_bvect_cart()) ; 
00530         
00531         // Matter part
00532         // -----------
00533         const Scalar& ee = et[i]->get_ener_euler() ;  // E
00534         const Scalar& pp = et[i]->get_press() ; // p
00535         const Scalar& psi4 = et[i]->get_psi4() ; // Psi^4
00536         Scalar rho = pow(psi4, double(2.5)) * (ee+pp) ; 
00537         rho.std_spectral_base() ;
00538 
00539         Vector jmom = rho * (et[i]->get_u_euler()) ; 
00540                 
00541         const Metric& gtilde = et[i]->get_gtilde() ; 
00542         const Metric_flat flat (mp.flat_met_cart()) ; 
00543         
00544         Vector vphi_cov = vphi.up_down(gtilde) ;
00545         
00546         Scalar integrand = contract(jmom, 0, vphi_cov, 0) ; 
00547               
00548         p_angu_mom->set(2) += integrand.integrale() ;
00549 
00550         // Extrinsic curvature part (0 if IWM)
00551         // -----------------------------------
00552         
00553         const Sym_tensor& aij = et[i]->get_tkij_auto() ;
00554         rho = pow(psi4, double(1.5)) ;  
00555         rho.std_spectral_base() ;
00556         
00557         // Construction of D_k \Phi^i
00558         Itbl type (2) ;
00559         type.set(0) = CON ;
00560         type.set(1) = COV ;
00561         
00562         Tensor dcov_vphi (mp, 2, type, mp.get_bvect_cart()) ;
00563         dcov_vphi.set(1,1) = 0. ;
00564         dcov_vphi.set(2,1) = ff ;
00565         dcov_vphi.set(3,1) = 0. ;
00566         dcov_vphi.set(2,2) = 0. ;
00567         dcov_vphi.set(3,2) = 0. ;
00568         dcov_vphi.set(3,3) = 0. ;
00569         dcov_vphi.set(1,2) = -ff ;
00570         dcov_vphi.set(1,3) = 0. ;
00571         dcov_vphi.set(2,3) = 0. ;
00572         dcov_vphi.inc_dzpuis(2) ;
00573         
00574         Connection gamijk (gtilde, flat) ;
00575         const Tensor& deltaijk = gamijk.get_delta() ;
00576         
00577         // Computation of \tilde D_i \tilde \Phi_j
00578         Sym_tensor kill_phi (mp, COV, mp.get_bvect_cart()) ;
00579         kill_phi = contract(gtilde.cov(), 1, dcov_vphi +
00580                     contract(deltaijk, 2, vphi, 0), 0) +
00581           contract(dcov_vphi + contract(deltaijk, 2, vphi, 0), 0,
00582                gtilde.cov(), 1) ; 
00583 
00584         integrand = rho * contract(aij, 0, 1, kill_phi, 0, 1) ; 
00585         
00586         p_angu_mom->set(2) += integrand.integrale() / (4*qpig) ;
00587         
00588         
00589     }  // End of the loop on the stars
00590 
00591     }   // End of the case where a new computation was necessary
00592   
00593     return *p_angu_mom ; 
00594   
00595 }
00596 
00597 
00598 
00599             //---------------------------------//
00600             //      Total energy           //
00601             //---------------------------------//
00602 
00603 double Binary::total_ener() const {
00604     /*
00605     if (p_total_ener == 0x0) {      // a new computation is requireed
00606     
00607     p_total_ener = new double ; 
00608         
00609         *p_total_ener = mass_adm() - star1.mass_b() - star2.mass_b() ; 
00610         
00611     }   // End of the case where a new computation was necessary
00612     
00613     */
00614     return *p_total_ener ; 
00615     
00616 }
00617 
00618 
00619             //---------------------------------//
00620             //   Error on the virial theorem   //
00621             //---------------------------------//
00622 
00623 double Binary::virial() const {
00624     
00625     if (p_virial == 0x0) {      // a new computation is requireed
00626     
00627     p_virial = new double ; 
00628         
00629         *p_virial = 1. - mass_kom() / mass_adm() ; 
00630         
00631     }
00632     
00633     return *p_virial ; 
00634     
00635 }

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