LORENE
et_rot_hydro.C
1 /*
2  * Function Etoile_rot::hydro_euler
3  *
4  * (see file etoile.h for documentation)
5  *
6  */
7 
8 /*
9  * Copyright (c) 2000-2001 Eric Gourgoulhon
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 
31 
32 /*
33  * $Id: et_rot_hydro.C,v 1.5 2016/12/05 16:17:54 j_novak Exp $
34  * $Log: et_rot_hydro.C,v $
35  * Revision 1.5 2016/12/05 16:17:54 j_novak
36  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
37  *
38  * Revision 1.4 2014/10/13 08:52:58 j_novak
39  * Lorene classes and functions now belong to the namespace Lorene.
40  *
41  * Revision 1.3 2014/10/06 15:13:09 j_novak
42  * Modified #include directives to use c++ syntax.
43  *
44  * Revision 1.2 2006/09/14 07:37:47 j_novak
45  * Removal of a test on u_euler.
46  *
47  * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
48  * LORENE
49  *
50  * Revision 2.6 2000/10/12 15:36:06 eric
51  * Le test sur u_euler est desormais effectue sur la totalite de u_euler
52  * (ie comprenant le shift).
53  *
54  * Revision 2.5 2000/10/06 15:08:13 eric
55  * Calcul 3-D de u_euler.
56  * uuu s'en deduit.
57  *
58  * Revision 2.4 2000/08/31 15:38:34 eric
59  * Appel du nouvel operateur Cmp::mult_rsint pour le calcul de uuu.
60  *
61  * Revision 2.3 2000/08/25 12:28:21 eric
62  * *** empty log message ***
63  *
64  * Revision 2.2 2000/08/18 14:01:49 eric
65  * *** empty log message ***
66  *
67  * Revision 2.1 2000/08/17 12:39:56 eric
68  * *** empty log message ***
69  *
70  * Revision 2.0 2000/07/21 16:31:00 eric
71  * *** empty log message ***
72  *
73  *
74  * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_rot_hydro.C,v 1.5 2016/12/05 16:17:54 j_novak Exp $
75  *
76  */
77 
78 // Headers C
79 #include <cstdlib>
80 
81 // Headers Lorene
82 #include "etoile.h"
83 #include "utilitaires.h"
84 
85 namespace Lorene {
87 
88  int nz = mp.get_mg()->get_nzone() ;
89  int nzm1 = nz - 1 ;
90 
91  // Computation of u_euler
92  // ----------------------
93 
94  const Coord& x = mp.x ;
95  const Coord& y = mp.y ;
96 
98 
99  u_euler.set(0) = - omega * y ; // Cartesian components of solid rotation
100  u_euler.set(1) = omega * x ;
101  u_euler.set(2) = 0 ;
102  u_euler.annule(nzm1) ;
103 
104  u_euler.set_triad( mp.get_bvect_cart() ) ; // Triad = Cartesian triad
105 
106  u_euler.set_std_base() ; // sets the standard bases for spectral expansions
107 
108  u_euler = ( u_euler - shift ) / nnn ;
109 
110  u_euler.set_std_base() ; // sets the standard bases for spectral expansions
111 
112  if ( (u_euler(0).get_etat() == ETATZERO) &&
113  (u_euler(1).get_etat() == ETATZERO) &&
114  (u_euler(2).get_etat() == ETATZERO) ) {
115 
116  u_euler = 0 ;
117  }
118 
119 
120  // Computation of uuu (norme of u_euler)
121  // ------------------
122 
123  // The scalar product is performed on the spherical components:
124 
125  Tenseur us = u_euler ;
126  us.change_triad( mp.get_bvect_spher() ) ;
127 
128  Cmp uuu2 = a_car() * ( us(0) * us(0) + us(1) * us(1) )
129  + b_car() * us(2) * us(2) ;
130 
131  uuu = sqrt( uuu2 ) ;
132 
133  if (uuu.get_etat() == ETATQCQ) {
134  ((uuu.set()).va).set_base( us(2).va.base ) ; // Same basis as
135  } // (Omega -N^phi) r sin(theta)
136 
137 
138  // Lorentz factor
139  // --------------
140 
141  Tenseur u2(mp) ;
142  u2 = unsurc2 * uuu2 ;
143 
144  Tenseur gam2 = 1 / (1 - u2) ;
145 
146  gam_euler = sqrt(gam2) ;
147 
148  gam_euler.set_std_base() ; // sets the standard spectral bases for
149  // a scalar field
150 
151  // Energy density E with respect to the Eulerian observer
152  //------------------------------------
153 
154  ener_euler = gam2 * ( ener + press ) - press ;
155 
157 
158  // Trace of the stress tensor with respect to the Eulerian observer
159  //------------------------------------
160 
161  s_euler = 3 * press + ( ener_euler + press ) * u2 ;
162 
163  s_euler.set_std_base() ;
164 
165  // The derived quantities are obsolete
166  // -----------------------------------
167 
168  del_deriv() ;
169 
170 
171 }
172 }
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
void set_triad(const Base_vect &new_triad)
Assigns a new vectorial basis (triad) of decomposition.
Definition: tenseur.C:690
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
void set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition: tenseur.C:1186
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
Definition: map.h:795
Lorene prototypes.
Definition: app_hor.h:67
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
Tenseur b_car
Square of the metric factor B.
Definition: etoile.h:1510
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
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: etoile_rot.C:344
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
void change_triad(const Base_vect &new_triad)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Definition: tenseur.C:684
Tenseur gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition: etoile.h:474
Map & mp
Mapping associated with the star.
Definition: etoile.h:432
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:465
int get_etat() const
Returns the logical state.
Definition: tenseur.h:710
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
Definition: et_rot_hydro.C:86
Active physical coordinates and mapping derivatives.
Definition: coord.h:90
Tenseur uuu
Norm of u_euler.
Definition: etoile.h:1521
double omega
Rotation angular velocity ([f_unit] )
Definition: etoile.h:1504
Tenseur a_car
Total conformal factor .
Definition: etoile.h:518
Tenseur ener
Total energy density in the fluid frame.
Definition: etoile.h:463
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:803
Coord y
y coordinate centered on the grid
Definition: map.h:739
Coord x
x coordinate centered on the grid
Definition: map.h:738
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
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition: tenseur.h:304