LORENE
star_rot_hydro.C
1/*
2 * Method Star_rot::hydro_euler
3 *
4 * (see file star_rot.h for documentation)
5 *
6 */
7
8/*
9 * Copyright (c) 2010 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: star_rot_hydro.C,v 1.4 2016/12/05 16:18:15 j_novak Exp $
34 * $Log: star_rot_hydro.C,v $
35 * Revision 1.4 2016/12/05 16:18:15 j_novak
36 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
37 *
38 * Revision 1.3 2014/10/13 08:53:39 j_novak
39 * Lorene classes and functions now belong to the namespace Lorene.
40 *
41 * Revision 1.2 2014/10/06 15:13:17 j_novak
42 * Modified #include directives to use c++ syntax.
43 *
44 * Revision 1.1 2010/01/25 18:15:52 e_gourgoulhon
45 * First version.
46 *
47 *
48 * $Header: /cvsroot/Lorene/C++/Source/Star/star_rot_hydro.C,v 1.4 2016/12/05 16:18:15 j_novak Exp $
49 *
50 */
51
52// Headers C
53#include <cstdlib>
54
55// Headers Lorene
56#include "star_rot.h"
57#include "utilitaires.h"
58
59namespace Lorene {
61
62 int nz = mp.get_mg()->get_nzone() ;
63 int nzm1 = nz - 1 ;
64
65 // Computation of u_euler
66 // ----------------------
67
68 const Coord& x = mp.x ;
69 const Coord& y = mp.y ;
70
71 u_euler.set_etat_qcq() ;
72
73 u_euler.set(1) = - omega * y ; // Cartesian components of solid rotation
74 u_euler.set(2) = omega * x ;
75 u_euler.set(3) = 0 ;
76 u_euler.annule_domain(nzm1) ;
77
78 u_euler.set_triad( mp.get_bvect_cart() ) ; // Triad = Cartesian triad
79
80 u_euler = ( u_euler + beta ) / nn ;
81
82 u_euler.std_spectral_base() ; // sets the standard bases for spectral expansions
83
84 if ( (u_euler(1).get_etat() == ETATZERO) &&
85 (u_euler(2).get_etat() == ETATZERO) &&
86 (u_euler(3).get_etat() == ETATZERO) ) {
87
88 u_euler.set_etat_zero() ;
89 }
90
91
92 // Computation of uuu (norme of u_euler)
93 // ------------------
94
95 // The scalar product is performed on the spherical components:
96
97 Vector us = u_euler ;
98 us.change_triad( mp.get_bvect_spher() ) ;
99
100 Scalar uuu2 = a_car * ( us(1)*us(1) + us(2)*us(2) ) + b_car * us(3)*us(3) ;
101
102 uuu = sqrt( uuu2 ) ;
103
104 if (uuu.get_etat() == ETATQCQ) {
105 // Same basis as (Omega -N^phi) r sin(theta)
106 (uuu.set_spectral_va()).set_base( us(3).get_spectral_va().get_base() ) ;
107 }
108
109 // Lorentz factor
110 // --------------
111
112 Scalar u2 = unsurc2 * uuu2 ;
113
114 Scalar gam2 = double(1) / (double(1) - u2) ;
115
116 gam_euler = sqrt(gam2) ;
117
118 gam_euler.std_spectral_base() ; // sets the standard spectral bases for
119 // a scalar field
120
121 // Energy density E with respect to the Eulerian observer
122 //--------------------------------------------------------
123
124 ener_euler = gam2 * ( ener + press ) - press ;
125
126 // Trace of the stress tensor with respect to the Eulerian observer
127 //------------------------------------
128
129 s_euler = 3 * press + ( ener_euler + press ) * u2 ;
130
131 // The derived quantities are obsolete
132 // -----------------------------------
133
134 del_deriv() ;
135
136}
137}
Active physical coordinates and mapping derivatives.
Definition coord.h:90
Tensor field of valence 0 (or component of a tensorial field).
Definition scalar.h:393
double unsurc2
: unsurc2=1 for a relativistic star, 0 for a Newtonian one.
Definition star_rot.h:111
Scalar b_car
Square of the metric factor B.
Definition star_rot.h:122
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
double omega
Rotation angular velocity ([f_unit] ).
Definition star_rot.h:113
Scalar uuu
Norm of u_euler.
Definition star_rot.h:133
Scalar a_car
Square of the metric factor A.
Definition star_rot.h:116
virtual void del_deriv() const
Deletes all the derived quantities.
Definition star_rot.C:310
Scalar ener
Total energy density in the fluid frame.
Definition star.h:193
Scalar nn
Lapse function N .
Definition star.h:225
Scalar ener_euler
Total energy density in the Eulerian frame.
Definition star.h:198
Scalar gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition star.h:204
Scalar s_euler
Trace of the stress scalar in the Eulerian frame.
Definition star.h:201
Scalar press
Fluid pressure.
Definition star.h:194
Vector u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition star.h:207
Map & mp
Mapping associated with the star.
Definition star.h:180
Vector beta
Shift vector.
Definition star.h:228
Tensor field of valence 1.
Definition vector.h:188
virtual void change_triad(const Base_vect &)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Cmp sqrt(const Cmp &)
Square root.
Definition cmp_math.C:223
Lorene prototypes.
Definition app_hor.h:67
Coord y
y coordinate centered on the grid
Definition map.h:739
Coord x
x coordinate centered on the grid
Definition map.h:738