LORENE
tslice_check_einstein.C
1/*
2 * Methods of class Time_slice to check Einstein equation solutions.
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
29
30/*
31 * $Id: tslice_check_einstein.C,v 1.11 2016/12/05 16:18:19 j_novak Exp $
32 * $Log: tslice_check_einstein.C,v $
33 * Revision 1.11 2016/12/05 16:18:19 j_novak
34 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
35 *
36 * Revision 1.10 2014/10/13 08:53:47 j_novak
37 * Lorene classes and functions now belong to the namespace Lorene.
38 *
39 * Revision 1.9 2014/10/06 15:13:22 j_novak
40 * Modified #include directives to use c++ syntax.
41 *
42 * Revision 1.8 2010/10/20 07:58:09 j_novak
43 * Better implementation of the explicit time-integration. Not fully-tested yet.
44 *
45 * Revision 1.7 2007/11/06 12:10:56 j_novak
46 * Corrected some mistakes.
47 *
48 * Revision 1.6 2007/11/06 11:53:34 j_novak
49 * Use of contravariant version of the 3+1 equations.
50 *
51 * Revision 1.5 2007/06/28 14:40:36 j_novak
52 * Dynamical check: the fields in the last domain are set to zero to avoid dzpuis problems
53 *
54 * Revision 1.4 2004/05/12 15:24:20 e_gourgoulhon
55 * Reorganized the #include 's, taking into account that
56 * time_slice.h contains now an #include "metric.h".
57 *
58 * Revision 1.3 2004/04/07 07:58:21 e_gourgoulhon
59 * Constructor as Minkowski slice: added call to std_spectral_base()
60 * after setting the lapse to 1.
61 *
62 * Revision 1.2 2004/04/05 12:38:45 j_novak
63 * Minor modifs to prevent some warnings.
64 *
65 * Revision 1.1 2004/04/05 11:54:20 j_novak
66 * First operational (but not tested!) version of checks of Eintein equations.
67 *
68 *
69 * $Header: /cvsroot/Lorene/C++/Source/Time_slice/tslice_check_einstein.C,v 1.11 2016/12/05 16:18:19 j_novak Exp $
70 *
71 */
72
73// C headers
74#include <cstdlib>
75#include <cassert>
76
77// Lorene headers
78#include "time_slice.h"
79#include "unites.h"
80
81namespace Lorene {
83 ostream& ost, bool verb) const {
84 using namespace Unites ;
85
86 bool vacuum = ( energy_density == 0x0 ) ;
87
88 Scalar field = trk() * trk() - contract( k_uu(), 0, 1, k_dd(), 0, 1 ) ;
89
90 field.dec_dzpuis() ; // dzpuis: 4 -> 3
91
92 field += gam().ricci_scal() ;
93
94 const Scalar* matter ;
95 if (vacuum)
96 matter = new Scalar (0*field) ;
97 else
98 matter = energy_density ;
99
100 if (verb) ost << endl ;
101 const char* comment = 0x0 ;
102 if (verb) comment = "Check of the Hamiltonian constraint" ;
103 Tbl resu = maxabs(field - (4*qpig) * (*matter), comment, ost,verb ) ;
104 if (verb) ost << endl ;
105
106 if (vacuum) delete matter ;
107
108 return resu ;
109
110}
111
113 ostream& ost, bool verb) const {
114 using namespace Unites ;
115
116 bool vacuum = ( momentum_density == 0x0 ) ;
117
118 Vector field = k_uu().divergence(gam()) - trk().derive_con(gam()) ;
119
120
121 const Vector* matter ;
122 if (vacuum)
123 matter = new Vector (0*field) ;
124 else {
125 assert (momentum_density->get_index_type(0) == CON) ;
126 matter = momentum_density ;
127 }
128
129 if (verb) ost << endl ;
130 const char* comment = 0x0 ;
131 if (verb) comment = "Check of the momentum constraint" ;
132 Tbl resu = maxabs(field - (2*qpig) * (*matter), comment, ost, verb ) ;
133
134 if (verb) ost << endl ;
135
136 if (vacuum) delete matter ;
137
138 return resu ;
139
140}
141
143 const Scalar* energy_density,
144 ostream& ost, bool verb) const {
145
146 using namespace Unites ;
147
148 bool vacuum = ( ( strain_tensor == 0x0 ) && ( energy_density == 0x0 ) ) ;
149
150 Sym_tensor dyn_lhs = (k_dd_evol.time_derive(jtime, scheme_order)).up_down(gam()) ;
151 int nz = dyn_lhs.get_mp().get_mg()->get_nzone() ;
152 dyn_lhs.annule_domain(nz-1) ;
153 dyn_lhs = dyn_lhs - k_dd().derive_lie(beta()).up_down(gam()) ;
154 dyn_lhs.annule_domain(nz-1) ;
155
156 const Sym_tensor* matter ;
157 if (vacuum)
158 matter = new Sym_tensor(0 * dyn_lhs) ;
159 else {
160 const Scalar* ener = 0x0 ;
161 const Sym_tensor* sij = 0x0 ;
162 bool new_e = false ;
163 bool new_s = false ;
164 if (energy_density == 0x0) {
165 ener = new Scalar ( 0 * nn() ) ;
166 new_e = true ;
167 }
168 else {
169 ener = energy_density ;
170 if (strain_tensor == 0x0) {
171 sij = new Sym_tensor(0 * dyn_lhs) ;
172 new_s = true ;
173 }
174 else
175 sij = strain_tensor ;
176 }
177 matter = new Sym_tensor( (sij->trace(gam()) - *ener)*gam().con()
178 - 2*(*sij) ) ;
179 if (new_e) delete ener ;
180 if (new_s) delete sij ;
181 }
182
183 Sym_tensor dyn_rhs = nn()*( (gam().ricci().up_down(gam()) + trk()*k_uu()
184 + qpig * (*matter))) ;
185 dyn_rhs.annule_domain(nz-1) ;
186 dyn_rhs = dyn_rhs - 2*nn()*contract(k_uu(), 1, k_dd().up(0, gam()), 1) ;
187 dyn_rhs.annule_domain(nz-1) ;
188 dyn_rhs = dyn_rhs - nn().derive_con(gam()).derive_con(gam()) ;
189 dyn_rhs.annule_domain(nz-1) ;
190
191 if (verb) ost << endl ;
192 const char* comment = 0x0 ;
193 if (verb) comment = "Check of the dynamical equations" ;
194 Tbl resu_tmp = maxabs(dyn_lhs - dyn_rhs, comment, ost, verb) ;
195
196 if (verb) ost << endl ;
197
198 if (vacuum) delete matter ;
199
200 Tbl resu(4, 4, resu_tmp.get_dim(0)) ;
201 resu.annule_hard() ; //## To initialize resu(0,0,...) which
202 // does not correspond to a valid component
203
204 for (int i=1; i<=3; i++) {
205 for (int j=1; j<=i; j++) {
206 Itbl idx(2) ;
207 idx.set(0) = i ;
208 idx.set(1) = j ;
209 int pos = dyn_lhs.position(idx) ;
210 assert (dyn_rhs.position(idx) == pos) ;
211
212 for (int lz=0; lz<resu_tmp.get_dim(0); lz++)
213 resu.set(i,j,lz) = resu_tmp(pos, lz) ;
214 }
215 }
216
217 return resu ;
218
219}
220}
Basic integer array class.
Definition itbl.h:122
int & set(int i)
Read/write of a particular element (index i ) (1D case).
Definition itbl.h:247
const Sym_tensor & ricci() const
Returns the Ricci tensor (given by the Connection p_connect ).
Definition metric.C:341
virtual const Scalar & ricci_scal() const
Returns the Ricci scalar.
Definition metric.C:353
Tensor field of valence 0 (or component of a tensorial field).
Definition scalar.h:393
const Vector & derive_con(const Metric &gam) const
Returns the "contravariant" derivative of *this with respect to some metric , by raising the index of...
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...
Class intended to describe valence-2 symmetric tensors.
Definition sym_tensor.h:226
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
const Vector & divergence(const Metric &) const
Returns the divergence of this with respect to a Metric .
Definition sym_tensor.C:352
Basic array class.
Definition tbl.h:161
void annule_hard()
Sets the Tbl to zero in a hard way.
Definition tbl.C:375
double & set(int i)
Read/write of a particular element (index i) (1D case).
Definition tbl.h:281
int get_dim(int i) const
Gives the i-th dimension (ie dim.dim[i]).
Definition tbl.h:403
virtual const Sym_tensor & k_dd() const
Extrinsic curvature tensor (covariant components ) at the current time step (jtime ).
int jtime
Time step index of the latest slice.
Definition time_slice.h:193
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:211
virtual const Scalar & trk() const
Trace K of the extrinsic curvature at the current time step (jtime ).
virtual const Scalar & nn() const
Lapse function N at the current time step (jtime ).
Tbl check_hamiltonian_constraint(const Scalar *energy_density=0x0, ostream &ost=cout, bool verb=true) const
Checks the level at which the hamiltonian constraint is verified.
int scheme_order
Order of the finite-differences scheme for the computation of time derivatives.
Definition time_slice.h:190
virtual const Vector & beta() const
shift vector at the current time step (jtime )
const Metric & gam() const
Induced metric at the current time step (jtime ).
virtual const Sym_tensor & k_uu() const
Extrinsic curvature tensor (contravariant components ) at the current time step (jtime ).
Tbl check_dynamical_equations(const Sym_tensor *strain_tensor=0x0, const Scalar *energy_density=0x0, ostream &ost=cout, bool verb=true) const
Checks the level at which the dynamical equations are verified.
Tbl check_momentum_constraint(const Vector *momentum_density=0x0, ostream &ost=cout, bool verb=true) const
Checks the level at which the momentum constraints are verified.
Tensor field of valence 1.
Definition vector.h:188
const Map & get_mp() const
Returns the mapping.
Definition tensor.h:874
int get_index_type(int i) const
Gives the type (covariant or contravariant) of the index number i .
Definition tensor.h:899
Tensor up_down(const Metric &gam) const
Computes a new tensor by raising or lowering all the indices of *this .
const Tensor & derive_con(const Metric &gam) const
Returns the "contravariant" derivative of this with respect to some metric , by raising the last inde...
Definition tensor.C:1023
void annule_domain(int l)
Sets the Tensor to zero in a given domain.
Definition tensor.C:675
virtual int position(const Itbl &ind) const
Returns the position in the array cmp of a component given by its indices.
Definition tensor_sym.C:248
Tensor trace(int ind1, int ind2) const
Trace on two different type indices.
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
Standard units of space, time and mass.