et_bin_bhns_extr_equil.C

00001 /*
00002  *  Method of class Et_bin_bhns_extr to compute an equilibrium configuration
00003  *  of a BH-NS binary system with an extreme mass ratio
00004  *
00005  *    (see file et_bin_bhns_extr.h for documentation).
00006  *
00007  */
00008 
00009 /*
00010  *   Copyright (c) 2004-2005 Keisuke Taniguchi
00011  *
00012  *   This file is part of LORENE.
00013  *
00014  *   LORENE is free software; you can redistribute it and/or modify
00015  *   it under the terms of the GNU General Public License version 2
00016  *   as published by the Free Software Foundation.
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 char et_bin_bhns_extr_equil_C[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_bhns_extr_equil.C,v 1.9 2005/02/28 23:11:05 k_taniguchi Exp $" ;
00030 
00031 /*
00032  * $Id: et_bin_bhns_extr_equil.C,v 1.9 2005/02/28 23:11:05 k_taniguchi Exp $
00033  * $Log: et_bin_bhns_extr_equil.C,v $
00034  * Revision 1.9  2005/02/28 23:11:05  k_taniguchi
00035  * Modification to include the case of the conformally flat background metric
00036  *
00037  * Revision 1.8  2005/01/25 17:33:19  k_taniguchi
00038  * Suppression of the filter for the source term of the shift vector.
00039  *
00040  * Revision 1.7  2005/01/03 18:01:12  k_taniguchi
00041  * Addition of the method to fix the position of the neutron star
00042  * in the coordinate system.
00043  *
00044  * Revision 1.6  2004/12/29 16:29:55  k_taniguchi
00045  * Suppression of "dzpius" for the shift vector.
00046  *
00047  * Revision 1.5  2004/12/22 18:26:53  k_taniguchi
00048  * Change an argument of poisson_vect_falloff.
00049  *
00050  * Revision 1.4  2004/12/06 17:59:50  k_taniguchi
00051  * Change the position of resize.
00052  *
00053  * Revision 1.3  2004/12/02 21:31:56  k_taniguchi
00054  * Set a filter for the shift vector.
00055  *
00056  * Revision 1.2  2004/12/02 15:05:36  k_taniguchi
00057  * Modification of the procedure for resize.
00058  *
00059  * Revision 1.1  2004/11/30 20:48:45  k_taniguchi
00060  * *** empty log message ***
00061  *
00062  *
00063  * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_bhns_extr_equil.C,v 1.9 2005/02/28 23:11:05 k_taniguchi Exp $
00064  *
00065  */
00066 
00067 // C headers
00068 #include <math.h>
00069 
00070 // Lorene headers
00071 #include "et_bin_bhns_extr.h"
00072 #include "etoile.h"
00073 #include "map.h"
00074 #include "coord.h"
00075 #include "param.h"
00076 #include "eos.h"
00077 #include "graphique.h"
00078 #include "utilitaires.h"
00079 #include "unites.h"
00080 
00081 void Et_bin_bhns_extr::equil_bhns_extr_ks(double ent_c, const double& mass,
00082                       const double& sepa, int mermax,
00083                       int mermax_poisson,
00084                       double relax_poisson,
00085                       int mermax_potvit,
00086                       double relax_potvit, int np_filter,
00087                       double thres_adapt,
00088                       Tbl& diff) {
00089 
00090     // Fundamental constants and units
00091     // -------------------------------
00092   using namespace Unites ;
00093 
00094     assert( kerrschild == true ) ;
00095 
00096     // Initializations
00097     // ---------------
00098 
00099     const Mg3d* mg = mp.get_mg() ;
00100     int nz = mg->get_nzone() ;      // total number of domains
00101 
00102     // The following is required to initialize mp_prev as a Map_et:
00103     Map_et& mp_et = dynamic_cast<Map_et&>(mp) ;
00104 
00105     // Domain and radial indices of points at the surface of the star:
00106     int l_b = nzet - 1 ;
00107     int i_b = mg->get_nr(l_b) - 1 ;
00108     int k_b ;
00109     int j_b ;
00110 
00111     // Value of the enthalpy defining the surface of the star
00112     double ent_b = 0 ;
00113 
00114     // Error indicators
00115     // ----------------
00116 
00117     double& diff_ent = diff.set(0) ;
00118     double& diff_vel_pot = diff.set(1) ;
00119     double& diff_logn = diff.set(2) ;
00120     double& diff_beta = diff.set(3) ;
00121     double& diff_shift_x = diff.set(4) ;
00122     double& diff_shift_y = diff.set(5) ;
00123     double& diff_shift_z = diff.set(6) ;
00124 
00125     // Parameters for the function Map_et::adapt
00126     // -----------------------------------------
00127 
00128     Param par_adapt ;
00129     int nitermax = 100 ;
00130     int niter ;
00131     int adapt_flag = 1 ;    //  1 = performs the full computation,
00132                     //  0 = performs only the rescaling by
00133                 //      the factor alpha_r
00134     int nz_search = nzet ;  // Number of domains for searching
00135                 //  the enthalpy isosurfaces
00136     double precis_secant = 1.e-14 ;
00137     double alpha_r ;
00138     double reg_map = 1. ; // 1 = regular mapping, 0 = contracting mapping
00139 
00140     Tbl ent_limit(nz) ;
00141 
00142     par_adapt.add_int(nitermax, 0) ; // maximum number of iterations to
00143                      // locate zeros by the secant method
00144     par_adapt.add_int(nzet, 1) ;    // number of domains where the adjustment
00145                     // to the isosurfaces of ent is to be
00146                     // performed
00147     par_adapt.add_int(nz_search, 2) ;   // number of domains to search for
00148                     // the enthalpy isosurface
00149     par_adapt.add_int(adapt_flag, 3) ; //  1 = performs the full computation,
00150                        //  0 = performs only the rescaling by
00151                        //      the factor alpha_r
00152     par_adapt.add_int(j_b, 4) ; //  theta index of the collocation point
00153                     //  (theta_*, phi_*)
00154     par_adapt.add_int(k_b, 5) ; //  theta index of the collocation point
00155                     //  (theta_*, phi_*)
00156     par_adapt.add_int_mod(niter, 0) ;  //  number of iterations actually
00157                        //  used in the secant method
00158     par_adapt.add_double(precis_secant, 0) ; // required absolute precision in
00159                          // the determination of zeros by
00160                          // the secant method
00161     par_adapt.add_double(reg_map, 1)    ;  // 1 = regular mapping,
00162                                            // 0 = contracting mapping
00163     par_adapt.add_double(alpha_r, 2) ;      // factor by which all the radial
00164                         // distances will be multiplied
00165     par_adapt.add_tbl(ent_limit, 0) ;   // array of values of the field ent
00166                         // to define the isosurfaces
00167 
00168     // Parameters for the function Map_et::poisson for logn_auto
00169     // ---------------------------------------------------------
00170 
00171     double precis_poisson = 1.e-16 ;
00172 
00173     Param par_poisson1 ;
00174 
00175     par_poisson1.add_int(mermax_poisson,  0) ; // maximum number of iterations
00176     par_poisson1.add_double(relax_poisson,  0) ; // relaxation parameter
00177     par_poisson1.add_double(precis_poisson, 1) ; // required precision
00178     par_poisson1.add_int_mod(niter, 0) ;  // number of iterations actually
00179                                           // used 
00180     par_poisson1.add_cmp_mod( ssjm1_logn ) ;
00181 
00182     // Parameters for the function Map_et::poisson for beta_auto
00183     // ---------------------------------------------------------
00184 
00185     Param par_poisson2 ;
00186 
00187     par_poisson2.add_int(mermax_poisson,  0) ; // maximum number of iterations
00188     par_poisson2.add_double(relax_poisson,  0) ; // relaxation parameter
00189     par_poisson2.add_double(precis_poisson, 1) ; // required precision
00190     par_poisson2.add_int_mod(niter, 0) ;  // number of iterations actually
00191                                           // used 
00192     par_poisson2.add_cmp_mod( ssjm1_beta ) ;
00193 
00194     // Parameters for the function Tenseur::poisson_vect
00195     // -------------------------------------------------
00196 
00197     Param par_poisson_vect ;
00198 
00199     par_poisson_vect.add_int(mermax_poisson, 0) ; // maximum number of
00200                                                   // iterations
00201     par_poisson_vect.add_double(relax_poisson, 0) ; // relaxation parameter
00202     par_poisson_vect.add_double(precis_poisson, 1) ; // required precision
00203     par_poisson_vect.add_cmp_mod( ssjm1_khi ) ;
00204     par_poisson_vect.add_tenseur_mod( ssjm1_wshift ) ;
00205     par_poisson_vect.add_int_mod(niter, 0) ;
00206 
00207     // External potential
00208     // See Eq (99) from Gourgoulhon et al. (2001)
00209     // -----------------------------------------
00210 
00211     Tenseur pot_ext = logn_comp + pot_centri + loggam ;
00212 
00213     Tenseur ent_jm1 = ent ;  // Enthalpy at previous step
00214 
00215     Tenseur source(mp) ;    // source term in the equation for logn_auto
00216                 // and beta_auto
00217 
00218     Tenseur source_shift(mp, 1, CON, ref_triad) ;  // source term in the
00219                            // equation for shift_auto
00220 
00221     //==========================================================//
00222     //                    Start of iteration                    //
00223     //==========================================================//
00224 
00225     for(int mer=0 ; mer<mermax ; mer++ ) {
00226 
00227         cout << "-----------------------------------------------" << endl ;
00228     cout << "step: " << mer << endl ;
00229     cout << "diff_ent = " << diff_ent << endl ;
00230 
00231     //------------------------------------------------------
00232     // Resolution of the elliptic equation for the velocity
00233     // scalar potential
00234     //------------------------------------------------------
00235 
00236     if (irrotational) {
00237         diff_vel_pot = velocity_pot_extr(mass, sepa, mermax_potvit,
00238                          precis_poisson, relax_potvit) ;
00239     }
00240 
00241     //-------------------------------------
00242     // Computation of the new radial scale
00243     //-------------------------------------
00244 
00245     // alpha_r (r = alpha_r r') is determined so that the enthalpy
00246     // takes the requested value ent_b at the stellar surface
00247 
00248     // Values at the center of the star:
00249     double logn_auto_c  = logn_auto()(0, 0, 0, 0) ;
00250     double pot_ext_c  = pot_ext()(0, 0, 0, 0) ;
00251 
00252     // Search for the reference point (theta_*, phi_*)
00253     // [notation of Bonazzola, Gourgoulhon & Marck PRD 58, 104020 (1998)]
00254     // at the surface of the star
00255     // ------------------------------------------------------------------
00256     double alpha_r2 = 0 ;
00257     for (int k=0; k<mg->get_np(l_b); k++) {
00258         for (int j=0; j<mg->get_nt(l_b); j++) {
00259 
00260             double pot_ext_b  = pot_ext()(l_b, k, j, i_b) ;
00261         double logn_auto_b  = logn_auto()(l_b, k, j, i_b) ;
00262 
00263         // See Eq (100) from Gourgoulhon et al. (2001)
00264         double alpha_r2_jk = ( ent_c - ent_b + pot_ext_c - pot_ext_b)
00265           / ( logn_auto_b - logn_auto_c ) ;
00266 
00267         if (alpha_r2_jk > alpha_r2) {
00268             alpha_r2 = alpha_r2_jk ;
00269             k_b = k ;
00270             j_b = j ;
00271         }
00272 
00273         }
00274     }
00275 
00276     alpha_r = sqrt(alpha_r2) ;
00277 
00278     cout << "k_b, j_b, alpha_r: " << k_b << "  " << j_b << "  " 
00279          << alpha_r << endl ;
00280 
00281     // New value of logn_auto
00282     // ----------------------
00283 
00284     logn_auto = alpha_r2 * logn_auto ;
00285     logn_auto_regu = alpha_r2 * logn_auto_regu ;
00286     logn_auto_c  = logn_auto()(0, 0, 0, 0) ;
00287 
00288     //------------------------------------------------------------
00289     // Change the values of the inner points of the second domain
00290     // by those of the outer points of the first domain
00291     //------------------------------------------------------------
00292 
00293     (logn_auto().va).smooth(nzet, (logn_auto.set()).va) ;
00294 
00295     //--------------------------------------------
00296     // First integral --> enthalpy in all space
00297     // See Eq (98) from Gourgoulhon et al. (2001)
00298     //--------------------------------------------
00299 
00300     ent = (ent_c + logn_auto_c + pot_ext_c) - logn_auto - pot_ext ;
00301 
00302     //----------------------------------------------------------
00303     // Change the enthalpy field to be set its maximum position
00304     // at the coordinate center
00305     //----------------------------------------------------------
00306 
00307     double dentdx = ent().dsdx()(0, 0, 0, 0) ;
00308     double dentdy = ent().dsdy()(0, 0, 0, 0) ;
00309 
00310     cout << "dH/dx|_center = " << dentdx << endl ;
00311     cout << "dH/dy|_center = " << dentdy << endl ;
00312 
00313     double dec_fact = 1. ;
00314 
00315     Tenseur func_in(mp) ;
00316     func_in.set_etat_qcq() ;
00317     func_in.set() = 1. - dec_fact * (dentdx/ent_c) * mp.x
00318       - dec_fact * (dentdy/ent_c) * mp.y ;
00319     func_in.set().annule(nzet, nz-1) ;
00320     func_in.set_std_base() ;
00321 
00322     Tenseur func_ex(mp) ;
00323     func_ex.set_etat_qcq() ;
00324     func_ex.set() = 1. ;
00325     func_ex.set().annule(0, nzet-1) ;
00326     func_ex.set_std_base() ;
00327 
00328     // New enthalpy field
00329     // ------------------
00330     ent.set() = ent() * (func_in() + func_ex()) ;
00331 
00332     (ent().va).smooth(nzet, (ent.set()).va) ;
00333 
00334     double dentdx_new = ent().dsdx()(0, 0, 0, 0) ;
00335     double dentdy_new = ent().dsdy()(0, 0, 0, 0) ;
00336     cout << "dH/dx|_new    = " << dentdx_new << endl ;
00337     cout << "dH/dy|_new    = " << dentdy_new << endl ;
00338 
00339     //-----------------------------------------------------
00340     // Adaptation of the mapping to the new enthalpy field
00341     //----------------------------------------------------
00342 
00343     // Shall the adaptation be performed (cusp) ?
00344     // ------------------------------------------
00345 
00346     double dent_eq = ent().dsdr().val_point(ray_eq_pi(),M_PI/2.,M_PI) ;
00347     double dent_pole = ent().dsdr().val_point(ray_pole(),0.,0.) ;
00348     double rap_dent = fabs( dent_eq / dent_pole ) ;
00349     cout << "| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ;
00350 
00351     if ( rap_dent < thres_adapt ) {
00352         adapt_flag = 0 ;    // No adaptation of the mapping
00353         cout << "******* FROZEN MAPPING  *********" << endl ;
00354     }
00355     else{
00356         adapt_flag = 1 ;    // The adaptation of the mapping is to be
00357                             // performed
00358     }
00359 
00360     ent_limit.set_etat_qcq() ;
00361     for (int l=0; l<nzet; l++) {    // loop on domains inside the star
00362         ent_limit.set(l) = ent()(l, k_b, j_b, i_b) ;
00363     }
00364 
00365     ent_limit.set(nzet-1) = ent_b ;
00366     Map_et mp_prev = mp_et ;
00367 
00368     mp.adapt(ent(), par_adapt) ;
00369 
00370     //----------------------------------------------------
00371     // Computation of the enthalpy at the new grid points
00372     //----------------------------------------------------
00373 
00374     mp_prev.homothetie(alpha_r) ;
00375 
00376     mp.reevaluate(&mp_prev, nzet+1, ent.set()) ;
00377 
00378     double fact_resize = 1. / alpha_r ;
00379     for (int l=nzet; l<nz-1; l++) {
00380         mp_et.resize(l, fact_resize) ;
00381     }
00382     mp_et.resize_extr(fact_resize) ;
00383 
00384     double ent_s_max = -1 ;
00385     int k_s_max = -1 ;
00386     int j_s_max = -1 ;
00387     for (int k=0; k<mg->get_np(l_b); k++) {
00388         for (int j=0; j<mg->get_nt(l_b); j++) {
00389             double xx = fabs( ent()(l_b, k, j, i_b) ) ;
00390         if (xx > ent_s_max) {
00391             ent_s_max = xx ;
00392             k_s_max = k ;
00393             j_s_max = j ;
00394         }
00395         }
00396     }
00397     cout << "Max. abs(enthalpy) at the boundary between domains nzet-1"
00398          << " and nzet : " << endl ;
00399     cout << "   " << ent_s_max << " reached for k = " << k_s_max
00400          << " and j = " << j_s_max << endl ;
00401 
00402     //-------------------
00403     // Equation of state
00404     //-------------------
00405 
00406     equation_of_state() ;   // computes new values for nbar (n), ener (e)
00407                             // and press (p) from the new ent (H)
00408 
00409     //----------------------------------------------------------
00410     // Matter source terms in the gravitational field equations
00411     //---------------------------------------------------------
00412 
00413     hydro_euler_extr(mass, sepa) ;  // computes new values for
00414                                     // ener_euler (E), s_euler (S)
00415                                     // and u_euler (U^i)
00416 
00417 
00418     //----------------------------------------------------
00419     // Auxiliary terms for the source of Poisson equation
00420     //----------------------------------------------------
00421 
00422     const Coord& xx = mp.x ;
00423     const Coord& yy = mp.y ;
00424     const Coord& zz = mp.z ;
00425 
00426     Tenseur r_bh(mp) ;
00427     r_bh.set_etat_qcq() ;
00428     r_bh.set() = pow( (xx+sepa)*(xx+sepa) + yy*yy + zz*zz, 0.5) ;
00429     r_bh.set_std_base() ;
00430 
00431     Tenseur xx_cov(mp, 1, COV, ref_triad) ;
00432     xx_cov.set_etat_qcq() ;
00433     xx_cov.set(0) = xx + sepa ;
00434     xx_cov.set(1) = yy ;
00435     xx_cov.set(2) = zz ;
00436     xx_cov.set_std_base() ;
00437 
00438     Tenseur xsr_cov(mp, 1, COV, ref_triad) ;
00439     xsr_cov = xx_cov / r_bh ;
00440     xsr_cov.set_std_base() ;
00441 
00442     Tenseur msr(mp) ;
00443     msr = ggrav * mass / r_bh ;
00444     msr.set_std_base() ;
00445 
00446     Tenseur lapse_bh(mp) ;
00447     lapse_bh = 1. / sqrt( 1.+2.*msr ) ;
00448     lapse_bh.set_std_base() ;
00449 
00450     Tenseur lapse_bh2(mp) ;  // lapse_bh * lapse_bh
00451     lapse_bh2 = 1. / (1.+2.*msr) ;
00452     lapse_bh2.set_std_base() ;
00453 
00454     Tenseur ldnu(mp) ;
00455     ldnu = flat_scalar_prod_desal(xsr_cov, d_logn_auto) ;
00456     ldnu.set_std_base() ;
00457 
00458     Tenseur ldbeta(mp) ;
00459     ldbeta = flat_scalar_prod_desal(xsr_cov, d_beta_auto) ;
00460     ldbeta.set_std_base() ;
00461 
00462     Tenseur lltkij(mp) ;
00463     lltkij.set_etat_qcq() ;
00464     lltkij.set() = 0 ;
00465     lltkij.set_std_base() ;
00466 
00467     for (int i=0; i<3; i++)
00468         for (int j=0; j<3; j++)
00469             lltkij.set() += xsr_cov(i) * xsr_cov(j) * tkij_auto(i, j) ;
00470 
00471     Tenseur lshift(mp) ;
00472     lshift = flat_scalar_prod_desal(xsr_cov, shift_auto) ;
00473     lshift.set_std_base() ;
00474 
00475     Tenseur d_ldnu(mp, 1, COV, ref_triad) ;
00476     d_ldnu = ldnu.gradient() ;        // (d/dx, d/dy, d/dz)
00477     d_ldnu.change_triad(ref_triad) ;  // --> (d/dX, d/dY, d/dZ)
00478 
00479     Tenseur ldldnu(mp) ;
00480     ldldnu = flat_scalar_prod_desal(xsr_cov, d_ldnu) ;
00481     ldldnu.set_std_base() ;
00482 
00483     Tenseur d_ldbeta(mp, 1, COV, ref_triad) ;
00484     d_ldbeta = ldbeta.gradient() ;      // (d/dx, d/dy, d/dz)
00485     d_ldbeta.change_triad(ref_triad) ;  // --> (d/dX, d/dY, d/dZ)
00486 
00487     Tenseur ldldbeta(mp) ;
00488     ldldbeta = flat_scalar_prod_desal(xsr_cov, d_ldbeta) ;
00489     ldldbeta.set_std_base() ;
00490 
00491     //------------------------------------------
00492     // Poisson equation for logn_auto (nu_auto)
00493     //------------------------------------------
00494 
00495     // Source
00496     // ------
00497 
00498     if (relativistic) {
00499 
00500         source = qpig * a_car % (ener_euler + s_euler) + akcar_auto
00501           - flat_scalar_prod_desal(d_logn_auto, d_beta_auto)
00502           + 2.*lapse_bh2 % msr % (ldnu % ldbeta + ldldnu)
00503           + lapse_bh2 % lapse_bh2 % msr % (2.*(ldnu + 4.*msr % ldnu)
00504                            - ldbeta) / r_bh
00505           - (4.*a_car % lapse_bh2 % lapse_bh2 % msr / 3. / nnn / r_bh)
00506           * (2.+3.*msr) * (3.+4.*msr) * lltkij
00507           + (2.*a_car % lapse_bh2 % lapse_bh2 % lapse_bh % msr
00508          / nnn / r_bh / r_bh) * (2.+10.*msr+9.*msr%msr) * lshift
00509           + (4.*pow(lapse_bh2, 3.) % msr % msr / 3. / r_bh / r_bh)
00510           % (2.*(a_car%lapse_bh2/nnn/nnn - 1.) * pow(2.+3.*msr, 2.)
00511          + (a_car - 1.) % pow(1.+3.*msr, 2.)
00512          - 3.*(a_car%lapse_bh/nnn - 1.)*(2.+10.*msr+9.*msr%msr)) ;
00513 
00514     }
00515     else {
00516         cout << "The computation of BH-NS binary systems"
00517          << " should be relativistic !!!" << endl ;
00518         abort() ;
00519     }
00520 
00521     source.set_std_base() ;
00522 
00523     // Resolution of the Poisson equation
00524     // ----------------------------------
00525 
00526     int k_falloff = 1 ;
00527 
00528     source().poisson_falloff(par_poisson1, logn_auto.set(), k_falloff) ;
00529 
00530     // Construct logn_auto_regu for et_bin_upmetr_extr.C
00531     // -------------------------------------------------
00532 
00533     logn_auto_regu = logn_auto ;
00534 
00535     // Check: has the Poisson equation been correctly solved ?
00536     // -----------------------------------------------------
00537 
00538     Tbl tdiff_logn = diffrel(logn_auto().laplacien(), source()) ;
00539 
00540     cout <<
00541       "Relative error in the resolution of the equation for logn_auto : "
00542          << endl ;
00543 
00544     for (int l=0; l<nz; l++) {
00545         cout << tdiff_logn(l) << "  " ;
00546     }
00547     cout << endl ;
00548     diff_logn = max(abs(tdiff_logn)) ;
00549 
00550     if (relativistic) {
00551 
00552         //--------------------------------
00553         // Poisson equation for beta_auto
00554         //--------------------------------
00555 
00556         // Source
00557         // ------
00558 
00559         source = qpig * a_car % s_euler + 0.75 * akcar_auto
00560           - 0.5 * flat_scalar_prod_desal(d_logn_auto, d_logn_auto)
00561           - 0.5 * flat_scalar_prod_desal(d_beta_auto, d_beta_auto)
00562           + lapse_bh2 % msr % (ldnu%ldnu + ldbeta%ldbeta + 2.*ldldbeta)
00563           + lapse_bh2 % lapse_bh2 % msr % (2.*(1.+4.*msr) * ldbeta
00564                            - ldnu) / r_bh
00565           - (a_car % lapse_bh2 % lapse_bh2 % msr / nnn / r_bh)
00566           * (2.+3.*msr) * (3.+4.*msr) * lltkij
00567           + (2.*a_car % lapse_bh2 % lapse_bh2 % lapse_bh % msr
00568          / nnn / r_bh / r_bh) * (2.+10.*msr+9.*msr%msr) * lshift
00569           + (2.*pow(lapse_bh2, 3.) % msr % msr / r_bh / r_bh)
00570           % ((a_car%lapse_bh2/nnn/nnn - 1.) * pow(2.+3.*msr, 2.)
00571          + (a_car - 1.) * pow(1.+3.*msr, 2.)
00572          - 2.*(a_car%lapse_bh/nnn - 1.)*(2.+10.*msr+9.*msr%msr)) ;
00573 
00574         source.set_std_base() ;
00575 
00576         // Resolution of the Poisson equation
00577         // ----------------------------------
00578 
00579         source().poisson_falloff(par_poisson2, beta_auto.set(),
00580                      k_falloff) ;
00581 
00582         // Check: has the Poisson equation been correctly solved ?
00583         // -----------------------------------------------------
00584 
00585         Tbl tdiff_beta = diffrel(beta_auto().laplacien(), source()) ;
00586 
00587         cout << "Relative error in the resolution of the equation for "
00588          << "beta_auto : " << endl ;
00589         for (int l=0; l<nz; l++) {
00590             cout << tdiff_beta(l) << "  " ;
00591         }
00592         cout << endl ;
00593         diff_beta = max(abs(tdiff_beta)) ;
00594 
00595         //----------------------------------------
00596         // Vector Poisson equation for shift_auto
00597         //----------------------------------------
00598 
00599         // Some auxiliary terms for the source
00600         // -----------------------------------
00601 
00602         Tenseur xx_con(mp, 1, CON, ref_triad) ;
00603         xx_con.set_etat_qcq() ;
00604         xx_con.set(0) = xx + sepa ;
00605         xx_con.set(1) = yy ;
00606         xx_con.set(2) = zz ;
00607         xx_con.set_std_base() ;
00608 
00609         Tenseur xsr_con(mp, 1, CON, ref_triad) ;
00610         xsr_con = xx_con / r_bh ;
00611         xsr_con.set_std_base() ;
00612 
00613         // Components of shift_auto with respect to the Cartesian triad
00614         //  (d/dx, d/dy, d/dz) of the mapping :
00615 
00616         Tenseur shift_auto_local = shift_auto ;
00617         shift_auto_local.change_triad( mp.get_bvect_cart() ) ;
00618 
00619         // Gradient (partial derivatives with respect to the Cartesian
00620         //           coordinates of the mapping)
00621         // dn(i, j) = D_i N^j
00622 
00623         Tenseur dn = shift_auto_local.gradient() ;
00624 
00625         // Return to the absolute reference frame
00626         dn.change_triad(ref_triad) ;
00627 
00628         // Trace of D_i N^j = divergence of N^j :
00629         Tenseur divn = contract(dn, 0, 1) ;
00630 
00631         // l^j D_j N^i
00632         Tenseur ldn_con = contract(xsr_con, 0, dn, 0) ;
00633 
00634         // D_j (l^k D_k N^i): dldn(j, i)
00635         Tenseur ldn_local = ldn_con ;
00636         ldn_local.change_triad( mp.get_bvect_cart() ) ;
00637         Tenseur dldn = ldn_local.gradient() ;
00638         dldn.change_triad(ref_triad) ;
00639 
00640         // l^j D_j (l^k D_k N^i)
00641         Tenseur ldldn = contract(xsr_con, 0, dldn, 0) ;
00642 
00643         // l_k D_j N^k
00644         Tenseur ldn_cov = contract(xsr_cov, 0, dn, 1) ;
00645 
00646         // l^j l_k D_j N^k
00647         Tenseur lldn_cov = contract(xsr_con, 0, ldn_cov, 0) ;
00648 
00649         // eta^{ij} l_k D_j N^k
00650         Tenseur eldn(mp, 1, CON, ref_triad) ;
00651         eldn.set_etat_qcq() ;
00652         eldn.set(0) = ldn_cov(0) ;
00653         eldn.set(1) = ldn_cov(1) ;
00654         eldn.set(2) = ldn_cov(2) ;
00655         eldn.set_std_base() ;
00656 
00657         // l^i D_j N^j
00658         Tenseur ldivn = xsr_con % divn ;
00659 
00660         // D_j (l^i D_k N^k): dldivn(j, i)
00661         Tenseur ldivn_local = ldivn ;
00662         ldivn_local.change_triad( mp.get_bvect_cart() ) ;
00663         Tenseur dldivn = ldivn_local.gradient() ;
00664         dldivn.change_triad(ref_triad) ;
00665 
00666         // l^j D_j (l^i D_k N^k)
00667         Tenseur ldldivn = contract(xsr_con, 0, dldivn, 0) ;
00668 
00669         // l_j N^j
00670         Tenseur ln = contract(xsr_cov, 0, shift_auto, 0) ;
00671 
00672         Tenseur vtmp =  6. * d_beta_auto - 8. * d_logn_auto ;
00673 
00674         Tenseur lvtmp = contract(xsr_con, 0, vtmp, 0) ;
00675 
00676         // eta^{ij} vtmp_j
00677         Tenseur evtmp(mp, 1, CON, ref_triad) ;
00678         evtmp.set_etat_qcq() ;
00679         evtmp.set(0) = vtmp(0) ;
00680         evtmp.set(1) = vtmp(1) ;
00681         evtmp.set(2) = vtmp(2) ;
00682         evtmp.set_std_base() ;
00683 
00684         // lapse_ns
00685         Tenseur lapse_ns(mp) ;
00686         lapse_ns = exp(logn_auto) ;
00687         lapse_ns.set_std_base() ;
00688 
00689         // Source
00690         // ------
00691 
00692         source_shift = (-4.*qpig) * nnn % a_car % (ener_euler + press)
00693           % u_euler
00694           + nnn % flat_scalar_prod_desal(tkij_auto, vtmp)
00695           - 2.*nnn % lapse_bh2 * msr / r_bh
00696           % flat_scalar_prod_desal(tkij_auto, xsr_cov)
00697           + 2.*lapse_bh2 * msr * (3.*ldldn + ldldivn) / 3.
00698           - lapse_bh2 * msr / r_bh
00699           * (4.*ldivn - lapse_bh2 % (3.*ldn_con + 8.*msr * ldn_con)
00700          - (eldn + 2.*lapse_bh2*(9.+11.*msr)*lldn_cov%xsr_con) / 3.)
00701           - 2.*lapse_bh2 % lapse_bh2 * msr / r_bh / r_bh
00702           * ( (4.+11.*msr) * shift_auto
00703           - lapse_bh2 * (12.+51.*msr+46.*msr*msr) * ln % xsr_con )
00704           / 3.
00705           + 8.*pow(lapse_bh2, 4.) * msr / r_bh / r_bh
00706           % (lapse_ns - 1.) * (2.+10.*msr+9.*msr*msr) * xsr_con / 3.
00707           + 2.*pow(lapse_bh2, 3.) * msr / r_bh * (2.+3.*msr)
00708           * ( (1.+2.*msr) * evtmp - (3.+2.*msr) * lvtmp * xsr_con) / 3. ;
00709 
00710         source_shift.set_std_base() ;
00711 
00712         // Resolution of the Poisson equation
00713         // ----------------------------------
00714 
00715         // Filter for the source of shift vector :
00716 
00717         for (int i=0; i<3; i++) {
00718             for (int l=0; l<nz; l++) {
00719             if (source_shift(i).get_etat() != ETATZERO)
00720             source_shift.set(i).filtre_phi(np_filter, l) ;
00721         }
00722         }
00723 
00724         // For Tenseur::poisson_vect, the triad must be the mapping
00725         // triad, not the reference one:
00726 
00727         source_shift.change_triad( mp.get_bvect_cart() ) ;
00728         /*
00729         for (int i=0; i<3; i++) {
00730             if(source_shift(i).dz_nonzero()) {
00731             assert( source_shift(i).get_dzpuis() == 4 ) ;
00732         }
00733         else {
00734             (source_shift.set(i)).set_dzpuis(4) ;
00735         }
00736         }
00737 
00738         source_shift.dec2_dzpuis() ;    // dzpuis 4 -> 2
00739         */
00740         double lambda_shift = double(1) / double(3) ;
00741 
00742         int* shift_falloff ;
00743         shift_falloff = new int[4] ;
00744         shift_falloff[0] = 1 ;
00745         shift_falloff[1] = 1 ;
00746         shift_falloff[2] = 2 ;
00747         shift_falloff[3] = 1 ;
00748 
00749         source_shift.poisson_vect_falloff(lambda_shift, par_poisson_vect,
00750                           shift_auto, w_shift,
00751                           khi_shift, shift_falloff) ;
00752 
00753         delete[] shift_falloff ;
00754 
00755         // Check: has the equation for shift_auto been correctly solved ?
00756         // --------------------------------------------------------------
00757 
00758         // Divergence of shift_auto :
00759         Tenseur divna = contract(shift_auto.gradient(), 0, 1) ;
00760         //      divna.dec2_dzpuis() ;    // dzpuis 2 -> 0
00761 
00762         // Grad(div) :
00763         Tenseur graddivn = divna.gradient() ;
00764         //      graddivn.inc2_dzpuis() ;    // dzpuis 2 -> 4
00765 
00766         // Full operator :
00767         Tenseur lap_shift(mp, 1, CON, mp.get_bvect_cart() ) ;
00768         lap_shift.set_etat_qcq() ;
00769         for (int i=0; i<3; i++) {
00770             lap_shift.set(i) = shift_auto(i).laplacien()
00771           + lambda_shift * graddivn(i) ;
00772         }
00773 
00774         Tbl tdiff_shift_x = diffrel(lap_shift(0), source_shift(0)) ;
00775         Tbl tdiff_shift_y = diffrel(lap_shift(1), source_shift(1)) ;
00776         Tbl tdiff_shift_z = diffrel(lap_shift(2), source_shift(2)) ;
00777 
00778         cout << "Relative error in the resolution of the equation for "
00779          << "shift_auto : " << endl ;
00780         cout << "x component : " ;
00781         for (int l=0; l<nz; l++) {
00782             cout << tdiff_shift_x(l) << "  " ;
00783         }
00784         cout << endl ;
00785         cout << "y component : " ;
00786         for (int l=0; l<nz; l++) {
00787             cout << tdiff_shift_y(l) << "  " ;
00788         }
00789         cout << endl ;
00790         cout << "z component : " ;
00791         for (int l=0; l<nz; l++) {
00792             cout << tdiff_shift_z(l) << "  " ;
00793         }
00794         cout << endl ;
00795 
00796         diff_shift_x = max(abs(tdiff_shift_x)) ;
00797         diff_shift_y = max(abs(tdiff_shift_y)) ;
00798         diff_shift_z = max(abs(tdiff_shift_z)) ;
00799 
00800         // Final result
00801         // ------------
00802         // The output of Tenseur::poisson_vect_falloff is on the mapping
00803         // triad, it should therefore be transformed to components on the
00804         // reference triad :
00805 
00806         shift_auto.change_triad( ref_triad ) ;
00807 
00808     }   // End of relativistic equations
00809 
00810     //------------------------------
00811     //  Relative change in enthalpy
00812     //------------------------------
00813 
00814     Tbl diff_ent_tbl = diffrel( ent(), ent_jm1() ) ;
00815     diff_ent = diff_ent_tbl(0) ;
00816     for (int l=1; l<nzet; l++) {
00817         diff_ent += diff_ent_tbl(l) ;
00818     }
00819     diff_ent /= nzet ;
00820 
00821     ent_jm1 = ent ;
00822 
00823     } // End of main loop
00824 
00825     //========================================================//
00826     //                    End of iteration                    //
00827     //========================================================//
00828 
00829 }
00830 
00831 
00832 void Et_bin_bhns_extr::equil_bhns_extr_cf(double ent_c, const double& mass,
00833                       const double& sepa, int mermax,
00834                       int mermax_poisson,
00835                       double relax_poisson,
00836                       int mermax_potvit,
00837                       double relax_potvit, int np_filter,
00838                       double thres_adapt,
00839                       Tbl& diff) {
00840 
00841     // Fundamental constants and units
00842     // -------------------------------
00843   using namespace Unites ;
00844 
00845     assert( kerrschild == false ) ;
00846 
00847     // Initializations
00848     // ---------------
00849 
00850     const Mg3d* mg = mp.get_mg() ;
00851     int nz = mg->get_nzone() ;      // total number of domains
00852 
00853     // The following is required to initialize mp_prev as a Map_et:
00854     Map_et& mp_et = dynamic_cast<Map_et&>(mp) ;
00855 
00856     // Domain and radial indices of points at the surface of the star:
00857     int l_b = nzet - 1 ;
00858     int i_b = mg->get_nr(l_b) - 1 ;
00859     int k_b ;
00860     int j_b ;
00861 
00862     // Value of the enthalpy defining the surface of the star
00863     double ent_b = 0 ;
00864 
00865     // Error indicators
00866     // ----------------
00867 
00868     double& diff_ent = diff.set(0) ;
00869     double& diff_vel_pot = diff.set(1) ;
00870     double& diff_logn = diff.set(2) ;
00871     double& diff_beta = diff.set(3) ;
00872     double& diff_shift_x = diff.set(4) ;
00873     double& diff_shift_y = diff.set(5) ;
00874     double& diff_shift_z = diff.set(6) ;
00875 
00876     // Parameters for the function Map_et::adapt
00877     // -----------------------------------------
00878 
00879     Param par_adapt ;
00880     int nitermax = 100 ;
00881     int niter ;
00882     int adapt_flag = 1 ;    //  1 = performs the full computation,
00883                     //  0 = performs only the rescaling by
00884                 //      the factor alpha_r
00885     int nz_search = nzet ;  // Number of domains for searching
00886                 //  the enthalpy isosurfaces
00887     double precis_secant = 1.e-14 ;
00888     double alpha_r ;
00889     double reg_map = 1. ; // 1 = regular mapping, 0 = contracting mapping
00890 
00891     Tbl ent_limit(nz) ;
00892 
00893     par_adapt.add_int(nitermax, 0) ; // maximum number of iterations to
00894                      // locate zeros by the secant method
00895     par_adapt.add_int(nzet, 1) ;    // number of domains where the adjustment
00896                     // to the isosurfaces of ent is to be
00897                     // performed
00898     par_adapt.add_int(nz_search, 2) ;   // number of domains to search for
00899                     // the enthalpy isosurface
00900     par_adapt.add_int(adapt_flag, 3) ; //  1 = performs the full computation,
00901                        //  0 = performs only the rescaling by
00902                        //      the factor alpha_r
00903     par_adapt.add_int(j_b, 4) ; //  theta index of the collocation point
00904                     //  (theta_*, phi_*)
00905     par_adapt.add_int(k_b, 5) ; //  theta index of the collocation point
00906                     //  (theta_*, phi_*)
00907     par_adapt.add_int_mod(niter, 0) ;  //  number of iterations actually
00908                        //  used in the secant method
00909     par_adapt.add_double(precis_secant, 0) ; // required absolute precision in
00910                          // the determination of zeros by
00911                          // the secant method
00912     par_adapt.add_double(reg_map, 1)    ;  // 1 = regular mapping,
00913                                            // 0 = contracting mapping
00914     par_adapt.add_double(alpha_r, 2) ;      // factor by which all the radial
00915                         // distances will be multiplied
00916     par_adapt.add_tbl(ent_limit, 0) ;   // array of values of the field ent
00917                         // to define the isosurfaces
00918 
00919     // Parameters for the function Map_et::poisson for logn_auto
00920     // ---------------------------------------------------------
00921 
00922     double precis_poisson = 1.e-16 ;
00923 
00924     Param par_poisson1 ;
00925 
00926     par_poisson1.add_int(mermax_poisson,  0) ; // maximum number of iterations
00927     par_poisson1.add_double(relax_poisson,  0) ; // relaxation parameter
00928     par_poisson1.add_double(precis_poisson, 1) ; // required precision
00929     par_poisson1.add_int_mod(niter, 0) ;  // number of iterations actually
00930                                           // used 
00931     par_poisson1.add_cmp_mod( ssjm1_logn ) ;
00932 
00933     // Parameters for the function Map_et::poisson for beta_auto
00934     // ---------------------------------------------------------
00935 
00936     Param par_poisson2 ;
00937 
00938     par_poisson2.add_int(mermax_poisson,  0) ; // maximum number of iterations
00939     par_poisson2.add_double(relax_poisson,  0) ; // relaxation parameter
00940     par_poisson2.add_double(precis_poisson, 1) ; // required precision
00941     par_poisson2.add_int_mod(niter, 0) ;  // number of iterations actually
00942                                           // used 
00943     par_poisson2.add_cmp_mod( ssjm1_beta ) ;
00944 
00945     // Parameters for the function Tenseur::poisson_vect
00946     // -------------------------------------------------
00947 
00948     Param par_poisson_vect ;
00949 
00950     par_poisson_vect.add_int(mermax_poisson, 0) ; // maximum number of
00951                                                   // iterations
00952     par_poisson_vect.add_double(relax_poisson, 0) ; // relaxation parameter
00953     par_poisson_vect.add_double(precis_poisson, 1) ; // required precision
00954     par_poisson_vect.add_cmp_mod( ssjm1_khi ) ;
00955     par_poisson_vect.add_tenseur_mod( ssjm1_wshift ) ;
00956     par_poisson_vect.add_int_mod(niter, 0) ;
00957 
00958     // External potential
00959     // See Eq (99) from Gourgoulhon et al. (2001)
00960     // -----------------------------------------
00961 
00962     Tenseur pot_ext = logn_comp + pot_centri + loggam ;
00963 
00964     Tenseur ent_jm1 = ent ;  // Enthalpy at previous step
00965 
00966     Tenseur source(mp) ;    // source term in the equation for logn_auto
00967                 // and beta_auto
00968 
00969     Tenseur source_shift(mp, 1, CON, ref_triad) ;  // source term in the
00970                            // equation for shift_auto
00971 
00972     //==========================================================//
00973     //                    Start of iteration                    //
00974     //==========================================================//
00975 
00976     for(int mer=0 ; mer<mermax ; mer++ ) {
00977 
00978         cout << "-----------------------------------------------" << endl ;
00979     cout << "step: " << mer << endl ;
00980     cout << "diff_ent = " << diff_ent << endl ;
00981 
00982     //------------------------------------------------------
00983     // Resolution of the elliptic equation for the velocity
00984     // scalar potential
00985     //------------------------------------------------------
00986 
00987     if (irrotational) {
00988         diff_vel_pot = velocity_pot_extr(mass, sepa, mermax_potvit,
00989                          precis_poisson, relax_potvit) ;
00990     }
00991 
00992     //-------------------------------------
00993     // Computation of the new radial scale
00994     //-------------------------------------
00995 
00996     // alpha_r (r = alpha_r r') is determined so that the enthalpy
00997     // takes the requested value ent_b at the stellar surface
00998 
00999     // Values at the center of the star:
01000     double logn_auto_c  = logn_auto()(0, 0, 0, 0) ;
01001     double pot_ext_c  = pot_ext()(0, 0, 0, 0) ;
01002 
01003     // Search for the reference point (theta_*, phi_*)
01004     // [notation of Bonazzola, Gourgoulhon & Marck PRD 58, 104020 (1998)]
01005     // at the surface of the star
01006     // ------------------------------------------------------------------
01007     double alpha_r2 = 0 ;
01008     for (int k=0; k<mg->get_np(l_b); k++) {
01009         for (int j=0; j<mg->get_nt(l_b); j++) {
01010 
01011             double pot_ext_b  = pot_ext()(l_b, k, j, i_b) ;
01012         double logn_auto_b  = logn_auto()(l_b, k, j, i_b) ;
01013 
01014         // See Eq (100) from Gourgoulhon et al. (2001)
01015         double alpha_r2_jk = ( ent_c - ent_b + pot_ext_c - pot_ext_b)
01016           / ( logn_auto_b - logn_auto_c ) ;
01017 
01018         if (alpha_r2_jk > alpha_r2) {
01019             alpha_r2 = alpha_r2_jk ;
01020             k_b = k ;
01021             j_b = j ;
01022         }
01023 
01024         }
01025     }
01026 
01027     alpha_r = sqrt(alpha_r2) ;
01028 
01029     cout << "k_b, j_b, alpha_r: " << k_b << "  " << j_b << "  " 
01030          << alpha_r << endl ;
01031 
01032     // New value of logn_auto
01033     // ----------------------
01034 
01035     logn_auto = alpha_r2 * logn_auto ;
01036     logn_auto_regu = alpha_r2 * logn_auto_regu ;
01037     logn_auto_c  = logn_auto()(0, 0, 0, 0) ;
01038 
01039     //------------------------------------------------------------
01040     // Change the values of the inner points of the second domain
01041     // by those of the outer points of the first domain
01042     //------------------------------------------------------------
01043 
01044     (logn_auto().va).smooth(nzet, (logn_auto.set()).va) ;
01045 
01046     //--------------------------------------------
01047     // First integral --> enthalpy in all space
01048     // See Eq (98) from Gourgoulhon et al. (2001)
01049     //--------------------------------------------
01050 
01051     ent = (ent_c + logn_auto_c + pot_ext_c) - logn_auto - pot_ext ;
01052 
01053     //---------------------------------------------------------
01054     // Change the enthalpy field to accelerate the convergence
01055     //---------------------------------------------------------
01056 
01057     double dentdx = ent().dsdx()(0, 0, 0, 0) ;
01058     double dentdy = ent().dsdy()(0, 0, 0, 0) ;
01059 
01060     cout << "dH/dx|_center = " << dentdx << endl ;
01061     cout << "dH/dy|_center = " << dentdy << endl ;
01062 
01063     double dec_fact = 1. ;
01064 
01065     Tenseur func_in(mp) ;
01066     func_in.set_etat_qcq() ;
01067     func_in.set() = 1. - dec_fact * (dentdx/ent_c) * mp.x ;
01068     func_in.set().annule(nzet, nz-1) ;
01069     func_in.set_std_base() ;
01070 
01071     Tenseur func_ex(mp) ;
01072     func_ex.set_etat_qcq() ;
01073     func_ex.set() = 1. ;
01074     func_ex.set().annule(0, nzet-1) ;
01075     func_ex.set_std_base() ;
01076 
01077     // New enthalpy field
01078     // ------------------
01079     ent.set() = ent() * (func_in() + func_ex()) ;
01080 
01081     (ent().va).smooth(nzet, (ent.set()).va) ;
01082 
01083     double dentdx_new = ent().dsdx()(0, 0, 0, 0) ;
01084 
01085     cout << "dH/dx|_new    = " << dentdx_new << endl ;
01086 
01087     //-----------------------------------------------------
01088     // Adaptation of the mapping to the new enthalpy field
01089     //----------------------------------------------------
01090 
01091     // Shall the adaptation be performed (cusp) ?
01092     // ------------------------------------------
01093 
01094     double dent_eq = ent().dsdr().val_point(ray_eq_pi(),M_PI/2.,M_PI) ;
01095     double dent_pole = ent().dsdr().val_point(ray_pole(),0.,0.) ;
01096     double rap_dent = fabs( dent_eq / dent_pole ) ;
01097     cout << "| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ;
01098 
01099     if ( rap_dent < thres_adapt ) {
01100         adapt_flag = 0 ;    // No adaptation of the mapping
01101         cout << "******* FROZEN MAPPING  *********" << endl ;
01102     }
01103     else{
01104         adapt_flag = 1 ;    // The adaptation of the mapping is to be
01105                             // performed
01106     }
01107 
01108     ent_limit.set_etat_qcq() ;
01109     for (int l=0; l<nzet; l++) {    // loop on domains inside the star
01110         ent_limit.set(l) = ent()(l, k_b, j_b, i_b) ;
01111     }
01112 
01113     ent_limit.set(nzet-1) = ent_b ;
01114     Map_et mp_prev = mp_et ;
01115 
01116     mp.adapt(ent(), par_adapt) ;
01117 
01118     //----------------------------------------------------
01119     // Computation of the enthalpy at the new grid points
01120     //----------------------------------------------------
01121 
01122     mp_prev.homothetie(alpha_r) ;
01123 
01124     mp.reevaluate_symy(&mp_prev, nzet+1, ent.set()) ;
01125 
01126     double fact_resize = 1. / alpha_r ;
01127     for (int l=nzet; l<nz-1; l++) {
01128         mp_et.resize(l, fact_resize) ;
01129     }
01130     mp_et.resize_extr(fact_resize) ;
01131 
01132     double ent_s_max = -1 ;
01133     int k_s_max = -1 ;
01134     int j_s_max = -1 ;
01135     for (int k=0; k<mg->get_np(l_b); k++) {
01136         for (int j=0; j<mg->get_nt(l_b); j++) {
01137             double xx = fabs( ent()(l_b, k, j, i_b) ) ;
01138         if (xx > ent_s_max) {
01139             ent_s_max = xx ;
01140             k_s_max = k ;
01141             j_s_max = j ;
01142         }
01143         }
01144     }
01145     cout << "Max. abs(enthalpy) at the boundary between domains nzet-1"
01146          << " and nzet : " << endl ;
01147     cout << "   " << ent_s_max << " reached for k = " << k_s_max
01148          << " and j = " << j_s_max << endl ;
01149 
01150     //-------------------
01151     // Equation of state
01152     //-------------------
01153 
01154     equation_of_state() ;   // computes new values for nbar (n), ener (e)
01155                             // and press (p) from the new ent (H)
01156 
01157     //----------------------------------------------------------
01158     // Matter source terms in the gravitational field equations
01159     //---------------------------------------------------------
01160 
01161     hydro_euler_extr(mass, sepa) ;  // computes new values for
01162                                     // ener_euler (E), s_euler (S)
01163                                     // and u_euler (U^i)
01164 
01165 
01166     //----------------------------------------------------
01167     // Auxiliary terms for the source of Poisson equation
01168     //----------------------------------------------------
01169 
01170     const Coord& xx = mp.x ;
01171     const Coord& yy = mp.y ;
01172     const Coord& zz = mp.z ;
01173 
01174     Tenseur r_bh(mp) ;
01175     r_bh.set_etat_qcq() ;
01176     r_bh.set() = pow( (xx+sepa)*(xx+sepa) + yy*yy + zz*zz, 0.5) ;
01177     r_bh.set_std_base() ;
01178 
01179     Tenseur xx_cov(mp, 1, COV, ref_triad) ;
01180     xx_cov.set_etat_qcq() ;
01181     xx_cov.set(0) = xx + sepa ;
01182     xx_cov.set(1) = yy ;
01183     xx_cov.set(2) = zz ;
01184     xx_cov.set_std_base() ;
01185 
01186     Tenseur msr(mp) ;
01187     msr = ggrav * mass / r_bh ;
01188     msr.set_std_base() ;
01189 
01190     Tenseur tmp(mp) ;
01191     tmp = 1. / ( 1. - 0.25*msr*msr ) ;
01192     tmp.set_std_base() ;
01193 
01194     Tenseur xdnu(mp) ;
01195     xdnu = flat_scalar_prod_desal(xx_cov, d_logn_auto) ;
01196     xdnu.set_std_base() ;
01197 
01198     Tenseur xdbeta(mp) ;
01199     xdbeta = flat_scalar_prod_desal(xx_cov, d_beta_auto) ;
01200     xdbeta.set_std_base() ;
01201 
01202     //------------------------------------------
01203     // Poisson equation for logn_auto (nu_auto)
01204     //------------------------------------------
01205 
01206     // Source
01207     // ------
01208 
01209     if (relativistic) {
01210 
01211         source = qpig * a_car % (ener_euler + s_euler) + akcar_auto
01212           - flat_scalar_prod_desal(d_logn_auto, d_beta_auto)
01213           - 0.5 * tmp % msr % msr % xdnu / r_bh / r_bh
01214           - tmp % msr % xdbeta / r_bh / r_bh ;
01215 
01216     }
01217     else {
01218         cout << "The computation of BH-NS binary systems"
01219          << " should be relativistic !!!" << endl ;
01220         abort() ;
01221     }
01222 
01223     source.set_std_base() ;
01224 
01225     // Resolution of the Poisson equation
01226     // ----------------------------------
01227 
01228     int k_falloff = 1 ;
01229 
01230     source().poisson_falloff(par_poisson1, logn_auto.set(), k_falloff) ;
01231 
01232     // Construct logn_auto_regu for et_bin_upmetr_extr.C
01233     // -------------------------------------------------
01234 
01235     logn_auto_regu = logn_auto ;
01236 
01237     // Check: has the Poisson equation been correctly solved ?
01238     // -----------------------------------------------------
01239 
01240     Tbl tdiff_logn = diffrel(logn_auto().laplacien(), source()) ;
01241 
01242     cout <<
01243       "Relative error in the resolution of the equation for logn_auto : "
01244          << endl ;
01245 
01246     for (int l=0; l<nz; l++) {
01247         cout << tdiff_logn(l) << "  " ;
01248     }
01249     cout << endl ;
01250     diff_logn = max(abs(tdiff_logn)) ;
01251 
01252     if (relativistic) {
01253 
01254         //--------------------------------
01255         // Poisson equation for beta_auto
01256         //--------------------------------
01257 
01258         // Source
01259         // ------
01260 
01261         source = qpig * a_car % s_euler + 0.75 * akcar_auto
01262           - 0.5 * flat_scalar_prod_desal(d_logn_auto, d_logn_auto)
01263           - 0.5 * flat_scalar_prod_desal(d_beta_auto, d_beta_auto)
01264           - tmp % msr % xdnu / r_bh / r_bh
01265           - 0.5 * tmp % msr %msr % xdbeta / r_bh / r_bh ;
01266 
01267         source.set_std_base() ;
01268 
01269         // Resolution of the Poisson equation
01270         // ----------------------------------
01271 
01272         source().poisson_falloff(par_poisson2, beta_auto.set(),
01273                      k_falloff) ;
01274 
01275         // Check: has the Poisson equation been correctly solved ?
01276         // -----------------------------------------------------
01277 
01278         Tbl tdiff_beta = diffrel(beta_auto().laplacien(), source()) ;
01279 
01280         cout << "Relative error in the resolution of the equation for "
01281          << "beta_auto : " << endl ;
01282         for (int l=0; l<nz; l++) {
01283             cout << tdiff_beta(l) << "  " ;
01284         }
01285         cout << endl ;
01286         diff_beta = max(abs(tdiff_beta)) ;
01287 
01288         //----------------------------------------
01289         // Vector Poisson equation for shift_auto
01290         //----------------------------------------
01291 
01292         // Some auxiliary terms for the source
01293         // -----------------------------------
01294 
01295         Tenseur bhtmp(mp, 1, COV, ref_triad) ;
01296         bhtmp.set_etat_qcq() ;
01297         bhtmp = tmp % msr % (3.*msr-8.) % xx_cov / r_bh / r_bh ;
01298         bhtmp.set_std_base() ;
01299 
01300         Tenseur vtmp =  6. * d_beta_auto - 8. * d_logn_auto ;
01301 
01302         // Source
01303         // ------
01304 
01305         source_shift = (-4.*qpig) * nnn % a_car % (ener_euler + press)
01306           % u_euler
01307           + nnn % flat_scalar_prod_desal(tkij_auto, vtmp+bhtmp) ;
01308 
01309         source_shift.set_std_base() ;
01310 
01311         // Resolution of the Poisson equation
01312         // ----------------------------------
01313 
01314         // Filter for the source of shift vector :
01315 
01316         for (int i=0; i<3; i++) {
01317             for (int l=0; l<nz; l++) {
01318             if (source_shift(i).get_etat() != ETATZERO)
01319             source_shift.set(i).filtre_phi(np_filter, l) ;
01320         }
01321         }
01322 
01323         // For Tenseur::poisson_vect, the triad must be the mapping
01324         // triad, not the reference one:
01325 
01326         source_shift.change_triad( mp.get_bvect_cart() ) ;
01327         /*
01328         for (int i=0; i<3; i++) {
01329             if(source_shift(i).dz_nonzero()) {
01330             assert( source_shift(i).get_dzpuis() == 4 ) ;
01331         }
01332         else {
01333             (source_shift.set(i)).set_dzpuis(4) ;
01334         }
01335         }
01336 
01337         source_shift.dec2_dzpuis() ;    // dzpuis 4 -> 2
01338         */
01339         double lambda_shift = double(1) / double(3) ;
01340 
01341         int* shift_falloff ;
01342         shift_falloff = new int[4] ;
01343         shift_falloff[0] = 1 ;
01344         shift_falloff[1] = 1 ;
01345         shift_falloff[2] = 2 ;
01346         shift_falloff[3] = 1 ;
01347 
01348         source_shift.poisson_vect_falloff(lambda_shift, par_poisson_vect,
01349                           shift_auto, w_shift,
01350                           khi_shift, shift_falloff) ;
01351 
01352         delete[] shift_falloff ;
01353 
01354         // Check: has the equation for shift_auto been correctly solved ?
01355         // --------------------------------------------------------------
01356 
01357         // Divergence of shift_auto :
01358         Tenseur divna = contract(shift_auto.gradient(), 0, 1) ;
01359         //      divna.dec2_dzpuis() ;    // dzpuis 2 -> 0
01360 
01361         // Grad(div) :
01362         Tenseur graddivn = divna.gradient() ;
01363         //      graddivn.inc2_dzpuis() ;    // dzpuis 2 -> 4
01364 
01365         // Full operator :
01366         Tenseur lap_shift(mp, 1, CON, mp.get_bvect_cart() ) ;
01367         lap_shift.set_etat_qcq() ;
01368         for (int i=0; i<3; i++) {
01369             lap_shift.set(i) = shift_auto(i).laplacien()
01370           + lambda_shift * graddivn(i) ;
01371         }
01372 
01373         Tbl tdiff_shift_x = diffrel(lap_shift(0), source_shift(0)) ;
01374         Tbl tdiff_shift_y = diffrel(lap_shift(1), source_shift(1)) ;
01375         Tbl tdiff_shift_z = diffrel(lap_shift(2), source_shift(2)) ;
01376 
01377         cout << "Relative error in the resolution of the equation for "
01378          << "shift_auto : " << endl ;
01379         cout << "x component : " ;
01380         for (int l=0; l<nz; l++) {
01381             cout << tdiff_shift_x(l) << "  " ;
01382         }
01383         cout << endl ;
01384         cout << "y component : " ;
01385         for (int l=0; l<nz; l++) {
01386             cout << tdiff_shift_y(l) << "  " ;
01387         }
01388         cout << endl ;
01389         cout << "z component : " ;
01390         for (int l=0; l<nz; l++) {
01391             cout << tdiff_shift_z(l) << "  " ;
01392         }
01393         cout << endl ;
01394 
01395         diff_shift_x = max(abs(tdiff_shift_x)) ;
01396         diff_shift_y = max(abs(tdiff_shift_y)) ;
01397         diff_shift_z = max(abs(tdiff_shift_z)) ;
01398 
01399         // Final result
01400         // ------------
01401         // The output of Tenseur::poisson_vect_falloff is on the mapping
01402         // triad, it should therefore be transformed to components on the
01403         // reference triad :
01404 
01405         shift_auto.change_triad( ref_triad ) ;
01406 
01407     }   // End of relativistic equations
01408 
01409     //------------------------------
01410     //  Relative change in enthalpy
01411     //------------------------------
01412 
01413     Tbl diff_ent_tbl = diffrel( ent(), ent_jm1() ) ;
01414     diff_ent = diff_ent_tbl(0) ;
01415     for (int l=1; l<nzet; l++) {
01416         diff_ent += diff_ent_tbl(l) ;
01417     }
01418     diff_ent /= nzet ;
01419 
01420     ent_jm1 = ent ;
01421 
01422     } // End of main loop
01423 
01424     //========================================================//
01425     //                    End of iteration                    //
01426     //========================================================//
01427 
01428 }
01429 

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