KADATH
metric_conf_factor.cpp
1 /*
2  Copyright 2018 Philippe Grandclement
3 
4  This file is part of Kadath.
5 
6  Kadath 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 3 of the License, or
9  (at your option) any later version.
10 
11  Kadath 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 Kadath. If not, see <http://www.gnu.org/licenses/>.
18 */
19 
20 #include "space.hpp"
21 #include "tensor.hpp"
22 #include "metric_nophi.hpp"
23 #include "term_eq.hpp"
24 #include "scalar.hpp"
25 #include "tensor_impl.hpp"
26 #include "system_of_eqs.hpp"
27 #include "metric_tensor.hpp"
28 #include "name_tools.hpp"
29 
30 namespace Kadath {
32  Metric(met.get_space()), p_met(&met), basis(met.get_basis()), fmet(met.get_space(), basis), conformal(ome), grad_conf (espace, COV, basis) {
33  type_tensor = met.get_type() ;
34 
35  // Compute the flat gradient of conformal factor (trick pass by Term_eq even if not needed)
36  for (int d=0 ; d<espace.get_nbr_domains() ; d++) {
37  Term_eq res (fmet.derive (COV, ' ', Term_eq (d, conformal))) ;
38  for (int i=1 ; i<=3 ; i++)
39  grad_conf.set(i).set_domain(d) = (*res.val_t)(i)(d) ;
40  }
41 }
42 
43 
45  Metric (so), p_met(so.p_met), basis(so.basis), fmet(so.fmet), conformal(so.conformal), grad_conf(so.grad_conf), place_syst(so.place_syst) {
46 }
47 
48 Metric_conf_factor::~Metric_conf_factor() {
49 }
50 
51 
52 int Metric_conf_factor::give_type(int dd) const {
53  return basis.get_basis(dd) ;
54 }
55 
56 void Metric_conf_factor::compute_cov (int dd) const {
57 
58 
59  int place = place_syst + (dd-syst->dom_min) ;
60  // Right storage : simple copy.
61  if (type_tensor==COV) {
62 
63  if (p_met_cov[dd]==0x0)
64  p_met_cov[dd] = new Term_eq(*syst->term[place]) ;
65  else
66  *p_met_cov[dd] = Term_eq(*syst->term[place]) ;
67  }
68  else {
69  //Need to work component by components...
70  bool doder = (syst->term[place]->der_t==0x0) ? false : true ;
71 
72  Term_eq** res = new Term_eq* [p_met->get_n_comp()] ;
73 
74  Scalar val (espace) ;
75  Scalar der (espace) ;
76 
77  Val_domain detval (espace.get_domain(dd)) ;
78  detval = (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->val_t)(2,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
79  + (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)
80  + (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
81  - (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,2)(dd)
82  - (*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(1,1)(dd)
83  - (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd) ;
84  Val_domain detder (espace.get_domain(dd)) ;
85  if (doder)
86  detder = (*syst->term[place]->der_t)(1,1)(dd)*(*syst->term[place]->val_t)(2,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
87  + (*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)
88  + (*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
89  - (*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,2)(dd)
90  - (*syst->term[place]->der_t)(2,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(1,1)(dd)
91  - (*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
92  + (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->der_t)(2,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
93  + (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(2,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)
94  + (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
95  - (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,2)(dd)
96  - (*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->der_t)(2,3)(dd)*(*syst->term[place]->val_t)(1,1)(dd)
97  - (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
98  + (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->val_t)(2,2)(dd)*(*syst->term[place]->der_t)(3,3)(dd)
99  + (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->der_t)(1,3)(dd)
100  + (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(2,3)(dd)
101  - (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(2,2)(dd)
102  - (*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->der_t)(1,1)(dd)
103  - (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(3,3)(dd) ;
104 
105  Val_domain cmpval(espace.get_domain(dd)) ;
106  Val_domain cmpder(espace.get_domain(dd)) ;
107 
108  // Compo 1 1
109  cmpval = (*syst->term[place]->val_t)(2,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
110  -(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd) ;
111  val.set_domain(dd) = cmpval / detval;
112  if (doder) {
113  cmpder = (*syst->term[place]->der_t)(2,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
114  + (*syst->term[place]->val_t)(2,2)(dd)*(*syst->term[place]->der_t)(3,3)(dd)
115  - (*syst->term[place]->der_t)(2,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
116  - (*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->der_t)(2,3)(dd) ;
117  der.set_domain(dd) = cmpder / detval - cmpval * detder / detval/detval;
118  res[0] = new Term_eq (dd, val, der) ;
119  }
120  else {
121  res[0] = new Term_eq (dd, val) ;
122  }
123 
124  // Compo 1 2
125  cmpval = (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
126  -(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd) ;
127  val.set_domain(dd) = cmpval / detval ;
128  if (doder) {
129  cmpder = (*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
130  + (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(2,3)(dd)
131  - (*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
132  - (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(3,3)(dd) ;
133  der.set_domain(dd) = cmpder / detval - cmpval * detder / detval/detval;
134  res[1] = new Term_eq (dd, val, der) ;
135  }
136  else {
137  res[1] = new Term_eq (dd, val) ;
138  }
139 
140  // Compo 1 3
141  cmpval = (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
142  -(*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,2)(dd) ;
143  val.set_domain(dd) = cmpval / detval ;
144  if (doder) {
145  cmpder = (*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
146  + (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(2,3)(dd)
147  - (*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,2)(dd)
148  - (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(2,2)(dd) ;
149  der.set_domain(dd) = cmpder / detval - cmpval * detder / detval/detval;
150  res[2] = new Term_eq (dd, val, der) ;
151  }
152  else {
153  res[2] = new Term_eq (dd, val) ;
154  }
155 
156  // Compo 2 2
157  cmpval = (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
158  -(*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd) ;
159  val.set_domain(dd) = cmpval / detval ;
160  if (doder) {
161  cmpder = (*syst->term[place]->der_t)(1,1)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
162  + (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->der_t)(3,3)(dd)
163  - (*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)
164  - (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(1,3)(dd) ;
165  der.set_domain(dd) = cmpder / detval - cmpval * detder / detval/detval;
166  res[3] = new Term_eq (dd, val, der) ;
167  }
168  else {
169  res[3] = new Term_eq (dd, val) ;
170  }
171 
172  // Compo 2 3
173  cmpval = (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,2)(dd)
174  -(*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->val_t)(2,3)(dd) ;
175  val.set_domain(dd) = cmpval / detval ;
176  if (doder) {
177  cmpder = (*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,2)(dd)
178  + (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(1,2)(dd)
179  - (*syst->term[place]->der_t)(1,1)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
180  - (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->der_t)(2,3)(dd) ;
181  der.set_domain(dd) = cmpder/ detval - cmpval * detder / detval/detval ;
182  res[4] = new Term_eq (dd, val, der) ;
183  }
184  else {
185  res[4] = new Term_eq (dd, val) ;
186  }
187 
188  // Compo 3 3
189  cmpval = (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->val_t)(2,2)(dd)
190  -(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(1,2)(dd) ;
191  val.set_domain(dd) = cmpval / detval ;
192  if (doder) {
193  cmpder = (*syst->term[place]->der_t)(1,1)(dd)*(*syst->term[place]->val_t)(2,2)(dd)
194  + (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->der_t)(2,2)(dd)
195  - (*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(1,2)(dd)
196  - (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(1,2)(dd) ;
197  der.set_domain(dd) = cmpder / detval - cmpval * detder / detval/detval;
198  res[5] = new Term_eq (dd, val, der) ;
199  }
200  else {
201  res[5] = new Term_eq (dd, val) ;
202  }
203 
204  // Value field :
205  Metric_tensor resval (espace, COV, basis) ;
206  resval.set(1,1) = res[0]->get_val_t() ;
207  resval.set(1,2) = res[1]->get_val_t() ;
208  resval.set(1,3) = res[2]->get_val_t() ;
209  resval.set(2,2) = res[3]->get_val_t() ;
210  resval.set(2,3) = res[4]->get_val_t() ;
211  resval.set(3,3) = res[5]->get_val_t() ;
212 
213  if (!doder) {
214  if (p_met_cov[dd]==0x0)
215  p_met_cov[dd] = new Term_eq(dd, resval) ;
216  else
217  *p_met_cov[dd] = Term_eq(dd, resval) ;
218  }
219  else {
220  // Der field :
221  Metric_tensor resder (espace, COV, basis) ;
222  resder.set(1,1) = res[0]->get_der_t() ;
223  resder.set(1,2) = res[1]->get_der_t() ;
224  resder.set(1,3) = res[2]->get_der_t() ;
225  resder.set(2,2) = res[3]->get_der_t() ;
226  resder.set(2,3) = res[4]->get_der_t() ;
227  resder.set(3,3) = res[5]->get_der_t() ;
228 
229  if (p_met_cov[dd]==0x0)
230  p_met_cov[dd] = new Term_eq(dd, resval, resder) ;
231  else
232  *p_met_cov[dd] = Term_eq(dd, resval, resder) ;
233  }
234  for (int i=0 ; i<p_met->get_n_comp() ; i++)
235  delete res [i] ;
236  delete [] res ;
237  }
238 }
239 
240 
241 void Metric_conf_factor::compute_con (int dd) const {
242 
243  int place = place_syst + (dd-syst->dom_min) ;
244  // Right storage : simple copy.
245  if (type_tensor==CON) {
246 
247  if (p_met_con[dd]==0x0)
248  p_met_con[dd] = new Term_eq(*syst->term[place]) ;
249  else
250  *p_met_con[dd] = Term_eq(*syst->term[place]) ;
251  }
252  else {
253  //Need to work component by components...
254  bool doder = (syst->term[place]->der_t==0x0) ? false : true ;
255  Term_eq** res = new Term_eq* [p_met->get_n_comp()] ;
256 
257  Scalar val (espace) ;
258  Scalar der (espace) ;
259 
260  Val_domain cmpval(espace.get_domain(dd)) ;
261  Val_domain cmpder(espace.get_domain(dd)) ;
262 
263  Val_domain detval (espace.get_domain(dd)) ;
264  detval = (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->val_t)(2,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
265  + (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)
266  + (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
267  - (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,2)(dd)
268  - (*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(1,1)(dd)
269  - (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd) ;
270  Val_domain detder (espace.get_domain(dd)) ;
271  if (doder)
272  detder = (*syst->term[place]->der_t)(1,1)(dd)*(*syst->term[place]->val_t)(2,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
273  + (*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)
274  + (*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
275  - (*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,2)(dd)
276  - (*syst->term[place]->der_t)(2,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(1,1)(dd)
277  - (*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
278  + (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->der_t)(2,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
279  + (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(2,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)
280  + (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
281  - (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,2)(dd)
282  - (*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->der_t)(2,3)(dd)*(*syst->term[place]->val_t)(1,1)(dd)
283  - (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
284  + (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->val_t)(2,2)(dd)*(*syst->term[place]->der_t)(3,3)(dd)
285  + (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->der_t)(1,3)(dd)
286  + (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(2,3)(dd)
287  - (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(2,2)(dd)
288  - (*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->der_t)(1,1)(dd)
289  - (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(3,3)(dd) ;
290 
291 
292  // Compo 1 1
293  cmpval = (*syst->term[place]->val_t)(2,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
294  -(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd) ;
295  val.set_domain(dd) = cmpval / detval;
296  if (doder) {
297  cmpder = (*syst->term[place]->der_t)(2,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
298  + (*syst->term[place]->val_t)(2,2)(dd)*(*syst->term[place]->der_t)(3,3)(dd)
299  - (*syst->term[place]->der_t)(2,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
300  - (*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->der_t)(2,3)(dd) ;
301  der.set_domain(dd) = cmpder / detval - cmpval * detder / detval/detval;
302  res[0] = new Term_eq (dd, val, der) ;
303  }
304  else {
305  res[0] = new Term_eq (dd, val) ;
306  }
307 
308  // Compo 1 2
309  cmpval = (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
310  -(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd) ;
311  val.set_domain(dd) = cmpval / detval;
312  if (doder) {
313  cmpder = (*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
314  + (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(2,3)(dd)
315  - (*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
316  - (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(3,3)(dd) ;
317  der.set_domain(dd) = cmpder / detval - cmpval * detder / detval/detval;
318  res[1] = new Term_eq (dd, val, der) ;
319  }
320  else {
321  res[1] = new Term_eq (dd, val) ;
322  }
323 
324  // Compo 1 3
325  cmpval = (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
326  -(*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,2)(dd) ;
327  val.set_domain(dd) = cmpval / detval;
328  if (doder) {
329  cmpder = (*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
330  + (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(2,3)(dd)
331  - (*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,2)(dd)
332  - (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(2,2)(dd) ;
333  der.set_domain(dd) = cmpder / detval - cmpval * detder / detval/detval;
334  res[2] = new Term_eq (dd, val, der) ;
335  }
336  else {
337  res[2] = new Term_eq (dd, val) ;
338  }
339 
340  // Compo 2 2
341  cmpval = (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
342  -(*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd) ;
343  val.set_domain(dd) = cmpval / detval;
344  if (doder) {
345  cmpder = (*syst->term[place]->der_t)(1,1)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
346  + (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->der_t)(3,3)(dd)
347  - (*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)
348  - (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(1,3)(dd) ;
349  der.set_domain(dd) = cmpder / detval - cmpval * detder / detval/detval;
350  res[3] = new Term_eq (dd, val, der) ;
351  }
352  else {
353  res[3] = new Term_eq (dd, val) ;
354  }
355 
356  // Compo 2 3
357  cmpval = (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,2)(dd)
358  -(*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->val_t)(2,3)(dd) ;
359  val.set_domain(dd) = cmpval / detval ;
360  if (doder) {
361  cmpder = (*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,2)(dd)
362  + (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(1,2)(dd)
363  - (*syst->term[place]->der_t)(1,1)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
364  - (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->der_t)(2,3)(dd) ;
365  der.set_domain(dd) = cmpder/ detval - cmpval * detder / detval/detval ;
366  res[4] = new Term_eq (dd, val, der) ;
367  }
368  else {
369  res[4] = new Term_eq (dd, val) ;
370  }
371 
372  // Compo 3 3
373  cmpval = (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->val_t)(2,2)(dd)
374  -(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(1,2)(dd) ;
375  val.set_domain(dd) = cmpval / detval;
376  if (doder) {
377  cmpder = (*syst->term[place]->der_t)(1,1)(dd)*(*syst->term[place]->val_t)(2,2)(dd)
378  + (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->der_t)(2,2)(dd)
379  - (*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(1,2)(dd)
380  - (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(1,2)(dd) ;
381  der.set_domain(dd) = cmpder / detval - cmpval * detder / detval/detval;
382  res[5] = new Term_eq (dd, val, der) ;
383  }
384  else {
385  res[5] = new Term_eq (dd, val) ;
386  }
387 
388  // Value field :
389  Metric_tensor resval (espace, CON, basis) ;
390  resval.set(1,1) = res[0]->get_val_t() ;
391  resval.set(1,2) = res[1]->get_val_t() ;
392  resval.set(1,3) = res[2]->get_val_t() ;
393  resval.set(2,2) = res[3]->get_val_t() ;
394  resval.set(2,3) = res[4]->get_val_t() ;
395  resval.set(3,3) = res[5]->get_val_t() ;
396 
397  if (!doder) {
398  if (p_met_con[dd]==0x0)
399  p_met_con[dd] = new Term_eq(dd, resval) ;
400  else
401  *p_met_con[dd] = Term_eq(dd, resval) ;
402  }
403  else {
404  // Der field :
405  Metric_tensor resder (espace, CON, basis) ;
406  resder.set(1,1) = res[0]->get_der_t() ;
407  resder.set(1,2) = res[1]->get_der_t() ;
408  resder.set(1,3) = res[2]->get_der_t() ;
409  resder.set(2,2) = res[3]->get_der_t() ;
410  resder.set(2,3) = res[4]->get_der_t() ;
411  resder.set(3,3) = res[5]->get_der_t() ;
412 
413  if (p_met_con[dd]==0x0)
414  p_met_con[dd] = new Term_eq(dd, resval, resder) ;
415  else
416  *p_met_con[dd] = Term_eq(dd, resval, resder) ;
417  }
418 
419  for (int i=0 ; i<p_met->get_n_comp() ; i++)
420  delete res [i] ;
421  delete [] res ;
422  }
423 }
424 
426  // Need both representation of the metric)
427  if (type_tensor == CON) {
428  if (p_met_con[dd]==0x0)
429  compute_con(dd) ;
430  if (p_met_cov[dd]==0x0)
431  compute_cov(dd) ;
432  }
433  else {
434  if (p_met_cov[dd]==0x0)
435  compute_cov(dd) ;
436  if (p_met_con[dd]==0x0)
437  compute_con(dd) ;
438  }
439 
440  Array<int> type_ind (3) ;
441  type_ind.set(0) = COV ; type_ind.set(1) = COV ; type_ind.set(2) = CON ;
442  Tensor res_val (espace, 3, type_ind, basis, 3) ;
443  Tensor res_der (espace, 3, type_ind, basis, 3) ;
444 
445  Term_eq flat_der (fmet.derive(COV, ' ', (*p_met_cov[dd]))) ;
446  bool doder = (flat_der.der_t==0x0) ? false : true ;
447 
448  Index pos (res_val) ;
449  do {
450 
451  Val_domain cmpval (espace.get_domain(dd)) ;
452  cmpval = 0 ;
453  for (int l=1 ; l<=espace.get_ndim() ; l++)
454  cmpval += conformal(dd) * (0.5*(*p_met_con[dd]->val_t)(pos(2)+1,l)(dd)*((*flat_der.val_t)(pos(0)+1, pos(1)+1, l)(dd) +
455  (*flat_der.val_t)(pos(1)+1, pos(0)+1, l)(dd) - (*flat_der.val_t)(l, pos(0)+1, pos(1)+1)(dd))) ;
456 
457  // Some corrections due to the conformal factor :
458  for (int l=1 ; l<=espace.get_ndim() ; l++)
459  cmpval += (*p_met_con[dd]->val_t)(pos(2)+1,l)(dd)*(*p_met_cov[dd]->val_t)(pos(0)+1,pos(1)+1)(dd)*grad_conf(l)(dd) ;
460  if (pos(1)==pos(2))
461  cmpval -= grad_conf(pos(0)+1)(dd) ;
462  if (pos(0)==pos(2))
463  cmpval -= grad_conf(pos(1)+1)(dd) ;
464 
465  res_val.set(pos).set_domain(dd) = cmpval ;
466 
467  if (doder) {
468  Val_domain cmpder (espace.get_domain(dd)) ;
469  cmpder = 0 ;
470  for (int l=1 ; l<=espace.get_ndim() ; l++)
471  cmpder += conformal(dd) * (0.5*(*p_met_con[dd]->der_t)(pos(2)+1,l)(dd)*((*flat_der.val_t)(pos(0)+1, pos(1)+1, l)(dd) +
472  (*flat_der.val_t)(pos(1)+1, pos(0)+1, l)(dd) - (*flat_der.val_t)(l, pos(0)+1, pos(1)+1)(dd))
473  + 0.5*(*p_met_con[dd]->val_t)(pos(2)+1,l)(dd)*((*flat_der.der_t)(pos(0)+1, pos(1)+1, l)(dd) +
474  (*flat_der.der_t)(pos(1)+1, pos(0)+1, l)(dd) - (*flat_der.der_t)(l, pos(0)+1, pos(1)+1)(dd))) ;
475 
476  // Some corrections due to the conformal factor :
477  for (int l=1 ; l<=espace.get_ndim() ; l++)
478  cmpder += ((*p_met_con[dd]->der_t)(pos(2)+1,l)(dd)*(*p_met_cov[dd]->val_t)(pos(0)+1,pos(1)+1)(dd)
479  +(*p_met_con[dd]->val_t)(pos(2)+1,l)(dd)*(*p_met_cov[dd]->der_t)(pos(0)+1,pos(1)+1)(dd)) *grad_conf(l)(dd) ;
480 
481 
482  res_der.set(pos).set_domain(dd) = cmpder ;
483  }
484  }
485  while (pos.inc()) ;
486 
487  if (!doder) {
488  if (p_christo[dd]==0x0)
489  p_christo[dd] = new Term_eq(dd, res_val) ;
490  else
491  *p_christo[dd] = Term_eq (dd, res_val) ;
492  }
493  else {
494  if (p_christo[dd]==0x0)
495  p_christo[dd] = new Term_eq(dd, res_val, res_der) ;
496  else
497  *p_christo[dd] = Term_eq(dd, res_val, res_der) ;
498  }
499 
500 }
501 
502 
504 
505  // Need christoffels
506  if (p_christo[dd]==0x0)
507  compute_christo(dd) ;
508 
509 
510  Array<int> indices (4) ;
511  indices.set(0) = CON ; indices.set(1) = COV ; indices.set(2) = COV ; indices.set(3) = COV ;
512  Tensor res_val (espace, 4, indices, basis, 3) ;
513  Tensor res_der (espace, 4, indices, basis, 3) ;
514 
515  Term_eq flat_der (fmet.derive(COV, ' ', (*p_christo[dd]))) ;
516  bool doder = (flat_der.der_t==0x0) ? false : true ;
517 
518  Index pos (res_val) ;
519  do {
520 
521  Val_domain cmpval (conformal(dd)*((*flat_der.val_t)(pos(2)+1, pos(1)+1,pos(3)+1,pos(0)+1)(dd)
522  - (*flat_der.val_t)(pos(3)+1, pos(1)+1, pos(2)+1, pos(0)+1)(dd))) ;
523 
524  // Corrections coming from the conformal part
525  cmpval -= (*p_christo[dd]->val_t)(pos(1)+1, pos(3)+1, pos(0)+1)(dd)*grad_conf(pos(2)+1)(dd) ;
526  cmpval += (*p_christo[dd]->val_t)(pos(1)+1, pos(2)+1, pos(0)+1)(dd)*grad_conf(pos(3)+1)(dd) ;
527 
528  for (int m=1 ; m<=espace.get_ndim() ; m++) {
529  cmpval += (*p_christo[dd]->val_t)(pos(2)+1,m, pos(0)+1)(dd)*(*p_christo[dd]->val_t)(pos(1)+1,pos(3)+1,m)(dd)
530  - (*p_christo[dd]->val_t)(pos(3)+1,m,pos(0)+1)(dd)*(*p_christo[dd]->val_t)(pos(1)+1,pos(2)+1,m)(dd) ;
531  }
532 
533 
534 
535  res_val.set(pos).set_domain(dd) = cmpval;
536  if (doder) {
537  Val_domain cmpder (conformal(dd)*((*flat_der.der_t)(pos(2)+1, pos(1)+1,pos(3)+1,pos(0)+1)(dd)
538  - (*flat_der.der_t)(pos(3)+1, pos(1)+1, pos(2)+1, pos(0)+1)(dd))) ;
539 
540  cmpder -= (*p_christo[dd]->der_t)(pos(1)+1, pos(3)+1, pos(0)+1)(dd)*grad_conf(pos(2)+1)(dd) ;
541  cmpder += (*p_christo[dd]->der_t)(pos(1)+1, pos(2)+1, pos(0)+1)(dd)*grad_conf(pos(3)+1)(dd) ;
542 
543 
544  for (int m=1 ; m<=espace.get_ndim() ; m++) {
545  cmpder += (*p_christo[dd]->der_t)(pos(2)+1,m, pos(0)+1)(dd)*(*p_christo[dd]->val_t)(pos(1)+1,pos(3)+1,m)(dd)
546  + (*p_christo[dd]->val_t)(pos(2)+1,m, pos(0)+1)(dd)*(*p_christo[dd]->der_t)(pos(1)+1,pos(3)+1,m)(dd)
547  - (*p_christo[dd]->der_t)(pos(3)+1,m,pos(0)+1)(dd)*(*p_christo[dd]->val_t)(pos(1)+1,pos(2)+1,m)(dd)
548  - (*p_christo[dd]->val_t)(pos(3)+1,m,pos(0)+1)(dd)*(*p_christo[dd]->der_t)(pos(1)+1,pos(2)+1,m)(dd);
549  }
550  res_der.set(pos).set_domain(dd) = cmpder ;
551  }
552  }
553  while (pos.inc()) ;
554 
555  if (!doder) {
556  if (p_riemann[dd]==0x0)
557  p_riemann[dd] = new Term_eq(dd, res_val) ;
558  else
559  *p_riemann[dd] = Term_eq (dd, res_val) ;
560  }
561  else {
562  if (p_riemann[dd]==0x0)
563  p_riemann[dd] = new Term_eq(dd, res_val, res_der) ;
564  else
565  *p_riemann[dd] = Term_eq(dd, res_val, res_der) ;
566  }
567 
568 }
569 
571  // Need christoffels
572  if (p_christo[dd]==0x0)
573  compute_christo(dd) ;
574 
575 
576  Array<int> indices (2) ;
577  indices.set(0) = COV ; indices.set(1) = COV ;
578  Tensor res_val (espace, 2, indices, basis, 3) ;
579  Tensor res_der (espace, 2, indices, basis, 3) ;
580 
581  Term_eq flat_der (fmet.derive(COV, ' ', (*p_christo[dd]))) ;
582  bool doder = (flat_der.der_t==0x0) ? false : true ;
583 
584  Index pos (res_val) ;
585  do {
586 
587  Val_domain cmpval (espace.get_domain(dd)) ;
588  cmpval = 0 ;
589 
590  for (int l=1 ; l<=espace.get_ndim() ; l++) {
591  cmpval += conformal(dd)*((*flat_der.val_t)(l, pos(1)+1,pos(0)+1,l)(dd)
592  - (*flat_der.val_t)(pos(1)+1, pos(0)+1,l ,l)(dd)) ;
593 
594  cmpval-=(*p_christo[dd]->val_t)(pos(0)+1, pos(1)+1, l)(dd)*grad_conf(l)(dd) ;
595  cmpval+=(*p_christo[dd]->val_t)(pos(0)+1, l, l)(dd)*grad_conf(pos(1)+1)(dd) ;
596 
597  for (int m=1 ; m<=espace.get_ndim() ; m++) {
598  cmpval += (*p_christo[dd]->val_t)(l,m, l)(dd)*(*p_christo[dd]->val_t)(pos(0)+1,pos(1)+1,m)(dd)
599  - (*p_christo[dd]->val_t)(pos(0)+1,l, m)(dd)*(*p_christo[dd]->val_t)(pos(1)+1,m,l)(dd) ;
600  }
601  }
602 
603  res_val.set(pos).set_domain(dd) = cmpval;
604  if (doder) {
605  Val_domain cmpder(espace.get_domain(dd)) ;
606  cmpder = 0 ;
607 
608  for (int l=1 ; l<=espace.get_ndim() ; l++) {
609  cmpder += conformal(dd)*((*flat_der.der_t)(l, pos(1)+1,pos(0)+1,l)(dd)
610  - (*flat_der.der_t)(pos(1)+1, pos(0)+1,l ,l)(dd)) ;
611 
612  cmpder-=(*p_christo[dd]->der_t)(pos(0)+1, pos(1)+1, l)(dd)*grad_conf(l)(dd) ;
613  cmpder+=(*p_christo[dd]->der_t)(pos(0)+1, l, l)(dd)*grad_conf(pos(1)+1)(dd) ;
614 
615 
616  for (int m=1 ; m<=espace.get_ndim() ; m++) {
617  cmpder += (*p_christo[dd]->der_t)(l,m, l)(dd)*(*p_christo[dd]->val_t)(pos(0)+1,pos(1)+1,m)(dd)
618  + (*p_christo[dd]->val_t)(l,m, l)(dd)*(*p_christo[dd]->der_t)(pos(0)+1,pos(1)+1,m)(dd)
619  - (*p_christo[dd]->der_t)(pos(0)+1,l, m)(dd)*(*p_christo[dd]->val_t)(pos(1)+1,m,l)(dd)
620  - (*p_christo[dd]->val_t)(pos(0)+1,l, m)(dd)*(*p_christo[dd]->der_t)(pos(1)+1,m,l)(dd);
621  }
622  }
623  res_der.set(pos).set_domain(dd) = cmpder ;
624  }
625  }
626  while (pos.inc()) ;
627 
628  if (!doder) {
629  if (p_ricci_tensor[dd]==0x0)
630  p_ricci_tensor[dd] = new Term_eq(dd, res_val) ;
631  else
632  *p_ricci_tensor[dd] = Term_eq (dd, res_val) ;
633  }
634  else {
635  if (p_ricci_tensor[dd]==0x0)
636  p_ricci_tensor[dd] = new Term_eq(dd, res_val, res_der) ;
637  else
638  *p_ricci_tensor[dd] = Term_eq(dd, res_val, res_der) ;
639  }
640 
641 
642  // Add flat metric part if needed
643  Term_eq ricci_flat (*fmet.give_ricci_tensor(dd)) ;
644  for (int i=1 ; i<=3 ; i++)
645  for (int j=i ; j<=3 ; j++) {
646  ricci_flat.set_val_t()->set(i,j).set_domain(dd) *= conformal(dd)*conformal(dd) ;
647  if (ricci_flat.der_t !=0x0)
648  ricci_flat.set_der_t()->set(i,j).set_domain(dd) *= conformal(dd)*conformal(dd) ;
649  }
650  *p_ricci_tensor[dd] = *p_ricci_tensor[dd] + ricci_flat ;
651 }
652 
653 Term_eq Metric_conf_factor::derive_flat (int type_der, char ind_der, const Term_eq& so) const {
654 
655  int dd = so.get_dom() ;
656 
657  if (p_met_con[dd]==0x0)
658  compute_con(dd) ;
659  if (p_met_cov[dd]==0x0)
660  compute_cov(dd) ;
661 
662 
663  // The partial derivative part
664  Term_eq res (fmet.derive_with_other (type_der, ind_der, so, this)) ;
665  return res ;
666 }
667 
668 
669 Term_eq Metric_conf_factor::derive (int type_der, char ind_der, const Term_eq& so) const {
670 
671  int dd = so.get_dom() ;
672 
673  if (p_christo[dd]==0x0)
674  compute_christo(dd) ;
675  // The partial derivative part
676  Term_eq res (conformal*fmet.derive_with_other (type_der, ind_der, so, this)) ;
677 
678  // Add the part containing the Christoffel :
679  //Must find a name for summation on Christofel :
680  bool found = false ;
681  int start = 97 ;
682  do {
683  bool same = false ;
684  if (ind_der==char(start))
685  same = true ;
686  for (int i=0 ; i<so.val_t->get_valence() ; i++)
687  if (so.val_t->get_name_ind()[i]==char(start))
688  same = true ;
689  if (!same)
690  found = true ;
691  else
692  start ++ ;
693  }
694  while ((!found) && (start<123)) ;
695  if (!found) {
696  cerr << "Trouble with indices in derive (you are not using tensors of order > 24, are you ?)" << endl ;
697  abort() ;
698  }
699  char name_sum = char(start) ;
700 
701  bool doder = ((so.der_t==0x0) || (p_christo[dd]->der_t==0x0)) ? false : true ;
702 
703  //loop on the components of so :
704  for (int cmp=0 ; cmp<so.val_t->get_valence() ; cmp ++) {
705 
706  int genre_indice = so.val_t->get_index_type(cmp) ;
707 
708  //Affecte names on the Christoffels :
709  p_christo[dd]->val_t->set_name_affected() ;
710  p_christo[dd]->val_t->set_name_ind(0, ind_der) ;
711  if (genre_indice==COV) {
712  p_christo[dd]->val_t->set_name_ind(1, so.val_t->get_name_ind()[cmp]) ;
713  p_christo[dd]->val_t->set_name_ind(2, name_sum) ;
714  }
715  else {
716  p_christo[dd]->val_t->set_name_ind(2, so.val_t->get_name_ind()[cmp]) ;
717  p_christo[dd]->val_t->set_name_ind(1, name_sum) ;
718  }
719 
720  if (doder) {
721  p_christo[dd]->der_t->set_name_affected() ;
722  p_christo[dd]->der_t->set_name_ind(0, ind_der) ;
723  if (genre_indice==COV) {
724  p_christo[dd]->der_t->set_name_ind(1, so.der_t->get_name_ind()[cmp]) ;
725  p_christo[dd]->der_t->set_name_ind(2, name_sum) ;
726  }
727  else {
728  p_christo[dd]->der_t->set_name_ind(2, so.der_t->get_name_ind()[cmp]) ;
729  p_christo[dd]->der_t->set_name_ind(1, name_sum) ;
730  }
731  }
732 
733  Term_eq auxi_christ (*p_christo[dd]) ;
734 
735  if (type_der==CON)
736  manipulate_ind(auxi_christ, 0) ;
737 
738  // Check if one inner summation is needed :
739  bool need_sum = false ;
740  char const * ind = p_christo[dd]->val_t->get_name_ind() ;
741  if ((ind[0]==ind[2]) || (ind[1]==ind[2]) || (ind[0]==ind[1]))
742  need_sum = true ;
743 
744  Term_eq* christ ;
745  if (need_sum)
746  if (!doder)
747  christ = new Term_eq (dd, auxi_christ.val_t->do_summation_one_dom(dd)) ;
748  else
749  christ = new Term_eq (dd, auxi_christ.val_t->do_summation_one_dom(dd),
750  auxi_christ.der_t->do_summation_one_dom(dd)) ;
751  else
752  christ = new Term_eq (auxi_christ) ;
753 
754  // Affecte names on the field :
755  Term_eq copie (so) ;
756  copie.val_t->set_name_affected() ;
757  for (int i=0 ; i<so.val_t->get_valence() ; i++)
758  if (i!=cmp)
759  copie.val_t->set_name_ind(i, so.val_t->get_name_ind()[i]) ;
760  else
761  copie.val_t->set_name_ind(i, name_sum) ;
762  if (doder) {
763  copie.der_t->set_name_affected() ;
764  for (int i=0 ; i<so.der_t->get_valence() ; i++)
765  if (i!=cmp)
766  copie.der_t->set_name_ind(i, so.der_t->get_name_ind()[i]) ;
767  else
768  copie.der_t->set_name_ind(i, name_sum) ;
769  }
770 
771  Term_eq part_christo ((*christ)*copie) ;
772  delete christ ;
773 
774  if (genre_indice==CON)
775  res = res + part_christo ;
776  else
777  res = res - part_christo ;
778  }
779  return res ;
780 }
781 
782 void Metric_conf_factor::set_system (System_of_eqs& ss, const char* name_met) {
783 
784  syst = &ss ;
785 
786  // Position in the system :
787  place_syst = ss.ndom*ss.nvar ;
788 
789  // unknown for the system (no name, the name is in the metric already)
790  ss.add_var (0x0, *p_met) ;
791 
792  if (ss.met!=0x0) {
793  cerr << "Metric already set for the system" << endl ;
794  abort() ;
795  }
796 
797  ss.met = this ;
798  ss.name_met = new char[LMAX] ;
799  trim_spaces (ss.name_met, name_met) ;
800 }
801 }
reference set(const Index &pos)
Read/write of an element.
Definition: array.hpp:186
int get_basis(int nd) const
Read only the basis in a given domain.
Definition: base_tensor.hpp:93
Class that gives the position inside a multi-dimensional Array.
Definition: index.hpp:38
bool inc(int increm, int var=0)
Increments the position of the Index.
Definition: index.hpp:99
Class to deal with a metric with a conformal decomposition (mainly used for AADS spacetimes) The true...
Definition: metric.hpp:562
virtual Term_eq derive(int, char, const Term_eq &) const
Computes the covariant derivative of a Term_eq (assumes Cartesian basis of decomposition).
Metric_flat fmet
Associated flat metric.
Definition: metric.hpp:567
virtual void compute_ricci_tensor(int) const
Computes the Ricci tensor, in a given Domain.
Vector grad_conf
flat gradient of the conformal factor
Definition: metric.hpp:569
Metric_conf_factor(Metric_tensor &, const Scalar &conf)
Constructor from a Metric_tensor and a conformal factor.
Metric_tensor * p_met
Pointer on the Metric_tensor describing the coformal metric.
Definition: metric.hpp:565
virtual void compute_christo(int) const
Computes the Christoffel symbols, in a given Domain.
virtual void set_system(System_of_eqs &syst, const char *name)
Associate the metric to a given system of equations.
int place_syst
Gives the location of the metric amongst the various unknowns of the associated System_of_eqs.
Definition: metric.hpp:570
virtual void compute_cov(int) const
Computes the covariant representation, in a given Domain.
virtual void compute_riemann(int) const
Computes the Riemann tensor, in a given Domain.
virtual Term_eq derive_flat(int, char, const Term_eq &) const
Computes the covariant flat derivative of a Term_eq.
virtual void compute_con(int) const
Computes the contravariant representation, in a given Domain.
Scalar conformal
The conformal factor (must be a purely radial function)
Definition: metric.hpp:568
virtual int give_type(int) const
Returns the type of tensorial basis of the covariant representation, in a given Domain.
const Base_tensor & basis
The tensorial basis used.
Definition: metric.hpp:566
Term_eq derive_with_other(int tder, char indder, const Term_eq &so, const Metric *othermet) const
Computes the flat covariant derivative.
virtual Term_eq derive(int, char, const Term_eq &) const
Computes the covariant derivative of a Term_eq (assumes Cartesian basis of decomposition).
Particular type of Tensor, dedicated to the desription of metrics.
Purely abstract class for metric handling.
Definition: metric.hpp:39
int type_tensor
States if one works in the CON or COV representation.
Definition: metric.hpp:86
MMPtr_array< Term_eq > p_ricci_tensor
Array of pointers on various Term_eq.
Definition: metric.hpp:70
const Term_eq * give_ricci_tensor(int dd) const
Gives the Ricci tensor, in a Domain.
Definition: metric.cpp:138
const Space & espace
The associated Space.
Definition: metric.hpp:42
MMPtr_array< Term_eq > p_met_cov
Array of pointers on various Term_eq.
Definition: metric.hpp:50
MMPtr_array< Term_eq > p_met_con
Array of pointers on various Term_eq.
Definition: metric.hpp:55
virtual void manipulate_ind(Term_eq &so, int ind) const
Uses the Metric to manipulate one of the index of a Term_eq (i.e.
Definition: metric.cpp:645
MMPtr_array< Term_eq > p_riemann
Array of pointers on various Term_eq.
Definition: metric.hpp:65
MMPtr_array< Term_eq > p_christo
Array of pointers on various Term_eq.
Definition: metric.hpp:60
const System_of_eqs * syst
Pointer of the system of equations where the metric is used (only one for now).
Definition: metric.hpp:44
The class Scalar does not really implements scalars in the mathematical sense but rather tensorial co...
Definition: scalar.hpp:67
Val_domain & set_domain(int)
Read/write of a particular Val_domain.
Definition: scalar.hpp:555
const Domain * get_domain(int i) const
returns a pointer on the domain.
Definition: space.hpp:1385
int get_ndim() const
Returns the number of dimensions.
Definition: space.hpp:1373
int get_nbr_domains() const
Returns the number of Domains.
Definition: space.hpp:1375
Class used to describe and solve a system of equations.
virtual void add_var(const char *name, double &var)
Addition of a variable (number case)
int ndom
Number of domains used.
char * name_met
Name by which the metric is recognized.
MMPtr_array< Term_eq > term
Pointers on the Term_eq corresponding to the unknown fields.
int dom_min
Smallest domain number.
Metric * met
Pointer on the associated Metric, if defined.
int nvar
Number of unknown fields.
Tensor handling.
Definition: tensor.hpp:149
void set_name_ind(int dd, char name)
Sets the name of one index ; the names must have been affected first.
void set_name_affected()
Affects the name of the indices.
Definition: tensor.hpp:435
Scalar & set(const Array< int > &ind)
Returns the value of a component (read/write version).
Definition: tensor_impl.hpp:91
char const * get_name_ind() const
Definition: tensor.hpp:424
int get_index_type(int i) const
Gives the type (covariant or contravariant) of a given index.
Definition: tensor.hpp:526
int get_n_comp() const
Returns the number of stored components.
Definition: tensor.hpp:514
int get_valence() const
Returns the valence.
Definition: tensor.hpp:509
Tensor do_summation_one_dom(int dd) const
Does the inner contraction of the Tensor in a given domain.
This class is intended to describe the manage objects appearing in the equations.
Definition: term_eq.hpp:62
Tensor * der_t
Pointer on the variation, if the Term_eq is a Tensor.
Definition: term_eq.hpp:69
void set_der_t(Tensor)
Sets the tensorial variation (only the values in the pertinent Domain are copied).
Definition: term_eq.cpp:171
int get_dom() const
Definition: term_eq.hpp:135
Tensor const & get_val_t() const
Definition: term_eq.hpp:355
Tensor const & get_der_t() const
Definition: term_eq.hpp:369
void set_val_t(Tensor)
Sets the tensorial value (only the values in the pertinent Domain are copied).
Definition: term_eq.cpp:155
Tensor * val_t
Pointer on the value, if the Term_eq is a Tensor.
Definition: term_eq.hpp:68
Class for storing the basis of decompositions of a field and its values on both the configuration and...
Definition: val_domain.hpp:69
Scalar & set(int)
Read/write access to a component.