LORENE
strot_diff_hydro.C
1 /*
2  * Function Star_rot_diff::hydro_euler
3  *
4  * (see file star_rot_diff.h for documentation)
5  *
6  */
7 
8 /*
9  * Copyright (c) 2025 Santiago Jaraba
10  *
11  * This file is part of LORENE.
12  *
13  * LORENE is free software; you can redistribute it and/or modify
14  * it under the terms of the GNU General Public License as published by
15  * the Free Software Foundation; either version 2 of the License, or
16  * (at your option) any later version.
17  *
18  * LORENE is distributed in the hope that it will be useful,
19  * but WITHOUT ANY WARRANTY; without even the implied warranty of
20  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21  * GNU General Public License for more details.
22  *
23  * You should have received a copy of the GNU General Public License
24  * along with LORENE; if not, write to the Free Software
25  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
26  *
27  */
28 
29 
30 // C headers
31 #include <cstdlib>
32 
33 // Lorene headers
34 #include "star_rot_diff.h"
35 #include "utilitaires.h"
36 
37 namespace Lorene {
39 
40  int nz = mp.get_mg()->get_nzone() ;
41  int nzm1 = nz - 1 ;
42 
43  // Computation of u_euler
44  // ----------------------
45 
46  Scalar x(mp) ;
47  Scalar y(mp) ;
48  x = mp.x ;
49  y = mp.y ;
50 
52 
53  // Cartesian components of differential rotation:
54 
55  u_euler.set(1) = - omega_field * y ;
56  u_euler.set(2) = omega_field * x ;
57  u_euler.set(3) = 0 ;
58  u_euler.annule_domain(nzm1) ;
59 
60  u_euler.set_triad( mp.get_bvect_cart() ) ; // Triad = Cartesian triad
61 
62  u_euler.std_spectral_base() ; // sets the standard bases for spectral expansions
63 
64  u_euler = ( u_euler + beta ) / nn ;
65 
66  u_euler.std_spectral_base() ; // sets the standard bases for spectral expansions
67 
68 //## Test
69  Vector utest(mp, CON, mp.get_bvect_spher()) ;
70  utest.set_etat_qcq() ;
71 
72  utest.set(1) = 0 ; // Spherical components of differential rotation
73  utest.set(2) = 0 ;
74  utest.set(3) = ( omega_field - nphi ) / nn;
75 
76  utest.set(3).annule_domain(nzm1) ;
77  utest.set(3).std_spectral_base() ;
78  utest.set(3).mult_rsint() ; // Multiplication by r sin(theta)
79 
80  utest.set_triad( mp.get_bvect_spher() ) ;
81 
82  utest.change_triad( mp.get_bvect_cart() ) ;
83 
84  for (int i=1; i<=3; i++) {
85  const Valeur& uu = u_euler(i).get_spectral_va() ;
86  Valeur& ut = utest.set(i).set_spectral_va() ;
87 
88  if (uu.get_etat() != ETATZERO) {
89  uu.coef() ;
90 
91  if (ut.get_etat() == ETATZERO) {
92  ut.set_etat_cf_qcq() ;
93  ut.c_cf = 0 ;
94  ut.c_cf->base = uu.c_cf->base ;
95  }
96  else {
97  ut.coef() ;
98  }
99 
100  Mtbl_cf diff = *(uu.c_cf) - *(ut.c_cf) ;
101  cout << "Star_rot_diff::hydro_euler: test u_euler(" << i << ") : "
102  << max( abs(diff) )(0) << endl ;
103 
104  }
105  }
106 //##
107 
108  if ( (u_euler(1).get_etat() == ETATZERO) &&
109  (u_euler(2).get_etat() == ETATZERO) &&
110  (u_euler(3).get_etat() == ETATZERO) ) {
111 
112  for (int i=1; i<=3; i++) u_euler.set(i) = 0 ;
113  }
114 
115 
116  // Computation of uuu (norme of u_euler)
117  // ------------------
118 
119  // The scalar product is performed on the spherical components:
120 
121  Vector us = u_euler ;
122  us.change_triad( mp.get_bvect_spher() ) ;
123 
124  Scalar uuu2 = a_car * ( us(1) * us(1) + us(2) * us(2) )
125  + b_car * us(3) * us(3) ;
126 
127  uuu = sqrt( uuu2 ) ;
128 
129  if (uuu.get_etat() == ETATQCQ) {
130  uuu.set_spectral_base( us(3).get_spectral_base() ) ; // Same basis as
131  } // (Omega -N^phi) r sin(theta)
132 
133  // Lorentz factor
134  // --------------
135 
136  Scalar u2(mp) ;
137  u2 = unsurc2 * uuu2 ;
138 
139  Scalar gam2 = 1 / (1 - u2) ;
140 
141  gam_euler = sqrt(gam2) ;
142 
143  gam_euler.std_spectral_base() ; // sets the standard spectral bases for
144  // a scalar field
145 
146  // Energy density E with respect to the Eulerian observer
147  //------------------------------------
148 
149  ener_euler = gam2 * ( ener + press ) - press ;
150 
152 
153  // Trace of the stress tensor with respect to the Eulerian observer
154  //------------------------------------
155 
156  s_euler = 3 * press + ( ener_euler + press ) * u2 ;
157 
159 
160  // The derived quantities are obsolete
161  // -----------------------------------
162 
163  del_deriv() ;
164 
165 
166 }
167 }
void annule_domain(int l)
Sets the Tensor to zero in a given domain.
Definition: tensor.C:676
virtual void set_etat_qcq()
Sets the logical state of all components to ETATQCQ (ordinary state).
Definition: tensor.C:491
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
Definition: valeur.h:312
Scalar a_car
Square of the metric factor A.
Definition: star_rot.h:116
Map & mp
Mapping associated with the star.
Definition: star.h:180
void set_etat_cf_qcq()
Sets the logical state to ETATQCQ (ordinary state) for values in the configuration space (Mtbl_cf c_c...
Definition: valeur.C:715
void coef() const
Computes the coeffcients of *this.
Definition: valeur_coef.C:151
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
Definition: map.h:810
Lorene prototypes.
Definition: app_hor.h:67
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:792
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:399
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: star_rot.C:310
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field. ...
Definition: scalar.C:790
Values and coefficients of a (real-value) function.
Definition: valeur.h:297
virtual void change_triad(const Base_vect &)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
int get_etat() const
Returns the logical state ETATNONDEF (undefined), ETATZERO (null) or ETATQCQ (ordinary).
Definition: scalar.h:566
Scalar s_euler
Trace of the stress scalar in the Eulerian frame.
Definition: star.h:201
Vector beta
Shift vector.
Definition: star.h:228
Tensor field of valence 1.
Definition: vector.h:188
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
virtual void std_spectral_base()
Sets the standard spectal bases of decomposition for each component.
Definition: vector.C:322
double unsurc2
: unsurc2=1 for a relativistic star, 0 for a Newtonian one.
Definition: star_rot.h:111
int get_etat() const
Returns the logical state.
Definition: valeur.h:760
Scalar nphi
Metric coefficient .
Definition: star_rot.h:125
Scalar gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition: star.h:204
Scalar ener
Total energy density in the fluid frame.
Definition: star.h:193
Scalar press
Fluid pressure.
Definition: star.h:194
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:467
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
Definition: cmp_math.C:438
Vector u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition: star.h:207
Scalar b_car
Square of the metric factor B.
Definition: star_rot.h:122
void set_spectral_base(const Base_val &)
Sets the spectral bases of the Valeur va
Definition: scalar.C:803
const Base_vect_cart & get_bvect_cart() const
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e.
Definition: map.h:818
Scalar & set(const Itbl &ind)
Returns the value of a component (read/write version).
Definition: tensor.C:664
Scalar nn
Lapse function N .
Definition: star.h:225
Coord y
y coordinate centered on the grid
Definition: map.h:754
Cmp abs(const Cmp &)
Absolute value.
Definition: cmp_math.C:413
Scalar uuu
Norm of u_euler.
Definition: star_rot.h:133
Coefficients storage for the multi-domain spectral method.
Definition: mtbl_cf.h:196
Coord x
x coordinate centered on the grid
Definition: map.h:753
void set_triad(const Base_vect &new_triad)
Assigns a new vectorial basis (triad) of decomposition.
Definition: tensor.C:529
Base_val base
Bases of the spectral expansions.
Definition: mtbl_cf.h:210
Scalar & set(int)
Read/write access to a component.
Definition: vector.C:302
Scalar ener_euler
Total energy density in the Eulerian frame.
Definition: star.h:198
Scalar omega_field
Field .
Definition: star_rot_diff.h:73
Tbl & set(int l)
Read/write of the value in a given domain (configuration space).
Definition: valeur.h:373