star_bin_equilibrium.C

00001 
00002 /*
00003  *   Method of class Star_bin to compute an equilibrium configuration
00004  *
00005  *  (see file star.h for documentation).
00006  */
00007 /*
00008  *   Copyright (c) 2004 Francois Limousin
00009  *
00010  *   This file is part of LORENE.
00011  *
00012  *   LORENE is free software; you can redistribute it and/or modify
00013  *   it under the terms of the GNU General Public License version 2
00014  *   as published by the Free Software Foundation.
00015  *
00016  *   LORENE is distributed in the hope that it will be useful,
00017  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00018  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00019  *   GNU General Public License for more details.
00020  *
00021  *   You should have received a copy of the GNU General Public License
00022  *   along with LORENE; if not, write to the Free Software
00023  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00024  *
00025  */
00026 
00027 char star_bin_equilibrium_C[] = "$Header: /cvsroot/Lorene/C++/Source/Star/star_bin_equilibrium.C,v 1.27 2006/08/01 14:26:01 f_limousin Exp $" ;
00028 
00029 /*
00030  * $Id: star_bin_equilibrium.C,v 1.27 2006/08/01 14:26:01 f_limousin Exp $
00031  * $Log: star_bin_equilibrium.C,v $
00032  * Revision 1.27  2006/08/01 14:26:01  f_limousin
00033  * Display
00034  *
00035  * Revision 1.26  2006/05/31 09:26:04  f_limousin
00036  * Modif. of the size of the different domains
00037  *
00038  * Revision 1.25  2006/04/11 14:24:44  f_limousin
00039  * New version of the code : improvement of the computation of some
00040  * critical sources, estimation of the dirac gauge, helical symmetry...
00041  *
00042  * Revision 1.24  2005/11/03 13:27:09  f_limousin
00043  * Final version for the letter.
00044  *
00045  * Revision 1.23  2005/09/14 12:48:02  f_limousin
00046  * Comment graphical outputs.
00047  *
00048  * Revision 1.22  2005/09/14 12:30:52  f_limousin
00049  * Saving of fields lnq and logn in class Star.
00050  *
00051  * Revision 1.21  2005/09/13 19:38:31  f_limousin
00052  * Reintroduction of the resolution of the equations in cartesian coordinates.
00053  *
00054  * Revision 1.20  2005/04/08 12:36:44  f_limousin
00055  * Just to avoid warnings...
00056  *
00057  * Revision 1.19  2005/02/24 16:27:21  f_limousin
00058  * Add mermax_poisson and relax_poisson in the parameters of the function.
00059  *
00060  * Revision 1.18  2005/02/24 16:04:13  f_limousin
00061  * Change the name of some variables (for instance dcov_logn --> dlogn).
00062  * Improve the resolution of the tensorial poisson equation for hh.
00063  *
00064  * Revision 1.17  2005/02/18 13:14:18  j_novak
00065  * Changing of malloc/free to new/delete + suppression of some unused variables
00066  * (trying to avoid compilation warnings).
00067  *
00068  * Revision 1.16  2005/02/17 17:32:53  f_limousin
00069  * Change the name of some quantities to be consistent with other classes
00070  * (for instance nnn is changed to nn, shift to beta, beta to lnq...)
00071  *
00072  * Revision 1.15  2005/02/11 18:13:47  f_limousin
00073  * Important modification : all the poisson equations for the metric
00074  * quantities are now solved on an affine mapping.
00075  *
00076  * Revision 1.14  2004/12/17 16:23:19  f_limousin
00077  * Modif. comments.
00078  *
00079  * Revision 1.13  2004/06/22 12:49:12  f_limousin
00080  * Change qq, qq_auto and qq_comp to beta, beta_auto and beta_comp.
00081  *
00082  * Revision 1.12  2004/05/27 12:41:00  p_grandclement
00083  * correction of some shadowed variables
00084  *
00085  * Revision 1.11  2004/05/25 14:18:00  f_limousin
00086  * Include filters
00087  *
00088  * Revision 1.10  2004/05/10 10:26:22  f_limousin
00089  * Minor changes to avoid warnings in the compilation of Lorene
00090  *
00091  * Revision 1.9  2004/04/08 16:32:48  f_limousin
00092  * The new variable is ln(Q) instead of Q=psi^2*N. It improves the
00093  * convergence of the code.
00094  *
00095  * Revision 1.8  2004/03/25 10:29:26  j_novak
00096  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
00097  *
00098  * Revision 1.7  2004/03/23 09:56:09  f_limousin
00099  * Many minor changes
00100  *
00101  * Revision 1.6  2004/02/27 21:16:32  e_gourgoulhon
00102  * Function contract_desal replaced by function contract with
00103  * argument desaliasing set to true.
00104  *
00105  * Revision 1.5  2004/02/27 09:51:51  f_limousin
00106  * Many minor changes.
00107  *
00108  * Revision 1.4  2004/02/21 17:05:13  e_gourgoulhon
00109  * Method Scalar::point renamed Scalar::val_grid_point.
00110  * Method Scalar::set_point renamed Scalar::set_grid_point.
00111  *
00112  * Revision 1.3  2004/01/20 15:17:48  f_limousin
00113  * First version
00114  *
00115  *
00116  * $Header: /cvsroot/Lorene/C++/Source/Star/star_bin_equilibrium.C,v 1.27 2006/08/01 14:26:01 f_limousin Exp $
00117  *
00118  */
00119 
00120 // C headers
00121 #include <math.h>
00122 
00123 // Lorene headers
00124 #include "cmp.h"
00125 #include "tenseur.h"
00126 #include "metrique.h"
00127 #include "star.h"
00128 #include "param.h"
00129 #include "graphique.h"
00130 #include "utilitaires.h"
00131 #include "tensor.h"
00132 #include "nbr_spx.h"
00133 #include "unites.h"
00134 
00135 
00136 void Star_bin::equilibrium(double ent_c, int mermax, int mermax_potvit, 
00137                int mermax_poisson, double relax_poisson, 
00138                double relax_potvit, double thres_adapt,
00139                Tbl& diff, double om) {
00140 
00141     // Fundamental constants and units
00142     // -------------------------------
00143     using namespace Unites ;
00144     
00145     // Initializations
00146     // ---------------
00147     
00148     const Mg3d* mg = mp.get_mg() ; 
00149     int nz = mg->get_nzone() ;      // total number of domains
00150     
00151     // The following is required to initialize mp_prev as a Map_et:
00152     Map_et& mp_et = dynamic_cast<Map_et&>(mp) ; 
00153     
00154     // Domain and radial indices of points at the surface of the star:
00155     int l_b = nzet - 1 ; 
00156     int i_b = mg->get_nr(l_b) - 1 ; 
00157     int k_b ;
00158     int j_b ; 
00159     
00160     // Value of the enthalpy defining the surface of the star
00161     double ent_b = 0 ; 
00162     
00163     // Error indicators
00164     // ----------------
00165     
00166     double& diff_ent = diff.set(0) ; 
00167     double& diff_vel_pot = diff.set(1) ; 
00168     double& diff_logn = diff.set(2) ; 
00169     double& diff_lnq = diff.set(3) ; 
00170     double& diff_beta_x = diff.set(4) ; 
00171     double& diff_beta_y = diff.set(5) ; 
00172     double& diff_beta_z = diff.set(6) ; 
00173     double& diff_h11 = diff.set(7) ; 
00174     double& diff_h21 = diff.set(8) ; 
00175     double& diff_h31 = diff.set(9) ; 
00176     double& diff_h22 = diff.set(10) ; 
00177     double& diff_h32 = diff.set(11) ; 
00178     double& diff_h33 = diff.set(12) ; 
00179 
00180 
00181 
00182     // Parameters for te function Map_et::adapt
00183     // -----------------------------------------
00184 
00185     Param par_adapt ;
00186     int nitermax = 100 ;
00187     int niter ; 
00188     int adapt_flag = 1 ;    //  1 = performs the full computation, 
00189     //  0 = performs only the rescaling by 
00190     //      the factor alpha_r
00191     //##    int nz_search = nzet + 1 ;  // Number of domains for searching the 
00192     // enthalpy
00193     int nz_search = nzet ;  // Number of domains for searching the enthalpy
00194                 //  isosurfaces
00195 
00196     double precis_secant = 1.e-14 ; 
00197     double alpha_r ; 
00198     double reg_map = 1. ; // 1 = regular mapping, 0 = contracting mapping
00199 
00200     Tbl ent_limit(nz) ; 
00201 
00202 
00203     par_adapt.add_int(nitermax, 0) ; // maximum number of iterations to 
00204     // locate zeros by the secant method
00205     par_adapt.add_int(nzet, 1) ;    // number of domains where the adjustment 
00206     // to the isosurfaces of ent is to be 
00207     // performed
00208     par_adapt.add_int(nz_search, 2) ;   // number of domains to search  
00209                                         // the enthalpy isosurface
00210     par_adapt.add_int(adapt_flag, 3) ; //  1 = performs the full computation, 
00211     //  0 = performs only the rescaling by 
00212     //      the factor alpha_r
00213     par_adapt.add_int(j_b, 4) ; //  theta index of the collocation point 
00214     //  (theta_*, phi_*)
00215     par_adapt.add_int(k_b, 5) ; //  theta index of the collocation point 
00216     //  (theta_*, phi_*)
00217 
00218     par_adapt.add_int_mod(niter, 0) ;  // number of iterations actually used in
00219     //  the secant method
00220     
00221     par_adapt.add_double(precis_secant, 0) ; // required absolute precision in 
00222     // the determination of zeros by 
00223     // the secant method
00224     par_adapt.add_double(reg_map, 1)    ;  // 1. = regular mapping, 
00225                                            // 0 = contracting mapping
00226     
00227     par_adapt.add_double(alpha_r, 2) ;      // factor by which all the radial 
00228                         // distances will be multiplied 
00229            
00230     par_adapt.add_tbl(ent_limit, 0) ;   // array of values of the field ent 
00231                         // to define the isosurfaces. 
00232  
00233 
00234     Cmp ssjm1logn (ssjm1_logn) ;
00235     Cmp ssjm1lnq (ssjm1_lnq) ;
00236     Cmp ssjm1h11 (ssjm1_h11) ;
00237     Cmp ssjm1h21 (ssjm1_h21) ;
00238     Cmp ssjm1h31 (ssjm1_h31) ;
00239     Cmp ssjm1h22 (ssjm1_h22) ;
00240     Cmp ssjm1h32 (ssjm1_h32) ;
00241     Cmp ssjm1h33 (ssjm1_h33) ;
00242  
00243 
00244     ssjm1logn.set_etat_qcq() ;
00245     ssjm1lnq.set_etat_qcq() ;
00246     ssjm1h11.set_etat_qcq() ;
00247     ssjm1h21.set_etat_qcq() ;
00248     ssjm1h31.set_etat_qcq() ;
00249     ssjm1h22.set_etat_qcq() ;
00250     ssjm1h32.set_etat_qcq() ;
00251     ssjm1h33.set_etat_qcq() ;
00252    
00253     
00254     double precis_poisson = 1.e-16 ;     
00255 
00256     // Parameters for the function Scalar::poisson for logn_auto
00257     // ---------------------------------------------------------------
00258  
00259     Param par_logn ;    
00260 
00261     par_logn.add_int(mermax_poisson,  0) ;  // maximum number of iterations
00262     par_logn.add_double(relax_poisson,  0) ; // relaxation parameter
00263     par_logn.add_double(precis_poisson, 1) ; // required precision
00264     par_logn.add_int_mod(niter, 0) ; // number of iterations actually used 
00265     par_logn.add_cmp_mod( ssjm1logn ) ; 
00266     
00267     // Parameters for the function Scalar::poisson for lnq_auto
00268     // ---------------------------------------------------------------
00269     
00270     Param par_lnq ; 
00271     
00272     par_lnq.add_int(mermax_poisson,  0) ;  // maximum number of iterations
00273     par_lnq.add_double(relax_poisson,  0) ; // relaxation parameter
00274     par_lnq.add_double(precis_poisson, 1) ; // required precision
00275     par_lnq.add_int_mod(niter, 0) ; // number of iterations actually used -
00276     par_lnq.add_cmp_mod( ssjm1lnq ) ; 
00277  
00278     // Parameters for the function Vector::poisson for beta method 2 
00279     // ---------------------------------------------------------------
00280     
00281     Param par_beta2 ; 
00282     
00283     par_beta2.add_int(mermax_poisson,  0) ;  // maximum number of iterations
00284     par_beta2.add_double(relax_poisson,  0) ; // relaxation parameter
00285     par_beta2.add_double(precis_poisson, 1) ; // required precision
00286     par_beta2.add_int_mod(niter, 0) ; // number of iterations actually used 
00287   
00288     Cmp ssjm1khi (ssjm1_khi) ;
00289     Tenseur ssjm1wbeta(mp, 1, CON, mp.get_bvect_cart()) ;
00290     ssjm1wbeta.set_etat_qcq() ;
00291     for (int i=0; i<3; i++) {
00292       ssjm1wbeta.set(i) = Cmp(ssjm1_wbeta(i+1)) ;
00293     }
00294     
00295     par_beta2.add_cmp_mod(ssjm1khi) ; 
00296     par_beta2.add_tenseur_mod(ssjm1wbeta) ; 
00297     
00298     // Parameters for the function Scalar::poisson for h11_auto
00299     // -------------------------------------------------------------
00300     
00301     Param par_h11 ; 
00302     
00303     par_h11.add_int(mermax_poisson,  0) ;  // maximum number of iterations
00304     par_h11.add_double(relax_poisson,  0) ; // relaxation parameter
00305     par_h11.add_double(precis_poisson, 1) ; // required precision
00306     par_h11.add_int_mod(niter, 0) ; // number of iterations actually used 
00307     par_h11.add_cmp_mod( ssjm1h11 ) ; 
00308     
00309     // Parameters for the function Scalar::poisson for h21_auto
00310     // -------------------------------------------------------------
00311     
00312     Param par_h21 ; 
00313     
00314     par_h21.add_int(mermax_poisson,  0) ;  // maximum number of iterations
00315     par_h21.add_double(relax_poisson,  0) ; // relaxation parameter
00316     par_h21.add_double(precis_poisson, 1) ; // required precision
00317     par_h21.add_int_mod(niter, 0) ; // number of iterations actually used 
00318     par_h21.add_cmp_mod( ssjm1h21 ) ; 
00319     
00320     // Parameters for the function Scalar::poisson for h31_auto
00321     // -------------------------------------------------------------
00322     
00323     Param par_h31 ; 
00324     
00325     par_h31.add_int(mermax_poisson,  0) ;  // maximum number of iterations
00326     par_h31.add_double(relax_poisson,  0) ; // relaxation parameter
00327     par_h31.add_double(precis_poisson, 1) ; // required precision
00328     par_h31.add_int_mod(niter, 0) ; // number of iterations actually used 
00329     par_h31.add_cmp_mod( ssjm1h31 ) ; 
00330     
00331     // Parameters for the function Scalar::poisson for h22_auto
00332     // -------------------------------------------------------------
00333     
00334     Param par_h22 ; 
00335     
00336     par_h22.add_int(mermax_poisson,  0) ;  // maximum number of iterations
00337     par_h22.add_double(relax_poisson,  0) ; // relaxation parameter
00338     par_h22.add_double(precis_poisson, 1) ; // required precision
00339     par_h22.add_int_mod(niter, 0) ; // number of iterations actually used 
00340     par_h22.add_cmp_mod( ssjm1h22 ) ; 
00341     
00342     // Parameters for the function Scalar::poisson for h32_auto
00343     // -------------------------------------------------------------
00344     
00345     Param par_h32 ; 
00346     
00347     par_h32.add_int(mermax_poisson,  0) ;  // maximum number of iterations
00348     par_h32.add_double(relax_poisson,  0) ; // relaxation parameter
00349     par_h32.add_double(precis_poisson, 1) ; // required precision
00350     par_h32.add_int_mod(niter, 0) ; // number of iterations actually used 
00351     par_h32.add_cmp_mod( ssjm1h32 ) ; 
00352     
00353     // Parameters for the function Scalar::poisson for h33_auto
00354     // -------------------------------------------------------------
00355     
00356     Param par_h33 ; 
00357     
00358     par_h33.add_int(mermax_poisson,  0) ;  // maximum number of iterations
00359     par_h33.add_double(relax_poisson,  0) ; // relaxation parameter
00360     par_h33.add_double(precis_poisson, 1) ; // required precision
00361     par_h33.add_int_mod(niter, 0) ; // number of iterations actually used 
00362     par_h33.add_cmp_mod( ssjm1h33 ) ; 
00363     
00364 
00365     // External potential
00366     // See Eq (99) from Gourgoulhon et al. (2001)
00367     // ------------------
00368     
00369     cout << "logn_comp" << norme(logn_comp) << endl ;
00370     cout << "pot_centri" << norme(pot_centri) << endl ;
00371     cout << "loggam" << norme(loggam) << endl ;
00372     Scalar pot_ext = logn_comp + pot_centri + loggam ;
00373     
00374     Scalar ent_jm1 = ent ;  // Enthalpy at previous step
00375     
00376     Scalar source_tot(mp) ; // source term in the equation for hij_auto, 
00377                             // logn_auto and beta_auto
00378                 
00379     Vector source_beta(mp, CON, mp.get_bvect_cart()) ;  // source term 
00380     // in the equation for beta_auto
00381 
00382 
00383 
00384     //=========================================================================
00385     //          Start of iteration
00386     //=========================================================================
00387 
00388     for(int mer=0 ; mer<mermax ; mer++ ) {
00389 
00390     cout << "-----------------------------------------------" << endl ;
00391     cout << "step: " << mer << endl ;
00392     cout << "diff_ent = " << diff_ent << endl ;    
00393 
00394     //-----------------------------------------------------
00395     // Resolution of the elliptic equation for the velocity
00396     // scalar potential
00397     //-----------------------------------------------------
00398 
00399     if (irrotational) {
00400         diff_vel_pot = velocity_potential(mermax_potvit, precis_poisson, 
00401                           relax_potvit) ; 
00402         
00403     }
00404 
00405     diff_vel_pot = 0. ; // to avoid the warning 
00406 
00407     //-----------------------------------------------------
00408     // Computation of the new radial scale
00409     //--------------------------------------------------
00410 
00411     // alpha_r (r = alpha_r r') is determined so that the enthalpy
00412     // takes the requested value ent_b at the stellar surface
00413     
00414     // Values at the center of the star:
00415     double logn_auto_c  = logn_auto.val_grid_point(0, 0, 0, 0) ; 
00416     double pot_ext_c  = pot_ext.val_grid_point(0, 0, 0, 0) ; 
00417 
00418     // Search for the reference point (theta_*, phi_*) [notation of
00419     //  Bonazzola, Gourgoulhon & Marck PRD 58, 104020 (1998)]
00420     //  at the surface of the star
00421     // ------------------------------------------------------------
00422     double alpha_r2 = 0 ; 
00423     for (int k=0; k<mg->get_np(l_b); k++) {
00424         for (int j=0; j<mg->get_nt(l_b); j++) {
00425         
00426         double pot_ext_b  = pot_ext.val_grid_point(l_b, k, j, i_b) ; 
00427         double logn_auto_b  = logn_auto.val_grid_point(l_b, k, j, i_b) ;
00428         // See Eq (100) from Gourgoulhon et al. (2001)
00429         double alpha_r2_jk = ( ent_c - ent_b + pot_ext_c - pot_ext_b) /
00430  
00431             ( logn_auto_b - logn_auto_c ) ;
00432           
00433         if (alpha_r2_jk > alpha_r2) {
00434             alpha_r2 = alpha_r2_jk ; 
00435             k_b = k ; 
00436             j_b = j ; 
00437         }
00438 
00439         }
00440     }
00441         
00442     alpha_r = sqrt(alpha_r2) ;
00443         
00444     cout << "k_b, j_b, alpha_r: " << k_b << "  " << j_b << "  " 
00445          <<  alpha_r << endl ;
00446 
00447     // New value of logn_auto 
00448     // ----------------------
00449 
00450     logn_auto = alpha_r2 * logn_auto ;
00451     logn_auto_c  = logn_auto.val_grid_point(0, 0, 0, 0) ;
00452 
00453     //------------------------------------------------------------
00454     // Change the values of the inner points of the second domain
00455     // by those of the outer points of the first domain
00456     //------------------------------------------------------------
00457 
00458     logn_auto.set_spectral_va().smooth(nzet, logn_auto.set_spectral_va()) ;
00459 
00460     //------------------------------------------
00461     // First integral   -->  enthalpy in all space
00462     // See Eq (98) from Gourgoulhon et al. (2001)
00463     //-------------------------------------------
00464 
00465     ent = (ent_c + logn_auto_c + pot_ext_c) - logn_auto - pot_ext ;
00466     cout.precision(8) ;
00467     cout << "pot" << norme(pot_ext) << endl ;
00468 
00469     (ent.set_spectral_va()).smooth(nzet, ent.set_spectral_va()) ;
00470 
00471     //----------------------------------------------------
00472     // Adaptation of the mapping to the new enthalpy field
00473     //----------------------------------------------------
00474     
00475     // Shall the adaptation be performed (cusp) ?
00476     // ------------------------------------------
00477     
00478     double dent_eq = ent.dsdr().val_point(ray_eq(),M_PI/2.,0.) ;
00479     double dent_pole = ent.dsdr().val_point(ray_pole(),0.,0.) ;
00480     double rap_dent = fabs( dent_eq / dent_pole ) ; 
00481     cout << "| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ; 
00482     
00483     if ( rap_dent < thres_adapt ) {
00484         adapt_flag = 0 ;    // No adaptation of the mapping 
00485         cout << "******* FROZEN MAPPING  *********" << endl ; 
00486     }
00487     else{
00488         adapt_flag = 1 ;    // The adaptation of the mapping is to be
00489         //  performed
00490     }
00491 
00492     ent_limit.set_etat_qcq() ; 
00493     for (int l=0; l<nzet; l++) {    // loop on domains inside the star
00494         ent_limit.set(l) = ent.val_grid_point(l, k_b, j_b, i_b) ; 
00495     }
00496     ent_limit.set(nzet-1) = ent_b  ; 
00497 
00498     Map_et mp_prev = mp_et ; 
00499 
00500     Cmp ent_cmp(ent) ;
00501     mp.adapt(ent_cmp, par_adapt) ; 
00502     ent = ent_cmp ;
00503 
00504     // Readjustment of the external boundary of domain l=nzet
00505     // to keep a fixed ratio with respect to star's surface
00506     
00507     if (nz>= 5) {
00508       double separation = 2. * fabs(mp.get_ori_x()) ;
00509       double ray_eqq = ray_eq() ;
00510       double ray_eqq_pi = ray_eq_pi() ;
00511       double new_rr_out_2 = (separation - ray_eqq) * 0.95 ; 
00512       double new_rr_out_3 = (separation + ray_eqq_pi) * 1.05 ; 
00513       
00514       double rr_in_1 = mp.val_r(1,-1., M_PI/2, 0.) ; 
00515       double rr_out_1 = mp.val_r(1, 1., M_PI/2, 0.) ; 
00516       double rr_out_2 = mp.val_r(2, 1., M_PI/2, 0.) ; 
00517       double rr_out_3 = mp.val_r(3, 1., M_PI/2, 0.) ; 
00518       
00519       mp.resize(1, 0.5*(new_rr_out_2 + rr_in_1) / rr_out_1) ; 
00520       mp.resize(2, new_rr_out_2 / rr_out_2) ; 
00521       mp.resize(3, new_rr_out_3 / rr_out_3) ;
00522 
00523       for (int dd=4; dd<=nz-2; dd++) {
00524         mp.resize(dd, new_rr_out_3 * pow(4., dd-3) / 
00525               mp.val_r(dd, 1., M_PI/2, 0.)) ;
00526       }
00527 
00528     }
00529     else {
00530       cout << "too small number of domains" << endl ;
00531     }
00532     
00533     //----------------------------------------------------
00534     // Computation of the enthalpy at the new grid points
00535     //----------------------------------------------------
00536     
00537     mp_prev.homothetie(alpha_r) ; 
00538     
00539     Cmp ent_cmp2 (ent) ;
00540     mp.reevaluate_symy(&mp_prev, nzet+1, ent_cmp2) ; 
00541     ent = ent_cmp2 ;
00542 
00543     double ent_s_max = -1 ; 
00544     int k_s_max = -1 ; 
00545     int j_s_max = -1 ; 
00546     for (int k=0; k<mg->get_np(l_b); k++) {
00547         for (int j=0; j<mg->get_nt(l_b); j++) {
00548         double xx = fabs( ent.val_grid_point(l_b, k, j, i_b) ) ;
00549         if (xx > ent_s_max) {
00550             ent_s_max = xx ; 
00551             k_s_max = k ; 
00552             j_s_max = j ; 
00553         }
00554         }
00555     }
00556     cout << "Max. abs(enthalpy) at the boundary between domains nzet-1"
00557          << " and nzet : " << endl ; 
00558     cout << "   " << ent_s_max << " reached for k = " << k_s_max <<
00559         " and j = " << j_s_max << endl ; 
00560 
00561     //----------------------------------------------------
00562     // Equation of state  
00563     //----------------------------------------------------
00564     
00565     equation_of_state() ;   // computes new values for nbar (n), ener (e) 
00566                 // and press (p) from the new ent (H)
00567     
00568     //---------------------------------------------------------
00569     // Matter source terms in the gravitational field equations 
00570     //---------------------------------------------------------
00571 
00572     hydro_euler() ;     // computes new values for ener_euler (E), 
00573                 // s_euler (S) and u_euler (U^i)
00574 
00575     
00576     // -------------------------------
00577     // AUXILIARY QUANTITIES
00578     // -------------------------------
00579 
00580     // Derivatives of N and logN
00581     //--------------------------
00582 
00583     const Vector dcov_logn_auto = logn_auto.derive_cov(flat) ;
00584     
00585     Tensor dcovdcov_logn_auto = (logn_auto.derive_cov(flat))
00586                                       .derive_cov(flat) ;
00587     dcovdcov_logn_auto.inc_dzpuis() ;
00588 
00589     // Derivatives of lnq, phi and Q
00590     //-------------------------------
00591 
00592     const Scalar phi (0.5 * (lnq - logn)) ;
00593     const Scalar phi_auto (0.5 * (lnq_auto - logn_auto)) ;
00594 
00595     const Vector dcov_phi_auto = phi_auto.derive_cov(flat) ;
00596     
00597     const Vector dcov_lnq = 2*dcov_phi + dcov_logn ;
00598     const Vector dcon_lnq = 2*dcon_phi + dcon_logn ;
00599     const Vector dcov_lnq_auto = lnq_auto.derive_cov(flat) ;
00600         Tensor dcovdcov_lnq_auto = dcov_lnq_auto.derive_cov(flat) ;
00601     dcovdcov_lnq_auto.inc_dzpuis() ;
00602 
00603     Scalar qq = exp(lnq) ;
00604     qq.std_spectral_base() ;
00605     const Vector& dcov_qq = qq.derive_cov(flat) ;
00606 
00607     const Scalar& divbeta_auto = beta_auto.divergence(flat) ;
00608     const Tensor& dcov_beta_auto = beta_auto.derive_cov(flat) ;
00609     Tensor dcovdcov_beta_auto = beta_auto.derive_cov(flat)
00610       .derive_cov(flat) ;
00611     dcovdcov_beta_auto.inc_dzpuis() ;
00612 
00613 
00614     // Derivatives of hij, gtilde... 
00615     //------------------------------
00616 
00617     Scalar psi2 (pow(psi4, 0.5)) ;
00618     psi2.std_spectral_base() ;
00619 
00620     const Tensor& dcov_hij = hij.derive_cov(flat) ;
00621     const Tensor& dcon_hij = hij.derive_con(flat) ;
00622     const Tensor& dcov_hij_auto = hij_auto.derive_cov(flat) ;
00623 
00624     const Sym_tensor& gtilde_cov = gtilde.cov() ;
00625     const Sym_tensor& gtilde_con = gtilde.con() ;
00626     const Tensor& dcov_gtilde = gtilde_cov.derive_cov(flat) ;
00627 
00628     Connection gamijk (gtilde, flat) ;
00629     const Tensor& deltaijk = gamijk.get_delta() ;
00630 
00631     // H^i and its derivatives ( = O in Dirac gauge)
00632     // ---------------------------------------------
00633 
00634     double lambda_dirac = 0. ;
00635 
00636     const Vector hdirac = lambda_dirac * hij.divergence(flat) ;
00637     const Vector hdirac_auto = lambda_dirac * hij_auto.divergence(flat) ;
00638 
00639     Tensor dcov_hdirac = lambda_dirac * hdirac.derive_cov(flat) ;
00640     dcov_hdirac.inc_dzpuis() ;
00641     Tensor dcov_hdirac_auto = lambda_dirac * hdirac_auto.derive_cov(flat) ;
00642     dcov_hdirac_auto.inc_dzpuis() ;
00643     Tensor dcon_hdirac_auto = lambda_dirac * hdirac_auto.derive_con(flat) ;
00644     dcon_hdirac_auto.inc_dzpuis() ;
00645 
00646 
00647     //--------------------------------------------------------
00648     // Poisson equation for logn_auto (nu_auto)
00649     //--------------------------------------------------------
00650   
00651     // Source 
00652     //--------
00653 
00654     int nr = mp.get_mg()->get_nr(0) ;
00655     int nt = mp.get_mg()->get_nt(0) ;
00656     int np = mp.get_mg()->get_np(0) ;
00657     
00658     Scalar source1(mp) ;
00659     Scalar source2(mp) ;
00660     Scalar source3(mp) ;
00661     Scalar source4(mp) ;
00662     Scalar source5(mp) ;
00663     Scalar source6(mp) ;
00664     Scalar source7(mp) ;
00665     Scalar source8(mp) ;
00666      
00667     source1 = qpig * psi4 % (ener_euler + s_euler) ; 
00668 
00669     source2 = psi4 % (kcar_auto + kcar_comp) ;
00670 
00671     source3 = - contract(dcov_logn_auto, 0, dcon_logn, 0, true) 
00672         - 2. * contract(contract(gtilde_con, 0, dcov_phi, 0), 
00673                 0, dcov_logn_auto, 0, true) ;
00674             
00675     source4 = - contract(hij, 0, 1, dcovdcov_logn_auto + 
00676                  dcov_logn_auto*dcov_logn, 0, 1) ;
00677 
00678     source5 = - contract(hdirac, 0, dcov_logn_auto, 0) ;
00679 
00680     source_tot = source1 + source2 + source3 + source4 + source5 ;
00681       
00682  
00683     cout << "moyenne de la source 1 pour logn_auto" << endl ;
00684     cout <<  norme(source1/(nr*nt*np)) << endl ;
00685     cout << "moyenne de la source 2 pour logn_auto" << endl ;
00686     cout <<  norme(source2/(nr*nt*np)) << endl ;
00687     cout << "moyenne de la source 3 pour logn_auto" << endl ;
00688     cout <<  norme(source3/(nr*nt*np)) << endl ;
00689     cout << "moyenne de la source 4 pour logn_auto" << endl ;
00690     cout <<  norme(source4/(nr*nt*np)) << endl ;
00691     cout << "moyenne de la source 5 pour logn_auto" << endl ;
00692     cout <<  norme(source5/(nr*nt*np)) << endl ;
00693     cout << "moyenne de la source pour logn_auto" << endl ;
00694     cout <<  norme(source_tot/(nr*nt*np)) << endl ;
00695 
00696     // Resolution of the Poisson equation 
00697     // ----------------------------------
00698 
00699     source_tot.poisson(par_logn, logn_auto) ; 
00700     ssjm1_logn = ssjm1logn ;
00701 
00702     cout << "logn_auto" << endl << norme(logn_auto/(nr*nt*np)) << endl ;
00703   
00704     // Check: has the Poisson equation been correctly solved ?
00705     // -----------------------------------------------------
00706 
00707     Tbl tdiff_logn = diffrel(logn_auto.laplacian(), source_tot) ;
00708     cout << 
00709         "Relative error in the resolution of the equation for logn_auto : "
00710          << endl ; 
00711     for (int l=0; l<nz; l++) {
00712         cout << tdiff_logn(l) << "  " ; 
00713     }
00714     cout << endl ;
00715     diff_logn = max(abs(tdiff_logn)) ; 
00716 
00717     //--------------------------------------------------------
00718     // Poisson equation for lnq_auto
00719     //--------------------------------------------------------
00720 
00721     // Source
00722     //--------
00723 
00724     source1 = qpig * psi4 % s_euler ;
00725 
00726     source2 = 0.75 * psi4 % (kcar_auto + kcar_comp) ;
00727 
00728     source3 = - contract(dcon_lnq, 0, dcov_lnq_auto, 0, true) ;
00729 
00730     source4 = 2. * contract(contract(gtilde_con, 0, dcov_phi, 0), 0, 
00731                 dcov_phi_auto + dcov_logn_auto, 0, true) ;
00732 
00733     source5 = 0.0625 * contract(gtilde_con, 0, 1, contract(
00734          dcov_hij_auto, 0, 1, dcov_gtilde, 0, 1), 0, 1) ;
00735         
00736     source6 = - 0.125 * contract(gtilde_con, 0, 1, contract(
00737                dcov_hij_auto, 0, 1, dcov_gtilde, 0, 2), 0, 1) ;
00738 
00739     source7 = - contract(hij, 0, 1, dcovdcov_lnq_auto + dcov_lnq_auto *
00740                  dcov_lnq, 0, 1) ;
00741 
00742     source8 = - 0.25 * contract(dcov_hdirac_auto, 0, 1) 
00743       - contract(hdirac, 0, dcov_lnq_auto, 0) ;
00744     
00745     source_tot = source1 + source2 + source3 + source4 + source5 + 
00746                  source6 + source7 + source8 ;
00747 
00748     
00749     cout << "moyenne de la source 1 pour lnq_auto" << endl ;
00750     cout <<  norme(source1/(nr*nt*np)) << endl ;
00751     cout << "moyenne de la source 2 pour lnq_auto" << endl ;
00752     cout <<  norme(source2/(nr*nt*np)) << endl ;
00753     cout << "moyenne de la source 3 pour lnq_auto" << endl ;
00754     cout <<  norme(source3/(nr*nt*np)) << endl ;
00755     cout << "moyenne de la source 4 pour lnq_auto" << endl ;
00756     cout <<  norme(source4/(nr*nt*np)) << endl ;
00757     cout << "moyenne de la source 5 pour lnq_auto" << endl ;
00758     cout <<  norme(source5/(nr*nt*np)) << endl ;
00759     cout << "moyenne de la source 6 pour lnq_auto" << endl ;
00760     cout <<  norme(source6/(nr*nt*np)) << endl ;
00761     cout << "moyenne de la source 7 pour lnq_auto" << endl ;
00762     cout <<  norme(source7/(nr*nt*np)) << endl ;
00763     cout << "moyenne de la source 8 pour lnq_auto" << endl ;
00764     cout <<  norme(source8/(nr*nt*np)) << endl ;
00765     cout << "moyenne de la source pour lnq_auto" << endl ;
00766     cout <<  norme(source_tot/(nr*nt*np)) << endl ;
00767     
00768 
00769     // Resolution of the Poisson equation 
00770     // ----------------------------------
00771 
00772     source_tot.poisson(par_lnq, lnq_auto) ; 
00773     ssjm1_lnq = ssjm1lnq ;
00774 
00775     cout << "lnq_auto" << endl << norme(lnq_auto/(nr*nt*np)) << endl ;
00776 
00777     // Check: has the Poisson equation been correctly solved 
00778     // -----------------------------------------------------
00779     
00780     Tbl tdiff_lnq = diffrel(lnq_auto.laplacian(), source_tot) ;
00781     cout << 
00782         "Relative error in the resolution of the equation for lnq : "
00783          << endl ; 
00784     for (int l=0; l<nz; l++) {
00785         cout << tdiff_lnq(l) << "  " ; 
00786     }
00787     cout << endl ;
00788     diff_lnq = max(abs(tdiff_lnq)) ; 
00789 
00790     //--------------------------------------------------------
00791     // Vector Poisson equation for beta_auto
00792     //--------------------------------------------------------
00793 
00794     // Source
00795     //--------
00796 
00797     Vector source1_beta(mp, CON, mp.get_bvect_cart()) ;
00798     Vector source2_beta(mp, CON, mp.get_bvect_cart()) ;
00799     Vector source3_beta(mp, CON, mp.get_bvect_cart()) ;
00800     Vector source4_beta(mp, CON, mp.get_bvect_cart()) ;
00801     Vector source5_beta(mp, CON, mp.get_bvect_cart()) ;
00802     Vector source6_beta(mp, CON, mp.get_bvect_cart()) ;
00803     Vector source7_beta(mp, CON, mp.get_bvect_cart()) ;
00804 
00805     source1_beta = (4.*qpig) * nn % psi4
00806                        %(ener_euler + press) * u_euler ;
00807 
00808     source2_beta = 2. * nn * contract(tkij_auto, 1, 
00809                         dcov_logn - 6 * dcov_phi, 0)  ;
00810 
00811     source3_beta = - 2. * nn * contract(tkij_auto, 0, 1, deltaijk, 
00812                           1, 2) ;
00813 
00814     source4_beta = - contract(hij, 0, 1, dcovdcov_beta_auto, 1, 2) ;
00815     
00816     source5_beta = - 0.3333333333333333 * contract(hij, 1, contract( 
00817                           dcovdcov_beta_auto, 0, 1), 0) ;
00818     
00819     source6_beta.set_etat_zero() ; //hdirac_auto.derive_lie(omdsdp) ;
00820     source6_beta.inc_dzpuis() ;
00821 
00822     source7_beta = contract(beta, 0, dcov_hdirac_auto, 1) ;
00823            + 2./3. * hdirac * divbeta_auto 
00824            - contract(hdirac, 0, dcov_beta_auto, 1) ;
00825 
00826     source_beta = source1_beta + source2_beta + source3_beta 
00827       + source4_beta + source5_beta + source6_beta + source7_beta ; 
00828     
00829 
00830     cout << "moyenne de la source 1 pour beta_auto" << endl ;
00831     cout <<  norme(source1_beta(1)/(nr*nt*np)) << endl ;
00832     cout <<  norme(source1_beta(2)/(nr*nt*np)) << endl ;
00833     cout <<  norme(source1_beta(3)/(nr*nt*np)) << endl ;
00834     cout << "moyenne de la source 2 pour beta_auto" << endl ;
00835     cout <<  norme(source2_beta(1)/(nr*nt*np)) << endl ;
00836     cout <<  norme(source2_beta(2)/(nr*nt*np)) << endl ;
00837     cout <<  norme(source2_beta(3)/(nr*nt*np)) << endl ;
00838     cout << "moyenne de la source 3 pour beta_auto" << endl ;
00839     cout <<  norme(source3_beta(1)/(nr*nt*np)) << endl ;
00840     cout <<  norme(source3_beta(2)/(nr*nt*np)) << endl ;
00841     cout <<  norme(source3_beta(3)/(nr*nt*np)) << endl ;
00842     cout << "moyenne de la source 4 pour beta_auto" << endl ;
00843     cout <<  norme(source4_beta(1)/(nr*nt*np)) << endl ;
00844     cout <<  norme(source4_beta(2)/(nr*nt*np)) << endl ;
00845     cout <<  norme(source4_beta(3)/(nr*nt*np)) << endl ;
00846     cout << "moyenne de la source 5 pour beta_auto" << endl ;
00847     cout <<  norme(source5_beta(1)/(nr*nt*np)) << endl ;
00848     cout <<  norme(source5_beta(2)/(nr*nt*np)) << endl ;
00849     cout <<  norme(source5_beta(3)/(nr*nt*np)) << endl ;
00850     cout << "moyenne de la source 6 pour beta_auto" << endl ;
00851     cout <<  norme(source6_beta(1)/(nr*nt*np)) << endl ;
00852     cout <<  norme(source6_beta(2)/(nr*nt*np)) << endl ;
00853     cout <<  norme(source6_beta(3)/(nr*nt*np)) << endl ;
00854     cout << "moyenne de la source 7 pour beta_auto" << endl ;
00855     cout <<  norme(source7_beta(1)/(nr*nt*np)) << endl ;
00856     cout <<  norme(source7_beta(2)/(nr*nt*np)) << endl ;
00857     cout <<  norme(source7_beta(3)/(nr*nt*np)) << endl ;
00858     cout << "moyenne de la source pour beta_auto" << endl ;
00859     cout <<  norme(source_beta(1)/(nr*nt*np)) << endl ;
00860     cout <<  norme(source_beta(2)/(nr*nt*np)) << endl ;
00861     cout <<  norme(source_beta(3)/(nr*nt*np)) << endl ;
00862 
00863     // Resolution of the Poisson equation 
00864     // ----------------------------------
00865 
00866     // Filter for the source of beta vector
00867     
00868     for (int i=1; i<=3; i++) {
00869       if (source_beta(i).get_etat() != ETATZERO)
00870         source_beta.set(i).filtre(4) ;
00871     }
00872     
00873     for (int i=1; i<=3; i++) {
00874       if(source_beta(i).dz_nonzero()) {
00875         assert( source_beta(i).get_dzpuis() == 4 ) ; 
00876       }
00877       else{
00878         (source_beta.set(i)).set_dzpuis(4) ; 
00879       }
00880     }
00881     
00882     double lambda = double(1) / double(3) ; 
00883 
00884     Tenseur source_p(mp, 1, CON, mp.get_bvect_cart() ) ;
00885     source_p.set_etat_qcq() ;
00886     for (int i=0; i<3; i++) {
00887       source_p.set(i) = Cmp(source_beta(i+1)) ;
00888     }
00889     
00890     Tenseur vect_auxi (mp, 1, CON, mp.get_bvect_cart()) ;
00891     vect_auxi.set_etat_qcq() ;
00892     for (int i=0; i<3 ;i++){
00893       vect_auxi.set(i) = 0. ;
00894     }
00895     Tenseur scal_auxi (mp) ;
00896     scal_auxi.set_etat_qcq() ;
00897     scal_auxi.set().annule_hard() ;
00898     scal_auxi.set_std_base() ;
00899     
00900     Tenseur resu_p(mp, 1, CON, mp.get_bvect_cart() ) ;
00901     resu_p.set_etat_qcq() ;
00902     for (int i=0; i<3 ;i++)
00903       resu_p.set(i).annule_hard() ;
00904     resu_p.set_std_base() ;
00905 
00906     //source_p.poisson_vect(lambda, par_beta2, resu_p, vect_auxi, 
00907     //            scal_auxi) ;
00908     
00909     source_p.poisson_vect_oohara(lambda, par_beta2, resu_p, scal_auxi) ;
00910 
00911     for (int i=1; i<=3; i++) 
00912       beta_auto.set(i) = resu_p(i-1) ;
00913 
00914     ssjm1_khi = ssjm1khi ;
00915     for (int i=0; i<3; i++){
00916         ssjm1_wbeta.set(i+1) = ssjm1wbeta(i) ;
00917     }
00918     
00919     cout << "beta_auto_x" << endl << norme(beta_auto(1)/(nr*nt*np)) 
00920          << endl ;
00921     cout << "beta_auto_y" << endl << norme(beta_auto(2)/(nr*nt*np)) 
00922          << endl ;
00923     cout << "beta_auto_z" << endl << norme(beta_auto(3)/(nr*nt*np)) 
00924          << endl ;
00925   
00926 
00927     // Check: has the equation for beta_auto been correctly solved ?
00928     // --------------------------------------------------------------
00929     
00930     Vector lap_beta = (beta_auto.derive_con(flat)).divergence(flat) 
00931         + lambda* beta_auto.divergence(flat).derive_con(flat) ;
00932     
00933     source_beta.dec_dzpuis() ;
00934     Tbl tdiff_beta_x = diffrel(lap_beta(1), source_beta(1)) ; 
00935     Tbl tdiff_beta_y = diffrel(lap_beta(2), source_beta(2)) ; 
00936     Tbl tdiff_beta_z = diffrel(lap_beta(3), source_beta(3)) ; 
00937 
00938     cout << 
00939         "Relative error in the resolution of the equation for beta_auto : "
00940          << endl ; 
00941     cout << "x component : " ;
00942     for (int l=0; l<nz; l++) {
00943         cout << tdiff_beta_x(l) << "  " ; 
00944     }
00945     cout << endl ;
00946     cout << "y component : " ;
00947     for (int l=0; l<nz; l++) {
00948         cout << tdiff_beta_y(l) << "  " ; 
00949     }
00950     cout << endl ;
00951     cout << "z component : " ;
00952     for (int l=0; l<nz; l++) {
00953         cout << tdiff_beta_z(l) << "  " ; 
00954     }
00955     cout << endl ;
00956     
00957     diff_beta_x = max(abs(tdiff_beta_x)) ; 
00958     diff_beta_y = max(abs(tdiff_beta_y)) ; 
00959     diff_beta_z = max(abs(tdiff_beta_z)) ; 
00960 
00961 
00962     if (!conf_flat) {
00963        
00964         //--------------------------------------------------------
00965         // Poisson equation for hij
00966         //--------------------------------------------------------
00967      
00968      
00969         // Declaration of all sources 
00970         //---------------------------
00971 
00972         Scalar source_tot_hij(mp) ;
00973         Tensor source_Sij(mp, 2, CON, mp.get_bvect_cart()) ;
00974         Tensor source_Rij(mp, 2, CON, mp.get_bvect_cart()) ;
00975         Tensor tens_temp(mp, 2, CON, mp.get_bvect_cart()) ;
00976 
00977         Tensor source_1 (mp, 2, CON, mp.get_bvect_cart()) ;
00978         Tensor source_2 (mp, 2, CON, mp.get_bvect_cart()) ;
00979         Tensor source_3a (mp, 2, CON, mp.get_bvect_cart()) ;
00980         Tensor source_3b (mp, 2, CON, mp.get_bvect_cart()) ;
00981         Tensor source_4 (mp, 2, CON, mp.get_bvect_cart()) ;
00982         Tensor source_5 (mp, 2, CON, mp.get_bvect_cart()) ;
00983         Tensor source_6 (mp, 2, CON, mp.get_bvect_cart()) ;
00984 
00985         // Sources
00986         //-----------
00987 
00988         source_1 = contract(dcon_hij, 1, dcov_lnq_auto, 0) ;
00989 
00990         source_2 = - contract(dcon_hij, 2, dcov_lnq_auto, 0) 
00991           - 2./3. * contract(hdirac, 0, dcov_lnq_auto, 0) * flat.con() ;
00992         
00993         // Lie derivative of A^{ij}
00994         // --------------------------
00995 
00996         Scalar decouple_logn = (logn_auto - 1.e-8)/(logn - 2.e-8) ;
00997 
00998         // Function exp(-(r-r_0)^2/sigma^2)
00999         // --------------------------------
01000         
01001         double r0 = mp.val_r(nz-2, 1, 0, 0) ;
01002         double sigma = 1.*r0 ;
01003         
01004         Scalar rr (mp) ;
01005         rr = mp.r ;
01006 
01007         Scalar ff (mp) ;
01008         ff = exp( -(rr - r0)*(rr - r0)/sigma/sigma ) ;
01009         for (int ii=0; ii<nz-1; ii++)
01010         ff.set_domain(ii) = 1. ;
01011         ff.set_outer_boundary(nz-1, 0) ;
01012         ff.std_spectral_base() ;
01013         
01014         //  ff.annule_domain(nz-1) ;
01015         //des_profile(ff, 0, 20, 0, 0) ;
01016 
01017         // Construction of Omega d/dphi
01018         // ----------------------------
01019         
01020         // Construction of D_k \Phi^i
01021         Itbl type (2) ;
01022         type.set(0) = CON ;
01023         type.set(1) = COV ;
01024 
01025         Tensor dcov_omdsdphi (mp, 2, type, mp.get_bvect_cart()) ;
01026         dcov_omdsdphi.set(1,1) = 0. ;
01027         dcov_omdsdphi.set(2,1) = om * ff ; 
01028         dcov_omdsdphi.set(3,1) = 0. ;
01029         dcov_omdsdphi.set(2,2) = 0. ;
01030         dcov_omdsdphi.set(3,2) = 0. ;
01031         dcov_omdsdphi.set(3,3) = 0. ;
01032         dcov_omdsdphi.set(1,2) = -om * ff ;
01033         dcov_omdsdphi.set(1,3) = 0. ;
01034         dcov_omdsdphi.set(2,3) = 0. ;
01035         dcov_omdsdphi.std_spectral_base() ;
01036 
01037         source_3a = contract(tkij_auto, 0, dcov_omdsdphi, 1) ;
01038         source_3a.inc_dzpuis() ;
01039 
01040         // Source 3b
01041         // ------------
01042 
01043         Vector omdsdp (mp, CON, mp.get_bvect_cart()) ;
01044         Scalar yya (mp) ;
01045         yya = mp.ya ;
01046         Scalar xxa (mp) ;
01047         xxa = mp.xa ;
01048         Scalar zza (mp) ;
01049         zza = mp.za ;
01050 
01051         if (fabs(mp.get_rot_phi()) < 1e-10){ 
01052         omdsdp.set(1) = - om * yya * ff ;
01053         omdsdp.set(2) = om * xxa * ff ;
01054         omdsdp.set(3).annule_hard() ;
01055         }
01056         else{
01057         omdsdp.set(1) = om * yya * ff ;
01058         omdsdp.set(2) = - om * xxa * ff ;
01059         omdsdp.set(3).annule_hard() ;
01060         }
01061 
01062         omdsdp.set(1).set_outer_boundary(nz-1, 0) ;
01063         omdsdp.set(2).set_outer_boundary(nz-1, 0) ;
01064         omdsdp.std_spectral_base() ;
01065 
01066         source_3b = - contract(omdsdp, 0, tkij_auto.derive_cov(flat), 2) ;
01067 
01068         // Source 4
01069         // ---------
01070 
01071         source_4 = - tkij_auto.derive_lie(beta) ;
01072         source_4.inc_dzpuis() ;
01073         source_4 += - 2./3. * beta.divergence(flat) * tkij_auto ;
01074 
01075         source_5 = dcon_hdirac_auto ;
01076         
01077         source_6 = - 2./3. * hdirac_auto.divergence(flat) * flat.con() ;
01078         source_6.inc_dzpuis() ;
01079 
01080         // Source terms for Sij
01081         //---------------------
01082         
01083         source_Sij = 8. * nn / psi4 * phi_auto.derive_con(gtilde) * 
01084           contract(gtilde_con, 0, dcov_phi, 0) ;
01085 
01086         source_Sij += 4. / psi4 * phi_auto.derive_con(gtilde) * 
01087           nn * contract(gtilde_con, 0, dcov_logn, 0) +
01088           4. / psi4 * nn * contract(gtilde_con, 0, dcov_logn, 0) *
01089           phi_auto.derive_con(gtilde) ;
01090 
01091         source_Sij += - nn / (3.*psi4) * gtilde_con * 
01092           ( 0.25 * contract(gtilde_con, 0, 1, contract(dcov_hij_auto, 0, 1,
01093                     dcov_gtilde, 0, 1), 0, 1)
01094            - 0.5 * contract(gtilde_con, 0, 1, contract(dcov_hij_auto, 0, 1,
01095                     dcov_gtilde, 0, 2), 0, 1)) ;
01096 
01097         source_Sij += - 8.*nn / (3.*psi4) * gtilde_con * 
01098      contract(dcov_phi_auto, 0, contract(gtilde_con, 0, dcov_phi, 0), 0) ;
01099 
01100         tens_temp = nn / (3.*psi4) * hdirac.divergence(flat)*hij_auto ;
01101         tens_temp.inc_dzpuis() ;
01102 
01103         source_Sij += tens_temp ;
01104        
01105         source_Sij += - 8./(3.*psi4) * contract(dcov_phi_auto, 0,
01106         nn*contract(gtilde_con, 0, dcov_logn, 0), 0) * gtilde_con ;
01107 
01108         source_Sij += 2.*nn* contract(gtilde_cov, 0, 1, tkij_auto *
01109                        (tkij_auto+tkij_comp), 1, 3) ;
01110 
01111         source_Sij += - 2. * qpig * nn * ( psi4 * stress_euler 
01112                   - 0.33333333333333333 * s_euler * gtilde_con ) ; 
01113         
01114         source_Sij += - 1./(psi4*psi2) * contract(gtilde_con, 1, 
01115               contract(gtilde_con, 1, qq*dcovdcov_lnq_auto + 
01116                    qq*dcov_lnq_auto*dcov_lnq, 0), 1) ;
01117 
01118         source_Sij += - 0.5/(psi4*psi2) * contract(contract(hij, 1,
01119                     dcov_hij_auto, 2), 1, dcov_qq, 0) -
01120           0.5/(psi4*psi2) * contract(contract(dcov_hij_auto, 2,
01121                           hij, 1), 1, dcov_qq, 0) ;
01122                     
01123         source_Sij += 0.5/(psi4*psi2) * contract(contract(hij, 0,
01124                     dcov_hij_auto, 2), 0, dcov_qq, 0) ;
01125 
01126         source_Sij += 1./(3.*psi4*psi2)*contract(gtilde_con, 0, 1,
01127                     qq*dcovdcov_lnq_auto + qq*dcov_lnq_auto*dcov_lnq, 0, 1)
01128                               *gtilde_con ;
01129 
01130         source_Sij += 1./(3.*psi4*psi2) * contract(hdirac, 0, 
01131                             dcov_qq, 0)*hij_auto ;
01132 
01133         // Source terms for Rij
01134         //---------------------
01135 
01136         source_Rij = contract(hij, 0, 1, dcov_hij_auto.derive_cov(flat), 
01137                   2, 3) ;
01138         source_Rij.inc_dzpuis() ;
01139 
01140 
01141         source_Rij += - contract(hij_auto, 1, dcov_hdirac, 1) -
01142           contract(dcov_hdirac, 1, hij_auto, 1) ;
01143         
01144         source_Rij += contract(hdirac, 0, dcov_hij_auto, 2) ;
01145 
01146         source_Rij += - contract(contract(dcov_hij_auto, 1, dcov_hij, 2),
01147                      1, 3) ;
01148 
01149         source_Rij += - contract(gtilde_cov, 0, 1, contract(contract(
01150             gtilde_con, 0, dcov_hij_auto, 2), 0, dcov_hij, 2), 1, 3) ;
01151 
01152         source_Rij += contract(contract(contract(contract(gtilde_cov, 0, 
01153          dcov_hij_auto, 1), 2, gtilde_con, 1), 0, dcov_hij, 1), 0, 3) +
01154           contract(contract(contract(contract(gtilde_cov, 0, 
01155          dcov_hij_auto, 1), 0, dcov_hij, 1), 0, 3), 0, gtilde_con, 1) ;
01156 
01157         source_Rij += 0.5 * contract(gtilde_con*gtilde_con, 1, 3, 
01158            contract(dcov_hij_auto, 0, 1, dcov_gtilde, 0, 1), 0, 1) ;
01159 
01160         source_Rij = source_Rij * 0.5 ;
01161 
01162         for(int i=1; i<=3; i++) 
01163         for(int j=1; j<=i; j++) {
01164 
01165             source_tot_hij = source_1(i,j) + source_1(j,i) 
01166             + source_2(i,j) + 2.*psi4/nn * (
01167                 source_4(i,j) - source_Sij(i,j)) 
01168               - 2.* source_Rij(i,j) +
01169               source_5(i,j) + source_5(j,i) + source_6(i,j) ;
01170             source_tot_hij.dec_dzpuis() ;
01171             
01172             source3 = 2.*psi4/nn * (source_3a(i,j) + source_3a(j,i) + 
01173                         source_3b(i,j)) ; 
01174     
01175             source_tot_hij = source_tot_hij + source3 ;
01176 
01177             //source_tot_hij.inc_dzpuis() ;
01178 
01179             cout << "source_mat" << endl 
01180              << norme((- 2. * qpig * nn * ( psi4 * stress_euler 
01181                    - 0.33333333333333333 * s_euler * gtilde_con ))
01182                   (i,j))/(nr*nt*np) << endl ;
01183             cout << "max source_mat" << endl 
01184              << max((- 2. * qpig * nn * ( psi4 * stress_euler 
01185                    - 0.33333333333333333 * s_euler * gtilde_con ))
01186                   (i,j)) << endl ;
01187         
01188             cout << "source1" << endl 
01189              << norme(source_1(i,j)/(nr*nt*np)) << endl ;
01190             cout << "max source1" << endl 
01191              << max(source_1(i,j)) << endl ;
01192             cout << "source2" << endl 
01193              << norme(source_2(i,j)/(nr*nt*np)) << endl ;
01194             cout << "max source2" << endl 
01195              << max(source_2(i,j)) << endl ;
01196             cout << "source3a" << endl 
01197              << norme(source_3a(i,j)/(nr*nt*np)) << endl ;
01198             cout << "max source3a" << endl 
01199              << max(source_3a(i,j)) << endl ;
01200                     cout << "source3b" << endl
01201                          << norme(source_3b(i,j)/(nr*nt*np)) << endl ;
01202                     cout << "max source3b" << endl
01203                          << max(source_3b(i,j)) << endl ;
01204             cout << "source4" << endl 
01205              << norme(source_4(i,j)/(nr*nt*np)) << endl ;
01206             cout << "max source4" << endl 
01207              << max(source_4(i,j)) << endl ;
01208             cout << "source5" << endl 
01209              << norme(source_5(i,j)/(nr*nt*np)) << endl ;
01210             cout << "max source5" << endl 
01211              << max(source_5(i,j)) << endl ;
01212             cout << "source6" << endl 
01213              << norme(source_6(i,j)/(nr*nt*np)) << endl ;
01214             cout << "max source6" << endl 
01215              << max(source_6(i,j)) << endl ;
01216             cout << "source_Rij" << endl 
01217              << norme(source_Rij(i,j)/(nr*nt*np)) << endl ;
01218             cout << "max source_Rij" << endl 
01219              << max(source_Rij(i,j)) << endl ;
01220             cout << "source_Sij" << endl 
01221              << norme(source_Sij(i,j)/(nr*nt*np)) << endl ;
01222             cout << "max source_Sij" << endl 
01223              << max(source_Sij(i,j)) << endl ;
01224             cout << "source_tot" << endl 
01225              << norme(source_tot_hij/(nr*nt*np)) << endl ;
01226             cout << "max source_tot" << endl 
01227              << max(source_tot_hij) << endl ;
01228                 
01229             // Resolution of the Poisson equations and
01230             // Check: has the Poisson equation been correctly solved ?
01231             // -----------------------------------------------------
01232 
01233             if(i==1 && j==1) {
01234          
01235             source_tot_hij.poisson(par_h11, hij_auto.set(1,1)) ; 
01236             
01237             Scalar laplac (hij_auto(1,1).laplacian()) ;
01238             laplac.dec_dzpuis() ;
01239             Tbl tdiff_h11 = diffrel(laplac, source_tot_hij) ;  
01240             cout << "Relative error in the resolution of the equation for "
01241                  << "h11_auto : " << endl ; 
01242             for (int l=0; l<nz; l++) {
01243                 cout << tdiff_h11(l) << "  " ; 
01244             }
01245             cout << endl ;
01246             diff_h11 = max(abs(tdiff_h11)) ; 
01247             }
01248                    
01249             if(i==2 && j==1) {
01250 
01251             source_tot_hij.poisson(par_h21, hij_auto.set(2,1)) ; 
01252         
01253             Scalar laplac (hij_auto(2,1).laplacian()) ;
01254             laplac.dec_dzpuis() ;
01255             Tbl tdiff_h21 = diffrel(laplac, source_tot_hij) ;  
01256     
01257             cout << 
01258                 "Relative error in the resolution of the equation for " 
01259                  << "h21_auto : "  << endl ; 
01260             for (int l=0; l<nz; l++) {
01261                 cout << tdiff_h21(l) << "  " ; 
01262             }
01263             cout << endl ;
01264             diff_h21 = max(abs(tdiff_h21)) ; 
01265             }
01266            
01267             if(i==3 && j==1) {
01268          
01269             source_tot_hij.poisson(par_h31, hij_auto.set(3,1)) ; 
01270 
01271             Scalar laplac (hij_auto(3,1).laplacian()) ;
01272             laplac.dec_dzpuis() ;
01273             Tbl tdiff_h31 = diffrel(laplac, source_tot_hij) ;  
01274 
01275             cout << 
01276                 "Relative error in the resolution of the equation for "
01277                  << "h31_auto : " << endl ; 
01278             for (int l=0; l<nz; l++) {
01279                 cout << tdiff_h31(l) << "  " ; 
01280             }
01281             cout << endl ;
01282             diff_h31 = max(abs(tdiff_h31)) ; 
01283             }
01284          
01285             if(i==2 && j==2) {
01286          
01287             source_tot_hij.poisson(par_h22, hij_auto.set(2,2)) ; 
01288 
01289             Scalar laplac (hij_auto(2,2).laplacian()) ;
01290             laplac.dec_dzpuis() ;
01291             Tbl tdiff_h22 = diffrel(laplac, source_tot_hij) ;  
01292 
01293             cout << 
01294                 "Relative error in the resolution of the equation for "
01295                  << "h22_auto : " << endl ; 
01296             for (int l=0; l<nz; l++) {
01297                 cout << tdiff_h22(l) << "  " ; 
01298             }
01299             cout << endl ;
01300             diff_h22 = max(abs(tdiff_h22)) ; 
01301             }
01302            
01303             if(i==3 && j==2) {
01304          
01305             source_tot_hij.poisson(par_h32, hij_auto.set(3,2)) ; 
01306     
01307             Scalar laplac (hij_auto(3,2).laplacian()) ;
01308             laplac.dec_dzpuis() ;
01309             Tbl tdiff_h32 = diffrel(laplac, source_tot_hij) ;  
01310 
01311             cout << 
01312                 "Relative error in the resolution of the equation for "
01313                  << "h32_auto : " << endl ; 
01314             for (int l=0; l<nz; l++) {
01315                 cout << tdiff_h32(l) << "  " ; 
01316             }
01317             cout << endl ;
01318             diff_h32 = max(abs(tdiff_h32)) ; 
01319             }
01320          
01321             if(i==3 && j==3) {
01322          
01323             source_tot_hij.poisson(par_h33, hij_auto.set(3,3)) ; 
01324 
01325             Scalar laplac (hij_auto(3,3).laplacian()) ;
01326             laplac.dec_dzpuis() ;
01327             Tbl tdiff_h33 = diffrel(laplac, source_tot_hij) ;  
01328 
01329             cout << 
01330                 "Relative error in the resolution of the equation for "
01331                  << "h33_auto : " << endl ;
01332             for (int l=0; l<nz; l++) {
01333                 cout << tdiff_h33(l) << "  " ;
01334             }
01335             cout << endl ;
01336             diff_h33 = max(abs(tdiff_h33)) ;
01337             }
01338 
01339         }
01340       
01341         cout << "Tenseur hij auto in cartesian coordinates" << endl ;
01342         for (int i=1; i<=3; i++)
01343         for (int j=1; j<=i; j++) {
01344             cout << "  Comp. " << i << " " << j << " :  " ;
01345             for (int l=0; l<nz; l++){
01346             cout << norme(hij_auto(i,j)/(nr*nt*np))(l) << " " ;
01347             }
01348             cout << endl ;
01349         }
01350         cout << endl ;
01351 
01352         ssjm1_h11 = ssjm1h11 ;
01353         ssjm1_h21 = ssjm1h21 ;
01354         ssjm1_h31 = ssjm1h31 ;
01355         ssjm1_h22 = ssjm1h22 ;
01356         ssjm1_h32 = ssjm1h32 ;
01357         ssjm1_h33 = ssjm1h33 ;
01358 
01359     }
01360 
01361     // End of relativistic equations    
01362        
01363            
01364     //-------------------------------------------------
01365     //  Relative change in enthalpy
01366     //-------------------------------------------------
01367 
01368     Tbl diff_ent_tbl = diffrel( ent, ent_jm1 ) ; 
01369     diff_ent = diff_ent_tbl(0) ; 
01370     for (int l=1; l<nzet; l++) {
01371         diff_ent += diff_ent_tbl(l) ; 
01372     }
01373     diff_ent /= nzet ; 
01374     
01375     
01376     ent_jm1 = ent ; 
01377 
01378 
01379         } // End of main loop
01380     
01381     //=========================================================================
01382     //          End of iteration
01383     //=========================================================================
01384 
01385 }  
01386 
01387 
01388 
01389 
01390 
01391     
01392     

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