et_bfrot_equilibre.C

00001 /*
00002  * Method of class Et_rot_bifluid to compute a static spherical configuration.
00003  *
00004  * (see file etoile.h for documentation).
00005  *
00006  */
00007 
00008 /*
00009  *   Copyright (c) 2001 Jerome Novak
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 as published by
00015  *   the Free Software Foundation; either version 2 of the License, or
00016  *   (at your option) any later version.
00017  *
00018  *   LORENE is distributed in the hope that it will be useful,
00019  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00020  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00021  *   GNU General Public License for more details.
00022  *
00023  *   You should have received a copy of the GNU General Public License
00024  *   along with LORENE; if not, write to the Free Software
00025  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00026  *
00027  */
00028 
00029 
00030 char et_bfrot_equilibre_C[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/et_bfrot_equilibre.C,v 1.17 2006/03/13 10:02:27 j_novak Exp $" ;
00031 
00032 /*
00033  * $Id: et_bfrot_equilibre.C,v 1.17 2006/03/13 10:02:27 j_novak Exp $
00034  * $Log: et_bfrot_equilibre.C,v $
00035  * Revision 1.17  2006/03/13 10:02:27  j_novak
00036  * Added things for triaxial perturbations.
00037  *
00038  * Revision 1.16  2004/09/01 10:56:05  r_prix
00039  * added option of converging baryon-mass to equilibrium_bi()
00040  *
00041  * Revision 1.15  2004/08/30 09:54:20  r_prix
00042  * experimental version of Kepler-limit finder for 2-fluid stars
00043  *
00044  * Revision 1.14  2004/03/25 10:29:03  j_novak
00045  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
00046  *
00047  * Revision 1.13  2003/12/11 12:43:35  r_prix
00048  * activated adaptive grid for 2-fluid star (taken from Etoile_rot)
00049  *
00050  * Revision 1.12  2003/12/04 14:28:26  r_prix
00051  * allow for the case of "slow-rot-style" EOS inversion, in which we need to adapt
00052  * the inner domain to n_outer=0 instead of mu_outer=0 ...
00053  * (this should only be used for comparison to analytic slow-rot solution!)
00054  *
00055  * Revision 1.11  2003/11/25 12:49:44  j_novak
00056  * Modified headers to compile on IRIX. Changed Mapping to be Map_af (speed
00057  * enhancement).
00058  *
00059  * Revision 1.10  2003/11/20 14:01:25  r_prix
00060  * changed member names to better conform to Lorene coding standards:
00061  * J_euler -> j_euler, EpS_euler -> enerps_euler, Delta_car -> delta_car
00062  *
00063  * Revision 1.9  2003/11/19 22:01:57  e_gourgoulhon
00064  * -- Relaxation on logn and dzeta performed only if mer >= 10.
00065  * -- err_grv2 is now evaluated also in the Newtonian case.
00066  *
00067  * Revision 1.8  2003/11/18 18:38:11  r_prix
00068  * use of new member EpS_euler: matter sources in equilibrium() and global quantities
00069  * no longer distinguish Newtonian/relativistic, as all terms should have the right limit...
00070  *
00071  * Revision 1.7  2003/11/17 13:49:43  r_prix
00072  * - moved superluminal check into hydro_euler()
00073  * - removed some warnings
00074  *
00075  * Revision 1.6  2003/11/13 12:14:35  r_prix
00076  * *) removed all use of etoile-type specific u_euler, press
00077  *   and use 3+1 components of Tmunu instead
00078  *
00079  * Revision 1.5  2002/10/16 14:36:35  j_novak
00080  * Reorganization of #include instructions of standard C++, in order to
00081  * use experimental version 3 of gcc.
00082  *
00083  * Revision 1.4  2002/04/05 09:09:36  j_novak
00084  * The inversion of the EOS for 2-fluids polytrope has been modified.
00085  * Some errors in the determination of the surface were corrected.
00086  *
00087  * Revision 1.3  2002/01/16 15:03:28  j_novak
00088  * *** empty log message ***
00089  *
00090  * Revision 1.2  2002/01/03 15:30:28  j_novak
00091  * Some comments modified.
00092  *
00093  * Revision 1.1.1.1  2001/11/20 15:19:28  e_gourgoulhon
00094  * LORENE
00095  *
00096  * Revision 1.1  2001/06/22  15:40:06  novak
00097  * Initial revision
00098  *
00099  *
00100  * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bfrot_equilibre.C,v 1.17 2006/03/13 10:02:27 j_novak Exp $
00101  *
00102  */
00103 
00104 // Headers C
00105 #include <math.h>
00106 
00107 // Headers Lorene
00108 #include "et_rot_bifluid.h"
00109 #include "param.h"
00110 #include "unites.h"     
00111 
00112 #include "graphique.h"
00113 #include "utilitaires.h"
00114 
00115 //-------------------------------------------------------------------------
00116 //-------------------------------------------------------------------------
00117 
00118 //                        Axial Equilibrium
00119 
00120 //-------------------------------------------------------------------------
00121 //-------------------------------------------------------------------------
00122 
00123 void Et_rot_bifluid::equilibrium_bi
00124 (double ent_c, double ent2_c, double omega0, double omega20, 
00125  const Tbl& ent_limit, const Tbl& ent2_limit, const Itbl& icontrol, 
00126  const Tbl& control, Tbl& diff,
00127  int mer_mass, double mbar1_wanted, double mbar2_wanted, double aexp_mass) 
00128 {
00129                  
00130   // Fundamental constants and units
00131   // -------------------------------
00132   using namespace Unites ;
00133 
00134   // For the display 
00135   // ---------------
00136   char display_bold[]="x[1m" ; display_bold[0] = 27 ;
00137   char display_normal[] = "x[0m" ; display_normal[0] = 27 ;
00138 
00139   // Grid parameters
00140   // ---------------
00141     
00142   const Mg3d* mg = mp.get_mg() ; 
00143   int nz = mg->get_nzone() ;        // total number of domains
00144 
00145   // The following is required to initialize mp_prev as a Map_af
00146   Map_et& mp_et = dynamic_cast<Map_et&>(mp) ;  // reference
00147     
00148   // Index of the point at phi=0, theta=pi/2 at the surface of the star:
00149   assert(mg->get_type_t() == SYM) ; 
00150   int l_b = nzet - 1 ; 
00151   int i_b = mg->get_nr(l_b) - 1 ; 
00152   int j_b = mg->get_nt(l_b) - 1 ; 
00153   int k_b = 0 ; 
00154     
00155   // Value of the enthalpies defining the surface of each fluid
00156   double ent1_b = ent_limit(nzet-1) ;
00157   double ent2_b = ent2_limit(nzet-1) ;
00158 
00159   // This value is chosen so that the grid contain both fluids
00160   //  double ent_b = (ent1_b > ent2_b ? ent1_b : ent2_b) ; 
00161     
00162   // Parameters to control the iteration
00163   // -----------------------------------
00164     
00165   int mer_max = icontrol(0) ; 
00166   int mer_rot = icontrol(1) ;
00167   int mer_change_omega = icontrol(2) ; 
00168   int mer_fix_omega = icontrol(3) ; 
00169   int mermax_poisson = icontrol(4) ; 
00170   int nzadapt = icontrol(5);        // number of domains for adaptive grid
00171   int kepler_fluid = icontrol(6);   // fluid-index for kepler-limit (0=none,3=both)
00172   int kepler_wait_steps = icontrol(7);
00173   int mer_triax = icontrol(8) ;
00174   
00175   int niter ;
00176 
00177   // Protections:
00178   if (mer_change_omega < mer_rot) {
00179     cout << "Et_rot_bifluid::equilibrium: mer_change_omega < mer_rot !" << endl ;
00180     cout << " mer_change_omega = " << mer_change_omega << endl ; 
00181     cout << " mer_rot = " << mer_rot << endl ; 
00182     abort() ; 
00183   }
00184   if (mer_fix_omega < mer_change_omega) {
00185     cout << "Et_rot_bifluid::equilibrium: mer_fix_omega < mer_change_omega !" 
00186      << endl ;
00187     cout << " mer_fix_omega = " << mer_fix_omega << endl ; 
00188     cout << " mer_change_omega = " << mer_change_omega << endl ; 
00189     abort() ; 
00190   }
00191 
00192   double precis = control(0) ; 
00193   double omega_ini = control(1) ;
00194   double omega2_ini = control(2) ;
00195   double relax = control(3) ;
00196   double relax_prev = double(1) - relax ;  
00197   double relax_poisson = control(4) ; 
00198   // some additional stuff for adaptive grid:
00199   double thres_adapt = control(5) ; 
00200   double precis_adapt = control(6) ; 
00201   double kepler_factor = control(7);
00202   if (kepler_factor <= 1.0)
00203     {
00204       cout << "ERROR: Kepler factor has to be greater than 1!!\n";
00205       abort();
00206     }
00207   double ampli_triax = control(8) ;
00208 
00209   // Error indicators
00210   // ----------------
00211     
00212   diff.set_etat_qcq() ; 
00213   double diff_ent ;
00214   double& diff_ent1 = diff.set(0) ; 
00215   double& diff_ent2 = diff.set(1) ; 
00216   double& diff_nuf = diff.set(2) ; 
00217   double& diff_nuq = diff.set(3) ; 
00218   //    double& diff_dzeta = diff.set(4) ; 
00219   //    double& diff_ggg = diff.set(5) ; 
00220   double& diff_shift_x = diff.set(6) ; 
00221   double& diff_shift_y = diff.set(7) ; 
00222   double& vit_triax = diff.set(8) ;  
00223 
00224   // Parameters for the function Map_et::adapt
00225   // -----------------------------------------
00226     
00227   Param par_adapt ; 
00228   int nitermax = 100 ;  
00229   int adapt_flag = 1 ;    //  1 = performs the full computation, 
00230                 //  0 = performs only the rescaling by 
00231                 //      the factor alpha_r
00232   int nz_search = nzet + 1 ;  // Number of domains for searching the enthalpy
00233                 //  isosurfaces
00234   double alpha_r ; 
00235   double reg_map = 1. ; // 1 = regular mapping, 0 = contracting mapping
00236 
00237   par_adapt.add_int(nitermax, 0) ; // maximum number of iterations to
00238                      // locate zeros by the secant method
00239   par_adapt.add_int(nzadapt, 1) ; // number of domains where the adjustment 
00240                     // to the isosurfaces of ent is to be 
00241                     // performed
00242   par_adapt.add_int(nz_search, 2) ; // number of domains to search for
00243                     // the enthalpy isosurface
00244   par_adapt.add_int(adapt_flag, 3) ; //  1 = performs the full computation, 
00245                        //  0 = performs only the rescaling by 
00246                        //      the factor alpha_r
00247   par_adapt.add_int(j_b, 4) ; //  theta index of the collocation point 
00248                     //  (theta_*, phi_*)
00249   par_adapt.add_int(k_b, 5) ; //  theta index of the collocation point 
00250                     //  (theta_*, phi_*)
00251 
00252   par_adapt.add_int_mod(niter, 0) ;  //  number of iterations actually used in 
00253                        //  the secant method
00254     
00255   par_adapt.add_double(precis_adapt, 0) ; // required absolute precision in 
00256                          // the determination of zeros by 
00257                          // the secant method
00258   par_adapt.add_double(reg_map, 1)  ;  // 1. = regular mapping, 0 = contracting mapping
00259     
00260   par_adapt.add_double(alpha_r, 2) ;    // factor by which all the radial 
00261                        // distances will be multiplied 
00262            
00263   par_adapt.add_tbl(ent_limit, 0) ; // array of values of the field ent 
00264                         // to define the isosurfaces.              
00265 
00266 
00267   // Parameters for the function Map_et::poisson for nuf
00268   // ----------------------------------------------------
00269 
00270   double precis_poisson = 1.e-16 ;     
00271 
00272   Param par_poisson_nuf ; 
00273   par_poisson_nuf.add_int(mermax_poisson,  0) ;  // maximum number of iterations
00274   par_poisson_nuf.add_double(relax_poisson,  0) ; // relaxation parameter
00275   par_poisson_nuf.add_double(precis_poisson, 1) ; // required precision
00276   par_poisson_nuf.add_int_mod(niter, 0) ;  //  number of iterations actually used 
00277   par_poisson_nuf.add_cmp_mod( ssjm1_nuf ) ; 
00278                        
00279   Param par_poisson_nuq ; 
00280   par_poisson_nuq.add_int(mermax_poisson,  0) ;  // maximum number of iterations
00281   par_poisson_nuq.add_double(relax_poisson,  0) ; // relaxation parameter
00282   par_poisson_nuq.add_double(precis_poisson, 1) ; // required precision
00283   par_poisson_nuq.add_int_mod(niter, 0) ;  //  number of iterations actually used 
00284   par_poisson_nuq.add_cmp_mod( ssjm1_nuq ) ; 
00285                        
00286   Param par_poisson_tggg ; 
00287   par_poisson_tggg.add_int(mermax_poisson,  0) ;  // maximum number of iterations
00288   par_poisson_tggg.add_double(relax_poisson,  0) ; // relaxation parameter
00289   par_poisson_tggg.add_double(precis_poisson, 1) ; // required precision
00290   par_poisson_tggg.add_int_mod(niter, 0) ;  //  number of iterations actually used 
00291   par_poisson_tggg.add_cmp_mod( ssjm1_tggg ) ; 
00292   double lambda_tggg ;
00293   par_poisson_tggg.add_double_mod( lambda_tggg ) ; 
00294     
00295   Param par_poisson_dzeta ; 
00296   double lbda_grv2 ;
00297   par_poisson_dzeta.add_double_mod( lbda_grv2 ) ; 
00298                        
00299   // Parameters for the function Tenseur::poisson_vect
00300   // -------------------------------------------------
00301 
00302   Param par_poisson_vect ; 
00303 
00304   par_poisson_vect.add_int(mermax_poisson,  0) ;  // maximum number of iterations
00305   par_poisson_vect.add_double(relax_poisson,  0) ; // relaxation parameter
00306   par_poisson_vect.add_double(precis_poisson, 1) ; // required precision
00307   par_poisson_vect.add_cmp_mod( ssjm1_khi ) ; 
00308   par_poisson_vect.add_tenseur_mod( ssjm1_wshift ) ; 
00309   par_poisson_vect.add_int_mod(niter, 0) ;   
00310 
00311                        
00312     // Initializations
00313     // ---------------
00314 
00315     // Initial angular velocities
00316   omega = 0 ; 
00317   omega2 = 0 ;
00318   
00319   double accrois_omega = (omega0 - omega_ini) / 
00320     double(mer_fix_omega - mer_change_omega) ; 
00321   double accrois_omega2 = (omega20 - omega2_ini) / 
00322     double(mer_fix_omega - mer_change_omega) ; 
00323 
00324 
00325   update_metric() ; // update of the metric coefficients
00326 
00327   equation_of_state() ; // update of the densities, pressure, etc...
00328     
00329   hydro_euler() ;   // update of the hydro quantities relative to the 
00330   //  Eulerian observer
00331 
00332     // Quantities at the previous step :    
00333 
00334   Map_et mp_prev = mp_et;
00335   Tenseur ent_prev = ent ;  
00336   Tenseur ent2_prev = ent2 ;
00337   Tenseur logn_prev = logn ;        
00338   Tenseur dzeta_prev = dzeta ;      
00339     
00340     // Creation of uninitialized tensors:
00341   Tenseur source_nuf(mp) ;    // source term in the equation for nuf
00342   Tenseur source_nuq(mp) ;    // source term in the equation for nuq
00343   Tenseur source_dzf(mp) ;  // matter source term in the eq. for dzeta
00344   Tenseur source_dzq(mp) ;  // quadratic source term in the eq. for dzeta
00345   Tenseur source_tggg(mp) ; // source term in the eq. for tggg
00346   Tenseur source_shift(mp, 1, CON, mp.get_bvect_cart()) ;  
00347   // source term for shift
00348   Tenseur mlngamma(mp) ;    // centrifugal potential
00349   Tenseur mlngamma2(mp) ;   // centrifugal potential
00350 
00351   Tenseur *outer_ent_p;     // pointer to the enthalpy field of the outer fluid
00352     
00353   // Preparations for the Poisson equations:
00354   // --------------------------------------
00355   if (nuf.get_etat() == ETATZERO) {
00356     nuf.set_etat_qcq() ; 
00357     nuf.set() = 0 ; 
00358   }
00359     
00360   if (relativistic) {
00361     if (nuq.get_etat() == ETATZERO) {
00362       nuq.set_etat_qcq() ; 
00363       nuq.set() = 0 ; 
00364     }
00365 
00366     if (tggg.get_etat() == ETATZERO) {
00367       tggg.set_etat_qcq() ; 
00368       tggg.set() = 0 ; 
00369     }
00370     
00371     if (dzeta.get_etat() == ETATZERO) {
00372       dzeta.set_etat_qcq() ; 
00373       dzeta.set() = 0 ; 
00374     }
00375   }
00376             
00377   ofstream fichconv("convergence.d") ;    // Output file for diff_ent
00378   fichconv << "#     diff_ent     GRV2        max_triax      vit_triax" << endl ; 
00379     
00380   ofstream fichfreq("frequency.d") ;    // Output file for  omega
00381   fichfreq << "#       f1 [Hz]     f2 [Hz]" << endl ; 
00382     
00383   ofstream fichevol("evolution.d") ;    // Output file for various quantities
00384   fichevol << "#           r_pole/r_eq       ent_c    ent2_c" << endl ; 
00385     
00386   diff_ent = 1 ; 
00387   double err_grv2 = 1 ; 
00388   double max_triax_prev = 0 ;       // Triaxial amplitude at previous step
00389     
00390     //=========================================================================
00391     //          Start of iteration
00392     //=========================================================================
00393 
00394   for(int mer=0 ; (diff_ent > precis) && (mer<mer_max) ; mer++ ) {
00395 
00396     cout << "-----------------------------------------------" << endl ;
00397     cout << "step: " << mer << endl ;
00398     cout << "diff_ent = " << display_bold << diff_ent << display_normal
00399      << endl ;    
00400     cout << "err_grv2 = " << err_grv2 << endl ;    
00401     fichconv << mer ;
00402     fichfreq << mer ;
00403     fichevol << mer ;
00404       
00405     if (mer >= mer_rot) {
00406     
00407       if (mer < mer_change_omega) {
00408     omega = omega_ini ; 
00409     omega2 = omega2_ini ;
00410       }
00411       else {
00412     if (mer <= mer_fix_omega) {
00413       omega = omega_ini + accrois_omega * 
00414         (mer - mer_change_omega) ;
00415       omega2 = omega2_ini + accrois_omega2 * 
00416         (mer - mer_change_omega) ;
00417     }
00418       }
00419     
00420     }
00421 
00422     //-----------------------------------------------
00423     //  Sources of the Poisson equations
00424     //-----------------------------------------------
00425     
00426     // Source for nu
00427     // -------------
00428     Tenseur beta = log(bbb) ; 
00429     beta.set_std_base() ; 
00430 
00431     // common source term for relativistic and Newtonian ! (enerps_euler has the right limit)
00432     source_nuf =  qpig * a_car * enerps_euler;
00433 
00434     if (relativistic) 
00435       source_nuq = ak_car - flat_scalar_prod(logn.gradient_spher(),logn.gradient_spher() + beta.gradient_spher()) ; 
00436     else 
00437       source_nuq = 0 ; 
00438 
00439     source_nuf.set_std_base() ;     
00440     source_nuq.set_std_base() ;     
00441 
00442     if (relativistic) {
00443       // Source for dzeta
00444       // ----------------
00445       source_dzf = 2 * qpig * a_car * sphph_euler;
00446       source_dzf.set_std_base() ; 
00447       
00448       source_dzq = 1.5 * ak_car - flat_scalar_prod(logn.gradient_spher(),logn.gradient_spher() ) ;      
00449       source_dzq.set_std_base() ;   
00450       
00451       // Source for tggg
00452       // ---------------
00453       
00454       source_tggg = 2 * qpig * nnn * a_car * bbb * (s_euler - sphph_euler);
00455       source_tggg.set_std_base() ; 
00456       
00457       (source_tggg.set()).mult_rsint() ; 
00458       
00459       
00460       // Source for shift
00461       // ----------------
00462       
00463       // Matter term 
00464       source_shift = (-4*qpig) * nnn * a_car * j_euler;
00465       
00466       // Quadratic terms:
00467       Tenseur vtmp =  3 * beta.gradient_spher() - logn.gradient_spher() ;
00468       vtmp.change_triad(mp.get_bvect_cart()) ; 
00469       
00470       Tenseur squad  = 2 * nnn * flat_scalar_prod(tkij, vtmp) ;     
00471       
00472       // The addition of matter terms and quadratic terms is performed
00473       //  component by component because u_euler is contravariant,
00474       //  while squad is covariant. 
00475       
00476       if (squad.get_etat() == ETATQCQ) {
00477     for (int i=0; i<3; i++) {
00478       source_shift.set(i) += squad(i) ; 
00479     }
00480       }
00481       
00482       source_shift.set_std_base() ;     
00483     }
00484     //----------------------------------------------
00485     // Resolution of the Poisson equation for nuf 
00486     //----------------------------------------------
00487     
00488     source_nuf().poisson(par_poisson_nuf, nuf.set()) ; 
00489     
00490 //     cout << "Test of the Poisson equation for nuf :" << endl ; 
00491 //     Tbl err = source_nuf().test_poisson(nuf(), cout, true) ; 
00492 //     diff_nuf = err(0, 0) ; 
00493     diff_nuf = 0 ;
00494     
00495     //---------------------------------------
00496     // Triaxial perturbation of nuf
00497     //---------------------------------------
00498     
00499     if (mer == mer_triax) {
00500     
00501     if ( mg->get_np(0) == 1 ) {
00502         cout << 
00503         "Et_rot_bifluid::equilibrium: np must be stricly greater than 1"
00504          << endl << " to set a triaxial perturbation !" << endl ; 
00505         abort() ; 
00506     }
00507     
00508     const Coord& phi = mp.phi ; 
00509     const Coord& sint = mp.sint ; 
00510     Cmp perturb(mp) ; 
00511     perturb = 1 + ampli_triax * sint*sint * cos(2*phi) ;
00512     nuf.set() = nuf() * perturb ;  
00513         
00514     nuf.set_std_base() ;    // set the bases for spectral expansions
00515                             //  to be the standard ones for a 
00516                             //  scalar field
00517 
00518     }
00519     
00520     // Monitoring of the triaxial perturbation
00521     // ---------------------------------------
00522     
00523     Valeur& va_nuf = nuf.set().va ; 
00524     va_nuf.coef() ;     // Computes the spectral coefficients
00525     double max_triax = 0 ; 
00526 
00527     if ( mg->get_np(0) > 1 ) {
00528 
00529     for (int l=0; l<nz; l++) {  // loop on the domains
00530         for (int j=0; j<mg->get_nt(l); j++) {
00531         for (int i=0; i<mg->get_nr(l); i++) {
00532         
00533             // Coefficient of cos(2 phi) : 
00534             double xcos2p = (*(va_nuf.c_cf))(l, 2, j, i) ; 
00535 
00536             // Coefficient of sin(2 phi) : 
00537             double xsin2p = (*(va_nuf.c_cf))(l, 3, j, i) ; 
00538             
00539             double xx = sqrt( xcos2p*xcos2p + xsin2p*xsin2p ) ; 
00540             
00541             max_triax = ( xx > max_triax ) ? xx : max_triax ;           
00542         }
00543         } 
00544     }
00545     }
00546     cout << "Triaxial part of nuf : " << max_triax << endl ; 
00547 
00548     if (relativistic) {
00549       
00550       //----------------------------------------------
00551       // Resolution of the Poisson equation for nuq 
00552       //----------------------------------------------
00553 
00554       source_nuq().poisson(par_poisson_nuq, nuq.set()) ; 
00555         
00556 //        cout << "Test of the Poisson equation for nuq :" << endl ; 
00557 //        err = source_nuq().test_poisson(nuq(), cout, true) ;
00558 //        diff_nuq = err(0, 0) ; 
00559       diff_nuq = 0 ;
00560     
00561       //---------------------------------------------------------
00562       // Resolution of the vector Poisson equation for the shift
00563       //---------------------------------------------------------
00564 
00565       if (source_shift.get_etat() != ETATZERO) {
00566 
00567     for (int i=0; i<3; i++) {
00568       if(source_shift(i).dz_nonzero()) {
00569         assert( source_shift(i).get_dzpuis() == 4 ) ; 
00570       }
00571       else{
00572         (source_shift.set(i)).set_dzpuis(4) ; 
00573       }
00574     }
00575 
00576       }
00577 
00578       double lambda_shift = double(1) / double(3) ; 
00579 
00580       if ( mg->get_np(0) == 1 ) {
00581     lambda_shift = 0 ; 
00582       }
00583     
00584       source_shift.poisson_vect(lambda_shift, par_poisson_vect, 
00585                 shift, w_shift, khi_shift) ;      
00586         
00587 //        cout << "Test of the Poisson equation for shift_x :" << endl ; 
00588 //        err = source_shift(0).test_poisson(shift(0), cout, true) ;
00589 //        diff_shift_x = err(0, 0) ; 
00590       diff_shift_x = 0 ; 
00591     
00592 //        cout << "Test of the Poisson equation for shift_y :" << endl ; 
00593 //        err = source_shift(1).test_poisson(shift(1), cout, true) ;
00594 //        diff_shift_y = err(0, 0) ; 
00595       diff_shift_y = 0 ;
00596         
00597       // Computation of tnphi and nphi from the Cartesian components
00598       //  of the shift
00599       // -----------------------------------------------------------
00600         
00601       fait_nphi() ; 
00602     
00603     }
00604 
00605 
00606     //----------------------------------------
00607     // Shall we search for the Kepler limit?
00608     //----------------------------------------
00609     bool kepler = false;
00610     bool too_fast = false;
00611 
00612     if ( (kepler_fluid > 0) && (mer > mer_fix_omega + kepler_wait_steps) )
00613       {
00614     if (kepler_fluid & 0x01)
00615       omega *= kepler_factor;
00616     if (kepler_fluid & 0x02)
00617       omega2 *= kepler_factor;
00618       }
00619 
00620 
00621     // ============================================================
00622     kepler = true;
00623     while (kepler)  
00624       { 
00625 
00626     // New computation of delta_car, gam_euler, enerps_euler etc...
00627     // ------------------------------------------------------
00628     hydro_euler() ; 
00629     
00630 
00631     //------------------------------------------------------
00632     //  First integral of motion 
00633     //------------------------------------------------------
00634     
00635     // Centrifugal potential : 
00636     if (relativistic) {
00637       mlngamma = - log( gam_euler ) ;
00638       mlngamma2 = - log( gam_euler2) ;
00639     }
00640     else {
00641       mlngamma = - 0.5 * uuu*uuu ;
00642       mlngamma2 = -0.5 * uuu2*uuu2 ;
00643     }
00644     
00645     // Central values of various potentials :
00646     double nuf_c = nuf()(0,0,0,0) ; 
00647     double nuq_c = nuq()(0,0,0,0) ; 
00648 
00649     // Scale factor to ensure that the enthalpy is equal to ent_b at 
00650     //  the equator for the "outer" fluid
00651     double alpha_r2 = 0;
00652     
00653     int j=j_b;
00654     
00655     // Boundary values of various potentials :
00656     double nuf_b  = nuf()(l_b, k_b, j, i_b) ; 
00657     double nuq_b  = nuq()(l_b, k_b, j, i_b) ; 
00658     double mlngamma_b  = mlngamma()(l_b, k_b, j, i_b) ; 
00659     double mlngamma2_b  = mlngamma2()(l_b, k_b, j, i_b) ; 
00660 
00661 
00662     // RP: "hack": adapt the radius correctly if using "slow-rot-style" EOS inversion
00663     // 
00664     if ( eos.identify() == 2 ) // only applies to Eos_bf_poly_newt
00665       {
00666         const Eos_bf_poly_newt &eos0 = dynamic_cast<const Eos_bf_poly_newt&>(eos);
00667         if (eos0.get_typeos() == 5)
00668           {
00669         double vn_b = uuu()(l_b, k_b, j, i_b);
00670         double vp_b = uuu2()(l_b, k_b, j, i_b);
00671         double D2_b = (vp_b - vn_b)*(vp_b - vn_b);
00672         double kdet = eos0.get_kap3() + eos0.get_beta()*D2_b;
00673         double kaps1 = kdet / ( eos0.get_kap2() - kdet );
00674         double kaps2 = kdet / ( eos0.get_kap1() - kdet );
00675         
00676         ent1_b = kaps1 * ( ent2_c - ent_c - mlngamma2_b + mlngamma_b );
00677         ent2_b = kaps2 * ( ent_c - ent2_c - mlngamma_b + mlngamma2_b );
00678         
00679         cout << "**********************************************************************\n";
00680         cout << "DEBUG: Rescaling domain for slow-rot-style EOS inversion \n";
00681         cout << "DEBUG: ent1_b = " << ent1_b << "; ent2_b = " << ent2_b << endl;
00682         cout << "**********************************************************************\n";
00683         
00684         adapt_flag = 0; // don't do adaptive-grid if using slow-rot-style inversion!
00685           }
00686       }
00687     
00688     double alpha1_r2 = ( ent_c - ent1_b - mlngamma_b + nuq_c - nuq_b) / ( nuf_b - nuf_c  ) ;
00689     double alpha2_r2 = ( ent2_c - ent2_b  - mlngamma2_b + nuq_c - nuq_b) / ( nuf_b - nuf_c  ) ;
00690     
00691     cout << "DEBUG: j= "<< j<<" ; alpha1 = " << alpha1_r2 <<" ; alpha2 = " << alpha2_r2 << endl;
00692 
00693     int outer_fluid = (alpha1_r2 > alpha2_r2) ? 1 : 2;  // index of 'outer' fluid (at equator!)
00694     
00695     outer_ent_p = (outer_fluid == 1) ? (&ent) : (&ent2);
00696     
00697     alpha_r2 = (outer_fluid == 1) ? alpha1_r2 : alpha2_r2 ;
00698     
00699     alpha_r = sqrt(alpha_r2);
00700     
00701     cout << "alpha_r = " << alpha_r << endl ; 
00702     
00703     // Readjustment of nu :
00704     // -------------------
00705     
00706     logn = alpha_r2 * nuf + nuq ;
00707     double nu_c =  logn()(0,0,0,0) ;
00708     
00709     
00710     // First integral   --> enthalpy in all space
00711     //-----------------
00712     
00713     ent = (ent_c + nu_c) - logn - mlngamma ;
00714     ent2 = (ent2_c + nu_c) - logn - mlngamma2 ;
00715 
00716 
00717     // now let's try to figure out if we have overstepped the Kepler-limit
00718     // (FIXME) we assume that the enthalpy of the _outer_ fluid being negative
00719     // inside the star 
00720     kepler = false;
00721     for (int l=0; l<nzet; l++) {
00722       int imax = mg->get_nr(l) - 1 ;
00723       if (l == l_b) imax-- ;    // The surface point is skipped
00724       for (int i=0; i<imax; i++) { 
00725         if ( (*outer_ent_p)()(l, 0, j_b, i) < 0. ) {
00726           kepler = true;
00727           cout << "(outer) ent < 0 for l, i : " << l << "  " << i 
00728            << "   ent = " << (*outer_ent_p)()(l, 0, j_b, i) << endl ;  
00729         } 
00730       }
00731     }
00732     
00733     if ( kepler ) 
00734       {
00735         cout << "**** KEPLERIAN VELOCITY REACHED ****" << endl ; 
00736         if (kepler_fluid & 0x01)
00737           omega /= kepler_factor ;    // Omega is decreased
00738         if (kepler_fluid & 0x02)
00739           omega2 /= kepler_factor;
00740 
00741         cout << "New rotation frequencies : " 
00742          << "Omega = " << omega/(2.*M_PI) * f_unit << " Hz; " 
00743          << "Omega2 = " << omega2/(2.*M_PI) * f_unit << endl ; 
00744 
00745         too_fast = true;
00746       }
00747 
00748       } /* while kepler */
00749 
00750     
00751     if ( too_fast ) 
00752       { // fact_omega is decreased for the next step 
00753     kepler_factor = sqrt( kepler_factor ) ; 
00754     cout << "**** New fact_omega : " << kepler_factor << endl ; 
00755       }
00756     // ============================================================
00757 
00758 
00759     // Cusp-check: shall the adaptation still be performed?
00760     // ------------------------------------------
00761     double dent_eq = (*outer_ent_p)().dsdr()(l_b, k_b, j_b, i_b) ; 
00762     double dent_pole = (*outer_ent_p)().dsdr()(l_b, k_b, 0, i_b) ;
00763     double rap_dent = fabs( dent_eq / dent_pole ) ; 
00764     cout << "| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ; 
00765     
00766     if ( rap_dent < thres_adapt ) {
00767       adapt_flag = 0 ;  // No adaptation of the mapping 
00768       cout << "******* FROZEN MAPPING  *********" << endl ; 
00769     }
00770     
00771     // Rescaling of the grid and adaption to (outer) enthalpy surface
00772     //---------------------------------------
00773     if (adapt_flag  && (nzadapt > 0) ) 
00774       {
00775     mp_prev = mp_et ; 
00776     
00777     mp.adapt( (*outer_ent_p)(), par_adapt) ; 
00778     
00779     mp_prev.homothetie(alpha_r) ;
00780     
00781     mp.reevaluate(&mp_prev, nzet+1, ent.set()) ; 
00782     mp.reevaluate(&mp_prev, nzet+1, ent2.set()) ; 
00783       }
00784     else
00785       mp.homothetie (alpha_r);
00786 
00787 
00788     //----------------------------------------------------
00789     // Equation of state  
00790     //----------------------------------------------------
00791     
00792     equation_of_state() ;   // computes new values for nbar1,2 , ener (e) 
00793                 // and press (p) from the new ent,ent2 
00794     
00795     //---------------------------------------------------------
00796     // Matter source terms in the gravitational field equations 
00797     //---------------------------------------------------------
00798 
00799     //## Computation of tnphi and nphi from the Cartesian components
00800     //  of the shift for the test in hydro_euler():
00801         
00802     fait_nphi() ; 
00803 
00804     hydro_euler() ;     // computes new values for ener_euler (E), 
00805                 // s_euler (S) and u_euler (U^i)
00806 
00807     if (relativistic) {
00808 
00809       //-------------------------------------------------------
00810       //    2-D Poisson equation for tggg
00811       //-------------------------------------------------------
00812 
00813       mp.poisson2d(source_tggg(), mp.cmp_zero(), par_poisson_tggg, tggg.set()) ; 
00814         
00815       //-------------------------------------------------------
00816       //    2-D Poisson equation for dzeta
00817       //-------------------------------------------------------
00818 
00819       mp.poisson2d(source_dzf(), source_dzq(), par_poisson_dzeta, dzeta.set()) ; 
00820         
00821       err_grv2 = lbda_grv2 - 1; 
00822       cout << "GRV2: " << err_grv2 << endl ; 
00823         
00824     }
00825     else {
00826       err_grv2 = grv2() ; 
00827     }
00828 
00829 
00830     //---------------------------------------
00831     // Computation of the metric coefficients (except for N^phi)
00832     //---------------------------------------
00833 
00834     // Relaxations on nu and dzeta :  
00835 
00836     if (mer >= 10) {
00837         logn = relax * logn + relax_prev * logn_prev ;
00838 
00839         dzeta = relax * dzeta + relax_prev * dzeta_prev ; 
00840     }
00841 
00842     // Update of the metric coefficients N, A, B and computation of K_ij :
00843 
00844     update_metric() ; 
00845     
00846     //-----------------------
00847     //  Informations display
00848     //-----------------------
00849 
00850     //    partial_display(cout) ; 
00851     fichfreq << "  " << omega / (2*M_PI) * f_unit ; 
00852     fichfreq << "  " << omega2 / (2*M_PI) * f_unit ; 
00853     fichevol << "  " << ray_pole() / ray_eq() ; 
00854     fichevol << "  " << ent_c ; 
00855     fichevol << "  " << ent2_c ; 
00856 
00857 
00858     //-----------------------------------------
00859     // Convergence towards given baryon masses  (if mer_mass > 0)
00860     //-----------------------------------------
00861     
00862     if ((mer_mass>0) && (mer > mer_mass)) {
00863       
00864       double xx, xprog, ax, fact; 
00865 
00866       // fluid 1
00867       xx = mass_b1() / mbar1_wanted - 1. ;
00868       cout << "Discrep. baryon mass1 <-> wanted bar. mass1 : " << xx << endl ; 
00869 
00870       xprog = ( mer > 2*mer_mass) ? 1. : double(mer - mer_mass)/double(mer_mass) ; 
00871       xx *= xprog ; 
00872       ax = 0.5 * ( 2. + xx ) / (1. + xx ) ; 
00873       fact = pow(ax, aexp_mass) ; 
00874       cout << "Fluid1:  xprog, xx, ax, fact : " << xprog << "  " << xx << "  " << ax << "  " << fact << endl ; 
00875       ent_c *= fact ; 
00876 
00877       // fluid 2
00878       xx = mass_b2() / mbar2_wanted - 1. ;
00879       cout << "Discrep. baryon mass2 <-> wanted bar. mass2 : " << xx << endl ; 
00880 
00881       xprog = ( mer > 2*mer_mass) ? 1. : double(mer - mer_mass)/double(mer_mass) ; 
00882       xx *= xprog ; 
00883       ax = 0.5 * ( 2. + xx ) / (1. + xx ) ; 
00884       fact = pow(ax, aexp_mass) ; 
00885       cout << "Fluid2: xprog, xx, ax, fact : " << xprog << "  " << xx << "  " << ax << "  " << fact << endl ; 
00886       ent2_c *= fact ; 
00887 
00888     } /* if mer > mer_mass */
00889     
00890 
00891 
00892 
00893 
00894     //-------------------------------------------------------------
00895     //  Relative change in enthalpies with respect to previous step 
00896     //-------------------------------------------------------------
00897 
00898     Tbl diff_ent_tbl = diffrel( ent(), ent_prev() ) ; 
00899     diff_ent1 = diff_ent_tbl(0) ; 
00900     for (int l=1; l<nzet; l++) {
00901       diff_ent1 += diff_ent_tbl(l) ; 
00902     }
00903     diff_ent1 /= nzet ; 
00904     diff_ent_tbl = diffrel( ent2(), ent2_prev() ) ;
00905     diff_ent2 = diff_ent_tbl(0) ; 
00906     for (int l=1; l<nzet; l++) {
00907       diff_ent2 += diff_ent_tbl(l) ; 
00908     }
00909     diff_ent2 /= nzet ;
00910     diff_ent = 0.5*(diff_ent1 + diff_ent2) ; 
00911     
00912     fichconv << "  " << log10( fabs(diff_ent) + 1.e-16 ) ;
00913     fichconv << "  " << log10( fabs(err_grv2) + 1.e-16 ) ;
00914     fichconv << "  " << log10( fabs(max_triax) + 1.e-16 ) ;
00915 
00916     vit_triax = 0 ;
00917     if ( (mer > mer_triax+1) && (max_triax_prev > 1e-13) ) {
00918     vit_triax = (max_triax - max_triax_prev) / max_triax_prev ; 
00919     }
00920     
00921     fichconv << "  " << vit_triax ;
00922     
00923     //------------------------------
00924     //  Recycling for the next step
00925     //------------------------------
00926     
00927     ent_prev = ent ; 
00928     ent2_prev = ent2 ; 
00929     logn_prev = logn ; 
00930     dzeta_prev = dzeta ; 
00931     max_triax_prev = max_triax ; 
00932     
00933     fichconv << endl ;
00934     fichfreq << endl ;
00935     fichevol << endl ;
00936     fichconv.flush() ; 
00937     fichfreq.flush() ; 
00938     fichevol.flush() ; 
00939 
00940   } // End of main loop
00941 
00942   //=========================================================================
00943   //            End of iteration
00944   //=========================================================================
00945 
00946   fichconv.close() ; 
00947   fichfreq.close() ; 
00948   fichevol.close() ; 
00949     
00950 
00951 }
00952 

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