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