LORENE
etoile_eqsph_falloff.C
1 /*
2  * Method of class Etoile to compute a static spherical configuration
3  * with the outer boundary condition at a finite radius
4  *
5  * (see file etoile.h for documentation).
6  *
7  */
8 
9 /*
10  * Copyright (c) 2004 Keisuke Taniguchi
11  *
12  * This file is part of LORENE.
13  *
14  * LORENE is free software; you can redistribute it and/or modify
15  * it under the terms of the GNU General Public License as published by
16  * the Free Software Foundation; either version 2 of the License, or
17  * (at your option) any later version.
18  *
19  * LORENE is distributed in the hope that it will be useful,
20  * but WITHOUT ANY WARRANTY; without even the implied warranty of
21  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
22  * GNU General Public License for more details.
23  *
24  * You should have received a copy of the GNU General Public License
25  * along with LORENE; if not, write to the Free Software
26  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
27  *
28  */
29 
30 
31 
32 
33 /*
34  * $Id: etoile_eqsph_falloff.C,v 1.3 2016/12/05 16:17:54 j_novak Exp $
35  * $Log: etoile_eqsph_falloff.C,v $
36  * Revision 1.3 2016/12/05 16:17:54 j_novak
37  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
38  *
39  * Revision 1.2 2014/10/13 08:52:58 j_novak
40  * Lorene classes and functions now belong to the namespace Lorene.
41  *
42  * Revision 1.1 2004/11/30 20:52:24 k_taniguchi
43  * *** empty log message ***
44  *
45  *
46  *
47  * $Header: /cvsroot/Lorene/C++/Source/Etoile/etoile_eqsph_falloff.C,v 1.3 2016/12/05 16:17:54 j_novak Exp $
48  *
49  */
50 
51 // Headers C
52 #include "math.h"
53 
54 // Headers Lorene
55 #include "etoile.h"
56 #include "param.h"
57 #include "unites.h"
58 
59 namespace Lorene {
60 void Etoile::equil_spher_falloff(double ent_c, double precis) {
61 
62  // Fundamental constants and units
63  // -------------------------------
64  using namespace Unites ;
65 
66  // Initializations
67  // ---------------
68 
69  const Mg3d* mg = mp.get_mg() ;
70  int nz = mg->get_nzone() ; // total number of domains
71 
72  // Index of the point at phi=0, theta=pi/2 at the surface of the star:
73  int l_b = nzet - 1 ;
74  int i_b = mg->get_nr(l_b) - 1 ;
75  int j_b = mg->get_nt(l_b) - 1 ;
76  int k_b = 0 ;
77 
78  // Value of the enthalpy defining the surface of the star
79  double ent_b = 0 ;
80 
81  // Initialization of the enthalpy field to the constant value ent_c :
82 
83  ent = ent_c ;
84  ent.annule(nzet, nz-1) ;
85 
86 
87  // Corresponding profiles of baryon density, energy density and pressure
88 
90 
91  // Initial metric
92  a_car = 1 ; // this value will remain unchanged in the Newtonian case
93  beta_auto = 0 ; // this value will remain unchanged in the Newtonian case
94 
95 
96  // Auxiliary quantities
97  // --------------------
98 
99  // Affine mapping for solving the Poisson equations
100  Map_af mpaff(mp);
101 
102  Param par_nul ; // Param (null) for Map_af::poisson.
103 
104  Tenseur ent_jm1(mp) ; // Enthalpy at previous step
105  ent_jm1 = 0 ;
106 
107  Tenseur source(mp) ;
108  Tenseur logn(mp) ;
109  Tenseur logn_quad(mp) ;
110  logn = 0 ;
111  logn_quad = 0 ;
112 
113  Cmp dlogn(mp) ;
114  Cmp dbeta(mp) ;
115 
116  double diff_ent = 1 ;
117  int mermax = 200 ; // Max number of iterations
118 
119  double alpha_r = 1 ;
120  int k_falloff = 1 ;
121 
122  //=========================================================================
123  // Start of iteration
124  //=========================================================================
125 
126  for(int mer=0 ; (diff_ent > precis) && (mer<mermax) ; mer++ ) {
127 
128  cout << "-----------------------------------------------" << endl ;
129  cout << "step: " << mer << endl ;
130  cout << "alpha_r: " << alpha_r << endl ;
131  cout << "diff_ent = " << diff_ent << endl ;
132 
133  //-----------------------------------------------------
134  // Resolution of Poisson equation for ln(N)
135  //-----------------------------------------------------
136 
137  // Matter part of ln(N)
138  // --------------------
139  if (relativistic) {
140  source = a_car * (ener + 3*press) ;
141  }
142  else {
143  source = nbar ;
144  }
145 
146  (source.set()).set_dzpuis(4) ;
147 
148  source.set_std_base() ; // Sets the standard spectral bases.
149 
151 
152  mpaff.poisson_falloff(source(), par_nul, logn_auto.set(), k_falloff) ;
153 
154  // NB: at this stage logn_auto is in machine units, not in c^2
155 
156  // Quadratic part of ln(N)
157  // -----------------------
158 
159  if (relativistic) {
160 
161  mpaff.dsdr(logn(), dlogn) ;
162  mpaff.dsdr(beta_auto(), dbeta) ;
163 
164  source = - dlogn * dbeta ;
165 
166  logn_quad.set_etat_qcq() ;
167 
168  mpaff.poisson_falloff(source(), par_nul, logn_quad.set(),
169  k_falloff) ;
170 
171  }
172 
173  //-----------------------------------------------------
174  // Computation of the new radial scale
175  //-----------------------------------------------------
176 
177  // alpha_r (r = alpha_r r') is determined so that the enthalpy
178  // takes the requested value ent_b at the stellar surface
179 
180  double nu_mat0_b = logn_auto()(l_b, k_b, j_b, i_b) ;
181  double nu_mat0_c = logn_auto()(0, 0, 0, 0) ;
182 
183  double nu_quad0_b = logn_quad()(l_b, k_b, j_b, i_b) ;
184  double nu_quad0_c = logn_quad()(0, 0, 0, 0) ;
185 
186  double alpha_r2 = ( ent_c - ent_b - nu_quad0_b + nu_quad0_c )
187  / ( qpig*(nu_mat0_b - nu_mat0_c) ) ;
188 
189  alpha_r = sqrt(alpha_r2) ;
190 
191  // New radial scale
192  mpaff.homothetie( alpha_r ) ;
193 
194  //--------------------
195  // First integral
196  //--------------------
197 
198  // Gravitation potential in units c^2 :
199  logn_auto = alpha_r2*qpig * logn_auto ;
200  logn = logn_auto + logn_quad ;
201 
202  // Enthalpy in all space
203  double logn_c = logn()(0, 0, 0, 0) ;
204  ent = ent_c - logn() + logn_c ;
205 
206  //---------------------
207  // Equation of state
208  //---------------------
209 
210  equation_of_state() ;
211 
212  if (relativistic) {
213 
214  //----------------------------
215  // Equation for beta = ln(AN)
216  //----------------------------
217 
218  mpaff.dsdr(logn(), dlogn) ;
219  mpaff.dsdr(beta_auto(), dbeta) ;
220 
221  source = 3 * qpig * a_car * press ;
222 
223  source = source()
224  - 0.5 * ( dlogn * dlogn + dbeta * dbeta ) ;
225 
226  source.set_std_base() ; // Sets the standard spectral bases.
227 
229 
230  mpaff.poisson_falloff(source(), par_nul, beta_auto.set(),
231  k_falloff) ;
232 
233 
234  // Metric coefficient A^2 update
235 
236  a_car = exp(2*(beta_auto - logn)) ;
237 
238  }
239 
240  // Relative difference with enthalpy at the previous step
241  // ------------------------------------------------------
242 
243  diff_ent = norme( diffrel(ent(), ent_jm1()) ) / nzet ;
244 
245  // Next step
246  // ---------
247 
248  ent_jm1 = ent ;
249 
250 
251  } // End of iteration loop
252 
253  //=========================================================================
254  // End of iteration
255  //=========================================================================
256 
257 
258  // The mapping is transfered to that of the star:
259  // ----------------------------------------------
260  mp = mpaff ;
261 
262 
263  // Sets value to all the Tenseur's of the star
264  // -------------------------------------------
265 
266  // ... hydro
267  ent.annule(nzet, nz-1) ; // enthalpy set to zero at the exterior of
268  // the star
269  ener_euler = ener ;
270  s_euler = 3 * press ;
271  gam_euler = 1 ;
272  u_euler = 0 ;
273 
274  // ... metric
275  nnn = exp( unsurc2 * logn ) ;
276  nnn.set_std_base() ;
277  shift = 0 ;
278  a_car.set_std_base() ;
279 
280  // Info printing
281  // -------------
282 
283  cout << endl
284  << "Characteristics of the star obtained by Etoile::equil_spher_falloff : "
285  << endl
286  << "-------------------------------------------------------------------"
287  << endl ;
288 
289  double ray = mp.val_r(l_b, 1., M_PI/2., 0) ;
290  cout << "Coordinate radius : " << ray / km << " km" << endl ;
291 
292  double rcirc = ray * sqrt( a_car()(l_b, k_b, j_b, i_b) ) ;
293 
294  double compact = qpig/(4.*M_PI) * mass_g() / rcirc ;
295 
296  cout << "Circumferential radius R : " << rcirc/km << " km" << endl ;
297  cout << "Baryon mass : " << mass_b()/msol << " Mo" << endl ;
298  cout << "Gravitational mass M : " << mass_g()/msol << " Mo" << endl ;
299  cout << "Compacity parameter GM/(c^2 R) : " << compact << endl ;
300 
301 
302  //-----------------
303  // Virial theorem
304  //-----------------
305 
306  //... Pressure term
307 
308  source = qpig * a_car * sqrt(a_car) * s_euler ;
309  source.set_std_base() ;
310  double vir_mat = source().integrale() ;
311 
312  //... Gravitational term
313 
314  Cmp tmp = beta_auto() - logn() ;
315 
316  source = - ( logn().dsdr() * logn().dsdr()
317  - 0.5 * tmp.dsdr() * tmp.dsdr() )
318  * sqrt(a_car()) ;
319 
320  source.set_std_base() ;
321 
322  double vir_grav = source().integrale() ;
323 
324  //... Relative error on the virial identity GRV3
325 
326  double grv3 = ( vir_mat + vir_grav ) / vir_mat ;
327 
328  cout << "Virial theorem GRV3 : " << endl ;
329  cout << " 3P term : " << vir_mat << endl ;
330  cout << " grav. term : " << vir_grav << endl ;
331  cout << " relative error : " << grv3 << endl ;
332 
333 }
334 }
const Cmp & dsdr() const
Returns of *this .
Definition: cmp_deriv.C:87
void annule(int l)
Sets the Tenseur to zero in a given domain.
Definition: tenseur.C:916
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition: cmp.h:446
Cmp exp(const Cmp &)
Exponential.
Definition: cmp_math.C:273
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
virtual void homothetie(double lambda)
Sets a new radial scale.
Definition: map_af.C:664
void set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition: tenseur.C:1186
virtual double mass_g() const
Gravitational mass.
Lorene prototypes.
Definition: app_hor.h:67
Standard units of space, time and mass.
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:777
Tenseur nnn
Total lapse function.
Definition: etoile.h:512
Tenseur s_euler
Trace of the stress tensor in the Eulerian frame.
Definition: etoile.h:471
double unsurc2
: unsurc2=1 for a relativistic star, 0 for a Newtonian one.
Definition: etoile.h:445
Tenseur press
Fluid pressure.
Definition: etoile.h:464
Tenseur shift
Total shift vector.
Definition: etoile.h:515
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition: cmp_math.C:507
virtual double val_r(int l, double xi, double theta, double pphi) const =0
Returns the value of the radial coordinate r for a given in a given domain.
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition: tenseur.C:840
Tenseur u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition: etoile.h:477
Tbl norme(const Cmp &)
Sums of the absolute values of all the values of the Cmp in each domain.
Definition: cmp_math.C:484
Tenseur nbar
Baryon density in the fluid frame.
Definition: etoile.h:462
Tenseur gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition: etoile.h:474
Parameter storage.
Definition: param.h:125
virtual double mass_b() const
Baryon mass.
Map & mp
Mapping associated with the star.
Definition: etoile.h:432
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:465
virtual void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Definition: etoile.C:569
int nzet
Number of domains of *mp occupied by the star.
Definition: etoile.h:435
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition: grilles.h:469
Tenseur a_car
Total conformal factor .
Definition: etoile.h:518
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one. ...
Definition: etoile.h:440
Multi-domain grid.
Definition: grilles.h:279
Tenseur ener
Total energy density in the fluid frame.
Definition: etoile.h:463
Affine radial mapping.
Definition: map.h:2042
Tenseur logn_auto
Total of the logarithm of the part of the lapse N generated principaly by the star.
Definition: etoile.h:487
Tenseur ent
Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case)
Definition: etoile.h:460
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition: grilles.h:474
virtual void equil_spher_falloff(double ent_c, double precis=1.e-14)
Computes a spherical static configuration with the outer boundary condition at a finite radius...
Tenseur beta_auto
Logarithm of the part of the product AN generated principaly by by the star.
Definition: etoile.h:509
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tenseur.C:652
Tenseur ener_euler
Total energy density in the Eulerian frame.
Definition: etoile.h:468
virtual void dsdr(const Cmp &ci, Cmp &resu) const
Computes of a Cmp.
Definition: map_af_deriv.C:282
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition: tenseur.h:304