LORENE
et_bin_nsbh_kinema.C
1 /*
2  * Copyright (c) 2000-2001 Eric Gourgoulhon
3  *
4  * This file is part of LORENE.
5  *
6  * LORENE is free software; you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License as published by
8  * the Free Software Foundation; either version 2 of the License, or
9  * (at your option) any later version.
10  *
11  * LORENE is distributed in the hope that it will be useful,
12  * but WITHOUT ANY WARRANTY; without even the implied warranty of
13  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  * GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License
17  * along with LORENE; if not, write to the Free Software
18  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19  *
20  */
21 
22 
23 char et_bin_kinema_nsbhC[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_nsbh_kinema.C,v 1.4 2014/10/13 08:52:56 j_novak Exp $" ;
24 
25 /*
26  * $Id: et_bin_nsbh_kinema.C,v 1.4 2014/10/13 08:52:56 j_novak Exp $
27  * $Log: et_bin_nsbh_kinema.C,v $
28  * Revision 1.4 2014/10/13 08:52:56 j_novak
29  * Lorene classes and functions now belong to the namespace Lorene.
30  *
31  * Revision 1.3 2008/08/19 06:42:00 j_novak
32  * Minor modifications to avoid warnings with gcc 4.3. Most of them concern
33  * cast-type operations, and constant strings that must be defined as const char*
34  *
35  * Revision 1.2 2005/10/18 13:12:33 p_grandclement
36  * update of the mixted binary codes
37  *
38  * Revision 1.1 2005/08/31 09:32:50 p_grandclement
39  * Forgot this file
40  *
41  *
42  * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_nsbh_kinema.C,v 1.4 2014/10/13 08:52:56 j_novak Exp $
43  *
44  */
45 
46 // Headers Lorene
47 #include "et_bin_nsbh.h"
48 #include "graphique.h"
49 
50 namespace Lorene {
51 void Et_bin_nsbh::kinematics(double omega, double) {
52 
53  int nz = mp.get_mg()->get_nzone() ;
54  int nzm1 = nz - 1 ;
55 
56  // --------------------
57  // Computation of B^i/N
58  // --------------------
59 
60  // 1/ Computation of - omega m^i
61 
62  const Coord& xa = mp.xa ;
63  const Coord& ya = mp.ya ;
64 
65  bsn.set_etat_qcq() ;
66 
67  bsn.set(0) = -omega * ya ;
68  bsn.set(1) = omega * xa ;
69  bsn.set(2) = 0 ;
70 
71  bsn.annule(nzm1, nzm1) ; // set to zero in the ZEC
72 
73  // 2/ Addition of shift and division by lapse
74  // See Eq (47) from Gourgoulhon et al. (2001)
75 
76  bsn = -( bsn + shift ) / nnn ;
77 
78  bsn.annule(nzm1, nzm1) ; // set to zero in the ZEC
79  bsn.set_std_base() ; // set the bases for spectral expansions
80 
81  //-------------------------
82  // Centrifugal potentatial
83  //-------------------------
84 
85  if (relativistic) {
86 
87  // Lorentz factor between the co-orbiting observer and the Eulerian one
88  // See Eq (23) from Gourgoulhon et al. (2001)
89 
90  Tenseur gam0 = 1 / sqrt( 1-sprod(bsn, bsn) ) ;
91  pot_centri = - log( gam0 ) ;
92  }
93  else {
94 
96 
97 
98  // See Eq (40) from Gourgoulhon et al. (2001)
99  pot_centri.set() = - 0.5 * omega * omega * (
100  xa *xa + ya * ya ) ;
101 
102  }
103 
104  pot_centri.annule(nzm1, nzm1) ; // set to zero in the external domain
105  pot_centri.set_std_base() ; // set the bases for spectral expansions
106 
107  // The derived quantities are obsolete
108  // -----------------------------------
109 
110  del_deriv() ;
111 
112 
113 }
114 }
Coord xa
Absolute x coordinate.
Definition: map.h:748
Cmp log(const Cmp &)
Neperian logarithm.
Definition: cmp_math.C:299
void annule(int l)
Sets the Tenseur to zero in a given domain.
Definition: tenseur.C:916
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
Tenseur pot_centri
Centrifugal potential.
Definition: etoile.h:956
Lorene prototypes.
Definition: app_hor.h:67
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:783
virtual void kinematics(double omega, double x_axe)
Computes the quantities bsn and pot_centri .
Tenseur nnn
Total lapse function.
Definition: etoile.h:512
Tenseur shift
Total shift vector.
Definition: etoile.h:515
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition: tenseur.C:840
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: etoile_bin.C:450
Tenseur bsn
3-vector shift, divided by N , of the rotating coordinates, .
Definition: etoile.h:953
Map & mp
Mapping associated with the star.
Definition: etoile.h:432
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:465
Active physical coordinates and mapping derivatives.
Definition: coord.h:90
Coord ya
Absolute y coordinate.
Definition: map.h:749
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one. ...
Definition: etoile.h:440
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tenseur.C:652
virtual Tenseur sprod(const Tenseur &t1, const Tenseur &t2) const
Performs the scalar product of two tensors by contracting the last index of t1 with the first index o...
Definition: etoile_bin.C:751
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition: tenseur.h:304