KADATH
metric_nophi_AADS.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_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), der_conf (conformal.der_r()) {
33  type_tensor = met.get_type() ;
34  for (int d=0 ; d<met.get_space().get_nbr_domains() ; d++)
35  if (basis.get_basis(d) != SPHERICAL_BASIS) {
36  cerr << "Metric_nophi_AADS only defined wrt spherical tensorial coordinates" << endl ;
37  abort() ;
38  }
39 
40 }
41 
43  Metric (so), p_met(so.p_met), basis(so.basis), fmet(so.fmet), conformal(so.conformal), der_conf(conformal.der_r()), place_syst(so.place_syst) {
44 }
45 
46 Metric_nophi_AADS::~Metric_nophi_AADS() {
47 }
48 
49 
50 int Metric_nophi_AADS::give_type(int dd) const {
51  return basis.get_basis(dd) ;
52 }
53 
54 void Metric_nophi_AADS::compute_cov (int dd) const {
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 
239 void Metric_nophi_AADS::compute_con (int dd) const {
240 
241  int place = place_syst + (dd-syst->dom_min) ;
242  // Right storage : simple copy.
243  if (type_tensor==CON) {
244 
245  if (p_met_con[dd]==0x0)
246  p_met_con[dd] = new Term_eq(*syst->term[place]) ;
247  else
248  *p_met_con[dd] = Term_eq(*syst->term[place]) ;
249  }
250  else {
251  //Need to work component by components...
252  bool doder = (syst->term[place]->der_t==0x0) ? false : true ;
253  Term_eq** res = new Term_eq* [p_met->get_n_comp()] ;
254 
255  Scalar val (espace) ;
256  Scalar der (espace) ;
257 
258  Val_domain cmpval(espace.get_domain(dd)) ;
259  Val_domain cmpder(espace.get_domain(dd)) ;
260 
261  Val_domain detval (espace.get_domain(dd)) ;
262  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)
263  + (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)
264  + (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
265  - (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,2)(dd)
266  - (*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(1,1)(dd)
267  - (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd) ;
268  Val_domain detder (espace.get_domain(dd)) ;
269  if (doder)
270  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)
271  + (*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)
272  + (*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
273  - (*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,2)(dd)
274  - (*syst->term[place]->der_t)(2,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(1,1)(dd)
275  - (*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
276  + (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->der_t)(2,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
277  + (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(2,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)
278  + (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
279  - (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,2)(dd)
280  - (*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->der_t)(2,3)(dd)*(*syst->term[place]->val_t)(1,1)(dd)
281  - (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
282  + (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->val_t)(2,2)(dd)*(*syst->term[place]->der_t)(3,3)(dd)
283  + (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->der_t)(1,3)(dd)
284  + (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(2,3)(dd)
285  - (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(2,2)(dd)
286  - (*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->der_t)(1,1)(dd)
287  - (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(3,3)(dd) ;
288 
289 
290  // Compo 1 1
291  cmpval = (*syst->term[place]->val_t)(2,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
292  -(*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd) ;
293  val.set_domain(dd) = cmpval / detval;
294  if (doder) {
295  cmpder = (*syst->term[place]->der_t)(2,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
296  + (*syst->term[place]->val_t)(2,2)(dd)*(*syst->term[place]->der_t)(3,3)(dd)
297  - (*syst->term[place]->der_t)(2,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
298  - (*syst->term[place]->val_t)(2,3)(dd)*(*syst->term[place]->der_t)(2,3)(dd) ;
299  der.set_domain(dd) = cmpder / detval - cmpval * detder / detval/detval;
300  res[0] = new Term_eq (dd, val, der) ;
301  }
302  else {
303  res[0] = new Term_eq (dd, val) ;
304  }
305 
306  // Compo 1 2
307  cmpval = (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
308  -(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd) ;
309  val.set_domain(dd) = cmpval / detval;
310  if (doder) {
311  cmpder = (*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
312  + (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(2,3)(dd)
313  - (*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
314  - (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(3,3)(dd) ;
315  der.set_domain(dd) = cmpder / detval - cmpval * detder / detval/detval;
316  res[1] = new Term_eq (dd, val, der) ;
317  }
318  else {
319  res[1] = new Term_eq (dd, val) ;
320  }
321 
322  // Compo 1 3
323  cmpval = (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
324  -(*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,2)(dd) ;
325  val.set_domain(dd) = cmpval / detval;
326  if (doder) {
327  cmpder = (*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
328  + (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(2,3)(dd)
329  - (*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(2,2)(dd)
330  - (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(2,2)(dd) ;
331  der.set_domain(dd) = cmpder / detval - cmpval * detder / detval/detval;
332  res[2] = new Term_eq (dd, val, der) ;
333  }
334  else {
335  res[2] = new Term_eq (dd, val) ;
336  }
337 
338  // Compo 2 2
339  cmpval = (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
340  -(*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd) ;
341  val.set_domain(dd) = cmpval / detval;
342  if (doder) {
343  cmpder = (*syst->term[place]->der_t)(1,1)(dd)*(*syst->term[place]->val_t)(3,3)(dd)
344  + (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->der_t)(3,3)(dd)
345  - (*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,3)(dd)
346  - (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(1,3)(dd) ;
347  der.set_domain(dd) = cmpder / detval - cmpval * detder / detval/detval;
348  res[3] = new Term_eq (dd, val, der) ;
349  }
350  else {
351  res[3] = new Term_eq (dd, val) ;
352  }
353 
354  // Compo 2 3
355  cmpval = (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,2)(dd)
356  -(*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->val_t)(2,3)(dd) ;
357  val.set_domain(dd) = cmpval / detval ;
358  if (doder) {
359  cmpder = (*syst->term[place]->der_t)(1,3)(dd)*(*syst->term[place]->val_t)(1,2)(dd)
360  + (*syst->term[place]->val_t)(1,3)(dd)*(*syst->term[place]->der_t)(1,2)(dd)
361  - (*syst->term[place]->der_t)(1,1)(dd)*(*syst->term[place]->val_t)(2,3)(dd)
362  - (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->der_t)(2,3)(dd) ;
363  der.set_domain(dd) = cmpder/ detval - cmpval * detder / detval/detval ;
364  res[4] = new Term_eq (dd, val, der) ;
365  }
366  else {
367  res[4] = new Term_eq (dd, val) ;
368  }
369 
370  // Compo 3 3
371  cmpval = (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->val_t)(2,2)(dd)
372  -(*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->val_t)(1,2)(dd) ;
373  val.set_domain(dd) = cmpval / detval;
374  if (doder) {
375  cmpder = (*syst->term[place]->der_t)(1,1)(dd)*(*syst->term[place]->val_t)(2,2)(dd)
376  + (*syst->term[place]->val_t)(1,1)(dd)*(*syst->term[place]->der_t)(2,2)(dd)
377  - (*syst->term[place]->der_t)(1,2)(dd)*(*syst->term[place]->val_t)(1,2)(dd)
378  - (*syst->term[place]->val_t)(1,2)(dd)*(*syst->term[place]->der_t)(1,2)(dd) ;
379  der.set_domain(dd) = cmpder / detval - cmpval * detder / detval/detval;
380  res[5] = new Term_eq (dd, val, der) ;
381  }
382  else {
383  res[5] = new Term_eq (dd, val) ;
384  }
385 
386  // Value field :
387  Metric_tensor resval (espace, CON, basis) ;
388  resval.set(1,1) = res[0]->get_val_t() ;
389  resval.set(1,2) = res[1]->get_val_t() ;
390  resval.set(1,3) = res[2]->get_val_t() ;
391  resval.set(2,2) = res[3]->get_val_t() ;
392  resval.set(2,3) = res[4]->get_val_t() ;
393  resval.set(3,3) = res[5]->get_val_t() ;
394 
395  if (!doder) {
396  if (p_met_con[dd]==0x0)
397  p_met_con[dd] = new Term_eq(dd, resval) ;
398  else
399  *p_met_con[dd] = Term_eq(dd, resval) ;
400  }
401  else {
402  // Der field :
403  Metric_tensor resder (espace, CON, basis) ;
404  resder.set(1,1) = res[0]->get_der_t() ;
405  resder.set(1,2) = res[1]->get_der_t() ;
406  resder.set(1,3) = res[2]->get_der_t() ;
407  resder.set(2,2) = res[3]->get_der_t() ;
408  resder.set(2,3) = res[4]->get_der_t() ;
409  resder.set(3,3) = res[5]->get_der_t() ;
410 
411  if (p_met_con[dd]==0x0)
412  p_met_con[dd] = new Term_eq(dd, resval, resder) ;
413  else
414  *p_met_con[dd] = Term_eq(dd, resval, resder) ;
415  }
416 
417  for (int i=0 ; i<p_met->get_n_comp() ; i++)
418  delete res [i] ;
419  delete [] res ;
420  }
421 }
422 
424  // Need both representation of the metric)
425  if (type_tensor == CON) {
426  if (p_met_con[dd]==0x0)
427  compute_con(dd) ;
428  if (p_met_cov[dd]==0x0)
429  compute_cov(dd) ;
430  }
431  else {
432  if (p_met_cov[dd]==0x0)
433  compute_cov(dd) ;
434  if (p_met_con[dd]==0x0)
435  compute_con(dd) ;
436  }
437 
438  Array<int> type_ind (3) ;
439  type_ind.set(0) = COV ; type_ind.set(1) = COV ; type_ind.set(2) = CON ;
440  Tensor res_val (espace, 3, type_ind, basis, 3) ;
441  Tensor res_der (espace, 3, type_ind, basis, 3) ;
442 
443  Term_eq flat_der (fmet.derive(COV, ' ', (*p_met_cov[dd]))) ;
444  bool doder = (flat_der.der_t==0x0) ? false : true ;
445 
446  Index pos (res_val) ;
447  do {
448 
449  Val_domain cmpval (espace.get_domain(dd)) ;
450  cmpval = 0 ;
451  for (int l=1 ; l<=espace.get_ndim() ; l++)
452  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) +
453  (*flat_der.val_t)(pos(1)+1, pos(0)+1, l)(dd) - (*flat_der.val_t)(l, pos(0)+1, pos(1)+1)(dd))) ;
454 
455  // Some corrections due to the conformal factor :
456  cmpval += (*p_met_con[dd]->val_t)(pos(2)+1,1)(dd)*(*p_met_cov[dd]->val_t)(pos(0)+1,pos(1)+1)(dd)*der_conf(dd) ;
457  if ((pos(1)==pos(2)) && (pos(0)==0))
458  cmpval -= der_conf(dd) ;
459  if ((pos(0)==pos(2)) && (pos(1)==0))
460  cmpval -= der_conf(dd) ;
461 
462  res_val.set(pos).set_domain(dd) = cmpval ;
463 
464  if (doder) {
465  Val_domain cmpder (espace.get_domain(dd)) ;
466  cmpder = 0 ;
467  for (int l=1 ; l<=espace.get_ndim() ; l++)
468  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) +
469  (*flat_der.val_t)(pos(1)+1, pos(0)+1, l)(dd) - (*flat_der.val_t)(l, pos(0)+1, pos(1)+1)(dd))
470  + 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) +
471  (*flat_der.der_t)(pos(1)+1, pos(0)+1, l)(dd) - (*flat_der.der_t)(l, pos(0)+1, pos(1)+1)(dd))) ;
472 
473  // Some corrections due to the conformal factor :
474  cmpder += ((*p_met_con[dd]->der_t)(pos(2)+1,1)(dd)*(*p_met_cov[dd]->val_t)(pos(0)+1,pos(1)+1)(dd)
475  +(*p_met_con[dd]->val_t)(pos(2)+1,1)(dd)*(*p_met_cov[dd]->der_t)(pos(0)+1,pos(1)+1)(dd)) *der_conf(dd) ;
476 
477 
478  res_der.set(pos).set_domain(dd) = cmpder ;
479  }
480  }
481  while (pos.inc()) ;
482 
483  if (!doder) {
484  if (p_christo[dd]==0x0)
485  p_christo[dd] = new Term_eq(dd, res_val) ;
486  else
487  *p_christo[dd] = Term_eq (dd, res_val) ;
488  }
489  else {
490  if (p_christo[dd]==0x0)
491  p_christo[dd] = new Term_eq(dd, res_val, res_der) ;
492  else
493  *p_christo[dd] = Term_eq(dd, res_val, res_der) ;
494  }
495 
496 }
497 
499 
500  // Need christoffels
501  if (p_christo[dd]==0x0)
502  compute_christo(dd) ;
503 
504 
505  Array<int> indices (4) ;
506  indices.set(0) = CON ; indices.set(1) = COV ; indices.set(2) = COV ; indices.set(3) = COV ;
507  Tensor res_val (espace, 4, indices, basis, 3) ;
508  Tensor res_der (espace, 4, indices, basis, 3) ;
509 
510  Term_eq flat_der (fmet.derive(COV, ' ', (*p_christo[dd]))) ;
511  bool doder = (flat_der.der_t==0x0) ? false : true ;
512 
513  Index pos (res_val) ;
514  do {
515 
516  Val_domain cmpval (conformal(dd)*((*flat_der.val_t)(pos(2)+1, pos(1)+1,pos(3)+1,pos(0)+1)(dd)
517  - (*flat_der.val_t)(pos(3)+1, pos(1)+1, pos(2)+1, pos(0)+1)(dd))) ;
518 
519  // Corrections coming from the conformal part
520  if (pos(2)==0)
521  cmpval -= (*p_christo[dd]->val_t)(pos(1)+1, pos(3)+1, pos(0)+1)(dd)*der_conf(dd) ;
522  if (pos(3)==0)
523  cmpval += (*p_christo[dd]->val_t)(pos(1)+1, pos(2)+1, pos(0)+1)(dd)*der_conf(dd) ;
524 
525  for (int m=1 ; m<=espace.get_ndim() ; m++) {
526  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)
527  - (*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) ;
528  }
529 
530 
531 
532  res_val.set(pos).set_domain(dd) = cmpval;
533  if (doder) {
534  Val_domain cmpder (conformal(dd)*((*flat_der.der_t)(pos(2)+1, pos(1)+1,pos(3)+1,pos(0)+1)(dd)
535  - (*flat_der.der_t)(pos(3)+1, pos(1)+1, pos(2)+1, pos(0)+1)(dd))) ;
536 
537  if (pos(2)==0)
538  cmpder -= (*p_christo[dd]->der_t)(pos(1)+1, pos(3)+1, pos(0)+1)(dd)*der_conf(dd) ;
539  if (pos(3)==0)
540  cmpder += (*p_christo[dd]->der_t)(pos(1)+1, pos(2)+1, pos(0)+1)(dd)*der_conf(dd) ;
541 
542 
543  for (int m=1 ; m<=espace.get_ndim() ; m++) {
544  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)
545  + (*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)
546  - (*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)
547  - (*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);
548  }
549  res_der.set(pos).set_domain(dd) = cmpder ;
550  }
551  }
552  while (pos.inc()) ;
553 
554  if (!doder) {
555  if (p_riemann[dd]==0x0)
556  p_riemann[dd] = new Term_eq(dd, res_val) ;
557  else
558  *p_riemann[dd] = Term_eq (dd, res_val) ;
559  }
560  else {
561  if (p_riemann[dd]==0x0)
562  p_riemann[dd] = new Term_eq(dd, res_val, res_der) ;
563  else
564  *p_riemann[dd] = Term_eq(dd, res_val, res_der) ;
565  }
566 
567 }
568 
570  // Need christoffels
571  if (p_christo[dd]==0x0)
572  compute_christo(dd) ;
573 
574 
575  Array<int> indices (2) ;
576  indices.set(0) = COV ; indices.set(1) = COV ;
577  Tensor res_val (espace, 2, indices, basis, 3) ;
578  Tensor res_der (espace, 2, indices, basis, 3) ;
579 
580  Term_eq flat_der (fmet.derive(COV, ' ', (*p_christo[dd]))) ;
581  bool doder = (flat_der.der_t==0x0) ? false : true ;
582 
583  Index pos (res_val) ;
584  do {
585 
586  Val_domain cmpval (espace.get_domain(dd)) ;
587  cmpval = 0 ;
588 
589  for (int l=1 ; l<=espace.get_ndim() ; l++) {
590  cmpval += conformal(dd)*((*flat_der.val_t)(l, pos(1)+1,pos(0)+1,l)(dd)
591  - (*flat_der.val_t)(pos(1)+1, pos(0)+1,l ,l)(dd)) ;
592  if (l==1)
593  cmpval-=(*p_christo[dd]->val_t)(pos(0)+1, pos(1)+1, l)(dd)*der_conf(dd) ;
594  if (pos(1)==0)
595  cmpval+=(*p_christo[dd]->val_t)(pos(0)+1, l, l)(dd)*der_conf(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  if (l==1)
613  cmpder-=(*p_christo[dd]->der_t)(pos(0)+1, pos(1)+1, l)(dd)*der_conf(dd) ;
614  if (pos(1)==0)
615  cmpder+=(*p_christo[dd]->der_t)(pos(0)+1, l, l)(dd)*der_conf(dd) ;
616 
617 
618  for (int m=1 ; m<=espace.get_ndim() ; m++) {
619  cmpder += (*p_christo[dd]->der_t)(l,m, l)(dd)*(*p_christo[dd]->val_t)(pos(0)+1,pos(1)+1,m)(dd)
620  + (*p_christo[dd]->val_t)(l,m, l)(dd)*(*p_christo[dd]->der_t)(pos(0)+1,pos(1)+1,m)(dd)
621  - (*p_christo[dd]->der_t)(pos(0)+1,l, m)(dd)*(*p_christo[dd]->val_t)(pos(1)+1,m,l)(dd)
622  - (*p_christo[dd]->val_t)(pos(0)+1,l, m)(dd)*(*p_christo[dd]->der_t)(pos(1)+1,m,l)(dd);
623  }
624  }
625  res_der.set(pos).set_domain(dd) = cmpder ;
626  }
627  }
628  while (pos.inc()) ;
629 
630  if (!doder) {
631  if (p_ricci_tensor[dd]==0x0)
632  p_ricci_tensor[dd] = new Term_eq(dd, res_val) ;
633  else
634  *p_ricci_tensor[dd] = Term_eq (dd, res_val) ;
635  }
636  else {
637  if (p_ricci_tensor[dd]==0x0)
638  p_ricci_tensor[dd] = new Term_eq(dd, res_val, res_der) ;
639  else
640  *p_ricci_tensor[dd] = Term_eq(dd, res_val, res_der) ;
641  }
642 
643 
644 }
645 
646 Term_eq Metric_nophi_AADS::derive_flat (int type_der, char ind_der, const Term_eq& so) const {
647 
648  int dd = so.get_dom() ;
649 
650  if (p_met_con[dd]==0x0)
651  compute_con(dd) ;
652  if (p_met_cov[dd]==0x0)
653  compute_cov(dd) ;
654 
655 
656  // The partial derivative part
657  Term_eq res (fmet.derive_with_other (type_der, ind_der, so, this)) ;
658  return res ;
659 }
660 
661 Term_eq Metric_nophi_AADS::derive (int type_der, char ind_der, const Term_eq& so) const {
662 
663  int dd = so.get_dom() ;
664 
665  if (p_christo[dd]==0x0)
666  compute_christo(dd) ;
667  // The partial derivative part
668  Term_eq res (conformal*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  bool doder = ((so.der_t==0x0) || (p_christo[dd]->der_t==0x0)) ? false : true ;
694 
695  //loop on the components of so :
696  for (int cmp=0 ; cmp<so.val_t->get_valence() ; cmp ++) {
697 
698  int genre_indice = so.val_t->get_index_type(cmp) ;
699 
700  //Affecte names on the Christoffels :
701  p_christo[dd]->val_t->set_name_affected() ;
702  p_christo[dd]->val_t->set_name_ind(0, ind_der) ;
703  if (genre_indice==COV) {
704  p_christo[dd]->val_t->set_name_ind(1, so.val_t->get_name_ind()[cmp]) ;
705  p_christo[dd]->val_t->set_name_ind(2, name_sum) ;
706  }
707  else {
708  p_christo[dd]->val_t->set_name_ind(2, so.val_t->get_name_ind()[cmp]) ;
709  p_christo[dd]->val_t->set_name_ind(1, name_sum) ;
710  }
711 
712  if (doder) {
713  p_christo[dd]->der_t->set_name_affected() ;
714  p_christo[dd]->der_t->set_name_ind(0, ind_der) ;
715  if (genre_indice==COV) {
716  p_christo[dd]->der_t->set_name_ind(1, so.der_t->get_name_ind()[cmp]) ;
717  p_christo[dd]->der_t->set_name_ind(2, name_sum) ;
718  }
719  else {
720  p_christo[dd]->der_t->set_name_ind(2, so.der_t->get_name_ind()[cmp]) ;
721  p_christo[dd]->der_t->set_name_ind(1, name_sum) ;
722  }
723  }
724 
725  Term_eq auxi_christ (*p_christo[dd]) ;
726 
727  if (type_der==CON)
728  manipulate_ind(auxi_christ, 0) ;
729 
730  // Check if one inner summation is needed :
731  bool need_sum = false ;
732  char const * ind = p_christo[dd]->val_t->get_name_ind() ;
733  if ((ind[0]==ind[2]) || (ind[1]==ind[2]) || (ind[0]==ind[1]))
734  need_sum = true ;
735 
736  Term_eq* christ ;
737  if (need_sum)
738  if (!doder)
739  christ = new Term_eq (dd, auxi_christ.val_t->do_summation_one_dom(dd)) ;
740  else
741  christ = new Term_eq (dd, auxi_christ.val_t->do_summation_one_dom(dd),
742  auxi_christ.der_t->do_summation_one_dom(dd)) ;
743  else
744  christ = new Term_eq (auxi_christ) ;
745 
746  // Affecte names on the field :
747  Term_eq copie (so) ;
748  copie.val_t->set_name_affected() ;
749  for (int i=0 ; i<so.val_t->get_valence() ; i++)
750  if (i!=cmp)
751  copie.val_t->set_name_ind(i, so.val_t->get_name_ind()[i]) ;
752  else
753  copie.val_t->set_name_ind(i, name_sum) ;
754  if (doder) {
755  copie.der_t->set_name_affected() ;
756  for (int i=0 ; i<so.der_t->get_valence() ; i++)
757  if (i!=cmp)
758  copie.der_t->set_name_ind(i, so.der_t->get_name_ind()[i]) ;
759  else
760  copie.der_t->set_name_ind(i, name_sum) ;
761  }
762 
763  Term_eq part_christo ((*christ)*copie) ;
764  delete christ ;
765 
766  if (genre_indice==CON)
767  res = res + part_christo ;
768  else
769  res = res - part_christo ;
770  }
771  return res ;
772 }
773 
774 void Metric_nophi_AADS::set_system (System_of_eqs& ss, const char* name_met) {
775 
776  syst = &ss ;
777 
778  // Position in the system :
779  place_syst = ss.ndom*ss.nvar ;
780 
781  // unknown for the system (no name, the name is in the metric already)
782  ss.add_var (0x0, *p_met) ;
783 
784  if (ss.met!=0x0) {
785  cerr << "Metric already set for the system" << endl ;
786  abort() ;
787  }
788 
789  ss.met = this ;
790  ss.name_met = new char[LMAX] ;
791  trim_spaces (ss.name_met, name_met) ;
792 }
793 
794 }
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
virtual Term_eq derive(int, char, const Term_eq &) const
Computes the covariant derivative of a Term_eq (assumes Cartesian basis of decomposition).
Term_eq derive_with_other(int tder, char indder, const Term_eq &so, const Metric *othermet) const
Computes the flat covariant derivative.
Class to deal with a metric independant of with a conformal decomposition (mainly used for AADS spac...
Metric_nophi_AADS(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.
Metric_flat_nophi fmet
Associated flat metric.
virtual void compute_christo(int) const
Computes the Christoffel symbols, in a given Domain.
virtual void compute_ricci_tensor(int) const
Computes the Ricci tensor, 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 Term_eq derive_flat(int, char, const Term_eq &) const
Computes the covariant flat derivative of a Term_eq.
Scalar der_conf
Radial derivative of the conformal factor.
virtual void set_system(System_of_eqs &syst, const char *name)
Associate the metric to a given system of equations.
const Base_tensor & basis
The tensorial basis used.
virtual void compute_con(int) const
Computes the contravariant representation, in a given Domain.
virtual void compute_cov(int) const
Computes the covariant representation, in a given Domain.
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_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.
Scalar conformal
The conformal factor (must be a purely radial function)
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_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.
const Space & get_space() const
Returns the Space.
Definition: tensor.hpp:499
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