strot_dirac_equilibrium.C

00001 /*
00002  *  Function Star_rot_Dirac::equilibrium
00003  *
00004  *    (see file star_rot_dirac.h for documentation).
00005  *
00006  */
00007 
00008 /*
00009  *   Copyright (c) 2005 Lap-Ming Lin & 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 version 2
00015  *   as published by the Free Software Foundation.
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 char strot_dirac_equilibrium_C[] = "$Header: /cvsroot/Lorene/C++/Source/Star/strot_dirac_equilibrium.C,v 1.11 2009/10/26 10:54:33 j_novak Exp $" ;
00029 
00030 /*
00031  * $Id: strot_dirac_equilibrium.C,v 1.11 2009/10/26 10:54:33 j_novak Exp $
00032  * $Log: strot_dirac_equilibrium.C,v $
00033  * Revision 1.11  2009/10/26 10:54:33  j_novak
00034  * Added the case of a NONSYM base in theta.
00035  *
00036  * Revision 1.10  2008/02/25 10:40:52  j_novak
00037  * Added the flag mer_hij to control the step from which the equation for h^{ij}
00038  * is being solved.
00039  *
00040  * Revision 1.9  2006/03/14 15:18:21  lm_lin
00041  *
00042  * Add convergence to a given baryon mass.
00043  *
00044  * Revision 1.8  2005/11/28 14:45:16  j_novak
00045  * Improved solution of the Poisson tensor equation in the case of a transverse
00046  * tensor.
00047  *
00048  * Revision 1.7  2005/09/16 14:04:49  j_novak
00049  * The equation for hij is now solved only for mer >  mer_fix_omega. It uses the
00050  * Poisson solver of the class Sym_tensor_trans.
00051  *
00052  * Revision 1.6  2005/04/20 14:26:29  j_novak
00053  * Removed some outputs.
00054  *
00055  * Revision 1.5  2005/04/05 09:24:05  j_novak
00056  * minor modifs
00057  *
00058  * Revision 1.4  2005/03/10 09:39:19  j_novak
00059  * The order of resolution has been changed in the iteration step.
00060  *
00061  * Revision 1.3  2005/02/17 17:30:09  f_limousin
00062  * Change the name of some quantities to be consistent with other classes
00063  * (for instance nnn is changed to nn, shift to beta, beta to lnq...)
00064  *
00065  * Revision 1.2  2005/02/09 13:36:42  lm_lin
00066  *
00067  * Calculate GRV2 during iterations.
00068  *
00069  * Revision 1.1  2005/01/31 08:51:48  j_novak
00070  * New files for rotating stars in Dirac gauge (still under developement).
00071  *
00072  *
00073  * $Header: /cvsroot/Lorene/C++/Source/Star/strot_dirac_equilibrium.C,v 1.11 2009/10/26 10:54:33 j_novak Exp $
00074  *
00075  */
00076 
00077 
00078 // C headers
00079 #include <math.h>
00080 #include <assert.h>
00081 
00082 // Lorene headers
00083 #include "star_rot_dirac.h"
00084 
00085 #include "utilitaires.h"
00086 #include "unites.h" 
00087 
00088 void Star_rot_Dirac::equilibrium(double ent_c, double omega0, 
00089              double fact_omega, int , const Tbl& ent_limit,
00090          const Itbl& icontrol, const Tbl& control,
00091          double mbar_wanted, double aexp_mass, Tbl& diff){  
00092     
00093 
00094   // Fundamental constants and units
00095   // --------------------------------
00096   using namespace Unites ;
00097 
00098   // For the display 
00099   // ---------------
00100   char display_bold[]="x[1m" ; display_bold[0] = 27 ;
00101   char display_normal[] = "x[0m" ; display_normal[0] = 27 ;
00102 
00103 
00104   // Grid parameters
00105   // ----------------
00106 
00107   const Mg3d* mg = mp.get_mg() ;
00108   int nz = mg->get_nzone() ;    // total number of domains
00109   int nzm1 = nz - 1 ;
00110 
00111   // Index of the point at phi=0, theta=pi/2 at the surface of the star:
00112   int type_t = mg->get_type_t() ; 
00113   assert( ( type_t == SYM) || (type_t == NONSYM) ) ; 
00114   int l_b = nzet - 1 ; 
00115   int i_b = mg->get_nr(l_b) - 1 ; 
00116   int j_b = (type_t == SYM ? mg->get_nt(l_b) - 1 : mg->get_nt(l_b)/2 ) ; 
00117   int k_b = 0 ;
00118 
00119   // Value of the enthalpy defining the surface of the star
00120   double ent_b = ent_limit(nzet-1) ;
00121 
00122   // Parameters to control the iteration
00123   // -----------------------------------
00124 
00125   int mer_max = icontrol(0) ;
00126   int mer_rot = icontrol(1) ;
00127   int mer_change_omega = icontrol(2) ; 
00128   int mer_fix_omega = icontrol(3) ; 
00129   int mer_mass = icontrol(4) ; 
00130   int delta_mer_kep = icontrol(5) ; 
00131   int mer_hij = icontrol(6) ;
00132   
00133   // Protections:
00134   if (mer_change_omega < mer_rot) {
00135     cout << "Star_rot_Dirac::equilibrium: mer_change_omega < mer_rot !" 
00136      << endl ;
00137     cout << " mer_change_omega = " << mer_change_omega << endl ; 
00138     cout << " mer_rot = " << mer_rot << endl ; 
00139     abort() ; 
00140   }
00141   if (mer_fix_omega < mer_change_omega) {
00142     cout << "Star_rot_Dirac::equilibrium: mer_fix_omega < mer_change_omega !" 
00143      << endl ;
00144     cout << " mer_fix_omega = " << mer_fix_omega << endl ; 
00145     cout << " mer_change_omega = " << mer_change_omega << endl ; 
00146     abort() ; 
00147   }
00148 
00149   // In order to converge to a given baryon mass, shall the central
00150   // enthalpy be varied or Omega ?
00151   bool change_ent = true ; 
00152   if (mer_mass < 0) {
00153     change_ent = false ; 
00154     mer_mass = abs(mer_mass) ;
00155   }
00156   
00157 
00158   double precis = control(0) ; 
00159   double omega_ini = control(1) ; 
00160   double relax = control(2) ;
00161   double relax_prev = double(1) - relax ;  
00162       
00163   // Error indicators
00164   // ----------------
00165 
00166   diff.annule_hard() ;
00167   double& diff_ent = diff.set(0) ; 
00168 
00169   double alpha_r = 1 ;
00170 
00171   // Initializations
00172   // ---------------
00173 
00174   // Initial angular velocities 
00175   omega = 0 ;
00176   
00177   double accrois_omega = (omega0 - omega_ini) / 
00178                 double(mer_fix_omega - mer_change_omega) ; 
00179  
00180 
00181   update_metric() ;    //update of the metric quantities 
00182 
00183   equation_of_state() ;  // update of the density, pressure,...etc
00184 
00185   hydro_euler() ; //update of the hydro quantities relative to the 
00186                   //  Eulerian observer
00187 
00188   // Quantities at the previous step :
00189   Scalar ent_prev = ent ;
00190   Scalar logn_prev = logn ;
00191   Scalar qqq_prev = qqq ;
00192   //  Vector beta_prev = beta ;
00193   // Sym_tensor_trans hh_prev = hh ;
00194 
00195   // Output files
00196   // -------------
00197 
00198   ofstream fichconv("convergence.d") ;    // Output file for diff_ent
00199   fichconv << "#     diff_ent    GRV2    max_triax   vit_triax" << endl ;   
00200 
00201   ofstream fichfreq("frequency.d") ;    // Output file for  omega
00202   fichfreq << "#       f [Hz]" << endl ; 
00203 
00204   ofstream fichevol("evolution.d") ;    // Output file for various quantities
00205   fichevol << 
00206     "#       |dH/dr_eq/dH/dr_pole|      r_pole/r_eq ent_c" 
00207        << endl ; 
00208 
00209   diff_ent = 1 ;
00210   double err_grv2 = 1 ;
00211 
00212 
00213 
00214 //=========================================================================
00215 //                   Start of iteration
00216 //=========================================================================
00217 
00218   for(int mer=0 ; (diff_ent > precis) && (mer<mer_max) ; mer++) {
00219 
00220  cout << "-----------------------------------------------" << endl ;
00221  cout << "step: " << mer << endl ;
00222  cout << "diff_ent = " << display_bold << diff_ent << display_normal
00223       << endl ;    
00224  cout << "err_grv2 = " << err_grv2 << endl ;    
00225  fichconv << mer ;
00226  fichfreq << mer ;
00227  fichevol << mer ;
00228 
00229 
00230  // switch on rotation 
00231  if (mer >= mer_rot) {
00232           
00233    if (mer < mer_change_omega) {
00234      omega = omega_ini ; 
00235    }
00236    else {
00237      if (mer <= mer_fix_omega) {
00238        omega = omega_ini + accrois_omega * 
00239      (mer - mer_change_omega) ;
00240      }
00241    }
00242 
00243 
00244  }
00245 
00246 
00247  //---------------------------------------------------//
00248  // Resolution of the Poisson equation for logn       //
00249  // Note: ln_f is due to the fluid part               //
00250  //       ln_q is due to the quadratic metric part    //
00251  //---------------------------------------------------//
00252 
00253  Scalar ln_f_new(mp) ;
00254  Scalar ln_q_new(mp) ;
00255 
00256  solve_logn_f( ln_f_new ) ;
00257  solve_logn_q( ln_q_new ) ;
00258 
00259  ln_f_new.std_spectral_base() ;
00260  ln_q_new.std_spectral_base() ;
00261 
00262 
00263  //--------------------------------------------------//
00264  // Resolution of the Poisson equation for shift     //
00265  //--------------------------------------------------//
00266 
00267  Vector beta_new(mp, CON, mp.get_bvect_spher()) ;
00268 
00269  solve_shift( beta_new ) ;
00270 
00271  //------------------------------------
00272  // Determination of the fluid velocity
00273  //------------------------------------
00274 
00275  if (mer > mer_fix_omega + delta_mer_kep) {
00276               
00277    omega *= fact_omega ;  // Increase of the angular velocity if 
00278  }              //  fact_omega != 1
00279  
00280  bool omega_trop_grand = false ;
00281  bool kepler = true ;
00282 
00283  while ( kepler ) {
00284 
00285    // Possible decrease of Omega to ensure a velocity < c
00286 
00287    bool superlum = true ;
00288 
00289    while ( superlum ){
00290    
00291      // New fluid velocity :
00292      //
00293 
00294      u_euler.set(1).set_etat_zero() ;
00295      u_euler.set(2).set_etat_zero() ;
00296 
00297      u_euler.set(3) = omega ;
00298      u_euler.set(3).annule(nzet,nzm1) ;   // nzet is defined in class Star
00299      u_euler.set(3).std_spectral_base() ;
00300      u_euler.set(3).mult_rsint() ;
00301      u_euler.set(3) += beta(3) ;
00302      u_euler.set(3).annule(nzet,nzm1) ;
00303      
00304      u_euler = u_euler / nn ; 
00305 
00306 
00307      // v2 (square of norm of u_euler)
00308      // -------------------------------
00309 
00310      v2 = contract(contract(gamma.cov(), 0, u_euler, 0), 0, u_euler, 0) ;
00311 
00312      // Is the new velocity larger than c in the equatorial plane ?
00313 
00314      superlum = false ;
00315      
00316      for (int l=0; l<nzet; l++) {
00317        for (int i=0; i<mg->get_nr(l); i++) {
00318       
00319      double u2 = v2.val_grid_point(l, 0, j_b, i) ; 
00320      if (u2 >= 1.) {     // superluminal velocity
00321        superlum = true ; 
00322        cout << "U > c  for l, i : " << l << "  " << i 
00323         << "   U = " << sqrt( u2 ) << endl ;  
00324      }
00325        }
00326      }
00327      if ( superlum ) {
00328        cout << "**** VELOCITY OF LIGHT REACHED ****" << endl ; 
00329        omega /= fact_omega ;    // Decrease of Omega
00330        cout << "New rotation frequency : " 
00331         << omega/(2.*M_PI) * f_unit <<  " Hz" << endl ; 
00332        omega_trop_grand = true ;  
00333      }
00334    }    // end of while ( superlum )
00335 
00336 
00337    // New computation of U (this time is not superluminal)
00338    // as well as of gam_euler, ener_euler,...etc
00339 
00340    hydro_euler() ;
00341 
00342 
00343 
00344  //--------------------------------//
00345  // First integral of motion       //
00346  //--------------------------------//
00347 
00348  Scalar mlngamma(mp) ;  //  -log( gam_euler )
00349 
00350  mlngamma = - log( gam_euler ) ;
00351 
00352  // Equatorial values of various potentials :
00353  double ln_f_b = ln_f_new.val_grid_point(l_b, k_b, j_b, i_b) ; 
00354  double ln_q_b = ln_q_new.val_grid_point(l_b, k_b, j_b, i_b) ;
00355  double mlngamma_b = mlngamma.val_grid_point(l_b, k_b, j_b, i_b) ;
00356 
00357  // Central values of various potentials :
00358  double ln_f_c = ln_f_new.val_grid_point(0,0,0,0) ;
00359  double ln_q_c = ln_q_new.val_grid_point(0,0,0,0) ;
00360  double mlngamma_c = 0 ;
00361 
00362  // Scale factor to ensure that the (log of) enthalpy is equal to 
00363  // ent_b at the equator
00364  double alpha_r2 = ( ent_c - ent_b + mlngamma_c - mlngamma_b
00365              + ln_q_c - ln_q_b) / ( ln_f_b - ln_f_c  ) ;
00366 
00367  alpha_r = sqrt(alpha_r2) ;
00368 
00369  cout << "alpha_r = " << alpha_r << endl ; 
00370 
00371  // Rescaling of the grid (no adaptation!)
00372  //---------------------------------------
00373  mp.homothetie(alpha_r) ;
00374 
00375  // Readjustment of logn :
00376  // -----------------------
00377 
00378  logn = alpha_r2 * ln_f_new + ln_q_new ;
00379 
00380  double logn_c = logn.val_grid_point(0,0,0,0) ;
00381 
00382  // First integral of motion -> (log of) enthalpy in all space
00383  // ----------------------------------------------------------
00384 
00385  ent = (ent_c + logn_c + mlngamma_c) - logn - mlngamma ;
00386 
00387 
00388  // --------------------------------------------------------------
00389  // Test: is the enthalpy negative somewhere in the equatorial plane
00390  //       inside the star?
00391  // --------------------------------------------------------
00392 
00393  kepler = false ; 
00394  for (int l=0; l<nzet; l++) {
00395    int imax = mg->get_nr(l) - 1 ;
00396    if (l == l_b) imax-- ;  // The surface point is skipped
00397    for (int i=0; i<imax; i++) { 
00398      if ( ent.val_grid_point(l, 0, j_b, i) < 0. ) {
00399        kepler = true ;
00400        cout << "ent < 0 for l, i : " << l << "  " << i 
00401         << "   ent = " << ent.val_grid_point(l, 0, j_b, i) << endl ;  
00402      } 
00403    }
00404  }
00405 
00406  if ( kepler ) {
00407    cout << "**** KEPLERIAN VELOCITY REACHED ****" << endl ; 
00408    omega /= fact_omega ;    // Omega is decreased
00409    cout << "New rotation frequency : " 
00410     << omega/(2.*M_PI) * f_unit << " Hz" << endl ; 
00411    omega_trop_grand = true ;  
00412  }
00413 
00414  } // End of while ( kepler )
00415 
00416  if ( omega_trop_grand ) {   // fact_omega is decreased for the
00417                              //  next step 
00418    fact_omega = sqrt( fact_omega ) ; 
00419    cout << "**** New fact_omega : " << fact_omega << endl ; 
00420  }
00421 
00422 
00423 //---------------------------------
00424  // Equation of state
00425  //---------------------------------
00426 
00427  equation_of_state() ;  // computes new values for nbar (n), ener (e),
00428                         // and press (p) from the new ent (H)
00429 
00430  hydro_euler() ;
00431 
00432  //---------------------------------------------//
00433  // Resolution of the Poisson equation for qqq  //
00434  //---------------------------------------------//
00435 
00436  Scalar q_new(mp) ;
00437 
00438  solve_qqq( q_new ) ;
00439 
00440  q_new.std_spectral_base() ;
00441 
00442  //----------------------------------------------//
00443  // Resolution of the Poisson equation for hh    //
00444  //----------------------------------------------//
00445 
00446  Sym_tensor_trans hij_new(mp, mp.get_bvect_spher(), flat) ;
00447 
00448  if (mer > mer_hij )
00449    solve_hij( hij_new ) ;
00450  else
00451    hij_new.set_etat_zero() ;
00452 
00453  hh = hij_new ; 
00454  qqq = q_new ;
00455  beta = beta_new ; 
00456 
00457  //---------------------------------------
00458  // Calculate error of the GRV2 identity 
00459  //---------------------------------------
00460 
00461  err_grv2 = grv2() ;
00462 
00463 
00464  //--------------------------------------
00465  // Relaxations on some quantities....?
00466  //
00467  // ** On logn and qqq?
00468  //--------------------------------------
00469 
00470  if (mer >= 10) {
00471    logn = relax * logn + relax_prev * logn_prev ;
00472 
00473    qqq = relax * qqq + relax_prev * qqq_prev ;
00474 
00475  }
00476 
00477  // Update of the metric quantities :
00478 
00479    update_metric() ;
00480 
00481  //---------------------------
00482  //  Informations display
00483  // More to come later......
00484  //---------------------------
00485 
00486  // partial_display(cout) ;    // What is partial_display(cout) ? 
00487  fichfreq << "  " << omega / (2*M_PI) * f_unit ; 
00488  fichevol << "  " << ent_c ; 
00489 
00490  
00491  //-----------------------------------------
00492  // Convergence towards a given baryon mass 
00493  //-----------------------------------------
00494 
00495  if (mer > mer_mass) {
00496         
00497    double xx ; 
00498    if (mbar_wanted > 0.) {
00499      xx = mass_b() / mbar_wanted - 1. ;
00500      cout << "Discrep. baryon mass <-> wanted bar. mass : " << xx 
00501       << endl ; 
00502    }
00503    else{
00504      xx = mass_g() / fabs(mbar_wanted) - 1. ;
00505      cout << "Discrep. grav. mass <-> wanted grav. mass : " << xx 
00506       << endl ; 
00507    }
00508    double xprog = ( mer > 2*mer_mass) ? 1. : 
00509      double(mer-mer_mass)/double(mer_mass) ; 
00510    xx *= xprog ; 
00511    double ax = .5 * ( 2. + xx ) / (1. + xx ) ; 
00512    double fact = pow(ax, aexp_mass) ; 
00513    cout << "  xprog, xx, ax, fact : " << xprog << "  " <<
00514      xx << "  " << ax << "  " << fact << endl ; 
00515    
00516    if ( change_ent ) {
00517      ent_c *= fact ; 
00518    }
00519    else {
00520      if (mer%4 == 0) omega *= fact ; 
00521    }
00522  }
00523     
00524 
00525  //-----------------------------------------------------------
00526  // Relative change in enthalpy with respect to previous step
00527  // ** Check:  Is diffrel(ent, ent_prev) ok?  
00528  //-----------------------------------------------------------
00529 
00530  Tbl diff_ent_tbl = diffrel( ent, ent_prev ) ;
00531  diff_ent = diff_ent_tbl(0) ;
00532  for (int l=1; l<nzet; l++) {
00533    diff_ent += diff_ent_tbl(l) ;
00534  }
00535  diff_ent /= nzet ;
00536 
00537  fichconv << " " << log10( fabs(diff_ent) + 1.e-16 ) ;
00538  fichconv << " " << log10( fabs(err_grv2) + 1.e-16 ) ;
00539 
00540  //------------------------------
00541  // Recycling for the next step
00542  //------------------------------
00543 
00544  ent_prev = ent ;
00545  logn_prev = logn ;
00546  qqq_prev = qqq ;
00547 
00548  fichconv << endl ;     
00549  fichfreq << endl ;     
00550  fichevol << endl ;
00551  fichconv.flush() ; 
00552  fichfreq.flush() ; 
00553  fichevol.flush() ; 
00554       
00555      
00556   } // End of main loop              
00557  
00558   //=================================================
00559   //             End of iteration
00560   //=================================================
00561 
00562   fichconv.close() ; 
00563   fichfreq.close() ; 
00564   fichevol.close() ; 
00565 
00566 
00567 }

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