star_rot_equil.C

00001 /*
00002  * Method Star_rot::equilibrium
00003  *
00004  * (see file star_h.h for documentation)
00005  *
00006  */
00007 
00008 /*
00009  *   Copyright (c) 2010 Eric Gourgoulhon
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 star_rot_equil_C[] = "$Header: /cvsroot/Lorene/C++/Source/Star/star_rot_equil.C,v 1.4 2010/01/26 16:49:53 e_gourgoulhon Exp $" ;
00031 
00032 /*
00033  * $Id: star_rot_equil.C,v 1.4 2010/01/26 16:49:53 e_gourgoulhon Exp $
00034  * $Log: star_rot_equil.C,v $
00035  * Revision 1.4  2010/01/26 16:49:53  e_gourgoulhon
00036  * Reformated some outputs to the screen.
00037  *
00038  * Revision 1.3  2010/01/26 10:49:43  e_gourgoulhon
00039  * First working version.
00040  *
00041  * Revision 1.2  2010/01/25 22:32:38  e_gourgoulhon
00042  * Start of real implementation
00043  *
00044  * Revision 1.1  2010/01/25 18:15:52  e_gourgoulhon
00045  * First version.
00046  *
00047  *
00048  * $Header: /cvsroot/Lorene/C++/Source/Star/star_rot_equil.C,v 1.4 2010/01/26 16:49:53 e_gourgoulhon Exp $
00049  *
00050  */
00051 
00052 // Headers C
00053 #include <math.h>
00054 
00055 // Headers Lorene
00056 #include "star_rot.h"
00057 #include "param.h"
00058 #include "tenseur.h"
00059 
00060 #include "graphique.h"
00061 #include "utilitaires.h"
00062 #include "unites.h"
00063 
00064 void Star_rot::equilibrium(double ent_c, double omega0, double fact_omega, 
00065                  int nzadapt, const Tbl& ent_limit, const Itbl& icontrol,
00066                  const Tbl& control, double mbar_wanted, 
00067                  double aexp_mass, Tbl& diff, Param*) {
00068                  
00069     // Fundamental constants and units
00070     // -------------------------------
00071 
00072     using namespace Unites ;        
00073     
00074     // For the display 
00075     // ---------------
00076     char display_bold[]="x[1m" ; display_bold[0] = 27 ;
00077     char display_normal[] = "x[0m" ; display_normal[0] = 27 ;
00078 
00079     // Grid parameters
00080     // ---------------
00081     
00082     const Mg3d* mg = mp.get_mg() ; 
00083     int nz = mg->get_nzone() ;      // total number of domains
00084     int nzm1 = nz - 1 ; 
00085     
00086     // The following is required to initialize mp_prev as a Map_et:
00087     Map_et& mp_et = dynamic_cast<Map_et&>(mp) ; 
00088         
00089     // Index of the point at phi=0, theta=pi/2 at the surface of the star:
00090     assert(mg->get_type_t() == SYM) ; 
00091     int l_b = nzet - 1 ; 
00092     int i_b = mg->get_nr(l_b) - 1 ; 
00093     int j_b = mg->get_nt(l_b) - 1 ; 
00094     int k_b = 0 ; 
00095     
00096     // Value of the enthalpy defining the surface of the star
00097     double ent_b = ent_limit(nzet-1) ;
00098     
00099     // Parameters to control the iteration
00100     // -----------------------------------
00101     
00102     int mer_max = icontrol(0) ; 
00103     int mer_rot = icontrol(1) ;
00104     int mer_change_omega = icontrol(2) ; 
00105     int mer_fix_omega = icontrol(3) ; 
00106     int mer_mass = icontrol(4) ; 
00107     int mermax_poisson = icontrol(5) ; 
00108     int mer_triax = icontrol(6) ; 
00109     int delta_mer_kep = icontrol(7) ; 
00110 
00111     // Protections:
00112     if (mer_change_omega < mer_rot) {
00113     cout << "Star_rot::equilibrium: mer_change_omega < mer_rot !" << endl ;
00114     cout << " mer_change_omega = " << mer_change_omega << endl ; 
00115     cout << " mer_rot = " << mer_rot << endl ; 
00116     abort() ; 
00117     }
00118     if (mer_fix_omega < mer_change_omega) {
00119     cout << "Star_rot::equilibrium: mer_fix_omega < mer_change_omega !" 
00120          << endl ;
00121     cout << " mer_fix_omega = " << mer_fix_omega << endl ; 
00122     cout << " mer_change_omega = " << mer_change_omega << endl ; 
00123     abort() ; 
00124     }
00125 
00126     // In order to converge to a given baryon mass, shall the central
00127     // enthalpy be varied or Omega ?
00128     bool change_ent = true ; 
00129     if (mer_mass < 0) {
00130     change_ent = false ; 
00131     mer_mass = abs(mer_mass) ;
00132     }
00133 
00134     double precis = control(0) ; 
00135     double omega_ini = control(1) ; 
00136     double relax = control(2) ;
00137     double relax_prev = double(1) - relax ;  
00138     double relax_poisson = control(3) ; 
00139     double thres_adapt = control(4) ; 
00140     double ampli_triax = control(5) ; 
00141     double precis_adapt = control(6) ; 
00142 
00143 
00144     // Error indicators
00145     // ----------------
00146     
00147     diff.set_etat_qcq() ; 
00148     double& diff_ent = diff.set(0) ; 
00149     double& diff_nuf = diff.set(1) ; 
00150     double& diff_nuq = diff.set(2) ; 
00151 //    double& diff_dzeta = diff.set(3) ; 
00152 //    double& diff_ggg = diff.set(4) ; 
00153     double& diff_shift_x = diff.set(5) ; 
00154     double& diff_shift_y = diff.set(6) ; 
00155     double& vit_triax = diff.set(7) ; 
00156     
00157     // Parameters for the function Map_et::adapt
00158     // -----------------------------------------
00159     
00160     Param par_adapt ; 
00161     int nitermax = 100 ;  
00162     int niter ; 
00163     int adapt_flag = 1 ;    //  1 = performs the full computation, 
00164                 //  0 = performs only the rescaling by 
00165                 //      the factor alpha_r
00166     int nz_search = nzet + 1 ;  // Number of domains for searching the enthalpy
00167                 //  isosurfaces
00168     double alpha_r ; 
00169     double reg_map = 1. ; // 1 = regular mapping, 0 = contracting mapping
00170 
00171     par_adapt.add_int(nitermax, 0) ; // maximum number of iterations to
00172                      // locate zeros by the secant method
00173     par_adapt.add_int(nzadapt, 1) ; // number of domains where the adjustment 
00174                     // to the isosurfaces of ent is to be 
00175                     // performed
00176     par_adapt.add_int(nz_search, 2) ;   // number of domains to search for
00177                     // the enthalpy isosurface
00178     par_adapt.add_int(adapt_flag, 3) ; //  1 = performs the full computation, 
00179                        //  0 = performs only the rescaling by 
00180                        //      the factor alpha_r
00181     par_adapt.add_int(j_b, 4) ; //  theta index of the collocation point 
00182                     //  (theta_*, phi_*)
00183     par_adapt.add_int(k_b, 5) ; //  theta index of the collocation point 
00184                     //  (theta_*, phi_*)
00185 
00186     par_adapt.add_int_mod(niter, 0) ;  //  number of iterations actually used in 
00187                        //  the secant method
00188     
00189     par_adapt.add_double(precis_adapt, 0) ; // required absolute precision in 
00190                          // the determination of zeros by 
00191                          // the secant method
00192     par_adapt.add_double(reg_map, 1)    ;  // 1. = regular mapping, 0 = contracting mapping
00193     
00194     par_adapt.add_double(alpha_r, 2) ;  // factor by which all the radial 
00195                        // distances will be multiplied 
00196            
00197     par_adapt.add_tbl(ent_limit, 0) ;   // array of values of the field ent 
00198                         // to define the isosurfaces.              
00199 
00200     // Parameters for the function Map_et::poisson for nuf
00201     // ----------------------------------------------------
00202 
00203     double precis_poisson = 1.e-16 ;     
00204 
00205     Cmp cssjm1_nuf(ssjm1_nuf) ; 
00206     Cmp cssjm1_nuq(ssjm1_nuq) ; 
00207     Cmp cssjm1_tggg(ssjm1_tggg) ; 
00208     Cmp cssjm1_khi(ssjm1_khi) ; 
00209     Tenseur cssjm1_wshift(mp, 1, CON, mp.get_bvect_cart() ) ;
00210     cssjm1_wshift.set_etat_qcq() ;
00211     for (int i=1; i<=3; i++) {
00212     cssjm1_wshift.set(i-1) = ssjm1_wshift(i) ; 
00213     }
00214 
00215     Tenseur cshift(mp, 1, CON, mp.get_bvect_cart() ) ;
00216     cshift.set_etat_qcq() ;
00217     for (int i=1; i<=3; i++) {
00218       cshift.set(i-1) = -beta(i) ; 
00219     }
00220 
00221     Tenseur cw_shift(mp, 1, CON, mp.get_bvect_cart() ) ;
00222     cw_shift.set_etat_qcq() ;
00223     for (int i=1; i<=3; i++) {
00224       cw_shift.set(i-1) = w_shift(i) ; 
00225     }
00226 
00227     Tenseur ckhi_shift(mp) ;
00228     ckhi_shift.set_etat_qcq() ;
00229     ckhi_shift.set() = khi_shift ; 
00230 
00231     Param par_poisson_nuf ; 
00232     par_poisson_nuf.add_int(mermax_poisson,  0) ;  // maximum number of iterations
00233     par_poisson_nuf.add_double(relax_poisson,  0) ; // relaxation parameter
00234     par_poisson_nuf.add_double(precis_poisson, 1) ; // required precision
00235     par_poisson_nuf.add_int_mod(niter, 0) ;  //  number of iterations actually used 
00236     par_poisson_nuf.add_cmp_mod( cssjm1_nuf ) ; 
00237                        
00238     Param par_poisson_nuq ; 
00239     par_poisson_nuq.add_int(mermax_poisson,  0) ;  // maximum number of iterations
00240     par_poisson_nuq.add_double(relax_poisson,  0) ; // relaxation parameter
00241     par_poisson_nuq.add_double(precis_poisson, 1) ; // required precision
00242     par_poisson_nuq.add_int_mod(niter, 0) ;  //  number of iterations actually used 
00243     par_poisson_nuq.add_cmp_mod( cssjm1_nuq ) ; 
00244                        
00245     Param par_poisson_tggg ; 
00246     par_poisson_tggg.add_int(mermax_poisson,  0) ;  // maximum number of iterations
00247     par_poisson_tggg.add_double(relax_poisson,  0) ; // relaxation parameter
00248     par_poisson_tggg.add_double(precis_poisson, 1) ; // required precision
00249     par_poisson_tggg.add_int_mod(niter, 0) ;  //  number of iterations actually used 
00250     par_poisson_tggg.add_cmp_mod( cssjm1_tggg ) ; 
00251     double lambda_tggg ;
00252     par_poisson_tggg.add_double_mod( lambda_tggg ) ; 
00253     
00254     Param par_poisson_dzeta ; 
00255     double lbda_grv2 ;
00256     par_poisson_dzeta.add_double_mod( lbda_grv2 ) ; 
00257                        
00258     // Parameters for the function Scalar::poisson_vect
00259     // -------------------------------------------------
00260 
00261     Param par_poisson_vect ; 
00262 
00263     par_poisson_vect.add_int(mermax_poisson,  0) ;  // maximum number of iterations
00264     par_poisson_vect.add_double(relax_poisson,  0) ; // relaxation parameter
00265     par_poisson_vect.add_double(precis_poisson, 1) ; // required precision
00266     par_poisson_vect.add_cmp_mod( cssjm1_khi ) ; 
00267     par_poisson_vect.add_tenseur_mod( cssjm1_wshift ) ; 
00268     par_poisson_vect.add_int_mod(niter, 0) ;   
00269 
00270                        
00271     // Initializations
00272     // ---------------
00273 
00274     // Initial angular velocity
00275     omega = 0 ; 
00276   
00277     double accrois_omega = (omega0 - omega_ini) / 
00278                 double(mer_fix_omega - mer_change_omega) ; 
00279 
00280 
00281     update_metric() ;   // update of the metric coefficients
00282 
00283     equation_of_state() ;   // update of the density, pressure, etc...
00284     
00285     hydro_euler() ; // update of the hydro quantities relative to the 
00286             //  Eulerian observer
00287 
00288     // Quantities at the previous step :    
00289     Map_et mp_prev = mp_et ; 
00290     Scalar ent_prev = ent ;     
00291     Scalar logn_prev = logn ;       
00292     Scalar dzeta_prev = dzeta ;     
00293     
00294     // Creation of uninitialized tensors:
00295     Scalar source_nuf(mp) ;    // source term in the equation for nuf
00296     Scalar source_nuq(mp) ;    // source term in the equation for nuq
00297     Scalar source_dzf(mp) ; // matter source term in the eq. for dzeta
00298     Scalar source_dzq(mp) ; // quadratic source term in the eq. for dzeta
00299     Scalar source_tggg(mp) ;    // source term in the eq. for tggg
00300     Vector source_shift(mp, CON, mp.get_bvect_cart()) ;  
00301                             // source term for shift
00302     Scalar mlngamma(mp) ;   // centrifugal potential
00303     
00304             
00305     ofstream fichconv("convergence.d") ;    // Output file for diff_ent
00306     fichconv << "#     diff_ent     GRV2       max_triax      vit_triax" << endl ; 
00307     
00308     ofstream fichfreq("frequency.d") ;    // Output file for  omega
00309     fichfreq << "#       f [Hz]" << endl ; 
00310     
00311     ofstream fichevol("evolution.d") ;    // Output file for various quantities
00312     fichevol << 
00313     "#       |dH/dr_eq/dH/dr_pole|      r_pole/r_eq ent_c" 
00314     << endl ; 
00315     
00316     diff_ent = 1 ; 
00317     double err_grv2 = 1 ; 
00318     double max_triax_prev = 0 ;     // Triaxial amplitude at previous step
00319 
00320     //=========================================================================
00321     //          Start of iteration
00322     //=========================================================================
00323 
00324     for(int mer=0 ; (diff_ent > precis) && (mer<mer_max) ; mer++ ) {
00325 
00326     cout << "-----------------------------------------------" << endl ;
00327     cout << "step: " << mer << endl ;
00328     cout << "diff_ent = " << display_bold << diff_ent << display_normal
00329          << endl ;    
00330     cout << "err_grv2 = " << err_grv2 << endl ;    
00331     fichconv << mer ;
00332     fichfreq << mer ;
00333     fichevol << mer ;
00334 
00335     if (mer >= mer_rot) {
00336     
00337         if (mer < mer_change_omega) {
00338         omega = omega_ini ; 
00339         }
00340         else {
00341         if (mer <= mer_fix_omega) {
00342             omega = omega_ini + accrois_omega * 
00343                   (mer - mer_change_omega) ;
00344         }
00345         }
00346 
00347     }
00348 
00349     //-----------------------------------------------
00350     //  Sources of the Poisson equations
00351     //-----------------------------------------------
00352     
00353     // Source for nu
00354     // -------------
00355     Scalar bet = log(bbb) ; 
00356     bet.std_spectral_base() ; 
00357 
00358     Vector d_logn = logn.derive_cov( mp.flat_met_spher() ) ; 
00359     Vector d_bet = bet.derive_cov( mp.flat_met_spher() ) ; 
00360 
00361     if (relativistic) {
00362         source_nuf =  qpig * a_car *( ener_euler + s_euler ) ; 
00363 
00364         source_nuq = ak_car - d_logn(1)*(d_logn(1)+d_bet(1))
00365         - d_logn(2)*(d_logn(2)+d_bet(2))
00366         - d_logn(3)*(d_logn(3)+d_bet(3)) ; 
00367     }
00368     else {
00369         source_nuf = qpig * nbar ; 
00370 
00371         source_nuq = 0 ; 
00372     }
00373     source_nuf.std_spectral_base() ;    
00374     source_nuq.std_spectral_base() ;    
00375 
00376     // Source for dzeta
00377     // ----------------
00378     source_dzf = 2 * qpig * a_car * (press + (ener_euler+press) * uuu*uuu ) ;
00379     source_dzf.std_spectral_base() ; 
00380   
00381     source_dzq = 1.5 * ak_car 
00382         - d_logn(1)*d_logn(1) - d_logn(2)*d_logn(2) - d_logn(3)*d_logn(3) ;     
00383     source_dzq.std_spectral_base() ;    
00384     
00385     // Source for tggg
00386     // ---------------
00387     
00388     source_tggg = 4 * qpig * nn * a_car * bbb * press ;
00389     source_tggg.std_spectral_base() ; 
00390     
00391     source_tggg.mult_rsint() ; 
00392     
00393 
00394     // Source for shift
00395     // ----------------
00396         
00397     // Matter term: 
00398     source_shift = (-4*qpig) * nn * a_car  * (ener_euler + press)
00399                 * u_euler ;
00400 
00401     // Quadratic terms:
00402     Vector vtmp =  6 * bet.derive_con( mp.flat_met_spher() ) 
00403         - 2 * logn.derive_con( mp.flat_met_spher() ) ;
00404     vtmp.change_triad(mp.get_bvect_cart()) ; 
00405 
00406     Vector squad  = nn * contract(tkij, 1, vtmp, 0) ; 
00407 
00408     source_shift = source_shift + squad.up(0, mp.flat_met_cart() ) ; 
00409 
00410     //----------------------------------------------
00411     // Resolution of the Poisson equation for nuf 
00412     //----------------------------------------------
00413 
00414     source_nuf.poisson(par_poisson_nuf, nuf) ; 
00415     
00416     cout << "Test of the Poisson equation for nuf :" << endl ; 
00417     Tbl err = source_nuf.test_poisson(nuf, cout, true) ; 
00418     diff_nuf = err(0, 0) ; 
00419 
00420     //---------------------------------------
00421     // Triaxial perturbation of nuf
00422     //---------------------------------------
00423     
00424     if (mer == mer_triax) {
00425     
00426         if ( mg->get_np(0) == 1 ) {
00427         cout << 
00428         "Star_rot::equilibrium: np must be stricly greater than 1"
00429         << endl << " to set a triaxial perturbation !" << endl ; 
00430         abort() ; 
00431         }
00432     
00433         const Coord& phi = mp.phi ; 
00434         const Coord& sint = mp.sint ; 
00435         Scalar perturb(mp) ; 
00436         perturb = 1 + ampli_triax * sint*sint * cos(2*phi) ;
00437         nuf = nuf * perturb ;  
00438         
00439         nuf.std_spectral_base() ;    // set the bases for spectral expansions
00440                     //  to be the standard ones for a 
00441                     //  scalar field
00442 
00443     }
00444     
00445     // Monitoring of the triaxial perturbation
00446     // ---------------------------------------
00447     
00448     const Valeur& va_nuf = nuf.get_spectral_va() ; 
00449     va_nuf.coef() ;     // Computes the spectral coefficients
00450     double max_triax = 0 ; 
00451 
00452     if ( mg->get_np(0) > 1 ) {
00453 
00454         for (int l=0; l<nz; l++) {  // loop on the domains
00455         for (int j=0; j<mg->get_nt(l); j++) {
00456             for (int i=0; i<mg->get_nr(l); i++) {
00457         
00458             // Coefficient of cos(2 phi) : 
00459             double xcos2p = (*(va_nuf.c_cf))(l, 2, j, i) ; 
00460 
00461             // Coefficient of sin(2 phi) : 
00462             double xsin2p = (*(va_nuf.c_cf))(l, 3, j, i) ; 
00463             
00464             double xx = sqrt( xcos2p*xcos2p + xsin2p*xsin2p ) ; 
00465             
00466             max_triax = ( xx > max_triax ) ? xx : max_triax ;           
00467             }
00468         } 
00469         }
00470         
00471     }
00472     
00473     cout << "Triaxial part of nuf : " << max_triax << endl ; 
00474 
00475     if (relativistic) {
00476         
00477         //----------------------------------------------
00478         // Resolution of the Poisson equation for nuq 
00479         //----------------------------------------------
00480 
00481         source_nuq.poisson(par_poisson_nuq, nuq) ; 
00482         
00483         cout << "Test of the Poisson equation for nuq :" << endl ; 
00484         err = source_nuq.test_poisson(nuq, cout, true) ;
00485         diff_nuq = err(0, 0) ; 
00486     
00487         //---------------------------------------------------------
00488         // Resolution of the vector Poisson equation for the shift
00489         //---------------------------------------------------------
00490 
00491 
00492         for (int i=1; i<=3; i++) {
00493           if(source_shift(i).get_etat() != ETATZERO) {
00494             if(source_shift(i).dz_nonzero()) {
00495             assert( source_shift(i).get_dzpuis() == 4 ) ; 
00496             }
00497             else{
00498             (source_shift.set(i)).set_dzpuis(4) ; 
00499             }
00500         }
00501         }
00502 
00503         double lambda_shift = double(1) / double(3) ; 
00504 
00505         if ( mg->get_np(0) == 1 ) {
00506         lambda_shift = 0 ; 
00507         }
00508     
00509         Tenseur csource_shift(mp, 1, CON, mp.get_bvect_cart() ) ;
00510         csource_shift.set_etat_qcq() ;
00511         for (int i=1; i<=3; i++) {
00512         csource_shift.set(i-1) = source_shift(i) ; 
00513         }
00514         csource_shift.set(2).set_etat_zero() ;  //## bizarre...
00515 
00516         csource_shift.poisson_vect(lambda_shift, par_poisson_vect, 
00517                           cshift, cw_shift, ckhi_shift) ;           
00518 
00519         for (int i=1; i<=3; i++) {
00520         beta.set(i) = - cshift(i-1) ;
00521         beta.set(i).set_dzpuis(0) ;     //## bizarre...
00522         w_shift.set(i) = cw_shift(i-1) ; 
00523         }
00524         khi_shift = ckhi_shift() ; 
00525 
00526         cout << "Test of the Poisson equation for shift_x :" << endl ; 
00527         err = source_shift(1).test_poisson(-beta(1), cout, true) ;
00528         diff_shift_x = err(0, 0) ; 
00529     
00530         cout << "Test of the Poisson equation for shift_y :" << endl ; 
00531         err = source_shift(2).test_poisson(-beta(2), cout, true) ;
00532         diff_shift_y = err(0, 0) ; 
00533         
00534         // Computation of tnphi and nphi from the Cartesian components
00535         //  of the shift
00536         // -----------------------------------------------------------
00537         
00538         fait_nphi() ; 
00539     
00540     }
00541 
00542     //-----------------------------------------
00543     // Determination of the fluid velociy U
00544     //-----------------------------------------
00545     
00546     if (mer > mer_fix_omega + delta_mer_kep) {
00547         
00548             omega *= fact_omega ;  // Increase of the angular velocity if 
00549     }              //  fact_omega != 1
00550     
00551     bool omega_trop_grand = false ; 
00552     bool kepler = true ; 
00553 
00554     while ( kepler ) {
00555     
00556         // Possible decrease of Omega to ensure a velocity < c 
00557     
00558         bool superlum = true ; 
00559     
00560         while ( superlum ) {
00561     
00562         // New fluid velocity U :
00563     
00564         Scalar tmp = omega - nphi ; 
00565         tmp.annule_domain(nzm1) ; 
00566         tmp.std_spectral_base() ;
00567     
00568         tmp.mult_rsint() ;      //  Multiplication by r sin(theta)
00569 
00570         uuu = bbb / nn * tmp ; 
00571     
00572         if (uuu.get_etat() == ETATQCQ) {
00573             // Same basis as (Omega -N^phi) r sin(theta) :
00574             (uuu.set_spectral_va()).set_base( tmp.get_spectral_va().get_base() ) ;   
00575         }
00576 
00577         // Is the new velocity larger than c in the equatorial plane ?
00578     
00579         superlum = false ; 
00580     
00581         for (int l=0; l<nzet; l++) {
00582             for (int i=0; i<mg->get_nr(l); i++) {
00583         
00584             double u1 = uuu.val_grid_point(l, 0, j_b, i) ; 
00585             if (u1 >= 1.) {     // superluminal velocity
00586                 superlum = true ; 
00587                 cout << "U > c  for l, i : " << l << "  " << i 
00588                  << "   U = " << u1 << endl ;  
00589             }
00590             }
00591         }
00592         if ( superlum ) {
00593             cout << "**** VELOCITY OF LIGHT REACHED ****" << endl ; 
00594             omega /= fact_omega ;    // Decrease of Omega
00595             cout << "New rotation frequency : " 
00596              << omega/(2.*M_PI) * f_unit <<  " Hz" << endl ; 
00597             omega_trop_grand = true ;  
00598         }
00599         }   // end of while ( superlum )
00600 
00601     
00602         // New computation of U (which this time is not superluminal)
00603         //  as well as of gam_euler, ener_euler, etc...
00604         // -----------------------------------
00605 
00606         hydro_euler() ; 
00607     
00608 
00609         //------------------------------------------------------
00610         //  First integral of motion 
00611         //------------------------------------------------------
00612     
00613         // Centrifugal potential : 
00614         if (relativistic) {
00615         mlngamma = - log( gam_euler ) ;
00616         }
00617         else {
00618         mlngamma = - 0.5 * uuu*uuu ; 
00619         }
00620     
00621         // Equatorial values of various potentials :
00622         double nuf_b  = nuf.val_grid_point(l_b, k_b, j_b, i_b) ; 
00623         double nuq_b  = nuq.val_grid_point(l_b, k_b, j_b, i_b) ; 
00624         double mlngamma_b  = mlngamma.val_grid_point(l_b, k_b, j_b, i_b) ; 
00625 
00626         // Central values of various potentials :
00627         double nuf_c = nuf.val_grid_point(0,0,0,0) ; 
00628         double nuq_c = nuq.val_grid_point(0,0,0,0) ; 
00629         double mlngamma_c = 0 ;
00630     
00631         // Scale factor to ensure that the enthalpy is equal to ent_b at 
00632         //  the equator
00633         double alpha_r2 = ( ent_c - ent_b + mlngamma_c - mlngamma_b
00634                 + nuq_c - nuq_b) / ( nuf_b - nuf_c  ) ;
00635         alpha_r = sqrt(alpha_r2) ;
00636         cout << "alpha_r = " << alpha_r << endl ; 
00637 
00638         // Readjustment of nu :
00639         // -------------------
00640 
00641         logn = alpha_r2 * nuf + nuq ;
00642         double nu_c =  logn.val_grid_point(0,0,0,0) ;
00643 
00644         // First integral   --> enthalpy in all space
00645         //-----------------
00646 
00647         ent = (ent_c + nu_c + mlngamma_c) - logn - mlngamma ;
00648 
00649         // Test: is the enthalpy negative somewhere in the equatorial plane
00650         //  inside the star ? If yes, this means that the Keplerian velocity
00651         //  has been overstep.
00652 
00653         kepler = false ; 
00654         for (int l=0; l<nzet; l++) {
00655         int imax = mg->get_nr(l) - 1 ;
00656         if (l == l_b) imax-- ;  // The surface point is skipped
00657         for (int i=0; i<imax; i++) { 
00658             if ( ent.val_grid_point(l, 0, j_b, i) < 0. ) {
00659             kepler = true ;
00660             cout << "ent < 0 for l, i : " << l << "  " << i 
00661                  << "   ent = " << ent.val_grid_point(l, 0, j_b, i) << endl ;  
00662             } 
00663         }
00664         }
00665     
00666         if ( kepler ) {
00667         cout << "**** KEPLERIAN VELOCITY REACHED ****" << endl ; 
00668         omega /= fact_omega ;    // Omega is decreased
00669         cout << "New rotation frequency : " 
00670              << omega/(2.*M_PI) * f_unit << " Hz" << endl ; 
00671         omega_trop_grand = true ;  
00672         }
00673 
00674     }   // End of while ( kepler )
00675     
00676     if ( omega_trop_grand ) {   // fact_omega is decreased for the
00677                     //  next step 
00678         fact_omega = sqrt( fact_omega ) ; 
00679         cout << "**** New fact_omega : " << fact_omega << endl ; 
00680     }
00681 
00682     //----------------------------------------------------
00683     // Adaptation of the mapping to the new enthalpy field
00684     //----------------------------------------------------
00685     
00686     // Shall the adaptation be performed (cusp) ?
00687     // ------------------------------------------
00688     
00689     double dent_eq = ent.dsdr().val_grid_point(l_b, k_b, j_b, i_b) ; 
00690     double dent_pole = ent.dsdr().val_grid_point(l_b, k_b, 0, i_b) ;
00691     double rap_dent = fabs( dent_eq / dent_pole ) ; 
00692     cout << "| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ; 
00693     
00694     if ( rap_dent < thres_adapt ) {
00695         adapt_flag = 0 ;    // No adaptation of the mapping 
00696         cout << "******* FROZEN MAPPING  *********" << endl ; 
00697     }
00698     else{
00699         adapt_flag = 1 ;    // The adaptation of the mapping is to be
00700                 //  performed
00701     }
00702 
00703     mp_prev = mp_et ; 
00704 
00705     Cmp cent(ent) ; 
00706     mp.adapt(cent, par_adapt) ; 
00707 
00708     //----------------------------------------------------
00709     // Computation of the enthalpy at the new grid points
00710     //----------------------------------------------------
00711     
00712     mp_prev.homothetie(alpha_r) ; 
00713     
00714     mp.reevaluate(&mp_prev, nzet+1, cent) ;
00715     ent = cent ; 
00716 
00717     //----------------------------------------------------
00718     // Equation of state  
00719     //----------------------------------------------------
00720     
00721     equation_of_state() ;   // computes new values for nbar (n), ener (e) 
00722                 // and press (p) from the new ent (H)
00723     
00724     //---------------------------------------------------------
00725     // Matter source terms in the gravitational field equations 
00726     //---------------------------------------------------------
00727 
00728     //## Computation of tnphi and nphi from the Cartesian components
00729     //  of the shift for the test in hydro_euler():
00730         
00731         fait_nphi() ; 
00732 
00733     hydro_euler() ;     // computes new values for ener_euler (E), 
00734                 // s_euler (S) and u_euler (U^i)
00735 
00736     if (relativistic) {
00737 
00738         //-------------------------------------------------------
00739         //  2-D Poisson equation for tggg
00740         //-------------------------------------------------------
00741 
00742         Cmp csource_tggg(source_tggg) ; 
00743         Cmp ctggg(tggg) ; 
00744         mp.poisson2d(csource_tggg, mp.cmp_zero(), par_poisson_tggg,
00745              ctggg) ; 
00746         tggg = ctggg ; 
00747 
00748         
00749         //-------------------------------------------------------
00750         //  2-D Poisson equation for dzeta
00751         //-------------------------------------------------------
00752 
00753         Cmp csource_dzf(source_dzf) ; 
00754         Cmp csource_dzq(source_dzq) ; 
00755         Cmp cdzeta(dzeta) ; 
00756         mp.poisson2d(csource_dzf, csource_dzq, par_poisson_dzeta,
00757              cdzeta) ;
00758         dzeta = cdzeta ; 
00759         
00760         err_grv2 = lbda_grv2 - 1; 
00761         cout << "GRV2: " << err_grv2 << endl ; 
00762         
00763     }
00764     else {
00765         err_grv2 = grv2() ; 
00766     }
00767 
00768 
00769     //---------------------------------------
00770     // Computation of the metric coefficients (except for N^phi)
00771     //---------------------------------------
00772 
00773     // Relaxations on nu and dzeta :  
00774 
00775     if (mer >= 10) {
00776         logn = relax * logn + relax_prev * logn_prev ;
00777 
00778         dzeta = relax * dzeta + relax_prev * dzeta_prev ; 
00779     }
00780 
00781     // Update of the metric coefficients N, A, B and computation of K_ij :
00782 
00783     update_metric() ; 
00784     
00785     //-----------------------
00786     //  Informations display
00787     //-----------------------
00788 
00789     partial_display(cout) ; 
00790     fichfreq << "  " << omega / (2*M_PI) * f_unit ; 
00791     fichevol << "  " << rap_dent ; 
00792     fichevol << "  " << ray_pole() / ray_eq() ; 
00793     fichevol << "  " << ent_c ; 
00794 
00795     //-----------------------------------------
00796     // Convergence towards a given baryon mass 
00797     //-----------------------------------------
00798 
00799     if (mer > mer_mass) {
00800         
00801         double xx ; 
00802         if (mbar_wanted > 0.) {
00803         xx = mass_b() / mbar_wanted - 1. ;
00804         cout << "Discrep. baryon mass <-> wanted bar. mass : " << xx 
00805              << endl ; 
00806         }
00807         else{
00808         xx = mass_g() / fabs(mbar_wanted) - 1. ;
00809         cout << "Discrep. grav. mass <-> wanted grav. mass : " << xx 
00810              << endl ; 
00811         }
00812         double xprog = ( mer > 2*mer_mass) ? 1. : 
00813                  double(mer-mer_mass)/double(mer_mass) ; 
00814         xx *= xprog ; 
00815         double ax = .5 * ( 2. + xx ) / (1. + xx ) ; 
00816         double fact = pow(ax, aexp_mass) ; 
00817         cout << "  xprog, xx, ax, fact : " << xprog << "  " <<
00818             xx << "  " << ax << "  " << fact << endl ; 
00819 
00820         if ( change_ent ) {
00821         ent_c *= fact ; 
00822         }
00823         else {
00824         if (mer%4 == 0) omega *= fact ; 
00825         }
00826     }
00827     
00828 
00829     //------------------------------------------------------------
00830     //  Relative change in enthalpy with respect to previous step 
00831     //------------------------------------------------------------
00832 
00833     Tbl diff_ent_tbl = diffrel( ent, ent_prev ) ; 
00834     diff_ent = diff_ent_tbl(0) ; 
00835     for (int l=1; l<nzet; l++) {
00836         diff_ent += diff_ent_tbl(l) ; 
00837     }
00838     diff_ent /= nzet ; 
00839     
00840     fichconv << "  " << log10( fabs(diff_ent) + 1.e-16 ) ;
00841     fichconv << "  " << log10( fabs(err_grv2) + 1.e-16 ) ;
00842     fichconv << "  " << log10( fabs(max_triax) + 1.e-16 ) ;
00843 
00844     vit_triax = 0 ;
00845     if ( (mer > mer_triax+1) && (max_triax_prev > 1e-13) ) {
00846         vit_triax = (max_triax - max_triax_prev) / max_triax_prev ; 
00847     }
00848 
00849     fichconv << "  " << vit_triax ;
00850     
00851     //------------------------------
00852     //  Recycling for the next step
00853     //------------------------------
00854     
00855     ent_prev = ent ; 
00856     logn_prev = logn ; 
00857     dzeta_prev = dzeta ; 
00858     max_triax_prev = max_triax ; 
00859     
00860     fichconv << endl ;
00861     fichfreq << endl ;
00862     fichevol << endl ;
00863     fichconv.flush() ; 
00864     fichfreq.flush() ; 
00865     fichevol.flush() ; 
00866     
00867     } // End of main loop
00868     
00869     //=========================================================================
00870     //          End of iteration
00871     //=========================================================================
00872 
00873     ssjm1_nuf = cssjm1_nuf ; 
00874     ssjm1_nuq = cssjm1_nuq ; 
00875     ssjm1_tggg = cssjm1_tggg ; 
00876     ssjm1_khi = cssjm1_khi ; 
00877     for (int i=1; i<=3; i++) {
00878     ssjm1_wshift.set(i) = cssjm1_wshift(i-1) ; 
00879     }
00880 
00881     fichconv.close() ; 
00882     fichfreq.close() ; 
00883     fichevol.close() ; 
00884     
00885 }

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