LORENE
star_bin_kinema_xcts.C
1 /*
2  * Method Star_bin_xcts::kinematics
3  * (see file star.h for documentation)
4  */
5 
6 /*
7  * Copyright (c) 2010 Michal Bejger
8  *
9  * This file is part of LORENE.
10  *
11  * LORENE is free software; you can redistribute it and/or modify
12  * it under the terms of the GNU General Public License version 2
13  * as published by the Free Software Foundation.
14  *
15  * LORENE is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18  * GNU General Public License for more details.
19  *
20  * You should have received a copy of the GNU General Public License
21  * along with LORENE; if not, write to the Free Software
22  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
23  *
24  */
25 
26 
27 
28 /*
29  * $Id: star_bin_kinema_xcts.C,v 1.7 2016/12/05 16:18:15 j_novak Exp $
30  * $Log: star_bin_kinema_xcts.C,v $
31  * Revision 1.7 2016/12/05 16:18:15 j_novak
32  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
33  *
34  * Revision 1.6 2014/10/13 08:53:38 j_novak
35  * Lorene classes and functions now belong to the namespace Lorene.
36  *
37  * Revision 1.5 2014/10/06 15:13:17 j_novak
38  * Modified #include directives to use c++ syntax.
39  *
40  * Revision 1.4 2010/12/09 10:43:53 m_bejger
41  * Small changes, annule --> annule_domain
42  *
43  * Revision 1.3 2010/10/26 19:57:02 m_bejger
44  * Various cleanups
45  *
46  * Revision 1.2 2010/06/04 19:57:52 m_bejger
47  * Added condition for mp.get_rot_phi
48  *
49  * Revision 1.1 2010/05/04 07:51:05 m_bejger
50  * Initial version
51  *
52  * $Header: /cvsroot/Lorene/C++/Source/Star/star_bin_kinema_xcts.C,v 1.7 2016/12/05 16:18:15 j_novak Exp $
53  *
54  */
55 
56 #include <cmath>
57 
58 // Headers Lorene
59 #include "star.h"
60 
61 namespace Lorene {
62 void Star_bin_xcts::kinematics(double omega, double x_axe) {
63 
64  int nzm1 = mp.get_mg()->get_nzone() - 1 ;
65 
66  // --------------------
67  // Computation of B^i/N
68  // --------------------
69 
70  // 1/ Computation of omega m^i
71 
72  const Coord& xa = mp.xa ;
73  const Coord& ya = mp.ya ;
74 
76 
77  if (fabs(mp.get_rot_phi()) < 1e-10) {
78 
79  bsn.set(1) = - omega * ya ;
80  bsn.set(2) = omega * (xa - x_axe) ;
81  bsn.set(3) = 0 ;
82 
83  } else {
84 
85  bsn.set(1) = omega * ya ;
86  bsn.set(2) = - omega * (xa - x_axe) ;
87  bsn.set(3) = 0 ;
88 
89  }
90 
91  bsn.annule_domain(nzm1) ; // set to zero in the ZEC
92 
93  // Addition of beta and division by lapse
94  // Eq. 47 from Gourgoulhon et al. (2001)
95  // New convention : l = Nn + B ==> B = \beta + \Omega d\phi
96  //---------------------------------------------------------
97 
98  bsn = ( bsn + beta ) / nn ;
99 
100  bsn.annule_domain(nzm1) ; // set to zero in the ZEC
102 
103  //----------------------
104  // Centrifugal potential
105  //----------------------
106 
107  // Lorentz factor between the co-orbiting observer
108  // and the Eulerian one
109  // Eq. 23 from Gourgoulhon et al. (2001)
110 
111  Sym_tensor gamma_cov (gamma.cov()) ;
112  gamma_cov.change_triad(mp.get_bvect_cart()) ;
113 
114  Scalar gam0 = 1 / sqrt(1 - contract(gamma_cov, 0, 1, bsn * bsn, 0, 1)) ;
115  pot_centri = - log( gam0 ) ;
116 
117  pot_centri.annule_domain(nzm1) ; // set to zero in the external domain
118  pot_centri.std_spectral_base() ; // set the bases for spectral expansions
119 
120  // The derived quantities are obsolete
121  // -----------------------------------
122 
123  del_deriv() ;
124 
125 }
126 }
Coord xa
Absolute x coordinate.
Definition: map.h:742
void annule_domain(int l)
Sets the Tensor to zero in a given domain.
Definition: tensor.C:675
Vector bsn
3-vector shift, divided by N, of the rotating coordinates, .
Definition: star.h:1126
Cmp log(const Cmp &)
Neperian logarithm.
Definition: cmp_math.C:299
Map & mp
Mapping associated with the star.
Definition: star.h:180
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
Metric gamma
3-metric
Definition: star.h:235
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
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:393
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
virtual void change_triad(const Base_vect &)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
double get_rot_phi() const
Returns the angle between the x –axis and X –axis.
Definition: map.h:787
Vector beta
Shift vector.
Definition: star.h:228
virtual void std_spectral_base()
Sets the standard spectal bases of decomposition for each component.
Definition: vector.C:322
Scalar pot_centri
Centrifugal potential.
Definition: star.h:1129
void kinematics(double omega, double x_axe)
Computes the quantities bsn and pot_centri.
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:465
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
Active physical coordinates and mapping derivatives.
Definition: coord.h:90
virtual const Sym_tensor & cov() const
Read-only access to the covariant representation.
Definition: metric.C:283
Coord ya
Absolute y coordinate.
Definition: map.h:743
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
Scalar nn
Lapse function N .
Definition: star.h:225
virtual void del_deriv() const
Deletes all the derived quantities.
virtual void change_triad(const Base_vect &new_triad)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Scalar & set(int)
Read/write access to a component.
Definition: vector.C:302
Class intended to describe valence-2 symmetric tensors.
Definition: sym_tensor.h:226