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