LORENE
tslice_dirac_max_setAB.C
1/*
2 * Methods of class Tslice_dirac_max dealing with the members potA and tildeB
3 *
4 * (see file time_slice.h for documentation).
5 *
6 */
7
8/*
9 * Copyright (c) 2007 Jerome Novak
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 version 2
15 * as published by the Free Software Foundation.
16 *
17 * LORENE is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
20 * GNU General Public License for more details.
21 *
22 * You should have received a copy of the GNU General Public License
23 * along with LORENE; if not, write to the Free Software
24 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
25 *
26 */
27
28
29
30/*
31 * $Id: tslice_dirac_max_setAB.C,v 1.12 2016/12/05 16:18:19 j_novak Exp $
32 * $Log: tslice_dirac_max_setAB.C,v $
33 * Revision 1.12 2016/12/05 16:18:19 j_novak
34 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
35 *
36 * Revision 1.11 2014/10/13 08:53:48 j_novak
37 * Lorene classes and functions now belong to the namespace Lorene.
38 *
39 * Revision 1.10 2014/10/06 15:13:22 j_novak
40 * Modified #include directives to use c++ syntax.
41 *
42 * Revision 1.9 2012/02/06 12:59:07 j_novak
43 * Correction of some errors.
44 *
45 * Revision 1.8 2011/07/22 13:21:02 j_novak
46 * Corrected an error on BC treatment.
47 *
48 * Revision 1.7 2010/10/20 07:58:10 j_novak
49 * Better implementation of the explicit time-integration. Not fully-tested yet.
50 *
51 * Revision 1.6 2008/12/04 18:22:49 j_novak
52 * Enhancement of the dzpuis treatment + various bug fixes.
53 *
54 * Revision 1.5 2008/12/02 15:02:22 j_novak
55 * Implementation of the new constrained formalism, following Cordero et al. 2009
56 * paper. The evolution eqs. are solved as a first-order system. Not tested yet!
57 *
58 * Revision 1.4 2007/06/05 07:38:37 j_novak
59 * Better treatment of dzpuis for A and tilde(B) potentials. Some errors in the bases manipulation have been also corrected.
60 *
61 * Revision 1.3 2007/05/24 12:10:41 j_novak
62 * Update of khi_evol and mu_evol.
63 *
64 * Revision 1.2 2007/04/25 15:21:01 j_novak
65 * Corrected an error in the initialization of tildeB in
66 * Tslice_dirac_max::initial_dat_cts. + New method for solve_hij_AB.
67 *
68 * Revision 1.1 2007/03/21 14:51:50 j_novak
69 * Introduction of potentials A and tilde(B) of h^{ij} into Tslice_dirac_max.
70 *
71 *
72 * $Header $
73 *
74 */
75
76// C headers
77#include <cassert>
78
79// Lorene headers
80#include "time_slice.h"
81#include "param.h"
82#include "unites.h"
83#include "proto.h"
84#include "graphique.h"
85
86namespace Lorene {
87void Tslice_dirac_max::set_AB_hh(const Scalar& A_in, const Scalar& B_in) {
88
89 A_hh_evol.update(A_in, jtime, the_time[jtime]) ;
90 B_hh_evol.update(B_in, jtime, the_time[jtime]) ;
91
92 // Reset of quantities depending on h^{ij}:
93 hh_evol.downdate(jtime) ;
94 trh_evol.downdate(jtime) ;
95 if (p_tgamma != 0x0) {
96 delete p_tgamma ;
97 p_tgamma = 0x0 ;
98 }
99 if (p_hdirac != 0x0) {
100 delete p_hdirac ;
101 p_hdirac = 0x0 ;
102 }
103 if (p_gamma != 0x0) {
104 delete p_gamma ;
105 p_gamma = 0x0 ;
106 }
107 gam_dd_evol.downdate(jtime) ;
108 gam_uu_evol.downdate(jtime) ;
109 adm_mass_evol.downdate(jtime) ;
110}
111
112void Tslice_dirac_max::hh_det_one(int j0, Param* par_bc, Param* par_mat) const {
113
114 assert (A_hh_evol.is_known(j0)) ; // The starting point
115 assert (B_hh_evol.is_known(j0)) ; // of the computation
116
117 const Map& mp = A_hh_evol[j0].get_mp() ;
118
119 // The representation of h^{ij} as an object of class Sym_tensor_trans :
120 Sym_tensor_trans hij(mp, *(ff.get_triad()), ff) ;
121 const Scalar* ptrace = 0x0 ;
122 if (trh_evol.is_known(j0-1)) ptrace = &trh_evol[j0-1] ;
123 hij.set_AtBtt_det_one(A_hh_evol[j0], B_hh_evol[j0], ptrace, par_bc, par_mat) ;
124
125 // Result set to trh_evol and hh_evol
126 // ----------------------------------
127 trh_evol.update(hij.the_trace(), j0, the_time[j0]) ;
128
129 // The longitudinal part of h^{ij}, which is zero by virtue of Dirac gauge :
130 Vector wzero(mp, CON, *(ff.get_triad())) ;
131 wzero.set_etat_zero() ;
132
133 // Temporary Sym_tensor with longitudinal part set to zero :
134 Sym_tensor hh_new(mp, CON, *(ff.get_triad())) ;
135 hh_new.set_longit_trans(wzero, hij) ;
136 hh_evol.update(hh_new, j0, the_time[j0]) ;
137
138 if (j0 == jtime) {
139 // Reset of quantities depending on h^{ij}:
140 if (p_tgamma != 0x0) {
141 delete p_tgamma ;
142 p_tgamma = 0x0 ;
143 }
144 if (p_hdirac != 0x0) {
145 delete p_hdirac ;
146 p_hdirac = 0x0 ;
147 }
148 if (p_gamma != 0x0) {
149 delete p_gamma ;
150 p_gamma = 0x0 ;
151 }
152 }
153 gam_dd_evol.downdate(j0) ;
154 gam_uu_evol.downdate(j0) ;
155 adm_mass_evol.downdate(j0) ;
156
157#ifndef NDEBUG
158 // Test
159 if (j0 == jtime) {
160 maxabs(tgam().determinant() - 1,
161 "Max. of absolute value of deviation from det tgam = 1") ;
162 }
163 else {
164 Metric tgam_j0( ff.con() + hh_evol[j0] ) ;
165 maxabs(tgam_j0.determinant() - 1,
166 "Max. of absolute value of deviation from det tgam = 1") ;
167 }
168#endif
169}
170
172 const {
173
174 const Map& mp = hijtt.get_mp() ;
175
176 // The representation of h^{ij} as an object of class Sym_tensor_trans :
177 Sym_tensor_trans hij(mp, *(ff.get_triad()), ff) ;
178 const Scalar* ptrace = 0x0 ;
179 if ( trh_evol.is_known( jtime - 1 ) ) ptrace = &trh_evol[jtime-1] ;
180 hij.set_tt_part_det_one(hijtt, ptrace, par_mat) ;
181
182 // Result set to trh_evol and hh_evol
183 // ----------------------------------
184 trh_evol.update(hij.the_trace(), jtime, the_time[jtime]) ;
185
186 // The longitudinal part of h^{ij}, which is zero by virtue of Dirac gauge :
187 Vector wzero(mp, CON, *(ff.get_triad())) ;
188 wzero.set_etat_zero() ;
189
190 // Temporary Sym_tensor with longitudinal part set to zero :
191 Sym_tensor hh_new(mp, CON, *(ff.get_triad())) ;
192 hh_new.set_longit_trans(wzero, hij) ;
193 hh_evol.update(hh_new, jtime, the_time[jtime]) ;
194
195 // Reset of quantities depending on h^{ij}:
196 if (p_tgamma != 0x0) {
197 delete p_tgamma ;
198 p_tgamma = 0x0 ;
199 }
200 if (p_hdirac != 0x0) {
201 delete p_hdirac ;
202 p_hdirac = 0x0 ;
203 }
204 if (p_gamma != 0x0) {
205 delete p_gamma ;
206 p_gamma = 0x0 ;
207 }
208 gam_dd_evol.downdate(jtime) ;
209 gam_uu_evol.downdate(jtime) ;
210 adm_mass_evol.downdate(jtime) ;
211
212#ifndef NDEBUG
213 // Test
214 maxabs(tgam().determinant() - 1,
215 "Max. of absolute value of deviation from det tgam = 1") ;
216#endif
217}
218 //----------------------------------------------//
219 // Equations for h^{ij} and \hat{A}^{ij} //
220 //----------------------------------------------//
221
222void Tslice_dirac_max::compute_sources( const Sym_tensor* p_strain_tens) const {
223
224 using namespace Unites ;
225
226 const Map& map = hh().get_mp() ;
227 const Base_vect& otriad = *hh().get_triad() ;
228 int nz = map.get_mg()->get_nzone() ;
229
230 Sym_tensor strain_tens(map, CON, otriad) ;
231 if (p_strain_tens != 0x0)
232 strain_tens = *(p_strain_tens) ;
233 else
234 strain_tens.set_etat_zero() ;
235
236 Sym_tensor aij = aa() ;
237 aij.annule_domain(nz-1) ;
238
239 const Sym_tensor& tgam_dd = tgam().cov() ; // {\tilde \gamma}_{ij}
240 const Sym_tensor& tgam_uu = tgam().con() ; // {\tilde \gamma}^{ij}
241 const Tensor_sym& dtgam = tgam_dd.derive_cov(ff) ;// D_k {\tilde \gamma}_{ij}
242 const Tensor_sym& dhh = hh().derive_cov(ff) ; // D_k h^{ij}
243 const Vector& dln_psi = ln_psi().derive_cov(ff) ; // D_i ln(Psi)
244 const Vector& tdln_psi_u = ln_psi().derive_con(tgam()) ; // tD^i ln(Psi)
245 Scalar log_N = log(nn()) ;
246 log_N.std_spectral_base() ;
247 const Vector& tdlnn_u = log_N.derive_con(tgam()) ; // tD^i ln(N)
248 const Scalar& div_beta = beta().divergence(ff) ; // D_k beta^k
249 Scalar qq = nn()*psi()*psi() ;
250 qq.annule_domain(nz-1) ;
251 const Vector& dqq = qq.derive_cov(ff) ; // D_i Q
252 Sym_tensor a_hat = hata() ;
253 a_hat.annule_domain(nz-1) ;
254 Scalar psi6 = psi4()*psi()*psi() ;
255 Sym_tensor sym_tmp(map, CON, otriad) ;
256 Scalar tmp(map) ;
257
258 //==================================
259 // Source for hij
260 //==================================
261
262 Sym_tensor source_hij = hh().derive_lie(beta()) + 2*(nn()/psi6 - 1.)*a_hat
263 - beta().ope_killing_conf(ff) + 0.6666666666666667*div_beta*hh() ;
264 source_hij.annule_domain(nz-1) ;
265 for (int i=1; i<=3; i++)
266 for (int j=i; j<=3; j++)
267 source_hij.set( i, j ).set_dzpuis(0) ;
268
269 tmp = 2.*A_hata_evol[jtime] + source_hij.compute_A(true) ;
270 tmp.set_spectral_va().ylm() ;
271 tmp.annule_domain(nz-1) ;
272 tmp.set_dzpuis(0) ;
273 source_A_hh_evol.update( tmp, jtime, the_time[jtime] ) ;
274
275 tmp = 2.*B_hata_evol[jtime] + source_hij.compute_tilde_B_tt(true) ;
276 tmp.set_spectral_va().ylm() ;
277 tmp.annule_domain(nz-1) ;
278 tmp.set_dzpuis(0) ;
279 source_B_hh_evol.update(tmp, jtime, the_time[jtime] ) ;
280
281 //==================================
282 // Source for \hat{A}^{ij}
283 //==================================
284
285 Sym_tensor source_aij = a_hat.derive_lie(beta())
286 + 1.666666666666667*a_hat*div_beta ;
287
288 // Quadratic part of the Ricci tensor of gam_tilde
289 // ------------------------------------------------
290
291 Sym_tensor ricci_star(map, CON, otriad) ;
292
293 ricci_star = contract(hh(), 0, 1, dhh.derive_cov(ff), 2, 3) ;
294
295 ricci_star.inc_dzpuis() ; // dzpuis : 3 --> 4
296
297 for (int i=1; i<=3; i++) {
298 for (int j=1; j<=i; j++) {
299 tmp = 0 ;
300 for (int k=1; k<=3; k++) {
301 for (int l=1; l<=3; l++) {
302 tmp += dhh(i,k,l) * dhh(j,l,k) ;
303 }
304 }
305 ricci_star.set(i,j) -= tmp ;
306 }
307 }
308
309 for (int i=1; i<=3; i++) {
310 for (int j=1; j<=i; j++) {
311 tmp = 0 ;
312 for (int k=1; k<=3; k++) {
313 for (int l=1; l<=3; l++) {
314 for (int m=1; m<=3; m++) {
315 for (int n=1; n<=3; n++) {
316
317 tmp += 0.5 * tgam_uu(i,k)* tgam_uu(j,l)
318 * dhh(m,n,k) * dtgam(m,n,l)
319 + tgam_dd(n,l) * dhh(m,n,k)
320 * (tgam_uu(i,k) * dhh(j,l,m) + tgam_uu(j,k) * dhh(i,l,m) )
321 - tgam_dd(k,l) *tgam_uu(m,n) * dhh(i,k,m) * dhh(j,l,n) ;
322 }
323 }
324 }
325 }
326 sym_tmp.set(i,j) = tmp ;
327 }
328 }
329 ricci_star += sym_tmp ; // a factor 1/2 is still missing, shall be put later
330
331 // Curvature scalar of conformal metric :
332 // -------------------------------------
333
334 Scalar tricci_scal =
335 0.25 * contract(tgam_uu, 0, 1,
336 contract(dhh, 0, 1, dtgam, 0, 1), 0, 1 )
337 - 0.5 * contract(tgam_uu, 0, 1,
338 contract(dhh, 0, 1, dtgam, 0, 2), 0, 1 ) ;
339
340 Scalar lap_A = A_hh_evol[jtime].laplacian(2) ;
341 Scalar tilde_lap_B(map) ;
342 tilde_laplacian( B_hh_evol[jtime], tilde_lap_B) ;
343 Sym_tensor_tt laplace_h(map, otriad, ff) ;
344 laplace_h.set_A_tildeB(lap_A, tilde_lap_B) ;
345 laplace_h.annule_domain(nz-1) ;
346
347 // sym_tmp.inc_dzpuis() ; // dzpuis : 3 --> 4
348
349 source_aij += (0.5*(qq - 1.))*laplace_h + qq*(0.5*ricci_star + 8.*tdln_psi_u * tdln_psi_u
350 + 4.*( tdln_psi_u * tdlnn_u + tdlnn_u * tdln_psi_u )
351 - 0.3333333333333333 * (tricci_scal + 8.*(contract(dln_psi, 0, tdln_psi_u, 0)
352 + contract(dln_psi, 0, tdlnn_u, 0) )
353 )*tgam_uu
354 ) ;
355
356 sym_tmp = contract(tgam_uu, 1, contract(tgam_uu, 1, dqq.derive_cov(ff), 0), 1) ;
357
358 for (int i=1; i<=3; i++) {
359 for (int j=1; j<=i; j++) {
360 tmp = 0 ;
361 for (int k=1; k<=3; k++) {
362 for (int l=1; l<=3; l++) {
363 tmp += ( tgam_uu(i,k)*dhh(l,j,k) + tgam_uu(k,j)*dhh(i,l,k)
364 - tgam_uu(k,l)*dhh(i,j,k) ) * dqq(l) ;
365 }
366 }
367 sym_tmp.set(i,j) += 0.5 * tmp ;
368 }
369 }
370
371 source_aij -= sym_tmp
372 - ( 0.3333333333333333*qq.derive_con(tgam()).divergence(tgam()) ) *tgam_uu ;
373
374 for (int i=1; i<=3; i++) {
375 for (int j=1; j<=i; j++) {
376 tmp = 0 ;
377 for (int k=1; k<=3; k++) {
378 for (int l=1; l<=3; l++) {
379 tmp += tgam_dd(k,l) * a_hat(i,k) * aij(j,l) ;
380 }
381 }
382 sym_tmp.set(i,j) = tmp ;
383 }
384 }
385
386 tmp = psi4() * strain_tens.trace(tgam()) ; // S = S_i^i
387
388 source_aij += (2.*nn())
389 * (
390 sym_tmp - qpig*psi6*( psi4()* strain_tens - (0.3333333333333333 * tmp) * tgam_uu )
391 ) ;
392
393 source_aij.annule_domain(nz-1) ;
394 for (int i=1; i<=3; i++)
395 for (int j=i; j<=3; j++)
396 source_aij.set(i,j).set_dzpuis(0) ;
397#ifndef NDEBUG
398 maxabs(source_aij, "source_aij tot") ;
399#endif
400
401 tmp = 0.5*lap_A + source_aij.compute_A(true) ;
402 tmp.annule_domain(nz-1) ;
403 tmp.set_dzpuis(0) ;
404 source_A_hata_evol.update( tmp, jtime, the_time[jtime] ) ;
405 tmp = 0.5*tilde_lap_B + source_aij.compute_tilde_B_tt(true) ;
406 tmp.annule_domain(nz-1) ;
407 tmp.set_dzpuis(0) ;
408 // Scalar dess = tmp - tilde_lap_B ;
409 // dess.set_spectral_va().ylm_i() ;
410 // des_profile(dess, 0, 8., 1, 1) ;
411 source_B_hata_evol.update( 0.5*tilde_lap_B, jtime, the_time[jtime] ) ;
412}
413
415
416 assert( source_A_hh_evol.is_known(jtime) ) ;
417 assert( source_B_hh_evol.is_known(jtime) ) ;
418 assert( source_A_hata_evol.is_known(jtime) ) ;
419 assert( source_B_hata_evol.is_known(jtime) ) ;
420
421 Scalar tmp_Ah = source_A_hh_evol[jtime] ;
422 Scalar tmp_Bh = source_B_hh_evol[jtime] ;
425
426 source_A_hh_evol.downdate(jtime) ;
427 source_B_hh_evol.downdate(jtime) ;
428 source_A_hata_evol.downdate(jtime) ;
429 source_B_hata_evol.downdate(jtime) ;
430
431 int jtime1 = jtime - depth + 1;
432 for (int j=jtime1; j <= jtime; j++) {
433 source_A_hh_evol.update( tmp_Ah, j, the_time[j] ) ;
434 source_B_hh_evol.update( tmp_Bh, j, the_time[j] ) ;
435 source_A_hata_evol.update( tmp_Aa, j, the_time[j] ) ;
436 source_B_hata_evol.update( tmp_Ba, j, the_time[j] ) ;
437 }
438}
439
440}
Vectorial bases (triads) with respect to which the tensorial components are defined.
Definition base_vect.h:105
Metric for tensor calculation.
Definition metric.h:90
virtual const Sym_tensor & con() const
Read-only access to the contravariant representation.
Definition metric.C:293
virtual const Sym_tensor & cov() const
Read-only access to the covariant representation.
Definition metric.C:283
virtual const Scalar & determinant() const
Returns the determinant.
Definition metric.C:395
Parameter storage.
Definition param.h:125
Tensor field of valence 0 (or component of a tensorial field).
Definition scalar.h:393
const Vector & derive_cov(const Metric &gam) const
Returns the gradient (1-form = covariant vector) of *this.
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
Valeur & set_spectral_va()
Returns va (read/write version).
Definition scalar.h:610
void set_dzpuis(int)
Modifies the dzpuis flag.
Definition scalar.C:814
const Vector & derive_con(const Metric &gam) const
Returns the "contravariant" derivative of *this with respect to some metric , by raising the index of...
Transverse symmetric tensors of rank 2.
Definition sym_tensor.h:611
void set_AtBtt_det_one(const Scalar &a_in, const Scalar &tbtt_in, const Scalar *h_prev=0x0, Param *par_bc=0x0, Param *par_mat=0x0, double precis=1.e-14, int it_max=100)
Assigns the derived member A and computes from its TT-part (see Sym_tensor::compute_tilde_B_tt() ).
const Scalar & the_trace() const
Returns the trace of the tensor with respect to metric *met_div.
void set_tt_part_det_one(const Sym_tensor_tt &hijtt, const Scalar *h_prev=0x0, Param *par_mat=0x0, double precis=1.e-14, int it_max=100)
Assignes the TT-part of the tensor.
Transverse and traceless symmetric tensors of rank 2.
Definition sym_tensor.h:933
void set_A_tildeB(const Scalar &a_in, const Scalar &tb_in, Param *par_bc=0x0, Param *par_mat=0x0)
Assigns the derived members A and .
Class intended to describe valence-2 symmetric tensors.
Definition sym_tensor.h:226
void set_longit_trans(const Vector &v, const Sym_tensor_trans &a)
Assigns the derived members p_longit_pot and p_transverse and updates the components accordingly.
const Scalar & compute_A(bool output_ylm=true, Param *par=0x0) const
Gives the field A (see member p_aaa ).
Scalar compute_tilde_B_tt(bool output_ylm=true, Param *par=0x0) const
Gives the field (see member p_tilde_b ) associated with the TT-part of the Sym_tensor .
Sym_tensor derive_lie(const Vector &v) const
Computes the Lie derivative of this with respect to some vector field v.
Definition sym_tensor.C:363
Symmetric tensors (with respect to two of their arguments).
Definition tensor.h:1050
Evolution_std< Sym_tensor > hh_evol
Values at successive time steps of the components .
Definition time_slice.h:533
virtual const Scalar & nn() const
Lapse function N at the current time step (jtime ).
virtual const Sym_tensor & hata() const
Conformal representation of the traceless part of the extrinsic curvature: .
const Scalar & ln_psi() const
Logarithm of at the current time step (jtime ).
virtual const Scalar & psi() const
Conformal factor relating the physical metric to the conformal one: .
Metric * p_tgamma
Pointer on the conformal metric at the current time step (jtime).
Definition time_slice.h:563
virtual Sym_tensor aa() const
Conformal representation of the traceless part of the extrinsic curvature: .
Vector * p_hdirac
Pointer on the vector (which vanishes in Dirac gauge), at the current time step (jtime).
Definition time_slice.h:574
const Scalar & psi4() const
Factor at the current time step (jtime ).
virtual const Metric & tgam() const
Conformal metric Returns the value at the current time step (jtime ).
Evolution_std< Scalar > A_hata_evol
Potential A associated with the symmetric tensor .
Definition time_slice.h:550
const Metric_flat & ff
Pointer on the flat metric with respect to which the conformal decomposition is performed.
Definition time_slice.h:510
Evolution_std< Scalar > B_hata_evol
Potential associated with the symmetric tensor .
Definition time_slice.h:555
int jtime
Time step index of the latest slice.
Definition time_slice.h:193
Metric * p_gamma
Pointer on the induced metric at the current time step (jtime).
Definition time_slice.h:242
Evolution_full< Tbl > adm_mass_evol
ADM mass at each time step, since the creation of the slice.
Definition time_slice.h:236
Evolution_std< Sym_tensor > gam_dd_evol
Values at successive time steps of the covariant components of the induced metric .
Definition time_slice.h:201
Evolution_std< Sym_tensor > gam_uu_evol
Values at successive time steps of the contravariant components of the induced metric .
Definition time_slice.h:206
virtual const Vector & beta() const
shift vector at the current time step (jtime )
int depth
Number of stored time slices.
Definition time_slice.h:182
Evolution_std< double > the_time
Time label of each slice.
Definition time_slice.h:196
virtual const Sym_tensor & hh(Param *par_bc=0x0, Param *par_mat=0x0) const
Deviation of the conformal metric from the flat metric : .
Evolution_std< Scalar > source_B_hh_evol
The potential of the source of equation for .
Definition time_slice.h:998
Evolution_std< Scalar > source_A_hh_evol
The A potential of the source of equation for .
Definition time_slice.h:992
Evolution_std< Scalar > A_hh_evol
The A potential of .
Definition time_slice.h:980
void hh_det_one(int j, Param *par_bc=0x0, Param *par_mat=0x0) const
Computes from the values of A and and using the condition , which fixes the trace of .
Evolution_std< Scalar > source_B_hata_evol
The potential of the source of equation for .
void compute_sources(const Sym_tensor *strain_tensor=0x0) const
Computes the sources source_A_XXX_evol and source_B_XXX_evol , for the solution of the evolution equa...
Evolution_std< Scalar > trh_evol
The trace, with respect to the flat metric ff , of .
void initialize_sources_copy() const
Copy the sources source_A_XXX_evol and source_B_XXX_evol to all time-steps.
virtual void set_AB_hh(const Scalar &A_in, const Scalar &B_in)
Sets the potentials A and of the TT part of (see the documentation of Sym_tensor for details).
Evolution_std< Scalar > B_hh_evol
The potential of .
Definition time_slice.h:986
Evolution_std< Scalar > source_A_hata_evol
The potential A of the source of equation for .
void ylm()
Computes the coefficients of *this.
Definition valeur_ylm.C:141
Tensor field of valence 1.
Definition vector.h:188
const Scalar & divergence(const Metric &) const
The divergence of this with respect to a Metric .
Definition vector.C:384
Sym_tensor ope_killing_conf(const Metric &gam) const
Computes the conformal Killing operator associated with a given metric.
Definition vector.C:465
Cmp log(const Cmp &)
Neperian logarithm.
Definition cmp_math.C:299
const Map & get_mp() const
Returns the mapping.
Definition tensor.h:874
const Tensor_sym & derive_cov(const Metric &gam) const
Returns the covariant derivative of this with respect to some metric .
virtual void inc_dzpuis(int inc=1)
Increases by inc units the value of dzpuis and changes accordingly the values in the compactified ext...
Definition tensor.C:825
const Base_vect * get_triad() const
Returns the vectorial basis (triad) on which the components are defined.
Definition tensor.h:879
void annule_domain(int l)
Sets the Tensor to zero in a given domain.
Definition tensor.C:675
const Tensor & derive_cov(const Metric &gam) const
Returns the covariant derivative of this with respect to some metric .
Definition tensor.C:1011
Tensor trace(int ind1, int ind2) const
Trace on two different type indices.
virtual void set_etat_zero()
Sets the logical state of all components to ETATZERO (zero state).
Definition tensor.C:506
Scalar & set(const Itbl &ind)
Returns the value of a component (read/write version).
Definition tensor.C:663
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
Tbl maxabs(const Tensor &aa, const char *comment=0x0, ostream &ost=cout, bool verb=true)
Maxima in each domain of the absolute values of the tensor components.
Lorene prototypes.
Definition app_hor.h:67
Map(const Mg3d &)
Constructor from a multi-domain 3D grid.
Definition map.C:142
Standard units of space, time and mass.