KADATH
domain_nucleus_symphi_for_metric.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 "spheric_symphi.hpp"
23 #include "term_eq.hpp"
24 #include "metric.hpp"
25 #include "scalar.hpp"
26 #include "tensor_impl.hpp"
27 namespace Kadath {
28 Term_eq Domain_nucleus_symphi::derive_flat_spher (int type_der, char ind_der, const Term_eq& so, const Metric* manipulator) const {
29 
30  assert ((type_der==COV) || (type_der==CON)) ;
31  int val_res = so.val_t->get_valence() + 1 ;
32 
33  bool donames = (ind_der==' ') ? false : true ;
34  bool need_sum = false ;
35 
36  Array<int> type_ind (val_res) ;
37  type_ind.set(0) = COV ;
38  for (int i=1 ; i<val_res ; i++)
39  type_ind.set(i) = so.val_t->get_index_type(i-1) ;
40 
41  if (donames) {
42  if (so.val_t->get_valence()>0)
43  assert (so.val_t->is_name_affected()) ;
44  // Need for summation ?
45  for (int i=1 ; i<val_res ; i++)
46  if (ind_der== so.val_t->get_name_ind()[i-1])
47  need_sum = true ;
48  }
49 
50  // Tensor for val
51  Base_tensor basis (so.val_t->get_space(), SPHERICAL_BASIS) ;
52  Tensor auxi_val (so.val_t->get_space(), val_res, type_ind, basis) ;
53 
54  if (donames) {
55  // Set the names of the indices :
56  auxi_val.set_name_affected() ;
57  auxi_val.set_name_ind(0, ind_der) ;
58  for (int i=1 ; i<val_res ; i++)
59  auxi_val.set_name_ind(i, so.val_t->get_name_ind()[i-1]) ;
60  }
61 
62  // Part derivative
63  //Loop on the components :
64  Index pos_auxi_bis(auxi_val) ;
65  Index pos_so_bis (*so.val_t) ;
66  do {
67  for (int i=0 ; i<val_res-1 ; i++)
68  pos_so_bis.set(i) = pos_auxi_bis(i+1) ;
69  switch (pos_auxi_bis(0)) {
70  case 0 :
71  // d/dr :
72  auxi_val.set(pos_auxi_bis).set_domain(num_dom) = (*so.val_t)(pos_so_bis)(num_dom).der_r() ;
73  break ;
74  case 1 :
75  // 1/r dtheta
76  auxi_val.set(pos_auxi_bis).set_domain(num_dom) = (*so.val_t)(pos_so_bis)(num_dom).der_var(2).div_r() ;
77  break ;
78  case 2 :
79  // 1/r sint d/dphi
80  auxi_val.set(pos_auxi_bis).set_domain(num_dom) = (*so.val_t)(pos_so_bis)(num_dom).der_var(3).div_r().div_sin_theta() ;
81  break ;
82  default :
83  cerr << "Bad indice in Domain_nucleus_symphi::derive_flat_spher" << endl ;
84  abort() ;
85  }
86  }
87  while (pos_auxi_bis.inc()) ;
88 
89  // Loop indice summation on connection symbols
90  for (int ind_sum=0 ; ind_sum<val_res-1 ; ind_sum++) {
91 
92  //Loop on the components :
93  Index pos_auxi(auxi_val) ;
94  Index pos_so (*so.val_t) ;
95 
96  do {
97  for (int i=0 ; i<val_res-1 ; i++)
98  pos_so.set(i) = pos_auxi(i+1) ;
99  // Different cases of the derivative index :
100  switch (pos_auxi(0)) {
101  case 0 :
102  // Dr nothing
103  break ;
104  case 1 :
105  // Dtheta
106  // Different cases of the source index
107  switch (pos_auxi(ind_sum+1)) {
108  case 0 :
109  //Dtheta S_r
110  pos_so.set(ind_sum) = 1 ;
111  auxi_val.set(pos_auxi).set_domain(num_dom) -= (*so.val_t)(pos_so)(num_dom).div_r() ;
112  break ;
113  case 1 :
114  // Dtheta S_theta
115  pos_so.set(ind_sum) = 0 ;
116  auxi_val.set(pos_auxi).set_domain(num_dom) += (*so.val_t)(pos_so)(num_dom).div_r() ;
117  break ;
118  case 2 :
119  //Dtheta S_phi
120  break ;
121  default :
122  cerr << "Bad indice in Domain_nucleus_symphi::derive_flat_spher" << endl ;
123  abort() ;
124  }
125  break ;
126  case 2 :
127  // Dphi
128  // Different cases of the source index
129  switch (pos_auxi(ind_sum+1)) {
130  case 0 :
131  //Dphi S_r
132  pos_so.set(ind_sum) = 2 ;
133  auxi_val.set(pos_auxi).set_domain(num_dom) -= (*so.val_t)(pos_so)(num_dom).div_r() ;
134  break ;
135  case 1 :
136  // Dphi S_theta
137  pos_so.set(ind_sum) = 2 ;
138  auxi_val.set(pos_auxi).set_domain(num_dom) -= (*so.val_t)(pos_so)(num_dom).div_r().mult_cos_theta().div_sin_theta() ;
139  break ;
140  case 2 :
141  //Dphi S_phi
142  pos_so.set(ind_sum) = 0 ;
143  auxi_val.set(pos_auxi).set_domain(num_dom) += (*so.val_t)(pos_so)(num_dom).div_r() ;
144  pos_so.set(ind_sum) = 1 ;
145  auxi_val.set(pos_auxi).set_domain(num_dom) += (*so.val_t)(pos_so)(num_dom).div_r().mult_cos_theta().div_sin_theta() ;
146  break ;
147  default :
148  cerr << "Bad indice in Domain_nucleus_symphi::derive_flat_spher" << endl ;
149  abort() ;
150  }
151  break ;
152  default :
153  cerr << "Bad indice in Domain_nucleus_symphi::derive_flat_spher" << endl ;
154  abort() ;
155  }
156  }
157  while (pos_auxi.inc()) ;
158  }
159 
160  if (so.der_t==0x0) {
161  // No need for derivative :
162  Term_eq auxi (num_dom, auxi_val) ;
163  // If derive contravariant : manipulate first indice :
164  if (type_der==CON)
165  manipulator->manipulate_ind (auxi, 0) ;
166 
167  if (!need_sum)
168  return auxi ;
169  else
171  }
172  else {
173  // Need to compute the derivative :
174  // Tensor for der
175  Tensor auxi_der (so.val_t->get_space(), val_res, type_ind, basis) ;
176 
177  if (donames) {
178  // Set the names of the indices :
179  auxi_der.set_name_affected() ;
180  auxi_der.set_name_ind(0, ind_der) ;
181  for (int i=1 ; i<val_res ; i++)
182  auxi_der.set_name_ind(i, so.der_t->get_name_ind()[i-1]) ;
183  }
184 
185  // Part derivative
186  //Loop on the components :
187  Index pos_auxi_der_bis(auxi_der) ;
188  do {
189  for (int i=0 ; i<val_res-1 ; i++)
190  pos_so_bis.set(i) = pos_auxi_der_bis(i+1) ;
191  switch (pos_auxi_der_bis(0)) {
192  case 0 :
193  // d/dr :
194  auxi_der.set(pos_auxi_der_bis).set_domain(num_dom) = (*so.der_t)(pos_so_bis)(num_dom).der_r() ;
195  break ;
196  case 1 :
197  // 1/r dtheta
198  auxi_der.set(pos_auxi_der_bis).set_domain(num_dom) = (*so.der_t)(pos_so_bis)(num_dom).der_var(2).div_r() ;
199  break ;
200  case 2 :
201  // 1/r sint d/dphi
202  auxi_der.set(pos_auxi_der_bis).set_domain(num_dom) = (*so.der_t)(pos_so_bis)(num_dom).der_var(3).div_r().div_sin_theta() ;
203  break ;
204  default :
205  cerr << "Bad indice in Domain_nucleus_symphi::derive_flat_spher" << endl ;
206  abort() ;
207  }
208  }
209  while (pos_auxi_der_bis.inc()) ;
210 
211  // Loop indice summation on connection symbols
212  for (int ind_sum=0 ; ind_sum<val_res-1 ; ind_sum++) {
213 
214  //Loop on the components :
215  Index pos_auxi_der(auxi_val) ;
216  Index pos_so (*so.val_t) ;
217 
218  do {
219  for (int i=0 ; i<val_res-1 ; i++)
220  pos_so.set(i) = pos_auxi_der(i+1) ;
221  // Different cases of the derivative index :
222  switch (pos_auxi_der(0)) {
223  case 0 :
224  // Dr nothing
225  break ;
226  case 1 :
227  // Dtheta
228  // Different cases of the source index
229  switch (pos_auxi_der(ind_sum+1)) {
230  case 0 :
231  //Dtheta S_r
232  pos_so.set(ind_sum) = 1 ;
233  auxi_der.set(pos_auxi_der).set_domain(num_dom) -= (*so.der_t)(pos_so)(num_dom).div_r() ;
234  break ;
235  case 1 :
236  // Dtheta S_theta
237  pos_so.set(ind_sum) = 0 ;
238  auxi_der.set(pos_auxi_der).set_domain(num_dom) += (*so.der_t)(pos_so)(num_dom).div_r() ;
239  break ;
240  case 2 :
241  //Dtheta S_phi
242  break ;
243  default :
244  cerr << "Bad indice in Domain_nucleus_symphi::derive_flat_spher" << endl ;
245  abort() ;
246  }
247  break ;
248  case 2 :
249  // Dphi
250  // Different cases of the source index
251  switch (pos_auxi_der(ind_sum+1)) {
252  case 0 :
253  //Dphi S_r
254  pos_so.set(ind_sum) = 2 ;
255  auxi_der.set(pos_auxi_der).set_domain(num_dom) -= (*so.der_t)(pos_so)(num_dom).div_r() ;
256  break ;
257  case 1 :
258  // Dphi S_theta
259  pos_so.set(ind_sum) = 2 ;
260  auxi_der.set(pos_auxi_der).set_domain(num_dom) -= (*so.der_t)(pos_so)(num_dom).div_r().mult_cos_theta().div_sin_theta() ;
261  break ;
262  case 2 :
263  //Dphi S_phi
264  pos_so.set(ind_sum) = 0 ;
265  auxi_der.set(pos_auxi_der).set_domain(num_dom) += (*so.der_t)(pos_so)(num_dom).div_r() ;
266  pos_so.set(ind_sum) = 1 ;
267  auxi_der.set(pos_auxi_der).set_domain(num_dom) += (*so.der_t)(pos_so)(num_dom).div_r().mult_cos_theta().div_sin_theta() ;
268  break ;
269  default :
270  cerr << "Bad indice in Domain_nucleus_symphi::derive_flat_spher" << endl ;
271  abort() ;
272  }
273  break ;
274  default :
275  cerr << "Bad indice in Domain_nucleus_symphi::derive_flat_spher" << endl ;
276  abort() ;
277  }
278  }
279  while (pos_auxi_der.inc()) ;
280  }
281 
282 
283  // Need for derivative :
284  Term_eq auxi (num_dom, auxi_val, auxi_der) ;
285  // If derive contravariant : manipulate first indice :
286  if (type_der==CON)
287  manipulator->manipulate_ind (auxi, 0) ;
288  if (!need_sum)
289  return auxi ;
290  else
293  }
294 }
295 
296 Term_eq Domain_nucleus_symphi::derive_flat_cart (int type_der, char ind_der, const Term_eq& so, const Metric* manipulator) const {
297 
298  bool doder = (so.der_t==0x0) ? false : true ;
299  if ((type_der==CON) && (manipulator->p_met_con[num_dom]->der_t==0x0))
300  doder = false ;
301 
302  assert ((type_der==COV) || (type_der==CON)) ;
303  int val_res = so.val_t->get_valence() + 1 ;
304 
305  bool doname = true ;
306  if (so.val_t->get_valence()>0)
307  if (!so.val_t->is_name_affected())
308  doname = false;
309 
310  Array<int> type_ind (val_res) ;
311  type_ind.set(0) = COV ;
312  for (int i=1 ; i<val_res ; i++)
313  type_ind.set(i) = so.val_t->get_index_type(i-1) ;
314 
315  // Need for summation ?
316  bool need_sum = false ;
317  if (doname)
318  for (int i=1 ; i<val_res ; i++)
319  if (ind_der== so.val_t->get_name_ind()[i-1])
320  need_sum = true ;
321 
322  // Tensor for val
323  Base_tensor basis (so.val_t->get_space(), CARTESIAN_BASIS) ;
324  Tensor auxi_val (so.val_t->get_space(), val_res, type_ind, basis) ;
325 
326 
327  // Set the names of the indices :
328  if (doname) {
329  auxi_val.set_name_affected() ;
330  auxi_val.set_name_ind(0, ind_der) ;
331  for (int i=1 ; i<val_res ; i++)
332  auxi_val.set_name_ind(i, so.val_t->get_name_ind()[i-1]) ;
333  }
334 
335  //Loop on the components :
336  Index pos_auxi(auxi_val) ;
337  Index pos_so (*so.val_t) ;
338  do {
339  for (int i=0 ; i<val_res-1 ; i++)
340  pos_so.set(i) = pos_auxi(i+1) ;
341  auxi_val.set(pos_auxi).set_domain(num_dom) = (*so.val_t)(pos_so)(num_dom).der_abs(pos_auxi(0)+1) ;
342  }
343  while (pos_auxi.inc()) ;
344 
345  if (!doder) {
346  // No need for derivative :
347  Term_eq auxi (num_dom, auxi_val) ;
348  // If derive contravariant : manipulate first indice :
349  if (type_der==CON)
350  manipulator->manipulate_ind (auxi, 0) ;
351 
352  if (!need_sum)
353  return auxi ;
354  else
356  }
357  else {
358  // Need to compute the derivative :
359  // Tensor for der
360  Tensor auxi_der (so.val_t->get_space(), val_res, type_ind, basis) ;
361  // Set the names of the indices :
362  auxi_der.set_name_affected() ;
363  auxi_der.set_name_ind(0, ind_der) ;
364  for (int i=1 ; i<val_res ; i++)
365  auxi_der.set_name_ind(i, so.der_t->get_name_ind()[i-1]) ;
366 
367  //Loop on the components :
368  Index pos_auxi_der(auxi_der) ;
369  do {
370  for (int i=0 ; i<val_res-1 ; i++)
371  pos_so.set(i) = pos_auxi_der(i+1) ;
372  auxi_der.set(pos_auxi_der).set_domain(num_dom) = (*so.der_t)(pos_so)(num_dom).der_abs(pos_auxi_der(0)+1) ;
373  }
374  while (pos_auxi_der.inc()) ;
375 
376  // Need for derivative :
377  Term_eq auxi (num_dom, auxi_val, auxi_der) ;
378 
379  // If derive contravariant : manipulate first indice :
380  if (type_der==CON)
381  manipulator->manipulate_ind (auxi, 0) ;
382 
383  if (!need_sum)
384  return auxi ;
385  else
388  }
389 
390 }}
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
virtual Term_eq derive_flat_cart(int, char, const Term_eq &, const Metric *) const
Computes the flat derivative of a Term_eq, in Cartesian coordinates.
virtual Term_eq derive_flat_spher(int, char, const Term_eq &, const Metric *) const
Computes the flat derivative of a Term_eq, in spherical orthonormal coordinates.
virtual Val_domain der_r(const Val_domain &) const
Compute the radial derivative of a scalar field.
virtual Val_domain div_r(const Val_domain &) const
Division by .
int num_dom
Number of the current domain (used by the Space)
Definition: space.hpp:63
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
Purely abstract class for metric handling.
Definition: metric.hpp:39
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
Val_domain & set_domain(int)
Read/write of a particular Val_domain.
Definition: scalar.hpp:555
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_valence() const
Returns the valence.
Definition: tensor.hpp:509
bool is_name_affected() const
Check whether the names of the indices have been affected.
Definition: tensor.hpp:429
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
Tensor * val_t
Pointer on the value, if the Term_eq is a Tensor.
Definition: term_eq.hpp:68
Val_domain div_r() const
Division by the radius.
Definition: val_domain.cpp:743
Val_domain div_sin_theta() const
Division by .
Val_domain mult_cos_theta() const
Multiplication by .