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