LORENE
tslice_dirac_max.C
1  /*
2  * Methods of class Tslice_dirac_max
3  *
4  * (see file time_slice.h for documentation).
5  *
6  */
7 
8 /*
9  * Copyright (c) 2004 Eric Gourgoulhon, Jose Luis Jaramillo & 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 char tslice_dirac_max_C[] = "$Header: /cvsroot/Lorene/C++/Source/Time_slice/tslice_dirac_max.C,v 1.26 2014/10/13 08:53:48 j_novak Exp $" ;
29 
30 /*
31  * $Id: tslice_dirac_max.C,v 1.26 2014/10/13 08:53:48 j_novak Exp $
32  * $Log: tslice_dirac_max.C,v $
33  * Revision 1.26 2014/10/13 08:53:48 j_novak
34  * Lorene classes and functions now belong to the namespace Lorene.
35  *
36  * Revision 1.25 2014/10/06 15:13:22 j_novak
37  * Modified #include directives to use c++ syntax.
38  *
39  * Revision 1.24 2008/12/04 18:22:49 j_novak
40  * Enhancement of the dzpuis treatment + various bug fixes.
41  *
42  * Revision 1.23 2008/12/02 15:02:22 j_novak
43  * Implementation of the new constrained formalism, following Cordero et al. 2009
44  * paper. The evolution eqs. are solved as a first-order system. Not tested yet!
45  *
46  * Revision 1.22 2007/11/06 14:47:07 j_novak
47  * New constructor from a rotating star in Dirac gauge (class Star_rot_Dirac).
48  * Evolution can take into account matter terms.
49  *
50  * Revision 1.21 2007/09/25 16:54:11 j_novak
51  * *** empty log message ***
52  *
53  * Revision 1.20 2007/09/25 16:52:15 j_novak
54  * *** empty log message ***
55  *
56  * Revision 1.19 2007/06/05 07:38:37 j_novak
57  * Better treatment of dzpuis for A and tilde(B) potentials. Some errors in the bases manipulation have been also corrected.
58  *
59  * Revision 1.18 2007/04/25 15:21:01 j_novak
60  * Corrected an error in the initialization of tildeB in
61  * Tslice_dirac_max::initial_dat_cts. + New method for solve_hij_AB.
62  *
63  * Revision 1.17 2007/03/21 14:51:50 j_novak
64  * Introduction of potentials A and tilde(B) of h^{ij} into Tslice_dirac_max.
65  *
66  * Revision 1.16 2004/12/28 14:21:48 j_novak
67  * Added the method Sym_tensor_trans::trace_from_det_one
68  *
69  * Revision 1.15 2004/07/08 12:29:01 j_novak
70  * use of new method Tensor::annule_extern_cn
71  *
72  * Revision 1.14 2004/06/30 08:02:40 j_novak
73  * Added filtering in l of khi_new and mu_new. ki_source is forced to go to
74  * zero at least as r^2.
75  *
76  * Revision 1.13 2004/06/17 06:59:41 e_gourgoulhon
77  * -- Method initial_data_cts: re-organized treatment of vanishing uu.
78  * -- Method hh_det_one: replaced the attenuation with tempo by a call
79  * to the new method Tensor::annule_extern_c2.
80  *
81  * Revision 1.12 2004/06/08 14:05:06 j_novak
82  * Added the attenuation of khi and mu in the last domain in ::det_one(). They are set to zero in the CED.
83  *
84  * Revision 1.11 2004/05/31 20:31:31 e_gourgoulhon
85  * -- Method hh_det_one takes now a time step as argument, to compute
86  * h^{ij} from khi and mu at some arbitrary time step and not only at
87  * the latest one.
88  * -- h^{ij} is no longer saved in binary files (method sauve);
89  * accordingly, the constructor from file calls the new version of
90  * hh_det_one to restore h^{ij}.
91  *
92  * Revision 1.10 2004/05/31 09:08:18 e_gourgoulhon
93  * Method sauve and constructor from binary file are now operational.
94  *
95  * Revision 1.9 2004/05/27 15:25:04 e_gourgoulhon
96  * Added constructors from binary file, as well as corresponding
97  * functions sauve and save.
98  *
99  * Revision 1.8 2004/05/17 19:54:10 e_gourgoulhon
100  * Method initial_data_cts: added arguments graph_device and method_poisson_vect.
101  *
102  * Revision 1.7 2004/05/12 15:24:20 e_gourgoulhon
103  * Reorganized the #include 's, taking into account that
104  * time_slice.h contains now an #include "metric.h".
105  *
106  * Revision 1.6 2004/05/10 09:16:32 e_gourgoulhon
107  * -- Method initial_data_cts: added a call to del_deriv() at the end.
108  * -- Methods set_trh and hh_det_one: added "adm_mass_evol.downdate(jtime)".
109  * -- Method trh() : the update is now performed via a call to hh_det_one().
110  *
111  * Revision 1.5 2004/05/06 15:23:55 e_gourgoulhon
112  * Added method initial_data_cts.
113  *
114  * Revision 1.4 2004/05/03 08:15:48 e_gourgoulhon
115  * Method hh_det_one(): added check at the end (deviation from det = 1).
116  *
117  * Revision 1.3 2004/04/08 16:44:19 e_gourgoulhon
118  * Added methods set_* and hh_det_one().
119  *
120  * Revision 1.2 2004/04/05 21:22:49 e_gourgoulhon
121  * Added constructor as standard time slice of Minkowski spacetime.
122  *
123  * Revision 1.1 2004/03/30 14:00:31 j_novak
124  * New class Tslide_dirac_max (first version).
125  *
126  *
127  * $Header: /cvsroot/Lorene/C++/Source/Time_slice/tslice_dirac_max.C,v 1.26 2014/10/13 08:53:48 j_novak Exp $
128  *
129  */
130 
131 // C headers
132 #include <cassert>
133 
134 // Lorene headers
135 #include "time_slice.h"
136 #include "utilitaires.h"
137 
138 
139 
140  //--------------//
141  // Constructors //
142  //--------------//
143 
144 
145 // Constructor from conformal decomposition
146 // ----------------------------------------
147 
148 namespace Lorene {
149 Tslice_dirac_max::Tslice_dirac_max(const Scalar& lapse_in, const Vector& shift_in,
150  const Metric_flat& ff_in, const Scalar& psi_in,
151  const Sym_tensor_trans& hh_in, const Sym_tensor& hata_in,
152  int depth_in)
153  : Time_slice_conf( lapse_in, shift_in, ff_in, psi_in, hh_in, hata_in,
154  0*lapse_in, depth_in),
155  A_hh_evol(depth_in), B_hh_evol(depth_in), source_A_hh_evol(depth_in),
156  source_B_hh_evol(depth_in), source_A_hata_evol(depth_in) ,
157  source_B_hata_evol(depth_in), trh_evol(hh_in.the_trace(), depth_in)
158 { }
159 
160 // Constructor as standard time slice of flat spacetime (Minkowski)
161 // ----------------------------------------------------------------
162 
164  const Metric_flat& ff_in, int depth_in)
165  : Time_slice_conf(mp, triad, ff_in, depth_in),
166  A_hh_evol(depth_in), B_hh_evol(depth_in),
167  source_A_hh_evol(depth_in), source_B_hh_evol(depth_in),
168  source_A_hata_evol(depth_in), source_B_hata_evol(depth_in),
169  trh_evol(depth_in) {
170 
171  double time_init = the_time[jtime] ;
172 
173  // All potentials identically zero:
174  Scalar tmp(mp) ;
175  tmp.set_etat_zero() ;
176  A_hh_evol.update(tmp, jtime, time_init) ;
177  B_hh_evol.update(tmp, jtime, time_init) ;
178 
179  source_A_hh_evol.update(tmp, jtime, time_init) ;
180  source_B_hh_evol.update(tmp, jtime, time_init) ;
181 
182  source_A_hata_evol.update(tmp, jtime, time_init) ;
183  source_B_hata_evol.update(tmp, jtime, time_init) ;
184 
185  // tr h identically zero:
186  trh_evol.update(tmp, jtime, time_init) ;
187 }
188 
189 
190 // Constructor from binary file
191 // ----------------------------
192 
194  const Metric_flat& ff_in, FILE* fich,
195  bool partial_read, int depth_in)
196  : Time_slice_conf(mp, triad, ff_in, fich, true, depth_in),
197  A_hh_evol(depth_in), B_hh_evol(depth_in),
198  source_A_hh_evol(depth_in), source_B_hh_evol(depth_in),
199  source_A_hata_evol(depth_in), source_B_hata_evol(depth_in),
200  trh_evol(depth_in) {
201 
202  if (partial_read) {
203  cout <<
204  "Constructor of Tslice_dirac_max from file: the case of partial reading\n"
205  << " is not ready yet !"
206  << endl ;
207  abort() ;
208  }
209 
210  // Reading of various fields
211  // -------------------------
212 
213  int jmin = jtime - depth + 1 ;
214  int indicator ;
215 
216  // h^{ij}
217  for (int j=jmin; j<=jtime; j++) {
218  fread_be(&indicator, sizeof(int), 1, fich) ;
219  if (indicator == 1) {
220  Sym_tensor hh_file(mp, triad, fich) ;
221  hh_evol.update(hh_file, j, the_time[j]) ;
222  }
223  }
224 
225  // A - hh
226  for (int j=jmin; j<=jtime; j++) {
227  fread_be(&indicator, sizeof(int), 1, fich) ;
228  if (indicator == 1) {
229  Scalar A_hh_file(mp, *(mp.get_mg()), fich) ;
230  A_hh_evol.update(A_hh_file, j, the_time[j]) ;
231  }
232  }
233 
234  // B - hh
235  for (int j=jmin; j<=jtime; j++) {
236  fread_be(&indicator, sizeof(int), 1, fich) ;
237  if (indicator == 1) {
238  Scalar B_hh_file(mp, *(mp.get_mg()), fich) ;
239  B_hh_evol.update(B_hh_file, j, the_time[j]) ;
240  }
241  }
242 }
243 
244 // Constructor from a rotating star
245 // --------------------------------
246 
247 Tslice_dirac_max::Tslice_dirac_max(const Star_rot_Dirac& star, double pdt, int depth_in)
248  : Time_slice_conf(star.get_nn(), star.get_beta(), star.get_mp().flat_met_spher(),
249  0.*star.get_nn(), star.get_hh(), 0.*star.get_aa(),
250  0.*star.get_nn(), depth_in),
251  A_hh_evol(depth_in), B_hh_evol(depth_in),
252  source_A_hh_evol(depth_in), source_B_hh_evol(depth_in),
253  source_A_hata_evol(depth_in), source_B_hata_evol(depth_in),
254  trh_evol(depth_in) {
255 
256  Scalar tmp = exp(star.get_ln_psi()) ;
257  tmp.std_spectral_base() ;
258  psi_evol.downdate(jtime) ;
259  psi_evol.update(tmp, jtime, the_time[jtime]) ;
260 
261  Sym_tensor tmp2 = psi4()*psi()*psi()*star.get_aa() ;
262  hata_evol.downdate(jtime) ;
263  hata_evol.update(tmp2, jtime, the_time[jtime]) ;
264 
265  A_hh() ;
266  B_hh() ;
267  A_hata() ;
268  B_hata() ;
269  compute_sources() ;
270 
271  // Update of various fields
272  // -------------------------
273  double ttime1 = the_time[jtime] ;
274  int jtime1 = jtime ;
275  for (int j=1; j < depth; j++) {
276  jtime1++ ;
277  ttime1 += pdt ;
278  psi_evol.update(psi_evol[jtime], jtime1, ttime1) ;
279  n_evol.update(n_evol[jtime], jtime1, ttime1) ;
280  beta_evol.update(beta_evol[jtime], jtime1, ttime1) ;
281  hh_evol.update(hh_evol[jtime], jtime1, ttime1) ;
282  trk_evol.update(trk_evol[jtime], jtime1, ttime1) ;
283  A_hh_evol.update(A_hh_evol[jtime], jtime1, ttime1) ;
284  B_hh_evol.update(B_hh_evol[jtime], jtime1, ttime1) ;
285  A_hata_evol.update(A_hata_evol[jtime], jtime1, ttime1) ;
286  B_hata_evol.update(B_hata_evol[jtime], jtime1, ttime1) ;
287  trh_evol.update(trh_evol[jtime], jtime1, ttime1) ;
288  k_dd_evol.update(k_dd_evol[jtime], jtime1, ttime1) ;
289  the_time.update(ttime1, jtime1, ttime1) ;
290  }
291  jtime += depth - 1 ;
292 
293  initialize_sources_copy() ;
294 }
295 
296 // Copy constructor
297 // ----------------
298 
300  : Time_slice_conf(tin),
301  A_hh_evol(tin.A_hh_evol),
302  B_hh_evol(tin.B_hh_evol),
303  source_A_hh_evol(tin.source_A_hh_evol),
304  source_B_hh_evol(tin.source_B_hh_evol),
305  source_A_hata_evol(tin.source_A_hata_evol),
306  source_B_hata_evol(tin.source_B_hata_evol),
307  trh_evol(tin.trh_evol) { }
308 
309 
310  //--------------//
311  // Destructor //
312  //--------------//
313 
315 
316 
317  //-----------------------//
318  // Mutators / assignment //
319  //-----------------------//
320 
321 void Tslice_dirac_max::operator=(const Tslice_dirac_max& tin) {
322 
324 
325  A_hh_evol = tin.A_hh_evol ;
326  B_hh_evol = tin.B_hh_evol ;
331  trh_evol = tin.trh_evol ;
332 
333 }
334 
335 
337 
338  Time_slice_conf::set_hh(hh_in) ;
339 
340  // Reset of quantities depending on h^{ij}:
341  A_hh_evol.downdate(jtime) ;
342  B_hh_evol.downdate(jtime) ;
343  source_A_hh_evol.downdate(jtime) ;
344  source_B_hh_evol.downdate(jtime) ;
345  source_A_hata_evol.downdate(jtime) ;
346  source_B_hata_evol.downdate(jtime) ;
347  trh_evol.downdate(jtime) ;
348 
349 }
350 
351 
353  const Scalar& trk_in, const Scalar& trk_point,
354  double pdt, double precis, int method_poisson_vect,
355  const char* graph_device, const Scalar* p_ener_dens,
356  const Vector* p_mom_dens, const Scalar* p_trace_stress) {
357 
358 
359  Time_slice_conf::initial_data_cts(uu, trk_in, trk_point, pdt, precis,
360  method_poisson_vect, graph_device,
361  p_ener_dens, p_mom_dens, p_trace_stress) ;
362 
363  int nz = trk_in.get_mp().get_mg()->get_nzone() ;
364  // Setting khi and mu for j <= jtime, taking into account u^{ij} = dh^{ij}/dt
365  //--------------------------------------------------------------------------
366  for (int j = jtime-depth+1 ; j <= jtime; j++) {
367 
368  // A and tildeB are computed from the value of hh
369  Scalar tmp = hh_evol[j].compute_A(true) ;
370  assert (tmp.get_etat() != ETATNONDEF) ;
371  if (tmp.get_etat() != ETATZERO) {
372  assert(tmp.get_mp().get_mg()->get_type_r(nz-1) == UNSURR) ;
373  tmp.annule_domain(nz-1) ;
374  }
375  tmp.set_dzpuis(0) ;
376  A_hh_evol.update(tmp, j, the_time[j]) ;
377 
378  tmp = hh_evol[jtime].compute_tilde_B_tt(true) ;
379  assert (tmp.get_etat() != ETATNONDEF) ;
380  if (tmp.get_etat() != ETATZERO) {
381  assert(tmp.get_mp().get_mg()->get_type_r(nz-1) == UNSURR) ;
382  tmp.annule_domain(nz-1) ;
383  }
384  tmp.set_dzpuis(0) ;
385  B_hh_evol.update(tmp, j, the_time[j]) ;
386  }
387 
388  cout << endl <<
389  "Tslice_dirac_max::initial_data_cts : variation of A and tilde(B) for J = "
390  << jtime << " :\n" ;
391  maxabs(A_hh_evol[jtime] - A_hh_evol[jtime-1], "A(h)^J - A(h)^{J-1}") ;
392 
393  maxabs(B_hh_evol[jtime] - B_hh_evol[jtime-1], "B(h)^J - B(h)^{J-1}") ;
394 
395  maxabs(A_hata_evol[jtime] - A_hata_evol[jtime-1], "A(hat{A})^J - A(hat{A})^{J-1}") ;
396 
397  maxabs(B_hata_evol[jtime] - B_hata_evol[jtime-1], "B(hat{A})^J - B(hat{A})^{J-1}") ;
398 
399  // Reset of derived quantities (at the new time step jtime)
400  // ---------------------------
401  del_deriv() ;
402 
403  compute_sources() ;
404  initialize_sources_copy() ; //## should be a call to the Runge-Kutta integrator
405 }
406 
407 
408 void Tslice_dirac_max::set_khi_mu(const Scalar& khi_in, const Scalar& mu_in) {
409 
410  const Map& mp = khi_in.get_mp() ;
411 
412  Sym_tensor_tt hh_tt(mp, mp.get_bvect_spher(), mp.flat_met_spher());
413  hh_tt.set_khi_mu(khi_in, mu_in, 2) ;
414 
415  Sym_tensor_trans hh_tmp(mp, mp.get_bvect_spher(), mp.flat_met_spher());
416  hh_tmp.trace_from_det_one(hh_tt) ;
417 
418  // Result set to trh_evol and hh_evol
419  // ----------------------------------
420  Scalar tmp = hh_tmp.the_trace() ;
421  tmp.dec_dzpuis(4) ;
422  trh_evol.update(tmp, jtime, the_time[jtime]) ;
423 
424  // The longitudinal part of h^{ij}, which is zero by virtue of Dirac gauge :
425  Vector wzero(mp, CON, *(ff.get_triad())) ;
426  wzero.set_etat_zero() ;
427 
428  // Temporary Sym_tensor with longitudinal part set to zero :
429  Sym_tensor hh_new(mp, CON, *(ff.get_triad())) ;
430 
431  hh_new.set_longit_trans(wzero, hh_tmp) ;
432 
433  hh_evol.update(hh_new, jtime, the_time[jtime]) ;
434 
435 }
436 
437 void Tslice_dirac_max::set_trh(const Scalar& trh_in) {
438 
439  trh_evol.update(trh_in, jtime, the_time[jtime]) ;
440  cout << "Tslice_dirac_max::set_trh : #### WARNING : \n"
441  << " this method does not check whether det(tilde gamma) = 1"
442  << endl ;
443 
444  // Reset of quantities depending on the trace:
445  hh_evol.downdate(jtime) ;
446  if (p_tgamma != 0x0) {
447  delete p_tgamma ;
448  p_tgamma = 0x0 ;
449  }
450  if (p_hdirac != 0x0) {
451  delete p_hdirac ;
452  p_hdirac = 0x0 ;
453  }
454  if (p_gamma != 0x0) {
455  delete p_gamma ;
456  p_gamma = 0x0 ;
457  }
458  source_A_hh_evol.downdate(jtime) ;
459  source_B_hh_evol.downdate(jtime) ;
460  source_A_hata_evol.downdate(jtime) ;
461  source_B_hata_evol.downdate(jtime) ;
462  gam_dd_evol.downdate(jtime) ;
463  gam_uu_evol.downdate(jtime) ;
464  adm_mass_evol.downdate(jtime) ;
465 
466 }
467 
468 
469  //----------------------------------------------------//
470  // Update of fields from base class Time_slice_conf //
471  //----------------------------------------------------//
472 
473 
474 const Sym_tensor& Tslice_dirac_max::hh(Param* par_bc, Param* par_mat) const {
475 
476  if (!( hh_evol.is_known(jtime) ) ) {
477 
478  assert (A_hh_evol.is_known(jtime)) ;
479  assert (B_hh_evol.is_known(jtime)) ;
480 
481  // Computation of h^{ij} to ensure det tgam_{ij} = det f_{ij} :
482 
483  hh_det_one(jtime, par_bc, par_mat) ;
484  }
485 
486  return hh_evol[jtime] ;
487 
488 }
489 
490 
492 
493  if( !(trk_evol.is_known(jtime)) ) {
494 
495  Scalar resu(ff.get_mp()) ;
496  resu.set_etat_zero() ;
497 
498  trk_evol.update(resu, jtime, the_time[jtime]) ;
499 
500  }
501 
502  return trk_evol[jtime] ;
503 
504 }
505 
506 
508 
509  if (p_hdirac == 0x0) {
510  p_hdirac = new Vector(ff.get_mp(), CON, ff.get_triad() ) ;
512  }
513 
514  return *p_hdirac ;
515 
516 }
517 
518 
519 
520 
521  //-----------------------------------//
522  // Update of fields from this class //
523  //-----------------------------------//
524 
526 
527  if (!( A_hh_evol.is_known(jtime) ) ) {
528  assert( hh_evol.is_known(jtime) ) ;
529 
530  A_hh_evol.update( hh_evol[jtime].compute_A(true), jtime, the_time[jtime] ) ;
531 
532  }
533  return A_hh_evol[jtime] ;
534 }
535 
537 
538  if (!( B_hh_evol.is_known(jtime) ) ) {
539  assert( hh_evol.is_known(jtime) ) ;
540 
541  B_hh_evol.update( hh_evol[jtime].compute_tilde_B_tt(true), jtime, the_time[jtime] ) ;
542 
543  }
544  return B_hh_evol[jtime] ;
545 }
546 
548 
549  if( !(trh_evol.is_known(jtime)) ) {
550 
551  // Computation of tr(h) to ensure det tgam_{ij} = det f_{ij} :
552  hh_det_one(jtime) ;
553 
554  }
555 
556  return trh_evol[jtime] ;
557 
558 }
559 
560 
561 
562  //------------------//
563  // output //
564  //------------------//
565 
566 ostream& Tslice_dirac_max::operator>>(ostream& flux) const {
567 
569 
570  flux << "Dirac gauge and maximal slicing" << '\n' ;
571 
572  if (A_hh_evol.is_known(jtime)) {
573  maxabs( A_hh_evol[jtime], "A_hh", flux) ;
574  }
575  if (B_hh_evol.is_known(jtime)) {
576  maxabs( B_hh_evol[jtime], "B_hh", flux) ;
577  }
578  if (trh_evol.is_known(jtime)) {
579  maxabs( trh_evol[jtime], "tr h", flux) ;
580  }
581 
582  return flux ;
583 
584 }
585 
586 
587 void Tslice_dirac_max::sauve(FILE* fich, bool partial_save) const {
588 
589  if (partial_save) {
590  cout <<
591  "Tslice_dirac_max::sauve : the partial_save case is not ready yet !"
592  << endl ;
593  abort() ;
594  }
595 
596  // Writing of quantities common to all derived classes of Time_slice_conf
597  // ----------------------------------------------------------------------
598 
599  Time_slice_conf::sauve(fich, true) ;
600 
601  // Writing of the other fields
602  // ---------------------------
603 
604  int jmin = jtime - depth + 1 ;
605 
606  // h^{ij}
607  assert( hh_evol.is_known(jtime) ) ;
608  for (int j=jmin; j<=jtime; j++) {
609  int indicator = (hh_evol.is_known(j)) ? 1 : 0 ;
610  fwrite_be(&indicator, sizeof(int), 1, fich) ;
611  if (indicator == 1) hh_evol[j].sauve(fich) ;
612  }
613 
614  // A_hh
615  A_hh() ; // forces the update at the current time step
616  for (int j=jmin; j<=jtime; j++) {
617  int indicator = (A_hh_evol.is_known(j)) ? 1 : 0 ;
618  fwrite_be(&indicator, sizeof(int), 1, fich) ;
619  if (indicator == 1) A_hh_evol[j].sauve(fich) ;
620  }
621 
622  // B_hh
623  B_hh() ; // forces the update at the current time step
624  for (int j=jmin; j<=jtime; j++) {
625  int indicator = (B_hh_evol.is_known(j)) ? 1 : 0 ;
626  fwrite_be(&indicator, sizeof(int), 1, fich) ;
627  if (indicator == 1) B_hh_evol[j].sauve(fich) ;
628  }
629 }
630 
631 
632 
633 
634 
635 
636 
637 }
Lorene::Time_slice_conf::del_deriv
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: time_slice_conf.C:347
Lorene::Time_slice_conf::hh_evol
Evolution_std< Sym_tensor > hh_evol
Values at successive time steps of the components .
Definition: time_slice.h:530
Lorene::Time_slice_conf::operator=
void operator=(const Time_slice_conf &)
Assignment to another Time_slice_conf.
Definition: time_slice_conf.C:376
Lorene::Tslice_dirac_max::set_trh
virtual void set_trh(const Scalar &trh_in)
Sets the trace, with respect to the flat metric ff , of .
Definition: tslice_dirac_max.C:437
Lorene::Evolution_std::update
virtual void update(const TyT &new_value, int j, double time_j)
Sets a new value at a given time step.
Definition: evolution_std.C:171
Lorene::Sym_tensor_trans::the_trace
const Scalar & the_trace() const
Returns the trace of the tensor with respect to metric *met_div.
Definition: sym_tensor_trans.C:270
Lorene::Sym_tensor::set_longit_trans
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.
Definition: sym_tensor_decomp.C:88
Lorene::Time_slice::depth
int depth
Number of stored time slices.
Definition: time_slice.h:179
Lorene::Sym_tensor
Class intended to describe valence-2 symmetric tensors.
Definition: sym_tensor.h:223
Lorene::Tslice_dirac_max::initial_data_cts
virtual void initial_data_cts(const Sym_tensor &uu, const Scalar &trk_in, const Scalar &trk_point, double pdt, double precis=1.e-12, int method_poisson_vect=6, const char *graph_device=0x0, const Scalar *ener_dens=0x0, const Vector *mom_dens=0x0, const Scalar *trace_stress=0x0)
Computes valid initial data by solving the constraint equations in the conformal thin-sandwich approa...
Definition: tslice_dirac_max.C:352
Lorene::Tslice_dirac_max::source_B_hata_evol
Evolution_std< Scalar > source_B_hata_evol
The potential of the source of equation for .
Definition: time_slice.h:1007
Lorene::Time_slice_conf::p_tgamma
Metric * p_tgamma
Pointer on the conformal metric at the current time step (jtime)
Definition: time_slice.h:560
Lorene::Tslice_dirac_max::set_hh
virtual void set_hh(const Sym_tensor &hh_in)
Sets the deviation of the conformal metric from the flat metric : .
Definition: tslice_dirac_max.C:336
Lorene::Time_slice::beta_evol
Evolution_std< Vector > beta_evol
Values at successive time steps of the shift vector .
Definition: time_slice.h:219
Lorene::Tslice_dirac_max::hdirac
virtual const Vector & hdirac() const
Vector which vanishes in Dirac gauge.
Definition: tslice_dirac_max.C:507
Lorene::Star_rot_Dirac
Class for relativistic rotating stars in Dirac gauge and maximal slicing.
Definition: star_rot_dirac.h:46
Lorene::Sym_tensor_tt
Transverse and traceless symmetric tensors of rank 2.
Definition: sym_tensor.h:938
Lorene::exp
Cmp exp(const Cmp &)
Exponential.
Definition: cmp_math.C:270
Lorene
Lorene prototypes.
Definition: app_hor.h:64
Lorene::Tslice_dirac_max
Spacelike time slice of a 3+1 spacetime with conformal decomposition in the maximal slicing and Dirac...
Definition: time_slice.h:968
Lorene::Time_slice::n_evol
Evolution_std< Scalar > n_evol
Values at successive time steps of the lapse function N.
Definition: time_slice.h:216
Lorene::Time_slice::adm_mass_evol
Evolution_full< Tbl > adm_mass_evol
ADM mass at each time step, since the creation of the slice.
Definition: time_slice.h:233
Lorene::Tensor::set_etat_zero
virtual void set_etat_zero()
Sets the logical state of all components to ETATZERO (zero state).
Definition: tensor.C:497
Lorene::Map::get_bvect_spher
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
Definition: map.h:783
Lorene::Map::get_mg
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:765
Lorene::Map::flat_met_spher
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:321
Lorene::Time_slice_conf::B_hata_evol
Evolution_std< Scalar > B_hata_evol
Potential associated with the symmetric tensor .
Definition: time_slice.h:552
Lorene::Mg3d::get_type_r
int get_type_r(int l) const
Returns the type of sampling in the radial direction in domain no.
Definition: grilles.h:474
Lorene::Metric_flat
Flat metric for tensor calculation.
Definition: metric.h:261
Lorene::Scalar
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:387
Lorene::Tensor::get_mp
const Map & get_mp() const
Returns the mapping.
Definition: tensor.h:861
Lorene::Tensor::annule_domain
void annule_domain(int l)
Sets the Tensor to zero in a given domain.
Definition: tensor.C:666
Lorene::Scalar::dec_dzpuis
virtual void dec_dzpuis(int dec=1)
Decreases by dec units the value of dzpuis and changes accordingly the values of the Scalar in the co...
Definition: scalar_r_manip.C:418
Lorene::Time_slice_conf
Spacelike time slice of a 3+1 spacetime with conformal decomposition.
Definition: time_slice.h:498
Lorene::Tslice_dirac_max::set_khi_mu
virtual void set_khi_mu(const Scalar &khi_in, const Scalar &mu_in)
Sets the potentials and of the TT part of (see the documentation of Sym_tensor_tt for details).
Definition: tslice_dirac_max.C:408
Lorene::Time_slice::gam_dd_evol
Evolution_std< Sym_tensor > gam_dd_evol
Values at successive time steps of the covariant components of the induced metric .
Definition: time_slice.h:198
Lorene::Tslice_dirac_max::source_B_hh_evol
Evolution_std< Scalar > source_B_hh_evol
The potential of the source of equation for .
Definition: time_slice.h:995
Lorene::Time_slice::k_dd_evol
Evolution_std< Sym_tensor > k_dd_evol
Values at successive time steps of the covariant components of the extrinsic curvature tensor .
Definition: time_slice.h:208
Lorene::Time_slice_conf::psi
virtual const Scalar & psi() const
Conformal factor relating the physical metric to the conformal one: .
Definition: time_slice_conf.C:693
Lorene::Sym_tensor_trans::trace_from_det_one
void trace_from_det_one(const Sym_tensor_tt &htt, double precis=1.e-14, int it_max=100)
Assigns the derived member p_tt and computes the trace so that *this + the flat metric has a determin...
Definition: sym_tensor_trans.C:315
Lorene::Time_slice_conf::set_hh
virtual void set_hh(const Sym_tensor &hh_in)
Sets the deviation of the conformal metric from the flat metric : .
Definition: time_slice_conf.C:489
Lorene::Tslice_dirac_max::source_A_hh_evol
Evolution_std< Scalar > source_A_hh_evol
The A potential of the source of equation for .
Definition: time_slice.h:989
Lorene::Tslice_dirac_max::source_A_hata_evol
Evolution_std< Scalar > source_A_hata_evol
The potential A of the source of equation for .
Definition: time_slice.h:1001
Lorene::Time_slice_conf::psi_evol
Evolution_std< Scalar > psi_evol
Values at successive time steps of the conformal factor relating the physical metric to the conform...
Definition: time_slice.h:517
Lorene::Mg3d::get_nzone
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:448
Lorene::Time_slice::operator>>
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator<<).
Definition: time_slice.C:411
Lorene::fread_be
int fread_be(int *aa, int size, int nb, FILE *fich)
Reads integer(s) from a binary file according to the big endian convention.
Definition: fread_be.C:69
Lorene::Star_rot_Dirac::get_aa
const Sym_tensor get_aa() const
Returns .
Definition: star_rot_dirac.h:229
Lorene::Time_slice::the_time
Evolution_std< double > the_time
Time label of each slice.
Definition: time_slice.h:193
Lorene::maxabs
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.
Definition: tensor_calculus_ext.C:606
Lorene::Scalar::set_etat_zero
virtual void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition: scalar.C:324
Lorene::Time_slice_conf::psi4
const Scalar & psi4() const
Factor at the current time step (jtime ).
Definition: time_slice_conf.C:707
Lorene::Tslice_dirac_max::sauve
virtual void sauve(FILE *fich, bool partial_save) const
Total or partial saves in a binary file.
Definition: tslice_dirac_max.C:587
Lorene::Scalar::set_dzpuis
void set_dzpuis(int)
Modifies the dzpuis flag.
Definition: scalar.C:808
Lorene::Vector
Tensor field of valence 1.
Definition: vector.h:188
Lorene::Tslice_dirac_max::~Tslice_dirac_max
virtual ~Tslice_dirac_max()
Destructor.
Definition: tslice_dirac_max.C:314
Lorene::Time_slice_conf::A_hata_evol
Evolution_std< Scalar > A_hata_evol
Potential A associated with the symmetric tensor .
Definition: time_slice.h:547
Lorene::Scalar::get_etat
int get_etat() const
Returns the logical state ETATNONDEF (undefined), ETATZERO (null) or ETATQCQ (ordinary).
Definition: scalar.h:554
Lorene::Tslice_dirac_max::B_hh
virtual const Scalar & B_hh() const
Returns the potential of .
Definition: tslice_dirac_max.C:536
Lorene::Tslice_dirac_max::compute_sources
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...
Definition: tslice_dirac_max_setAB.C:219
Lorene::Time_slice_conf::B_hata
virtual const Scalar & B_hata() const
Returns the potential of .
Definition: time_slice_conf.C:678
Lorene::Time_slice::trk_evol
Evolution_std< Scalar > trk_evol
Values at successive time steps of the trace K of the extrinsic curvature.
Definition: time_slice.h:224
Lorene::Time_slice::gam_uu_evol
Evolution_std< Sym_tensor > gam_uu_evol
Values at successive time steps of the contravariant components of the induced metric .
Definition: time_slice.h:203
Lorene::Time_slice_conf::p_hdirac
Vector * p_hdirac
Pointer on the vector (which vanishes in Dirac gauge), at the current time step (jtime).
Definition: time_slice.h:571
Lorene::Scalar::std_spectral_base
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field.
Definition: scalar.C:784
Lorene::Time_slice::jtime
int jtime
Time step index of the latest slice.
Definition: time_slice.h:190
Lorene::Time_slice_conf::A_hata
virtual const Scalar & A_hata() const
Returns the potential A of .
Definition: time_slice_conf.C:664
Lorene::Param
Parameter storage.
Definition: param.h:125
Lorene::Tslice_dirac_max::trh
virtual const Scalar & trh() const
Computes the trace h, with respect to the flat metric ff , of .
Definition: tslice_dirac_max.C:547
Lorene::Sym_tensor_tt::set_khi_mu
void set_khi_mu(const Scalar &khi_i, const Scalar &mu_i, int dzp=0, Param *par1=0x0, Param *par2=0x0, Param *par3=0x0)
Sets the component , as well as the angular potential (see member p_khi and p_mu ).
Definition: sym_tensor_tt_etamu.C:288
Lorene::fwrite_be
int fwrite_be(const int *aa, int size, int nb, FILE *fich)
Writes integer(s) into a binary file according to the big endian convention.
Definition: fwrite_be.C:70
Lorene::Star_rot_Dirac::get_ln_psi
const Scalar & get_ln_psi() const
Returns .
Definition: star_rot_dirac.h:199
Lorene::Sym_tensor_trans
Transverse symmetric tensors of rank 2.
Definition: sym_tensor.h:608
Lorene::Tslice_dirac_max::Tslice_dirac_max
Tslice_dirac_max(const Scalar &lapse_in, const Vector &shift_in, const Metric_flat &ff_in, const Scalar &psi_in, const Sym_tensor_trans &hh_in, const Sym_tensor &hata_in, int depth_in=3)
Constructor from conformal decomposition.
Definition: tslice_dirac_max.C:149
Lorene::Time_slice::p_gamma
Metric * p_gamma
Pointer on the induced metric at the current time step (jtime)
Definition: time_slice.h:239
Lorene::Tslice_dirac_max::B_hh_evol
Evolution_std< Scalar > B_hh_evol
The potential of .
Definition: time_slice.h:983
Lorene::Tslice_dirac_max::A_hh_evol
Evolution_std< Scalar > A_hh_evol
The A potential of .
Definition: time_slice.h:977
Lorene::Tslice_dirac_max::hh
virtual const Sym_tensor & hh(Param *par_bc=0x0, Param *par_mat=0x0) const
Deviation of the conformal metric from the flat metric : .
Definition: tslice_dirac_max.C:474
Lorene::Tslice_dirac_max::A_hh
virtual const Scalar & A_hh() const
Returns the potential A of .
Definition: tslice_dirac_max.C:525
Lorene::Map
Base class for coordinate mappings.
Definition: map.h:670
Lorene::Tslice_dirac_max::hh_det_one
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 .
Definition: tslice_dirac_max_setAB.C:109
Lorene::Tslice_dirac_max::trk
virtual const Scalar & trk() const
Trace K of the extrinsic curvature at the current time step (jtime ).
Definition: tslice_dirac_max.C:491
Lorene::Time_slice_conf::hata_evol
Evolution_std< Sym_tensor > hata_evol
Values at successive time steps of the components .
Definition: time_slice.h:542
Lorene::Base_vect
Vectorial bases (triads) with respect to which the tensorial components are defined.
Definition: base_vect.h:105
Lorene::Time_slice::sauve
virtual void sauve(FILE *fich, bool partial_save) const
Total or partial saves in a binary file.
Definition: time_slice.C:507