KADATH
domain_polar_inner_adapted_term_eq.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 "headcpp.hpp"
21 #include "utilities.hpp"
22 #include "adapted_polar.hpp"
23 #include "point.hpp"
24 #include "array_math.hpp"
25 #include "val_domain.hpp"
26 #include "scalar.hpp"
27 #include "tensor_impl.hpp"
28 #include "vector.hpp"
29 namespace Kadath {
31  return derive_r(so) ;
32 }
33 
35 
36  assert (so.get_dom() == num_dom) ;
37 
38  Tensor val_res(*so.val_t, false) ;
39 
40  for (int cmp = 0 ; cmp<so.val_t->get_n_comp() ; cmp++) {
41  Array<int>ind (val_res.indices(cmp)) ;
42  val_res.set(ind).set_domain(num_dom) = (*so.val_t)(ind)(num_dom).der_var(1) / (*der_rad_term_eq->val_t)()(num_dom) ;
43  }
44 
45  bool doder = ((so.der_t==0x0) || (der_rad_term_eq->der_t==0x0)) ? false : true ;
46  if (doder) {
47  Tensor der_res (*so.der_t, false) ;
48  for (int cmp = 0 ; cmp<so.val_t->get_n_comp() ; cmp++) {
49  Array<int>ind (val_res.indices(cmp)) ;
50  der_res.set(ind).set_domain(num_dom) = ((*so.der_t)(ind)(num_dom).der_var(1)*(*der_rad_term_eq->val_t)()(num_dom)
51  - (*so.val_t)(ind)(num_dom).der_var(1)*(*der_rad_term_eq->der_t)()(num_dom)) /
53  }
54  return Term_eq (num_dom, val_res, der_res) ;
55  }
56  else
57  return Term_eq (num_dom, val_res) ;
58 }
59 
61 
62  assert (so.get_dom() == num_dom) ;
63  Term_eq* dtprime ;
64  Tensor val_res(*so.val_t, false) ;
65 
66  for (int cmp = 0 ; cmp<so.val_t->get_n_comp() ; cmp++) {
67  Array<int>ind (val_res.indices(cmp)) ;
68  val_res.set(ind).set_domain(num_dom) = (*so.val_t)(ind)(num_dom).der_var(2) ;
69  }
70  bool doder = (so.der_t==0x0) ? false : true ;
71 
72  if (doder) {
73  Tensor der_res (*so.der_t, false) ;
74  for (int cmp = 0 ; cmp<so.val_t->get_n_comp() ; cmp++) {
75  Array<int>ind (val_res.indices(cmp)) ;
76  der_res.set(ind).set_domain(num_dom) = (*so.der_t)(ind)(num_dom).der_var(2) ;
77  }
78  dtprime = new Term_eq (num_dom, val_res, der_res) ;
79  }
80  else
81  dtprime = new Term_eq (num_dom, val_res) ;
82 
83 
84  Term_eq res (*dtprime - (*dt_rad_term_eq) * derive_r(so)) ;
85  delete dtprime ;
86  return res ;
87 
88 }
89 
91 
92  assert (so.get_dom() == num_dom) ;
93 
94  int valso = so.get_val_t().get_valence() ;
95  Array<int> type_ind (valso+1) ;
96  type_ind.set(0) = COV ;
97  for (int i=0 ; i<valso ; i++)
98  type_ind.set(i+1) = so.get_val_t().get_index_type(i) ;
99 
100  Base_tensor basis (sp) ;
101  basis.set_basis(num_dom) = SPHERICAL_BASIS ;
102  Tensor res (sp, valso+1 , type_ind, basis) ;
103 
104  Term_eq comp_r (derive_r(so)) ;
105  Term_eq comp_t (derive_t(so)/(*rad_term_eq)) ;
106 
107  // Loop on cmp :
108  for (int nc=0 ; nc<so.get_val_t().get_n_comp() ; nc++) {
109 
110  Array<int> ind (so.get_val_t().indices(nc)) ;
111  Array<int> indtarget (valso+1) ;
112  for (int i=0 ; i<valso ; i++)
113  indtarget.set(i+1) = ind(i) ;
114 
115  // R comp :
116  indtarget.set(0) = 1 ;
117  res.set(indtarget).set_domain(num_dom) = comp_r.get_val_t()(ind)(num_dom);
118  // theta comp :
119  indtarget.set(0) = 2 ;
120  res.set(indtarget).set_domain(num_dom) = comp_t.get_val_t()(ind)(num_dom) ;
121  }
122  bool doder = ((so.der_t==0x0) || (inner_radius_term_eq->der_t==0x0)) ? false : true ;
123 
124  if (!doder) {
125  return Term_eq (num_dom, res) ;
126  }
127  else {
128 
129  Tensor resder (sp, valso+1 , type_ind, basis) ;
130  // Loop on cmp :
131  for (int nc=0 ; nc<so.get_val_t().get_n_comp() ; nc++) {
132 
133  Array<int> ind (so.get_val_t().indices(nc)) ;
134  Array<int> indtarget (valso+1) ;
135  for (int i=0 ; i<valso ; i++)
136  indtarget.set(i+1) = ind(i) ;
137 
138  // R comp :
139  indtarget.set(0) = 1 ;
140  resder.set(indtarget).set_domain(num_dom) = comp_r.get_der_t()(ind)(num_dom);
141  // theta comp :
142  indtarget.set(0) = 2 ;
143  resder.set(indtarget).set_domain(num_dom) = comp_t.get_der_t()(ind)(num_dom) ;
144  }
145  return Term_eq (num_dom, res, resder) ;
146  }
147 }
148 
150 
152 
153  Scalar val_norme (sp) ;
154  val_norme.set_domain(num_dom) = sqrt((*grad.val_t)(1)(num_dom)*(*grad.val_t)(1)(num_dom)
155  + (*grad.val_t)(2)(num_dom)*(*grad.val_t)(2)(num_dom)) ;
156  bool doder = ((grad.der_t==0x0) || (inner_radius_term_eq->der_t==0x0)) ? false : true ;
157  if (doder) {
158  Scalar der_norme (sp) ;
159  der_norme.set_domain(num_dom) = ((*grad.der_t)(1)(num_dom)*(*grad.val_t)(1)(num_dom)
160  + (*grad.der_t)(2)(num_dom)*(*grad.val_t)(2)(num_dom)) / val_norme(num_dom) ;
161  Term_eq norme (num_dom, val_norme, der_norme) ;
162  if (normal_spher==0x0)
163  normal_spher = new Term_eq (grad / norme) ;
164  else
165  *normal_spher = Term_eq(grad/norme) ;
166  }
167  else {
168  Term_eq norme (num_dom, val_norme) ;
169  if (normal_spher==0x0)
170  normal_spher = new Term_eq (grad / norme) ;
171  else
172  *normal_spher = Term_eq(grad/norme) ;
173  }
174 }
175 
177 
178 
179  do_normal_spher() ;
180 
181  Base_tensor basis (sp) ;
182  basis.set_basis(num_dom) = CARTESIAN_BASIS ;
183  Vector val (sp, CON, basis) ;
184 
187 
188  bool doder = (normal_spher->der_t==0x0) ? false : true ;
189  if (doder) {
190  Vector der (sp, CON, basis) ;
191 
194 
195  if (normal_cart==0x0)
196  normal_cart = new Term_eq (num_dom, val, der) ;
197  else
198  *normal_cart = Term_eq(num_dom, val,der) ;
199  }
200  else {
201  if (normal_cart==0x0)
202  normal_cart = new Term_eq (num_dom, val) ;
203  else
204  *normal_cart = Term_eq(num_dom, val) ;
205  }
206 }
207 
209 
210 
211  switch (bound) {
212  case INNER_BC : {
213  // Deformed surface
214  if (normal_spher == 0x0)
215  do_normal_spher() ;
216 
217  int valso = so.get_val_t().get_valence() ;
218 
219  Term_eq grad (flat_grad_spher(so)) ;
220 
221  Tensor res (so.get_val_t(), false) ;
222 
223  // Loop on cmp :
224  for (int nc=0 ; nc<so.get_val_t().get_n_comp() ; nc++) {
225 
226  Array<int> ind (so.get_val_t().indices(nc)) ;
227  Array<int> indgrad (valso+1) ;
228  for (int i=0 ; i<valso ; i++)
229  indgrad.set(i+1) = ind(i) ;
230 
231  indgrad.set(0) = 1 ;
232  res.set(ind).set_domain(num_dom) = (*grad.val_t)(indgrad)(num_dom) * (*normal_spher->val_t)(1)(num_dom) ;
233  indgrad.set(0) = 2 ;
234  res.set(ind).set_domain(num_dom) += (*grad.val_t)(indgrad)(num_dom) * (*normal_spher->val_t)(2)(num_dom) ;
235  }
236 
237  bool doder = ((grad.der_t == 0x0) || (normal_spher->der_t == 0x0)) ? false : true ;
238  if (doder) {
239 
240  Tensor der (so.get_val_t(), false) ;
241 
242  // Loop on cmp :
243  for (int nc=0 ; nc<so.get_val_t().get_n_comp() ; nc++) {
244 
245  Array<int> ind (so.get_val_t().indices(nc)) ;
246  Array<int> indgrad (valso+1) ;
247  for (int i=0 ; i<valso ; i++)
248  indgrad.set(i+1) = ind(i) ;
249 
250  indgrad.set(0) = 1 ;
251  der.set(ind).set_domain(num_dom) = (*grad.der_t)(indgrad)(num_dom) * (*normal_spher->val_t)(1)(num_dom)
252  + (*grad.val_t)(indgrad)(num_dom) * (*normal_spher->der_t)(1)(num_dom) ;
253  indgrad.set(0) = 2 ;
254  der.set(ind).set_domain(num_dom) += (*grad.der_t)(indgrad)(num_dom) * (*normal_spher->val_t)(2)(num_dom)
255  + (*grad.val_t)(indgrad)(num_dom) * (*normal_spher->der_t)(2)(num_dom) ;
256  }
257 
258  return Term_eq (num_dom, res, der) ;
259  }
260  else
261  return Term_eq (num_dom, res) ;
262  }
263  case OUTER_BC : {
264  return derive_r(so) ;
265  }
266  default :
267  cerr << "Unknown boundary in Domain_polar_shell_inner_adapted::der_normal" << endl ;
268  abort() ;
269  }
270 }
271 
272 
274 
275  if (so.get_val_t().get_valence() != 0) {
276  cerr << "Domain_polar_shell_inner_adapted::lap_term_eq only defined for scalars" << endl ;
277  abort() ;
278  }
279 
280  // Angular part :
281 
283 
284  Term_eq dert (derive_t(so));
285  Term_eq p1 (derive_t(dert)) ;
287 
288  Term_eq dr(derive_r(so)) ;
289 
290  Term_eq res (derive_r(dr) + 2 *dr / (*rad_term_eq) + (p1+p2)/(*rad_term_eq)/(*rad_term_eq));
291  return res ;
292 }
293 
295  if (mm!=0) {
296  cerr << "Domain_polar_shell_inner_adapted::lap2_term_eq not defined for m != 0 (for now)" << endl ;
297  abort() ;
298  }
299 
300  if (so.get_val_t().get_valence() != 0) {
301  cerr << "Domain_polar_shell_inner_adapted::lap2_term_eq only defined for scalars" << endl ;
302  abort() ;
303  }
304 
305  // Angular part :
306  Term_eq dert (derive_t(so));
307  Term_eq p1 (derive_t(dert)) ;
308 
309  Term_eq dr(derive_r(so)) ;
310 
311  Term_eq res (derive_r(dr) + dr / (*rad_term_eq) + p1/(*rad_term_eq)/(*rad_term_eq));
312 
313 
314  return res ;
315 }
316 
318 
319 
320  return so * (*rad_term_eq) ;
321 }
322 
324 
325 
326  return so / (*rad_term_eq) ;
327 }
328 
330 
331  Tensor der (*so->val_t, false) ;
332  for (int cmp=0 ; cmp<so->val_t->get_n_comp() ; cmp ++) {
333 
334 
335 
336  Val_domain derr ((*(*so->val_t).cmp[cmp])(num_dom).der_var(1) / (*der_rad_term_eq->val_t)()(num_dom)) ;
337 
338  if (!derr.check_if_zero()) {
339  Val_domain res(this) ;
340  res.allocate_conf() ;
341  Index pos (nbr_points) ;
342  do {
343  res.set(pos) = (*(*inner_radius_term_eq->der_t).cmp[0])(num_dom)(pos) / 2. * (1-((*coloc[0])(pos(0)))) * derr(pos) ;
344  }
345  while (pos.inc()) ;
346  res.set_base() = (*(*so->val_t).cmp[cmp])(num_dom).get_base() ;
347  der.cmp[cmp]->set_domain(num_dom) = res ;
348  }
349  else
350  der.cmp[cmp]->set_domain(num_dom).set_zero() ;
351  }
352 
353  so->set_der_t (der) ;
354 }
355 
357  int dom = so.get_dom() ;
358  assert (dom == num_dom) ;
359 
360  int valence = so.val_t->get_valence() ;
361 
362  Array<int> type_ind (valence+1) ;
363  type_ind.set(0) = COV ;
364  for (int i=1 ; i<valence+1 ; i++)
365  type_ind.set(i) = so.val_t->get_index_type(i-1) ;
366 
367  Base_tensor basis (so.get_val_t().get_space()) ;
368  basis.set_basis(dom) = SPHERICAL_BASIS ;
369 
370  Term_eq comp_r (derive_r(so)) ;
371  Term_eq comp_t (derive_t(so) / (*rad_term_eq)) ;
372 
373  Tensor val_res (so.get_val_t().get_space(), valence+1, type_ind, basis) ;
374  {
375  Index pos_so (*so.val_t) ;
376  Index pos_res (val_res) ;
377  do {
378  for (int i=1 ; i<valence+1 ; i++)
379  pos_res.set(i) = pos_so(i-1) ;
380  // R part
381  pos_res.set(0) = 0 ;
382  val_res.set(pos_res).set_domain(num_dom) = (*comp_r.val_t)(pos_so)(num_dom) ;
383  // Theta part
384  pos_res.set(0) = 1 ;
385  val_res.set(pos_res).set_domain(num_dom) = (*comp_t.val_t)(pos_so)(num_dom) ;
386  }
387  while (pos_so.inc()) ;
388 }
389 
390  bool doder = ((so.der_t==0x0) || (inner_radius_term_eq->der_t==0x0)) ? false : true ;
391  if (doder) {
392  Tensor der_res (so.get_val_t().get_space(), valence+1, type_ind, basis) ;
393  Index pos_so (*so.der_t) ;
394  Index pos_res (val_res) ;
395  do {
396  for (int i=1 ; i<valence+1 ; i++)
397  pos_res.set(i) = pos_so(i-1) ;
398  // R part
399  pos_res.set(0) = 0 ;
400  der_res.set(pos_res).set_domain(num_dom) = (*comp_r.der_t)(pos_so)(num_dom) ;
401  // Theta part
402  pos_res.set(0) = 1 ;
403  der_res.set(pos_res).set_domain(num_dom) = (*comp_t.der_t)(pos_so)(num_dom) ;
404  }
405  while (pos_so.inc()) ;
406 
407  return Term_eq (num_dom, val_res, der_res) ;
408  }
409  else
410  return Term_eq (num_dom, val_res) ;
411 }
412 
413 
415  int dom = so.get_dom() ;
416  assert (dom == num_dom) ;
417 
418  int valence = so.val_t->get_valence() ;
419 
420  Array<int> type_ind (valence+1) ;
421  type_ind.set(0) = COV ;
422  for (int i=1 ; i<valence+1 ; i++)
423  type_ind.set(i) = so.val_t->get_index_type(i-1) ;
424 
425  Base_tensor basis (so.get_val_t().get_space()) ;
426  basis.set_basis(dom) = CARTESIAN_BASIS ;
427 
428  Term_eq comp_r (derive_r(so)) ;
429  Term_eq comp_t (derive_t(so) / (*rad_term_eq)) ;
430 
431  Tensor val_res (so.get_val_t().get_space(), valence+1, type_ind, basis) ;
432  {
433  Index pos_so (*so.val_t) ;
434  Index pos_res (val_res) ;
435  do {
436  for (int i=1 ; i<valence+1 ; i++)
437  pos_res.set(i) = pos_so(i-1) ;
438 
439  // X part
440  pos_res.set(0) = 0 ;
441  val_res.set(pos_res).set_domain(num_dom) = mult_sin_theta((*comp_r.val_t)(pos_so)(num_dom)) + mult_cos_theta((*comp_t.val_t)(pos_so)(num_dom)) ;
442  // Z part
443  pos_res.set(0) = 1 ;
444  val_res.set(pos_res).set_domain(num_dom) =
445  mult_cos_theta((*comp_r.val_t)(pos_so)(num_dom)) - mult_sin_theta((*comp_t.val_t)(pos_so)(num_dom)) ;
446  }
447  while (pos_so.inc()) ;
448 }
449 
450  bool doder = ((so.der_t==0x0) || (rad_term_eq->der_t==0x0)) ? false : true ;
451  if (doder) {
452  Tensor der_res (so.get_val_t().get_space(), valence+1, type_ind, basis) ;
453  Index pos_so (*so.der_t) ;
454  Index pos_res (val_res) ;
455  do {
456  for (int i=1 ; i<valence+1 ; i++)
457  pos_res.set(i) = pos_so(i-1) ;
458 
459  Val_domain auxi (mult_sin_theta((*comp_r.der_t)(pos_so)(num_dom)) + mult_cos_theta((*comp_t.der_t)(pos_so)(num_dom))) ;
460  // X part
461  pos_res.set(0) = 0 ;
462  der_res.set(pos_res).set_domain(num_dom) = mult_sin_theta((*comp_r.der_t)(pos_so)(num_dom)) + mult_cos_theta((*comp_t.der_t)(pos_so)(num_dom)) ;
463  // Z part
464  pos_res.set(0) = 1 ;
465  der_res.set(pos_res).set_domain(num_dom) =
466  mult_cos_theta((*comp_r.der_t)(pos_so)(num_dom)) - mult_sin_theta((*comp_t.der_t)(pos_so)(num_dom)) ;
467  }
468  while (pos_so.inc()) ;
469 
470  return Term_eq (num_dom, val_res, der_res) ;
471  }
472  else
473  return Term_eq (num_dom, val_res) ;
474 }
475 
477  return partial_cart(so) ;
478 }
479 
480 const Term_eq* Domain_polar_shell_inner_adapted::give_normal (int bound, int tipe) const {
481  assert (bound==INNER_BC) ;
482  switch (tipe) {
483  case CARTESIAN_BASIS :
484  if (normal_cart==0x0)
485  do_normal_cart() ;
486  return normal_cart ;
487 
488  case SPHERICAL_BASIS :
489  if (normal_spher==0x0)
490  do_normal_spher() ;
491  return normal_spher ;
492 
493  default:
494  cerr << "Unknown type of tensorial basis in Domain_polar_shell_inner_adapted::give_normal" << endl ;
495  abort() ;
496  }
497 }
498 
499 
500 
502  cerr << "Domain_polar_shell_inner_adapted::integ_term_eq not implemented yet" << endl ;
503  abort() ;
504 }
505 }
506 
reference set(const Index &pos)
Read/write of an element.
Definition: array.hpp:186
Describes the tensorial basis used by the various tensors.
Definition: base_tensor.hpp:49
int & set_basis(int nd)
Read/write the basis in a given domain.
Definition: base_tensor.hpp:91
virtual Term_eq lap2_term_eq(const Term_eq &, int) const
Returns the flat 2d-Laplacian of Term_eq, for a given harmonic.
virtual Term_eq der_normal_term_eq(const Term_eq &, int) const
Returns the normal derivative of a Term_eq.
Term_eq * der_rad_term_eq
Pointer on the Term_eq containing the .
void do_normal_cart() const
Computes the normal wrt the inner boundary, in Cartesian coordinates.
virtual void update_term_eq(Term_eq *) const
Update the value of a field, after the shape of the Domain has been changed by the system.
Term_eq * dt_rad_term_eq
Pointer on the Term_eq containing the .
Term_eq * rad_term_eq
Pointer on the Term_eq containing the radius.
Term_eq * inner_radius_term_eq
Pointer on the inner boundary , as a Term_eq.
Term_eq derive_t(const Term_eq &so) const
Computes .
const Space & sp
The corresponding Space ; required for updating fields whene the mapping changes.
virtual const Term_eq * give_normal(int, int) const
Returns the vector normal to a surface.
Term_eq * normal_cart
Pointer on the Term_eq containing the normal vector to the inner boundary, in Cartesian coordinates.
virtual Val_domain mult_sin_theta(const Val_domain &) const
Multiplication by .
virtual Term_eq dr_term_eq(const Term_eq &) const
Radial derivative of a Term_eq.
Term_eq derive_r(const Term_eq &so) const
Computes .
virtual Term_eq lap_term_eq(const Term_eq &, int) const
Returns the flat Laplacian of Term_eq, for a given harmonic.
virtual Term_eq div_r_term_eq(const Term_eq &) const
Division by of a Term_eq.
Term_eq flat_grad_spher(const Term_eq &so) const
Computes the flat gradient of a field, in orthonormal spherical coordinates.
virtual Term_eq partial_spher(const Term_eq &) const
Computes the part of the gradient containing the partial derivative of the field, in spherical orthon...
virtual Term_eq partial_cart(const Term_eq &) const
Computes the part of the gradient containing the partial derivative of the field, in Cartesian coordi...
virtual Term_eq mult_r_term_eq(const Term_eq &) const
Multiplication by of a Term_eq.
Term_eq * normal_spher
Pointer on the Term_eq containing the normal vector to the inner boundary, in orthonormal spherical c...
virtual Term_eq integ_term_eq(const Term_eq &, int) const
Surface integral of a Term_eq.
virtual Term_eq grad_term_eq(const Term_eq &) const
Gradient of Term_eq.
virtual Val_domain mult_cos_theta(const Val_domain &) const
Multiplication by .
void do_normal_spher() const
Computes the normal wrt the inner boundary, in orthonormal spherical coordinates.
Term_eq do_comp_by_comp(const Term_eq &so, Val_domain(Domain::*pfunc)(const Val_domain &) const) const
Function used to apply the same operation to all the components of a tensor, in the current domain.
Definition: domain.cpp:283
int num_dom
Number of the current domain (used by the Space)
Definition: space.hpp:63
virtual Val_domain div_sin_theta(const Val_domain &) const
Division by .
Definition: domain.cpp:1248
Dim_array nbr_points
Number of colocation points.
Definition: space.hpp:65
Memory_mapped_array< Array< double > * > coloc
Colocation points in each dimension (stored in ndim 1d- arrays)
Definition: space.hpp:75
virtual Val_domain mult_cos_theta(const Val_domain &) const
Multiplication by .
Definition: domain.cpp:1224
Class that gives the position inside a multi-dimensional Array.
Definition: index.hpp:38
int & set(int i)
Read/write of the position in a given dimension.
Definition: index.hpp:72
bool inc(int increm, int var=0)
Increments the position of the Index.
Definition: index.hpp:99
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
Tensor handling.
Definition: tensor.hpp:149
Memory_mapped_array< Scalar * > cmp
Array of size n_comp of pointers onto the components.
Definition: tensor.hpp:179
Scalar & set(const Array< int > &ind)
Returns the value of a component (read/write version).
Definition: tensor_impl.hpp:91
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
virtual Array< int > indices(int pos) const
Gives the values of the indices corresponding to a location in the array used for storage of the comp...
Definition: tensor.hpp:484
int get_valence() const
Returns the valence.
Definition: tensor.hpp:509
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
void set_der_t(Tensor)
Sets the tensorial variation (only the values in the pertinent Domain are copied).
Definition: term_eq.cpp:171
int get_dom() const
Definition: term_eq.hpp:135
Tensor const & get_val_t() const
Definition: term_eq.hpp:355
Tensor const & get_der_t() const
Definition: term_eq.hpp:369
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
bool check_if_zero() const
Check whether the logical state is zero or not.
Definition: val_domain.hpp:142
double & set(const Index &pos)
Read/write the value of the field in the configuration space.
Definition: val_domain.cpp:171
Base_spectral & set_base()
Sets the basis of decomposition.
Definition: val_domain.hpp:126
void allocate_conf()
Allocates the values in the configuration space and destroys the values in the coefficients space.
Definition: val_domain.cpp:209
A class derived from Tensor to deal specificaly with objects of valence 1 (and so also 1-forms).
Definition: vector.hpp:41
Scalar & set(int)
Read/write access to a component.