et_bfrot_global.C

00001 /*
00002  *   Copyright (c) 2001 Jerome Novak
00003  *
00004  *   This file is part of LORENE.
00005  *
00006  *   LORENE is free software; you can redistribute it and/or modify
00007  *   it under the terms of the GNU General Public License as published by
00008  *   the Free Software Foundation; either version 2 of the License, or
00009  *   (at your option) any later version.
00010  *
00011  *   LORENE is distributed in the hope that it will be useful,
00012  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00013  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00014  *   GNU General Public License for more details.
00015  *
00016  *   You should have received a copy of the GNU General Public License
00017  *   along with LORENE; if not, write to the Free Software
00018  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00019  *
00020  */
00021 
00022 
00023 char et_bfrot_global_C[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/et_bfrot_global.C,v 1.14 2004/09/03 13:55:07 j_novak Exp $" ;
00024 
00025 /*
00026  * $Id: et_bfrot_global.C,v 1.14 2004/09/03 13:55:07 j_novak Exp $
00027  * $Log: et_bfrot_global.C,v $
00028  * Revision 1.14  2004/09/03 13:55:07  j_novak
00029  * Use of enerps_euler instead of ener_euler in the calculation of mom_angu()
00030  *
00031  * Revision 1.13  2004/03/25 10:29:03  j_novak
00032  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
00033  *
00034  * Revision 1.12  2003/11/20 14:01:26  r_prix
00035  * changed member names to better conform to Lorene coding standards:
00036  * J_euler -> j_euler, EpS_euler -> enerps_euler, Delta_car -> delta_car
00037  *
00038  * Revision 1.11  2003/11/18 18:38:11  r_prix
00039  * use of new member EpS_euler: matter sources in equilibrium() and global quantities
00040  * no longer distinguish Newtonian/relativistic, as all terms should have the right limit...
00041  *
00042  * Revision 1.10  2003/11/13 12:13:15  r_prix
00043  * *) removed all uses of Etoile-type specific quantities: u_euler, press
00044  *    and exclusively use ener_euler, s_euler, sphph_euler, J_euler instead
00045  * *) this also fixes a bug in previous expression for mass_g() in 2-fluid case
00046  *
00047  * NOTE: some current Newtonian expressions are still completely broken..
00048  *
00049  * Revision 1.9  2003/09/17 08:27:50  j_novak
00050  * New methods: mass_b1() and mass_b2().
00051  *
00052  * Revision 1.8  2003/02/07 10:47:43  j_novak
00053  * The possibility of having gamma5 xor gamma6 =0 has been introduced for
00054  * tests. The corresponding parameter files have been added.
00055  *
00056  * Revision 1.7  2002/10/16 14:36:36  j_novak
00057  * Reorganization of #include instructions of standard C++, in order to
00058  * use experimental version 3 of gcc.
00059  *
00060  * Revision 1.6  2002/10/14 14:20:08  j_novak
00061  * Error corrected for angu_mom()
00062  *
00063  * Revision 1.5  2002/05/10 09:26:52  j_novak
00064  * Added new class Et_rot_mag for magnetized rotating neutron stars (under development)
00065  *
00066  * Revision 1.4  2002/04/05 09:09:36  j_novak
00067  * The inversion of the EOS for 2-fluids polytrope has been modified.
00068  * Some errors in the determination of the surface were corrected.
00069  *
00070  * Revision 1.3  2002/01/08 14:43:53  j_novak
00071  * better determination of surfaces for 2-fluid stars
00072  *
00073  * Revision 1.2  2002/01/03 15:30:28  j_novak
00074  * Some comments modified.
00075  *
00076  * Revision 1.1.1.1  2001/11/20 15:19:28  e_gourgoulhon
00077  * LORENE
00078  *
00079  * Revision 1.4  2001/08/31  15:07:12  novak
00080  * Retour a la version 1.2, sans la routine prolonge_c1
00081  *
00082  * Revision 1.2  2001/08/28 15:32:17  novak
00083  * The surface is now defined by the baryonic density
00084  *
00085  * Revision 1.1  2001/06/22 15:39:46  novak
00086  * Initial revision
00087  *
00088  *
00089  * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bfrot_global.C,v 1.14 2004/09/03 13:55:07 j_novak Exp $
00090  *
00091  */
00092 
00093 
00094 // Headers C
00095 #include <stdlib.h>
00096 #include <math.h>
00097 
00098 // Headers Lorene
00099 #include "et_rot_bifluid.h"
00100 #include "unites.h"
00101 
00102 //--------------------------//
00103 //  Baryon mass     //
00104 //--------------------------//
00105 
00106 double Et_rot_bifluid::mass_b1() const {
00107 
00108   if (p_mass_b1 == 0x0) {    // a new computation is required
00109     
00110     assert(nbar.get_etat() == ETATQCQ);
00111     Cmp dens1 = a_car() * bbb() * gam_euler() * eos.get_m1()* nbar();
00112     dens1.std_base_scal() ; 
00113 
00114     p_mass_b1 = new double( dens1.integrale() ) ;
00115 
00116   }
00117     
00118   return *p_mass_b1 ; 
00119 
00120 } 
00121 
00122 double Et_rot_bifluid::mass_b2() const {
00123 
00124   if (p_mass_b2 == 0x0) {    // a new computation is required
00125       assert(nbar2.get_etat() == ETATQCQ);
00126     
00127       Cmp dens2 = a_car() * bbb() * gam_euler2() * eos.get_m2() * nbar2();
00128       dens2.std_base_scal() ; 
00129 
00130       p_mass_b2 = new double( dens2.integrale() ) ;
00131 
00132   }
00133     
00134   return *p_mass_b2 ; 
00135 
00136 } 
00137 
00138 double Et_rot_bifluid::mass_b() const {
00139 
00140   if (p_mass_b == 0x0) 
00141     p_mass_b = new double(mass_b1() + mass_b2() ) ;
00142 
00143   return *p_mass_b;
00144 
00145 } 
00146 
00147 
00148 //----------------------------//
00149 //  Gravitational mass    //
00150 //----------------------------//
00151 
00152 double Et_rot_bifluid::mass_g() const {
00153 
00154   if (p_mass_g == 0x0) {    // a new computation is required
00155     
00156       Tenseur vtmp (j_euler);
00157       vtmp.change_triad( mp.get_bvect_spher() );  // change to spherical triad
00158 
00159       Tenseur tJphi (mp);
00160       tJphi = vtmp(2);  // get the phi tetrad-component: J^ph r sin(th)
00161 
00162       // relativistic AND Newtonian: enerps_euler has right limit and tnphi->0
00163       Tenseur source = nnn * enerps_euler + 2 * b_car * tnphi * tJphi;
00164       source = a_car * bbb * source ;
00165       
00166       source.set_std_base() ;
00167 
00168       p_mass_g = new double( source().integrale() ) ;
00169       
00170   }
00171     
00172   return *p_mass_g ; 
00173 
00174 } 
00175         
00176 //----------------------------//
00177 //  Angular momentum      //
00178 //----------------------------//
00179 
00180 double Et_rot_bifluid::angu_mom() const {
00181 
00182   if (p_angu_mom == 0x0) {    // a new computation is required
00183     
00184     Cmp dens(mp) ;
00185 
00186     // this should work for both relativistic AND Newtonian, 
00187     // provided j_euler has the right limit...
00188     Tenseur tmp = j_euler;
00189     tmp.change_triad( mp.get_bvect_spher() ) ;  
00190     dens = tmp(2) ;
00191     dens.mult_rsint() ;
00192 
00193     dens = a_car() * b_car()* bbb() * dens ; 
00194 
00195     dens.std_base_scal() ; 
00196       
00197     p_angu_mom = new double( dens.integrale() ) ;
00198       
00199   }
00200     
00201   return *p_angu_mom ; 
00202 
00203 }
00204 
00205 
00206 //----------------------------//
00207 //       GRV2         //
00208 //----------------------------//
00209 
00210 double Et_rot_bifluid::grv2() const {
00211 
00212   using namespace Unites ;
00213 
00214   if (p_grv2 == 0x0) {    // a new computation is required
00215     
00216     Tenseur sou_m(mp);
00217     // should work for both relativistic AND Newtonian, provided sphph_euler has right limit..
00218     sou_m =  2 * qpig * a_car * sphph_euler ;
00219 
00220     Tenseur sou_q =  1.5 * ak_car - flat_scalar_prod(logn.gradient_spher(),logn.gradient_spher() ) ;    
00221 
00222     p_grv2 = new double( double(1) - lambda_grv2(sou_m(), sou_q()) ) ;  
00223 
00224   }
00225     
00226   return *p_grv2 ; 
00227 
00228 }
00229 
00230 
00231 //----------------------------//
00232 //       GRV3         //
00233 //----------------------------//
00234 
00235 double Et_rot_bifluid::grv3(ostream* ost) const {
00236 
00237   using namespace Unites ;
00238 
00239   if (p_grv3 == 0x0) {    // a new computation is required
00240 
00241     Tenseur source(mp) ; 
00242     
00243     // Gravitational term [cf. Eq. (43) of Gourgoulhon & Bonazzola
00244     // ------------------       Class. Quantum Grav. 11, 443 (1994)]
00245 
00246     if (relativistic) {
00247       Tenseur alpha = dzeta - logn ; 
00248       Tenseur beta = log( bbb ) ; 
00249       beta.set_std_base() ; 
00250         
00251       source = 0.75 * ak_car - flat_scalar_prod(logn.gradient_spher(), logn.gradient_spher() )
00252     + 0.5 * flat_scalar_prod(alpha.gradient_spher(), beta.gradient_spher() ) ; 
00253         
00254       Cmp aa = alpha() - 0.5 * beta() ; 
00255       Cmp daadt = aa.srdsdt() ; // 1/r d/dth
00256         
00257       // What follows is valid only for a mapping of class Map_radial : 
00258       const Map_radial* mpr = dynamic_cast<const Map_radial*>(&mp) ; 
00259       if (mpr == 0x0) {
00260     cout << "Et_rot_bifluid::grv3: the mapping does not belong"
00261          << " to the class Map_radial !" << endl ; 
00262     abort() ; 
00263       }
00264         
00265       // Computation of 1/tan(theta) * 1/r daa/dtheta
00266       if (daadt.get_etat() == ETATQCQ) {
00267     Valeur& vdaadt = daadt.va ; 
00268     vdaadt = vdaadt.ssint() ;   // division by sin(theta)
00269     vdaadt = vdaadt.mult_ct() ; // multiplication by cos(theta)
00270       }
00271         
00272       Cmp temp = aa.dsdr() + daadt ; 
00273       temp = ( bbb() - a_car()/bbb() ) * temp ; 
00274       temp.std_base_scal() ; 
00275         
00276       // Division by r 
00277       Valeur& vtemp = temp.va ; 
00278       vtemp = vtemp.sx() ;    // division by xi in the nucleus
00279       // Id in the shells
00280       // division by xi-1 in the ZEC
00281       vtemp = (mpr->xsr) * vtemp ; // multiplication by xi/r in the nucleus
00282       //          by 1/r in the shells
00283       //          by r(xi-1) in the ZEC
00284 
00285       // In the ZEC, a multiplication by r has been performed instead
00286       //   of the division:             
00287       temp.set_dzpuis( temp.get_dzpuis() + 2 ) ;  
00288         
00289       source = bbb() * source() + 0.5 * temp ; 
00290 
00291     }
00292     else{
00293       source = - 0.5 * flat_scalar_prod(logn.gradient_spher(),logn.gradient_spher() ) ; 
00294     }
00295     
00296     source.set_std_base() ; 
00297 
00298     double int_grav = source().integrale() ; 
00299 
00300     // Matter term
00301     // -----------
00302 
00303     // should work for relativistic AND Newtonian, provided s_euler has the right limit...
00304     source  = qpig * a_car * bbb * s_euler ;
00305 
00306     source.set_std_base() ; 
00307 
00308     double int_mat = source().integrale() ; 
00309 
00310     // Virial error
00311     // ------------
00312     if (ost != 0x0) {
00313       *ost << "Et_rot_bifluid::grv3 : gravitational term : " << int_grav 
00314        << endl ;
00315       *ost << "Et_rot_bifluid::grv3 : matter term :        " << int_mat 
00316        << endl ;
00317     }
00318 
00319     p_grv3 = new double( (int_grav + int_mat) / int_mat ) ;      
00320 
00321   }
00322     
00323   return *p_grv3 ; 
00324 
00325 }
00326 
00327 
00328 //----------------------------//
00329 //       R_circ       //
00330 //----------------------------//
00331 
00332 double Et_rot_bifluid::r_circ2() const {
00333 
00334   if (p_r_circ2 == 0x0) {    // a new computation is required
00335     
00336     // Index of the point at phi=0, theta=pi/2 at the surface of the star:
00337     const Mg3d* mg = mp.get_mg() ; 
00338     assert(mg->get_type_t() == SYM) ; 
00339     int l_b = nzet - 1 ; 
00340     int i_b = mg->get_nr(l_b) - 1 ; 
00341     int j_b = mg->get_nt(l_b) - 1 ; 
00342     int k_b = 0 ; 
00343     
00344     p_r_circ2 = new double( bbb()(l_b, k_b, j_b, i_b) * ray_eq2() ) ; 
00345 
00346   }
00347     
00348   return *p_r_circ2 ; 
00349 
00350 }
00351 
00352 
00353 //----------------------------//
00354 //     Flattening         //
00355 //----------------------------//
00356 
00357 double Et_rot_bifluid::aplat2() const {
00358 
00359   if (p_aplat2 == 0x0) {    // a new computation is required
00360     
00361     p_aplat2 = new double( ray_pole2() / ray_eq2() ) ;   
00362 
00363   }
00364     
00365   return *p_aplat2 ; 
00366 
00367 }
00368 
00369 
00370 
00371 //----------------------------//
00372 //     Quadrupole moment      //
00373 //----------------------------//
00374 
00375 double Et_rot_bifluid::mom_quad() const {
00376 
00377   using namespace Unites ;
00378 
00379   if (p_mom_quad == 0x0) {    // a new computation is required
00380     
00381     // Source for of the Poisson equation for nu
00382     // -----------------------------------------
00383 
00384     Tenseur source(mp) ; 
00385     
00386     if (relativistic) {
00387       Tenseur beta = log(bbb) ; 
00388       beta.set_std_base() ; 
00389       source =  qpig * a_car * enerps_euler 
00390     + ak_car - flat_scalar_prod(logn.gradient_spher(), 
00391                     logn.gradient_spher() + beta.gradient_spher()) ; 
00392     }
00393     else {
00394       source = qpig * (nbar + nbar2); 
00395     }
00396     source.set_std_base() ;     
00397 
00398     // Multiplication by -r^2 P_2(cos(theta))
00399     //  [cf Eq.(7) of Salgado et al. Astron. Astrophys. 291, 155 (1994) ]
00400     // ------------------------------------------------------------------
00401     
00402     // Multiplication by r^2 : 
00403     // ----------------------
00404     Cmp& csource = source.set() ; 
00405     csource.mult_r() ; 
00406     csource.mult_r() ; 
00407     if (csource.check_dzpuis(2)) {
00408       csource.inc2_dzpuis() ; 
00409     }
00410         
00411     // Muliplication by cos^2(theta) :
00412     // -----------------------------
00413     Cmp temp = csource ; 
00414     
00415     // What follows is valid only for a mapping of class Map_radial : 
00416     assert( dynamic_cast<const Map_radial*>(&mp) != 0x0 ) ; 
00417         
00418     if (temp.get_etat() == ETATQCQ) {
00419       Valeur& vtemp = temp.va ; 
00420       vtemp = vtemp.mult_ct() ; // multiplication by cos(theta)
00421       vtemp = vtemp.mult_ct() ; // multiplication by cos(theta)
00422     }
00423     
00424     // Muliplication by -P_2(cos(theta)) :
00425     // ----------------------------------
00426     source = 0.5 * source() - 1.5 * temp ; 
00427     
00428     // Final result
00429     // ------------
00430 
00431     p_mom_quad = new double( source().integrale() / qpig ) ;     
00432 
00433   }
00434     
00435   return *p_mom_quad ; 
00436 
00437 }
00438 
00439 
00440 //--------------------------//
00441 //  Stellar surface     //
00442 //--------------------------//
00443 
00444 const Itbl& Et_rot_bifluid::l_surf() const {
00445 
00446   if (p_l_surf == 0x0) {    // a new computation is required
00447     
00448     assert(p_xi_surf == 0x0) ;  // consistency check
00449     
00450     int np = mp.get_mg()->get_np(0) ;   
00451     int nt = mp.get_mg()->get_nt(0) ;   
00452     
00453     p_l_surf = new Itbl(np, nt) ;
00454     p_xi_surf = new Tbl(np, nt) ;
00455     
00456     double nb0 = 0 ;    // definition of the surface
00457     double precis = 1.e-15 ; 
00458     int nitermax = 100 ; 
00459     int niter ; 
00460           
00461     // Cmp defining the surface of the star (via the density fields)
00462     // 
00463     Cmp surf(mp) ;
00464     surf = -0.2*nbar()(0,0,0,0) ;
00465     surf.annule(0, nzet-1) ;
00466     surf += nbar() ; ;
00467     surf = prolonge_c1(surf, nzet) ;
00468     
00469     (surf.va).equipot(nb0, nzet, precis, nitermax, niter, *p_l_surf, 
00470             *p_xi_surf) ; 
00471     
00472   }
00473    
00474   return *p_l_surf ; 
00475     
00476 }
00477 const Itbl& Et_rot_bifluid::l_surf2() const {
00478 
00479   if (p_l_surf2 == 0x0) {    // a new computation is required
00480     
00481     assert(p_xi_surf2 == 0x0) ;  // consistency check
00482     
00483     int np = mp.get_mg()->get_np(0) ;   
00484     int nt = mp.get_mg()->get_nt(0) ;   
00485     
00486     p_l_surf2 = new Itbl(np, nt) ;
00487     p_xi_surf2 = new Tbl(np, nt) ;
00488     
00489     double nb0 = 0 ;    // definition of the surface
00490     double precis = 1.e-15 ; 
00491     int nitermax = 100 ; 
00492     int niter ; 
00493     
00494     // Cmp defining the surface of the star (via the density fields)
00495     // 
00496     Cmp surf2(mp) ;
00497     surf2 = -0.2*nbar2()(0,0,0,0) ;
00498     surf2.annule(0, nzet-1) ;
00499     surf2 += nbar2() ; ;
00500     surf2 = prolonge_c1(surf2, nzet) ;
00501     
00502     (surf2.va).equipot(nb0, nzet, precis, nitermax, niter, *p_l_surf2, 
00503             *p_xi_surf2) ; 
00504     
00505   }
00506    
00507   return *p_l_surf2 ; 
00508     
00509 }
00510 
00511 const Tbl& Et_rot_bifluid::xi_surf2() const {
00512 
00513   if (p_xi_surf2 == 0x0) {    // a new computation is required
00514     
00515     assert(p_l_surf2 == 0x0) ;  // consistency check
00516     
00517     l_surf2() ;  // the computation is done by l_surf2()
00518     
00519   }
00520    
00521   return *p_xi_surf2 ; 
00522     
00523 }
00524 
00525 //--------------------------//
00526 //  Coordinate radii    //
00527 //--------------------------//
00528 
00529 double Et_rot_bifluid::ray_eq2() const {
00530 
00531   if (p_ray_eq2 == 0x0) {    // a new computation is required
00532     
00533     const Mg3d& mg = *(mp.get_mg()) ;
00534     
00535     int type_t = mg.get_type_t() ; 
00536     int nt = mg.get_nt(0) ;     
00537     
00538     if ( type_t == SYM ) {
00539       assert( ( mg.get_type_p() == SYM) || (mg.get_type_p() == NONSYM) ) ; 
00540       int k = 0 ; 
00541       int j = nt-1 ; 
00542       int l = l_surf2()(k, j) ; 
00543       double xi = xi_surf2()(k, j) ; 
00544       double theta = M_PI / 2 ; 
00545       double phi = 0 ; 
00546         
00547       p_ray_eq2 = new double( mp.val_r(l, xi, theta, phi) ) ;
00548 
00549     }
00550     else {
00551       cout << "Et_rot_bifluid::ray_eq2 : the case type_t = " << type_t
00552        << " is not contemplated yet !" << endl ;
00553       abort() ; 
00554     }
00555 
00556   }
00557     
00558   return *p_ray_eq2 ; 
00559 
00560 } 
00561 
00562 
00563 double Et_rot_bifluid::ray_eq2_pis2() const {
00564 
00565   if (p_ray_eq2_pis2 == 0x0) {    // a new computation is required
00566     
00567     const Mg3d& mg = *(mp.get_mg()) ;
00568     
00569     int type_t = mg.get_type_t() ; 
00570     int type_p = mg.get_type_p() ; 
00571     int nt = mg.get_nt(0) ;     
00572     int np = mg.get_np(0) ;     
00573     
00574     if ( type_t == SYM ) {
00575       
00576       int j = nt-1 ; 
00577       double theta = M_PI / 2 ; 
00578       double phi = M_PI / 2 ;
00579       
00580       switch (type_p) {
00581     
00582       case SYM : {
00583     int k = np / 2  ; 
00584     int l = l_surf2()(k, j) ; 
00585     double xi = xi_surf2()(k, j) ; 
00586     p_ray_eq2_pis2 = new double( mp.val_r(l, xi, theta, phi) ) ;
00587     break ; 
00588       }
00589       
00590       case NONSYM : {
00591     assert( np % 4 == 0 ) ; 
00592     int k = np / 4  ; 
00593     int l = l_surf2()(k, j) ; 
00594     double xi = xi_surf2()(k, j) ; 
00595     p_ray_eq2_pis2 = new double( mp.val_r(l, xi, theta, phi) ) ;
00596     break ; 
00597       }
00598         
00599       default : {
00600     cout << "Et_rot_bifluid::ray_eq2_pis2 : the case type_p = " 
00601          << type_p << " is not contemplated yet !" << endl ;
00602     abort() ; 
00603       }
00604       } 
00605       
00606     }
00607     else {
00608       cout << "Et_rot_bifluid::ray_eq2_pis2 : the case type_t = " << type_t
00609        << " is not contemplated yet !" << endl ;
00610       abort() ; 
00611     }
00612     
00613   }
00614     
00615   return *p_ray_eq2_pis2 ; 
00616 
00617 } 
00618 
00619 
00620 double Et_rot_bifluid::ray_eq2_pi() const {
00621 
00622   if (p_ray_eq2_pi == 0x0) {    // a new computation is required
00623     
00624     const Mg3d& mg = *(mp.get_mg()) ;
00625     
00626     int type_t = mg.get_type_t() ; 
00627     int type_p = mg.get_type_p() ; 
00628     int nt = mg.get_nt(0) ;     
00629     int np = mg.get_np(0) ;     
00630     
00631     if ( type_t == SYM ) {
00632       
00633       switch (type_p) {
00634     
00635       case SYM : {
00636     p_ray_eq2_pi = new double( ray_eq2() ) ;
00637     break ; 
00638       }     
00639       
00640       case NONSYM : {
00641     int k = np / 2  ; 
00642     int j = nt-1 ; 
00643     int l = l_surf2()(k, j) ; 
00644     double xi = xi_surf2()(k, j) ; 
00645     double theta = M_PI / 2 ; 
00646     double phi = M_PI ; 
00647     
00648     p_ray_eq2_pi = new double( mp.val_r(l, xi, theta, phi) ) ;
00649     break ;
00650       }
00651       
00652       default : {
00653     
00654     cout << "Et_rot_bifluid::ray_eq2_pi : the case type_t = " << type_t
00655          << " and type_p = " << type_p << endl ; 
00656     cout << " is not contemplated yet !" << endl ;
00657     abort() ; 
00658     break ; 
00659       }
00660       }
00661     }
00662     
00663   }
00664   
00665   return *p_ray_eq2_pi ; 
00666   
00667 } 
00668 
00669 double Et_rot_bifluid::ray_pole2() const {
00670 
00671   if (p_ray_pole2 == 0x0) {    // a new computation is required
00672     
00673     assert( ((mp.get_mg())->get_type_t() == SYM) 
00674         || ((mp.get_mg())->get_type_t() == NONSYM) ) ;  
00675     
00676     int k = 0 ; 
00677     int j = 0 ; 
00678     int l = l_surf2()(k, j) ; 
00679     double xi = xi_surf2()(k, j) ; 
00680     double theta = 0 ; 
00681     double phi = 0 ; 
00682         
00683     p_ray_pole2 = new double( mp.val_r(l, xi, theta, phi) ) ;
00684 
00685   }
00686     
00687   return *p_ray_pole2 ; 
00688 
00689 } 
00690 
00691 
00692 

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