LORENE
bhole_equations_bin.C
1/*
2 * Copyright (c) 2001 Philippe Grandclement
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
24
25/*
26 * $Id: bhole_equations_bin.C,v 1.6 2016/12/05 16:17:45 j_novak Exp $
27 * $Log: bhole_equations_bin.C,v $
28 * Revision 1.6 2016/12/05 16:17:45 j_novak
29 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
30 *
31 * Revision 1.5 2014/10/13 08:52:40 j_novak
32 * Lorene classes and functions now belong to the namespace Lorene.
33 *
34 * Revision 1.4 2014/10/06 15:12:58 j_novak
35 * Modified #include directives to use c++ syntax.
36 *
37 * Revision 1.3 2006/04/27 09:12:32 p_grandclement
38 * First try at irrotational black holes
39 *
40 * Revision 1.2 2002/10/16 14:36:33 j_novak
41 * Reorganization of #include instructions of standard C++, in order to
42 * use experimental version 3 of gcc.
43 *
44 * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
45 * LORENE
46 *
47 * Revision 2.6 2001/05/07 09:11:56 phil
48 * *** empty log message ***
49 *
50 * Revision 2.5 2001/04/30 09:30:50 phil
51 * filtre pour poisson vectoriel
52 *
53 * Revision 2.4 2001/04/26 12:04:34 phil
54 * *** empty log message ***
55 *
56 * Revision 2.3 2001/04/25 15:08:48 phil
57 * vire fait_tkij
58 *
59 * Revision 2.2 2001/04/02 12:16:11 phil
60 * *** empty log message ***
61 *
62 * Revision 2.1 2001/03/22 10:41:36 phil
63 * modification prototypage
64 *
65 * Revision 2.0 2001/03/01 08:18:04 phil
66 * *** empty log message ***
67 *
68 *
69 * $Header: /cvsroot/Lorene/C++/Source/Bhole_binaire/bhole_equations_bin.C,v 1.6 2016/12/05 16:17:45 j_novak Exp $
70 *
71 */
72
73//standard
74#include <cstdlib>
75#include <cmath>
76
77// Lorene
78#include "nbr_spx.h"
79#include "tenseur.h"
80#include "bhole.h"
81#include "proto.h"
82#include "utilitaires.h"
83#include "graphique.h"
84
85// Resolution pour le lapse
86namespace Lorene {
87void Bhole_binaire::solve_lapse (double precision, double relax) {
88
89 assert ((relax >0) && (relax<=1)) ;
90
91 cout << "-----------------------------------------------" << endl ;
92 cout << "Resolution LAPSE" << endl ;
93
94 Tenseur lapse_un_old (hole1.n_auto) ;
95 Tenseur lapse_deux_old (hole2.n_auto) ;
96
97 Tenseur auxi_un (flat_scalar_prod(hole1.tkij_tot, hole1.tkij_auto)) ;
98 Cmp kk_un (auxi_un(0, 0)+auxi_un(1, 1)+auxi_un(2, 2)) ;
99
100 Tenseur auxi_deux (flat_scalar_prod(hole2.tkij_tot, hole2.tkij_auto)) ;
101 Cmp kk_deux (auxi_deux(0, 0)+auxi_deux(1, 1)+auxi_deux(2, 2)) ;
102
103 // Les sources
104
105 Cmp source_un
106(-2*flat_scalar_prod(hole1.grad_n_tot, hole1.psi_auto.gradient())()/hole1.psi_tot()
107 +hole1.n_tot()*pow(hole1.psi_tot(), 4.)*kk_un) ;
108 source_un.std_base_scal() ;
109
110 Cmp source_deux
111(-2*flat_scalar_prod(hole2.grad_n_tot, hole2.psi_auto.gradient())()/hole2.psi_tot()
112 +hole2.n_tot()*pow(hole2.psi_tot(), 4.)*kk_deux) ;
113 source_deux.std_base_scal() ;
114
115 //On resout
116 hole1.n_auto.set() = hole1.n_auto() - 1./2. ;
117 hole2.n_auto.set() = hole2.n_auto() - 1./2. ;
118
119 dirichlet_binaire (source_un, source_deux, -1., -1., hole1.n_auto.set(),
120 hole2.n_auto.set(), 0, precision) ;
121
122 hole1.n_auto.set() = hole1.n_auto() + 1./2. ;
123 hole2.n_auto.set() = hole2.n_auto() + 1./2. ;
124
125 hole1.n_auto.set().raccord(1) ;
126 hole2.n_auto.set().raccord(1) ;
127
128 // La relaxation :
129 hole1.n_auto.set() = relax*hole1.n_auto() + (1-relax)*lapse_un_old() ;
130 hole2.n_auto.set() = relax*hole2.n_auto() + (1-relax)*lapse_deux_old() ;
131
132 hole1.fait_n_comp (hole2) ;
133 hole2.fait_n_comp (hole1) ;
134}
135
136//Resolution sur Psi
137void Bhole_binaire::solve_psi (double precision, double relax) {
138
139 assert ((relax>0) && (relax<=1)) ;
140
141 cout << "-----------------------------------------------" << endl ;
142 cout << "Resolution PSI" << endl ;
143
144 Tenseur psi_un_old (hole1.psi_auto) ;
145 Tenseur psi_deux_old (hole2.psi_auto) ;
146
147 Tenseur auxi_un (flat_scalar_prod(hole1.tkij_tot, hole1.tkij_auto)) ;
148 Cmp kk_un (auxi_un(0, 0)+auxi_un(1, 1)+auxi_un(2, 2)) ;
149
150 Tenseur auxi_deux (flat_scalar_prod(hole2.tkij_tot, hole2.tkij_auto)) ;
151 Cmp kk_deux (auxi_deux(0, 0)+auxi_deux(1, 1)+auxi_deux(2, 2)) ;
152
153 // Les sources
154 Cmp source_un (-pow(hole1.psi_tot(), 5.)*kk_un/8.) ;
155 source_un.std_base_scal() ;
156
157 Cmp source_deux (-pow(hole2.psi_tot(), 5.)*kk_deux/8.) ;
158 source_deux.std_base_scal() ;
159
160 // Les valeurs limites :
161 int np_un = hole1.mp.get_mg()->get_np(1) ;
162 int nt_un = hole1.mp.get_mg()->get_nt(1) ;
163 Valeur lim_un (hole1.mp.get_mg()->get_angu()) ;
164 lim_un = 1 ;
165 for (int k=0 ; k<np_un ; k++)
166 for (int j=0 ; j<nt_un ; j++)
167 lim_un.set(0, k, j, 0) = -0.5/hole1.rayon*hole1.psi_tot()(1, k, j, 0) ;
168 lim_un.std_base_scal() ;
169
170 int np_deux = hole2.mp.get_mg()->get_np(1) ;
171 int nt_deux = hole2.mp.get_mg()->get_nt(1) ;
172 Valeur lim_deux (hole2.mp.get_mg()->get_angu()) ;
173 lim_deux = 1 ;
174 for (int k=0 ; k<np_deux ; k++)
175 for (int j=0 ; j<nt_deux ; j++)
176 lim_deux.set(0, k, j, 0) = -0.5/hole2.rayon*hole2.psi_tot()(1, k, j, 0) ;
177 lim_deux.std_base_scal() ;
178
179 //On resout
180 neumann_binaire (source_un, source_deux, lim_un, lim_deux,
181 hole1.psi_auto.set(), hole2.psi_auto.set(), 0, precision) ;
182
183 hole1.psi_auto.set() = hole1.psi_auto() + 1./2. ;
184 hole2.psi_auto.set() = hole2.psi_auto() + 1./2. ;
185
186 hole1.psi_auto.set().raccord(1) ;
187 hole2.psi_auto.set().raccord(1) ;
188
189 // La relaxation :
190 hole1.psi_auto.set() = relax*hole1.psi_auto() + (1-relax)*psi_un_old() ;
191 hole2.psi_auto.set() = relax*hole2.psi_auto() + (1-relax)*psi_deux_old() ;
192
193 hole1.fait_psi_comp (hole2) ;
194 hole2.fait_psi_comp (hole1) ;
195}
196
197
198// Resolution sur shift a omega bloque.
199void Bhole_binaire::solve_shift (double precision, double relax) {
200
201 cout << "------------------------------------------------" << endl ;
202 cout << "Resolution shift : Omega = " << omega << endl ;
203
204 // On determine les sources
205 Tenseur source_un (
206 -6*flat_scalar_prod (hole1.taij_tot, hole1.psi_auto.gradient())/hole1.psi_tot
207 +2*flat_scalar_prod (hole1.tkij_tot, hole1.n_auto.gradient())) ;
208 if (source_un.get_etat() == ETATZERO) {
209 source_un.set_etat_qcq() ;
210 for (int i=0 ; i<3 ; i++)
211 source_un.set(i).set_etat_zero() ;
212 }
213 source_un.set_std_base() ;
214
215 Tenseur source_deux (
216 -6*flat_scalar_prod (hole2.taij_tot, hole2.psi_auto.gradient())/hole2.psi_tot
217 +2*flat_scalar_prod (hole2.tkij_tot, hole2.n_auto.gradient())) ;
218 if (source_deux.get_etat() == ETATZERO) {
219 source_deux.set_etat_qcq() ;
220 for (int i=0 ; i<3 ; i++)
221 source_deux.set(i).set_etat_zero() ;
222 }
223 source_deux.set_std_base() ;
224
225 // On filtre les hautes frequences.
226 for (int i=0 ; i<3 ; i++) {
227 if (source_un(i).get_etat() != ETATZERO)
228 source_un.set(i).filtre(4) ;
229 if (source_deux(i).get_etat() != ETATZERO)
230 source_deux.set(i).filtre(4) ;
231 }
232
233 // Les alignemenents pour le signe des CL.
234 double orientation_un = hole1.mp.get_rot_phi() ;
235 assert ((orientation_un==0) || (orientation_un == M_PI)) ;
236
237 double orientation_deux = hole2.mp.get_rot_phi() ;
238 assert ((orientation_deux==0) || (orientation_deux == M_PI)) ;
239
240 int aligne_un = (orientation_un == 0) ? 1 : -1 ;
241 int aligne_deux = (orientation_deux == 0) ? 1 : -1 ;
242
243 // On regarde si toutes les composantes sont nulles :
244 int ind1 = 0 ;
245 int ind2 = 0 ;
246 for (int i=0 ; i<3 ; i++) {
247 if (source_un(i).get_etat() == ETATQCQ)
248 ind1 = 1 ;
249 if (source_deux(i).get_etat() == ETATQCQ)
250 ind2 = 1 ;
251 }
252
253 if (ind1==0)
254 source_un.set_etat_zero() ;
255 if (ind2==0)
256 source_deux.set_etat_zero() ;
257
258 // On determine les Cl en fonction de omega :
259 int np_un = hole1.mp.get_mg()->get_np (1) ;
260 int nt_un = hole1.mp.get_mg()->get_nt (1) ;
261
262 Mtbl xa_mtbl_un (source_un.get_mp()->get_mg()) ;
263 xa_mtbl_un.set_etat_qcq() ;
264 Mtbl ya_mtbl_un (source_un.get_mp()->get_mg()) ;
265 ya_mtbl_un.set_etat_qcq() ;
266
267 xa_mtbl_un = source_un.get_mp()->xa ;
268 ya_mtbl_un = source_un.get_mp()->ya ;
269
270 Mtbl x_mtbl_un (source_un.get_mp()->get_mg()) ;
271 x_mtbl_un.set_etat_qcq() ;
272 Mtbl y_mtbl_un (source_un.get_mp()->get_mg()) ;
273 y_mtbl_un.set_etat_qcq() ;
274
275 x_mtbl_un = source_un.get_mp()->x ;
276 y_mtbl_un = source_un.get_mp()->y ;
277
278 // Les bases
279 Base_val** bases_un = hole1.mp.get_mg()->std_base_vect_cart() ;
280 Base_val** bases_deux = hole2.mp.get_mg()->std_base_vect_cart() ;
281
282 Valeur lim_x_un (*hole1.mp.get_mg()->get_angu()) ;
283 lim_x_un = 1 ; // Juste pour affecter dans espace des configs ;
284 lim_x_un.set_etat_c_qcq() ;
285 for (int k=0 ; k<np_un ; k++)
286 for (int j=0 ; j<nt_un ; j++)
287 lim_x_un.set(0, k, j, 0) = aligne_un*omega*ya_mtbl_un(0, 0, 0, 0) + aligne_un*hole1.omega_local*y_mtbl_un(1,k,j,0) ;
288 lim_x_un.base = *bases_un[0] ;
289
290 Valeur lim_y_un (*hole1.mp.get_mg()->get_angu()) ;
291 lim_y_un = 1 ; // Juste pour affecter dans espace des configs ;
292 lim_y_un.set_etat_c_qcq() ;
293 for (int k=0 ; k<np_un ; k++)
294 for (int j=0 ; j<nt_un ; j++)
295 lim_y_un.set(0, k, j, 0) = -aligne_un*omega*xa_mtbl_un(0, 0, 0, 0) - aligne_un*hole1.omega_local*x_mtbl_un(1,k,j,0) ;;
296 lim_y_un.base = *bases_un[1] ;
297
298 Valeur lim_z_un (*hole1.mp.get_mg()->get_angu()) ;
299 lim_z_un = 1 ;
300 for (int k=0 ; k<np_un ; k++)
301 for (int j=0 ; j<nt_un ; j++)
302 lim_z_un.set(0, k, j, 0) = 0 ;
303 lim_z_un.base = *bases_un[2] ;
304
305 // On determine les Cl en fonction de omega :
306 int np_deux = hole2.mp.get_mg()->get_np (1) ;
307 int nt_deux = hole2.mp.get_mg()->get_nt (1) ;
308
309 Mtbl xa_mtbl_deux (source_deux.get_mp()->get_mg()) ;
310 xa_mtbl_deux.set_etat_qcq() ;
311 Mtbl ya_mtbl_deux (source_deux.get_mp()->get_mg()) ;
312 ya_mtbl_deux.set_etat_qcq() ;
313
314 xa_mtbl_deux = source_deux.get_mp()->xa ;
315 ya_mtbl_deux = source_deux.get_mp()->ya ;
316
317 Mtbl x_mtbl_deux (source_deux.get_mp()->get_mg()) ;
318 x_mtbl_deux.set_etat_qcq() ;
319 Mtbl y_mtbl_deux (source_deux.get_mp()->get_mg()) ;
320 y_mtbl_deux.set_etat_qcq() ;
321
322 x_mtbl_deux = source_deux.get_mp()->x ;
323 y_mtbl_deux = source_deux.get_mp()->y ;
324
325 Valeur lim_x_deux (*hole2.mp.get_mg()->get_angu()) ;
326 lim_x_deux = 1 ; // Juste pour affecter dans espace des configs ;
327 lim_x_deux.set_etat_c_qcq() ;
328 for (int k=0 ; k<np_deux ; k++)
329 for (int j=0 ; j<nt_deux ; j++)
330 lim_x_deux.set(0, k, j, 0) = aligne_deux*omega*ya_mtbl_deux(0, 0, 0, 0) + aligne_deux*hole2.omega_local*y_mtbl_deux(1,k,j,0) ;
331 lim_x_deux.base = *bases_deux[0] ;
332
333 Valeur lim_y_deux (*hole2.mp.get_mg()->get_angu()) ;
334 lim_y_deux = 1 ; // Juste pour affecter dans espace des configs ;
335 lim_y_deux.set_etat_c_qcq() ;
336 for (int k=0 ; k<np_deux ; k++)
337 for (int j=0 ; j<nt_deux ; j++)
338 lim_y_deux.set(0, k, j, 0) = -aligne_deux*omega*xa_mtbl_deux(0, 0, 0, 0) - aligne_deux*hole2.omega_local*x_mtbl_deux(1,k,j,0) ;
339 lim_y_deux.base = *bases_deux[1] ;
340
341 Valeur lim_z_deux (*hole2.mp.get_mg()->get_angu()) ;
342 lim_z_deux = 1 ;
343 for (int k=0 ; k<np_deux ; k++)
344 for (int j=0 ; j<nt_deux ; j++)
345 lim_z_deux.set(0, k, j, 0) = 0 ;
346 lim_z_deux.base = *bases_deux[2] ;
347
348 for (int i=0 ; i<3 ; i++) {
349 delete bases_un[i] ;
350 delete bases_deux[i] ;
351 }
352 delete [] bases_un ;
353 delete [] bases_deux ;
354
355 // On resout le truc :
356 Tenseur shift_un_old (hole1.shift_auto) ;
357 Tenseur shift_deux_old (hole2.shift_auto) ;
358
359 poisson_vect_binaire (1./3., source_un, source_deux,
360 lim_x_un, lim_y_un, lim_z_un,
361 lim_x_deux, lim_y_deux, lim_z_deux,
362 hole1.shift_auto, hole2.shift_auto, 0, precision) ;
363
364 for (int i=0 ; i<3 ; i++) {
365 hole1.shift_auto.set(i).raccord(1) ;
366 hole2.shift_auto.set(i).raccord(1) ;
367 }
368
369 // On regularise les shift.
370 hole1.shift_auto = relax*hole1.shift_auto +
371 (1-relax)*shift_un_old ;
372 hole2.shift_auto = relax*hole2.shift_auto +
373 (1-relax)*shift_deux_old ;
374
375 double diff_un = regle (hole1.shift_auto, hole2.shift_auto, omega, hole1.omega_local) ;
376 double diff_deux = regle (hole2.shift_auto, hole1.shift_auto, omega, hole2.omega_local) ;
377 hole1.regul = diff_un ;
378 hole2.regul = diff_deux ;
379
380 cout << "Difference relatives : " << diff_un << " " << diff_deux << endl ;
381}
382}
Bases of the spectral expansions.
Definition base_val.h:325
double omega
Position of the axis of rotation.
Definition bhole.h:769
Bhole hole1
Black hole one.
Definition bhole.h:762
void solve_shift(double precis, double relax)
Solves the equation for the shift, using the Oohara-Nakarmure scheme~:
Bhole hole2
Black hole two.
Definition bhole.h:763
void solve_psi(double precis, double relax)
Solves the equation for the conformal factor~:
void solve_lapse(double precis, double relax)
Solves the equation for the lapse~:
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition cmp.h:446
void std_base_scal()
Sets the spectral bases of the Valeur va to the standard ones for a scalar.
Definition cmp.C:647
void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition cmp.C:292
void filtre(int n)
Sets the n lasts coefficients in r to 0 in the external domain.
Definition cmp_manip.C:77
Multi-domain array.
Definition mtbl.h:118
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition mtbl.C:302
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition tenseur.h:304
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition tenseur.C:830
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition tenseur.C:642
const Map * get_mp() const
Returns pointer on the mapping.
Definition tenseur.h:702
void set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition tenseur.C:1176
void set_etat_zero()
Sets the logical state to ETATZERO (zero state).
Definition tenseur.C:651
int get_etat() const
Returns the logical state.
Definition tenseur.h:710
Values and coefficients of a (real-value) function.
Definition valeur.h:297
void set_etat_c_qcq()
Sets the logical state to ETATQCQ (ordinary state) for values in the configuration space (Mtbl c ).
Definition valeur.C:704
Tbl & set(int l)
Read/write of the value in a given domain (configuration space).
Definition valeur.h:373
Base_val base
Bases on which the spectral expansion is performed.
Definition valeur.h:315
void std_base_scal()
Sets the bases for spectral expansions (member base ) to the standard ones for a scalar.
Definition valeur.C:827
Cmp pow(const Cmp &, int)
Power .
Definition cmp_math.C:351
Tenseur flat_scalar_prod(const Tenseur &t1, const Tenseur &t2)
Scalar product of two Tenseur when the metric is : performs the contraction of the last index of t1 w...
Lorene prototypes.
Definition app_hor.h:67