LORENE
star_bhns_vel_pot.C
1 /*
2  * Method of class Star_bhns to compute the velocity scalar potential $\psi$
3  * by solving the continuity equation.
4  *
5  * (see file star_bhns.h for documentation).
6  *
7  */
8 
9 /*
10  * Copyright (c) 2005-2007 Keisuke Taniguchi
11  *
12  * This file is part of LORENE.
13  *
14  * LORENE is free software; you can redistribute it and/or modify
15  * it under the terms of the GNU General Public License version 2
16  * as published by the Free Software Foundation.
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  * $Id: star_bhns_vel_pot.C,v 1.5 2018/11/16 14:34:37 j_novak Exp $
33  * $Log: star_bhns_vel_pot.C,v $
34  * Revision 1.5 2018/11/16 14:34:37 j_novak
35  * Changed minor points to avoid some compilation warnings.
36  *
37  * Revision 1.4 2016/12/05 16:18:16 j_novak
38  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
39  *
40  * Revision 1.3 2014/10/13 08:53:41 j_novak
41  * Lorene classes and functions now belong to the namespace Lorene.
42  *
43  * Revision 1.2 2008/05/15 19:20:29 k_taniguchi
44  * Change of some parameters.
45  *
46  * Revision 1.1 2007/06/22 01:33:14 k_taniguchi
47  * *** empty log message ***
48  *
49  *
50  * $Header: /cvsroot/Lorene/C++/Source/Star_bhns/star_bhns_vel_pot.C,v 1.5 2018/11/16 14:34:37 j_novak Exp $
51  *
52  */
53 
54 // C++ headers
55 //#include <>
56 
57 // C headers
58 //#include <>
59 
60 // Lorene headers
61 #include "star_bhns.h"
62 #include "eos.h"
63 #include "param.h"
64 #include "cmp.h"
65 #include "tenseur.h"
66 #include "utilitaires.h"
67 #include "unites.h"
68 
69 // Local prototype
70 namespace Lorene {
71 Cmp raccord_c1(const Cmp& uu, int l1) ;
72 
73 double Star_bhns::velo_pot_bhns(const double&, const double&,
74  bool kerrschild, int mermax, double precis,
75  double relax) {
76 
77  // Fundamental constants and units
78  // -------------------------------
79  using namespace Unites ;
80 
81  int nzm1 = mp.get_mg()->get_nzone() - 1 ;
82 
83  //--------------------------------
84  // Specific relativistic enthalpy ---> hhh
85  //--------------------------------
86 
87  Scalar hhh = exp(ent) ;
88  hhh.std_spectral_base() ;
89 
90  //---------------------------------------------------
91  // Computation of W^i = psi4 h Gamma_n B^i/lapse_tot
92  // See Eq (62) from Gourgoulhon et al. (2001)
93  //---------------------------------------------------
94 
95  Vector www = hhh * gam_euler * bsn * psi4 ;
96 
97  www.change_triad( mp.get_bvect_cart() ) ; // components on the mapping
98  // Cartesian basis
99 
100  //-------------------------------------------------
101  // Constant value of W^i at the center of the star
102  //-------------------------------------------------
103 
104  Vector v_orb(mp, CON, mp.get_bvect_cart()) ;
105  v_orb.set_etat_qcq() ;
106 
107  for (int i=1; i<=3; i++) {
108  v_orb.set(i) = www(i).val_grid_point(0, 0, 0, 0) ;
109  }
110 
111  v_orb.set_triad( *(www.get_triad()) ) ;
112  v_orb.std_spectral_base() ;
113 
114  //-------------------------------------------------
115  // Source and coefficients a,b for poisson_compact (independent from psi0)
116  //-------------------------------------------------
117 
118  Scalar dndh_log = eos.der_nbar_ent(ent, nzet) ;
119 
120  // In order to avoid any division by zero in the computation of zeta_h
121  // the value of dndh_log is set to 1 in the external domains:
122  for (int l=nzet; l<=nzm1; l++) {
123  dndh_log.set_domain(l) = 1. ;
124  }
125 
126  Scalar zeta_h( ent / dndh_log ) ;
127  zeta_h.std_spectral_base() ;
128 
129  double erreur ;
130 
131  Scalar source(mp) ;
132  Vector bb(mp, CON, mp.get_bvect_spher()) ;
133  Metric_flat flat_spher( mp.flat_met_spher() ) ;
134 
135  if (kerrschild) {
136 
137  cout << "!!!!! WARNING: Not yet available !!!!!" << endl ;
138  abort() ;
139 
140  } // End of Kerr-Schild
141  else { // Isotropic coordinates with the maximal slicing
142 
143  Scalar lnlappsi(mp) ;
144  lnlappsi = log( lapconf_tot * confo_tot ) ;
145  lnlappsi.std_spectral_base() ;
146 
147  bb = (1. - zeta_h) * ent.derive_con(flat_spher)
148  + zeta_h * lnlappsi.derive_con(flat_spher) ;
149 
150  Vector dentmb(mp, COV, mp.get_bvect_cart()) ;
151  dentmb.set_etat_qcq() ;
152  for (int i=1; i<=3; i++) {
153  dentmb.set(i) = ent.deriv(i)
155  - (d_confo_auto(i) + d_confo_comp(i)) / confo_tot ;
156  }
157  dentmb.std_spectral_base() ;
158 
159  source = contract(www - v_orb, 0, ent.derive_cov(flat), 0)
160  +zeta_h*(contract(v_orb, 0, dentmb, 0)
162  ) ;
163 
164  source.annule(nzet, nzm1) ;
165 
166  } // End of Isotropic coordinates
167 
168 
169  //----------------------------------------------------
170  // Resolution by means of Map_radial::poisson_compact
171  //----------------------------------------------------
172 
173  Param par ;
174  int niter ;
175  par.add_int(mermax) ;
176  par.add_double(precis, 0) ;
177  par.add_double(relax, 1) ;
178  par.add_int_mod(niter) ;
179 
180  if (psi0.get_etat() == ETATZERO) {
181  psi0 = 0. ;
182  }
183 
184  Cmp source_cmp(source) ;
185  Cmp zeta_h_cmp(zeta_h) ;
186  Cmp psi0_cmp(psi0) ;
187 
188  Tenseur bb_cmp(mp, 1, CON, mp.get_bvect_spher()) ;
189  bb_cmp.set_etat_qcq() ;
190  Cmp bb_cmp1(bb(1)) ;
191  Cmp bb_cmp2(bb(2)) ;
192  Cmp bb_cmp3(bb(3)) ;
193  bb_cmp.set(0) = bb_cmp1 ;
194  bb_cmp.set(1) = bb_cmp2 ;
195  bb_cmp.set(2) = bb_cmp3 ;
196 
197  source_cmp.va.ylm() ;
198 
199  mp.poisson_compact(source_cmp, zeta_h_cmp, bb_cmp, par, psi0_cmp) ;
200 
201  psi0 = psi0_cmp ;
202 
203  //-----------------------
204  // Check of the solution
205  //-----------------------
206 
207  Scalar bb_dpsi0 = contract(bb, 0, psi0.derive_cov(flat_spher), 0) ;
208 
209  Scalar oper = zeta_h * psi0.laplacian() + bb_dpsi0 ;
210 
211  source.set_spectral_va().ylm_i() ;
212 
213  erreur = diffrel(oper, source)(0) ;
214 
215  cout << "Check of the resolution of the continuity equation : "
216  << endl ;
217  cout << " norme(source) : " << norme(source)(0)
218  << " diff oper/source : " << erreur << endl ;
219 
220  //--------------------------
221  // Computation of grad(psi)
222  //--------------------------
223 
224  for (int i=1; i<=3; i++)
225  d_psi.set(i) = (psi0.derive_cov(flat))(i) + v_orb(i) ;
226 
227 
228  // C^1 continuation of d_psi outside the star
229  // (to ensure a smooth enthalpy field accross the stellar surface)
230  // ----------------------------------------------------------------
231 
232  d_psi.annule(nzet, nzm1) ;
233 
234  for (int i=1; i<=3; i++) {
235  Cmp d_psi_i( d_psi(i) ) ;
236  d_psi_i = raccord_c1(d_psi_i, nzet) ;
237  d_psi.set(i) = d_psi_i ;
238  }
239 
240  return erreur ;
241 
242 }
243 }
Scalar psi4
Fourth power of the total conformal factor.
Definition: star_bhns.h:176
Vector d_lapconf_auto
Derivative of the lapconf function generated by the star .
Definition: star_bhns.h:130
Scalar lapconf_tot
Total lapconf function.
Definition: star_bhns.h:119
Cmp log(const Cmp &)
Neperian logarithm.
Definition: cmp_math.C:299
virtual void set_etat_qcq()
Sets the logical state of all components to ETATQCQ (ordinary state).
Definition: tensor.C:490
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition: cmp.h:446
Cmp exp(const Cmp &)
Exponential.
Definition: cmp_math.C:273
void add_int(const int &n, int position=0)
Adds the address of a new int to the list.
Definition: param.C:249
Map & mp
Mapping associated with the star.
Definition: star.h:180
void ylm_i()
Inverse of ylm()
Definition: valeur_ylm_i.C:134
Cmp der_nbar_ent(const Cmp &ent, int nzet, int l_min=0, Param *par=0x0) const
Computes the logarithmic derivative from the log-enthalpy and extra parameters.
Definition: eos.C:431
virtual void annule(int l_min, int l_max)
Sets the Scalar to zero in several domains.
Definition: scalar.C:397
Vector bsn
3-vector shift, divided by N , of the rotating coordinates, .
Definition: star_bhns.h:99
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 Eos & eos
Equation of state of the stellar matter.
Definition: star.h:185
Standard units of space, time and mass.
void ylm()
Computes the coefficients of *this.
Definition: valeur_ylm.C:141
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:777
const Scalar & laplacian(int ced_mult_r=4) const
Returns the Laplacian of *this.
Definition: scalar_deriv.C:436
Flat metric for tensor calculation.
Definition: metric.h:261
Vector d_lapconf_comp
Derivative of the lapconf function generated by the companion black hole.
Definition: star_bhns.h:135
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:393
virtual void poisson_compact(const Cmp &source, const Cmp &aa, const Tenseur &bb, const Param &par, Cmp &psi) const =0
Resolution of the elliptic equation in the case where the stellar interior is covered by a single do...
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.
int get_etat() const
Returns the logical state ETATNONDEF (undefined), ETATZERO (null) or ETATQCQ (ordinary).
Definition: scalar.h:560
Scalar ent
Log-enthalpy.
Definition: star.h:190
const Vector & derive_con(const Metric &gam) const
Returns the "contravariant" derivative of *this with respect to some metric , by raising the index of...
Definition: scalar_deriv.C:402
Tensor field of valence 1.
Definition: vector.h:188
Tbl & set_domain(int l)
Read/write of the value in a given domain.
Definition: scalar.h:621
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition: cmp_math.C:507
Vector d_psi
Gradient of (in the irrotational case) (Spherical components with respect to the mapping of the star...
Definition: star_bhns.h:82
int nzet
Number of domains of *mp occupied by the star.
Definition: star.h:183
Metric_flat flat
Flat metric defined on the mapping (Spherical components with respect to the mapping of the star )...
Definition: star_bhns.h:192
Scalar confo_tot
Total conformal factor.
Definition: star_bhns.h:163
Scalar gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition: star.h:204
Tbl norme(const Cmp &)
Sums of the absolute values of all the values of the Cmp in each domain.
Definition: cmp_math.C:484
const Base_vect * get_triad() const
Returns the vectorial basis (triad) on which the components are defined.
Definition: tensor.h:879
Vector d_confo_comp
Derivative of the conformal factor generated by the companion black hole.
Definition: star_bhns.h:173
Parameter storage.
Definition: param.h:125
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:465
const Scalar & deriv(int i) const
Returns of *this , where .
Definition: scalar_deriv.C:359
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
Tbl & set(int l)
Read/write of the value in a given domain.
Definition: cmp.h:724
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
Vector d_confo_auto
Derivative of the conformal factor generated by the star .
Definition: star_bhns.h:168
void add_double(const double &x, int position=0)
Adds the the address of a new double to the list.
Definition: param.C:318
void set_triad(const Base_vect &new_triad)
Assigns a new vectorial basis (triad) of decomposition.
Definition: tensor.C:528
Valeur & set_spectral_va()
Returns va (read/write version)
Definition: scalar.h:610
Scalar & set(int)
Read/write access to a component.
Definition: vector.C:302
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tenseur.C:652
const Vector & derive_cov(const Metric &gam) const
Returns the gradient (1-form = covariant vector) of *this
Definition: scalar_deriv.C:390
virtual void annule(int l_min, int l_max)
Sets the Tensor to zero in several domains.
Definition: tensor.C:680
Scalar psi0
Scalar potential of the non-translational part of fluid 4-velocity (in the irrotational case) ...
Definition: star_bhns.h:77
Valeur va
The numerical value of the Cmp.
Definition: cmp.h:464
const Metric_flat & flat_met_spher() const
Returns the flat metric associated with the spherical coordinates and with components expressed in th...
Definition: map.C:324
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition: tenseur.h:304
void add_int_mod(int &n, int position=0)
Adds the address of a new modifiable int to the list.
Definition: param.C:388
double velo_pot_bhns(const double &mass_bh, const double &sepa, bool kerrschild, int mermax, double precis, double relax)
Computes the non-translational part of the velocity scalar potential by solving the continuity equat...