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