KADATH
domain_shell_outer_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.hpp"
23 #include "point.hpp"
24 #include "array_math.hpp"
25 #include "scalar.hpp"
26 #include "tensor_impl.hpp"
27 #include "vector.hpp"
28 
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 
47  if (doder) {
48  Tensor der_res (*so.der_t, false) ;
49  for (int cmp = 0 ; cmp<so.val_t->get_n_comp() ; cmp++) {
50  Array<int>ind (val_res.indices(cmp)) ;
51  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)
52  - (*so.val_t)(ind)(num_dom).der_var(1)*(*der_rad_term_eq->der_t)()(num_dom)) /
54  }
55  return Term_eq (num_dom, val_res, der_res) ;
56  }
57  else
58  return Term_eq (num_dom, val_res) ;
59 }
60 
62 
63  assert (so.get_dom() == num_dom) ;
64  Term_eq* dtprime ;
65  Tensor val_res(*so.val_t, false) ;
66 
67  for (int cmp = 0 ; cmp<so.val_t->get_n_comp() ; cmp++) {
68  Array<int>ind (val_res.indices(cmp)) ;
69  val_res.set(ind).set_domain(num_dom) = (*so.val_t)(ind)(num_dom).der_var(2) ;
70  }
71  bool doder = (so.der_t==0x0) ? false : true ;
72 
73  if (doder) {
74  Tensor der_res (*so.der_t, false) ;
75  for (int cmp = 0 ; cmp<so.val_t->get_n_comp() ; cmp++) {
76  Array<int>ind (val_res.indices(cmp)) ;
77  der_res.set(ind).set_domain(num_dom) = (*so.der_t)(ind)(num_dom).der_var(2) ;
78  }
79  dtprime = new Term_eq (num_dom, val_res, der_res) ;
80  }
81  else
82  dtprime = new Term_eq (num_dom, val_res) ;
83 
84 
85  Term_eq res (*dtprime - (*dt_rad_term_eq) * derive_r(so)) ;
86  delete dtprime ;
87  return res ;
88 
89 }
90 
92 
93  Term_eq* dpprime ;
94  Tensor val_res(*so.val_t, false) ;
95 
96  for (int cmp = 0 ; cmp<so.val_t->get_n_comp() ; cmp++) {
97  Array<int>ind (val_res.indices(cmp)) ;
98  val_res.set(ind).set_domain(num_dom) = (*so.val_t)(ind)(num_dom).der_var(3) ;
99  }
100  bool doder = (so.der_t==0x0) ? false : true ;
101 
102  if (doder) {
103  Tensor der_res (*so.der_t, false) ;
104  for (int cmp = 0 ; cmp<so.val_t->get_n_comp() ; cmp++) {
105  Array<int>ind (val_res.indices(cmp)) ;
106  der_res.set(ind).set_domain(num_dom) = (*so.der_t)(ind)(num_dom).der_var(3) ;
107  }
108  dpprime = new Term_eq (num_dom, val_res, der_res) ;
109  }
110  else
111  dpprime = new Term_eq (num_dom, val_res) ;
112 
113 
114  Term_eq res (*dpprime - (*dp_rad_term_eq) * derive_r(so)) ;
115  delete dpprime ;
116  return res ;
117 
118 }
119 
120 
122 
123  assert (so.get_dom() == num_dom) ;
124 
125  int valso = so.get_val_t().get_valence() ;
126  Array<int> type_ind (valso+1) ;
127  type_ind.set(0) = COV ;
128  for (int i=0 ; i<valso ; i++)
129  type_ind.set(i+1) = so.get_val_t().get_index_type(i) ;
130 
131  Base_tensor basis (sp) ;
132  basis.set_basis(num_dom) = SPHERICAL_BASIS ;
133  Tensor res (sp, valso+1 , type_ind, basis) ;
134 
135  Term_eq comp_r (derive_r(so)) ;
136  Term_eq comp_t (derive_t(so)/(*rad_term_eq)) ;
138 
139 
140  // Loop on cmp :
141  for (int nc=0 ; nc<so.get_val_t().get_n_comp() ; nc++) {
142 
143  Array<int> ind (so.get_val_t().indices(nc)) ;
144  Array<int> indtarget (valso+1) ;
145  for (int i=0 ; i<valso ; i++)
146  indtarget.set(i+1) = ind(i) ;
147 
148  // R comp :
149  indtarget.set(0) = 1 ;
150  res.set(indtarget).set_domain(num_dom) = comp_r.get_val_t()(ind)(num_dom);
151  // theta comp :
152  indtarget.set(0) = 2 ;
153  res.set(indtarget).set_domain(num_dom) = comp_t.get_val_t()(ind)(num_dom) ;
154  // Phi comp :
155  indtarget.set(0) = 3 ;
156  res.set(indtarget).set_domain(num_dom) = comp_p.get_val_t()(ind)(num_dom) ;
157  }
158  bool doder = ((so.der_t==0x0) || (outer_radius_term_eq->der_t==0x0)) ? false : true ;
159 
160  if (!doder) {
161  return Term_eq (num_dom, res) ;
162  }
163  else {
164 
165  Tensor resder (sp, valso+1 , type_ind, basis) ;
166  // Loop on cmp :
167  for (int nc=0 ; nc<so.get_val_t().get_n_comp() ; nc++) {
168 
169  Array<int> ind (so.get_val_t().indices(nc)) ;
170  Array<int> indtarget (valso+1) ;
171  for (int i=0 ; i<valso ; i++)
172  indtarget.set(i+1) = ind(i) ;
173 
174  // R comp :
175  indtarget.set(0) = 1 ;
176  resder.set(indtarget).set_domain(num_dom) = comp_r.get_der_t()(ind)(num_dom);
177  // theta comp :
178  indtarget.set(0) = 2 ;
179  resder.set(indtarget).set_domain(num_dom) = comp_t.get_der_t()(ind)(num_dom) ;
180  // Phi comp :
181  indtarget.set(0) = 3 ;
182  resder.set(indtarget).set_domain(num_dom) = comp_p.get_der_t()(ind)(num_dom) ;
183  }
184  return Term_eq (num_dom, res, resder) ;
185  }
186 }
187 
189 
191 
192  Scalar val_norme (sp) ;
193  val_norme.set_domain(num_dom) = sqrt((*grad.val_t)(1)(num_dom)*(*grad.val_t)(1)(num_dom)
194  + (*grad.val_t)(2)(num_dom)*(*grad.val_t)(2)(num_dom)
195  + (*grad.val_t)(3)(num_dom)*(*grad.val_t)(3)(num_dom)) ;
196  bool doder = ((grad.der_t==0x0) || (outer_radius_term_eq->der_t==0x0)) ? false : true ;
197  if (doder) {
198  Scalar der_norme (sp) ;
199  der_norme.set_domain(num_dom) = ((*grad.der_t)(1)(num_dom)*(*grad.val_t)(1)(num_dom)
200  + (*grad.der_t)(2)(num_dom)*(*grad.val_t)(2)(num_dom)
201  + (*grad.der_t)(3)(num_dom)*(*grad.val_t)(3)(num_dom)) / val_norme(num_dom) ;
202  Term_eq norme (num_dom, val_norme, der_norme) ;
203  if (normal_spher==0x0)
204  normal_spher = new Term_eq (grad / norme) ;
205  else
206  *normal_spher = Term_eq (grad/norme) ;
207  }
208  else {
209  Term_eq norme (num_dom, val_norme) ;
210  if (normal_spher==0x0)
211  normal_spher = new Term_eq (grad / norme) ;
212  else
213  *normal_spher = Term_eq (grad/norme) ;
214  }
215 }
216 
217 
218 
220 
221 
222  do_normal_spher() ;
223 
224  Base_tensor basis (sp) ;
225  basis.set_basis(num_dom) = CARTESIAN_BASIS ;
226  Vector val (sp, CON, basis) ;
227 
233 
234  bool doder = (normal_spher->der_t==0x0) ? false : true ;
235  if (doder) {
236  Vector der (sp, CON, basis) ;
237 
243 
244  if (normal_cart==0x0)
245  normal_cart = new Term_eq (num_dom, val, der) ;
246  else
247  *normal_cart = Term_eq(num_dom, val,der) ;
248  }
249  else {
250  if (normal_cart==0x0)
251  normal_cart = new Term_eq (num_dom, val) ;
252  else
253  *normal_cart = Term_eq(num_dom, val) ;
254  }
255 }
256 
258 
259 
260  switch (bound) {
261  case OUTER_BC : {
262  // Deformed surface
263  if (normal_spher == 0x0)
264  do_normal_spher() ;
265 
266  int valso = so.get_val_t().get_valence() ;
267 
268  Term_eq grad (flat_grad_spher(so)) ;
269 
270  Tensor res (so.get_val_t(), false) ;
271 
272  // Loop on cmp :
273  for (int nc=0 ; nc<so.get_val_t().get_n_comp() ; nc++) {
274 
275  Array<int> ind (so.get_val_t().indices(nc)) ;
276  Array<int> indgrad (valso+1) ;
277  for (int i=0 ; i<valso ; i++)
278  indgrad.set(i+1) = ind(i) ;
279 
280  indgrad.set(0) = 1 ;
281  res.set(ind).set_domain(num_dom) = (*grad.val_t)(indgrad)(num_dom) * (*normal_spher->val_t)(1)(num_dom) ;
282  indgrad.set(0) = 2 ;
283  res.set(ind).set_domain(num_dom) += (*grad.val_t)(indgrad)(num_dom) * (*normal_spher->val_t)(2)(num_dom) ;
284  indgrad.set(0) = 3 ;
285  res.set(ind).set_domain(num_dom) += (*grad.val_t)(indgrad)(num_dom) * (*normal_spher->val_t)(3)(num_dom) ;
286  }
287 
288  bool doder = ((grad.der_t == 0x0) || (normal_spher->der_t == 0x0)) ? false : true ;
289  if (doder) {
290 
291  Tensor der (so.get_val_t(), false) ;
292 
293  // Loop on cmp :
294  for (int nc=0 ; nc<so.get_val_t().get_n_comp() ; nc++) {
295 
296  Array<int> ind (so.get_val_t().indices(nc)) ;
297  Array<int> indgrad (valso+1) ;
298  for (int i=0 ; i<valso ; i++)
299  indgrad.set(i+1) = ind(i) ;
300 
301  indgrad.set(0) = 1 ;
302  der.set(ind).set_domain(num_dom) = (*grad.der_t)(indgrad)(num_dom) * (*normal_spher->val_t)(1)(num_dom)
303  + (*grad.val_t)(indgrad)(num_dom) * (*normal_spher->der_t)(1)(num_dom) ;
304  indgrad.set(0) = 2 ;
305  der.set(ind).set_domain(num_dom) += (*grad.der_t)(indgrad)(num_dom) * (*normal_spher->val_t)(2)(num_dom)
306  + (*grad.val_t)(indgrad)(num_dom) * (*normal_spher->der_t)(2)(num_dom) ;
307  indgrad.set(0) = 3 ;
308  der.set(ind).set_domain(num_dom) += (*grad.der_t)(indgrad)(num_dom) * (*normal_spher->val_t)(3)(num_dom)
309  + (*grad.val_t)(indgrad)(num_dom) * (*normal_spher->der_t)(3)(num_dom) ;
310  }
311 
312 
313  return Term_eq (num_dom, res, der) ;
314  }
315  else
316  return Term_eq (num_dom, res) ;
317  }
318  case INNER_BC : {
319  return derive_r(so) ;
320  }
321  default :
322  cerr << "Unknown boundary in Domain_shell_outer_adapted::der_normal" << endl ;
323  abort() ;
324  }
325 }
326 
328 #ifndef REMOVE_ALL_CHECKS
329  if (mm!=0) {
330  cerr << "Domain_shell_outer_adapted::lap_term_eq not defined for m != 0 (for now)" << endl ;
331  abort() ;
332  }
333 
334  if (so.get_val_t().get_valence() != 0) {
335  cerr << "Domain_shell_outer_adapted::lap_term_eq only defined for scalars" << endl ;
336  abort() ;
337  }
338 #endif
339 
340  // Angular part :
341  Term_eq dert (derive_t(so));
342  Term_eq p1 (derive_t(dert)) ;
344  Term_eq der2p (derive_p(derive_p(so))) ;
346 
347  Term_eq dr(derive_r(so)) ;
348 
349  Term_eq res (derive_r(dr) + 2 *dr / (*rad_term_eq) + (p1+p2+p3)/(*rad_term_eq)/(*rad_term_eq));
350 
351 
352  return res ;
353 }
354 
356 
357 
358  return so * (*rad_term_eq) ;
359 }
360 
362 
363  Tensor der (*so->val_t, false) ;
364  for (int cmp=0 ; cmp<so->val_t->get_n_comp() ; cmp ++) {
365 
366 
367 
368  Val_domain derr ((*(*so->val_t).cmp[cmp])(num_dom).der_var(1) / (*der_rad_term_eq->val_t)()(num_dom)) ;
369 
370  if (!derr.check_if_zero()) {
371  Val_domain res(this) ;
372  res.allocate_conf() ;
373  Index pos (nbr_points) ;
374  do {
375  res.set(pos) = (*(*outer_radius_term_eq->der_t).cmp[0])(num_dom)(pos) / 2. * (1+((*coloc[0])(pos(0)))) * derr(pos) ;
376  }
377  while (pos.inc()) ;
378  res.set_base() = (*(*so->val_t).cmp[cmp])(num_dom).get_base() ;
379  der.cmp[cmp]->set_domain(num_dom) = res ;
380  }
381  else
382  der.cmp[cmp]->set_domain(num_dom).set_zero() ;
383  }
384 
385  so->set_der_t (der) ;
386 }
387 
389  int dom = so.get_dom() ;
390  assert (dom == num_dom) ;
391 
392  int valence = so.val_t->get_valence() ;
393 
394  Array<int> type_ind (valence+1) ;
395  type_ind.set(0) = COV ;
396  for (int i=1 ; i<valence+1 ; i++)
397  type_ind.set(i) = so.val_t->get_index_type(i-1) ;
398 
399  Base_tensor basis (so.get_val_t().get_space()) ;
400  basis.set_basis(dom) = SPHERICAL_BASIS ;
401 
402  Term_eq comp_r (derive_r(so)) ;
403  Term_eq comp_t (derive_t(so) / (*rad_term_eq)) ;
405 
406  Tensor val_res (so.get_val_t().get_space(), valence+1, type_ind, basis) ;
407  {
408  Index pos_so (*so.val_t) ;
409  Index pos_res (val_res) ;
410  do {
411  for (int i=1 ; i<valence+1 ; i++)
412  pos_res.set(i) = pos_so(i-1) ;
413  // R part
414  pos_res.set(0) = 0 ;
415  val_res.set(pos_res).set_domain(num_dom) = (*comp_r.val_t)(pos_so)(num_dom) ;
416  // Theta part
417  pos_res.set(0) = 1 ;
418  val_res.set(pos_res).set_domain(num_dom) = (*comp_t.val_t)(pos_so)(num_dom) ;
419  // Phi part
420  pos_res.set(0) = 2 ;
421  val_res.set(pos_res).set_domain(num_dom) = (*comp_p.val_t)(pos_so)(num_dom) ;
422  }
423  while (pos_so.inc()) ;
424 }
425 
426  bool doder = ((so.der_t==0x0) || (outer_radius_term_eq->der_t==0x0)) ? false : true ;
427  if (doder) {
428  Tensor der_res (so.get_val_t().get_space(), valence+1, type_ind, basis) ;
429  Index pos_so (*so.der_t) ;
430  Index pos_res (val_res) ;
431  do {
432  for (int i=1 ; i<valence+1 ; i++)
433  pos_res.set(i) = pos_so(i-1) ;
434  // R part
435  pos_res.set(0) = 0 ;
436  der_res.set(pos_res).set_domain(num_dom) = (*comp_r.der_t)(pos_so)(num_dom) ;
437  // Theta part
438  pos_res.set(0) = 1 ;
439  der_res.set(pos_res).set_domain(num_dom) = (*comp_t.der_t)(pos_so)(num_dom) ;
440  // Phi part
441  pos_res.set(0) = 2 ;
442  der_res.set(pos_res).set_domain(num_dom) = (*comp_p.der_t)(pos_so)(num_dom) ;
443  }
444  while (pos_so.inc()) ;
445 
446  return Term_eq (num_dom, val_res, der_res) ;
447  }
448  else
449  return Term_eq (num_dom, val_res) ;
450 }
451 
452 
454 
455  int dom = so.get_dom() ;
456  assert (dom == num_dom) ;
457 
458  int valence = so.val_t->get_valence() ;
459  int val_res = so.val_t->get_valence() + 1 ;
460 
461  Array<int> type_ind (val_res) ;
462  type_ind.set(0) = COV ;
463  for (int i=1 ; i<val_res ; i++)
464  type_ind.set(i) = so.val_t->get_index_type(i-1) ;
465 
466  Base_tensor basis (so.get_val_t().get_space()) ;
467  basis.set_basis(dom) = SPHERICAL_BASIS ;
468 
469  Tensor auxi_val (so.get_val_t().get_space(), val_res, type_ind, basis) ;
470  for (int cmp=0 ; cmp<auxi_val.get_n_comp() ; cmp ++)
471  auxi_val.set(auxi_val.indices(cmp)) = 0 ;
472 
473  for (int ind_sum=0 ; ind_sum<valence ; ind_sum++) {
474 
475  //Loop on the components :
476  Index pos_auxi(auxi_val) ;
477  Index pos_so (*so.val_t) ;
478 
479  do {
480  for (int i=0 ; i<valence ; i++)
481  pos_so.set(i) = pos_auxi(i+1) ;
482  // Different cases of the derivative index :
483  switch (pos_auxi(0)) {
484  case 0 :
485  // Dr nothing
486  break ;
487  case 1 :
488  // Dtheta
489  // Different cases of the source index
490  switch (pos_auxi(ind_sum+1)) {
491  case 0 :
492  //Dtheta S_r
493  pos_so.set(ind_sum) = 1 ;
494  auxi_val.set(pos_auxi).set_domain(dom) -= (*so.val_t)(pos_so)(dom).div_r() ;
495  break ;
496  case 1 :
497  // Dtheta S_theta
498  pos_so.set(ind_sum) = 0 ;
499  auxi_val.set(pos_auxi).set_domain(dom) += (*so.val_t)(pos_so)(dom).div_r() ;
500  break ;
501  case 2 :
502  //Dtheta S_phi
503  break ;
504  default :
505  cerr << "Bad indice in Domain_shell_outer_adapted::connection_spher" << endl ;
506  abort() ;
507  }
508  break ;
509  case 2 :
510  // Dphi
511  // Different cases of the source index
512  switch (pos_auxi(ind_sum+1)) {
513  case 0 :
514  //Dphi S_r
515  pos_so.set(ind_sum) = 2 ;
516  auxi_val.set(pos_auxi).set_domain(dom) -= (*so.val_t)(pos_so)(dom).div_r() ;
517  break ;
518  case 1 :
519  // Dphi S_theta
520  pos_so.set(ind_sum) = 2 ;
521  auxi_val.set(pos_auxi).set_domain(dom) -= (*so.val_t)(pos_so)(dom).div_r().mult_cos_theta().div_sin_theta() ;
522  break ;
523  case 2 :
524  //Dphi S_phi
525  pos_so.set(ind_sum) = 0 ;
526  auxi_val.set(pos_auxi).set_domain(dom) += (*so.val_t)(pos_so)(dom).div_r() ;
527  pos_so.set(ind_sum) = 1 ;
528  auxi_val.set(pos_auxi).set_domain(dom) += (*so.val_t)(pos_so)(dom).div_r().mult_cos_theta().div_sin_theta() ;
529  break ;
530  default :
531  cerr << "Bad indice in Domain_shell_outer_adapted::connection_spher" << endl ;
532  abort() ;
533  }
534  break ;
535  default :
536  cerr << "Bad indice in Domain_shell_outer_adapted::connection_spher" << endl ;
537  abort() ;
538  }
539  }
540  while (pos_auxi.inc()) ;
541  }
542 
543 
544  if ((so.der_t==0x0) || (outer_radius_term_eq->der_t==0x0)) {
545  // No need for derivative :
546  return Term_eq (dom, auxi_val) ;
547  }
548  else {
549 
550  // Need to compute the derivative :
551  // Tensor for der
552  Tensor auxi_der (so.get_val_t().get_space(), val_res, type_ind, basis) ;
553  for (int cmp=0 ; cmp<auxi_der.get_n_comp() ; cmp ++)
554  auxi_der.set(auxi_der.indices(cmp)) = 0 ;
555 
556  // Loop indice summation on connection symbols
557  for (int ind_sum=0 ; ind_sum<valence ; ind_sum++) {
558 
559  //Loop on the components :
560  Index pos_auxi_der(auxi_der) ;
561  Index pos_so (*so.der_t) ;
562 
563  do {
564  for (int i=0 ; i<valence ; i++)
565  pos_so.set(i) = pos_auxi_der(i+1) ;
566  // Different cases of the derivative index :
567  switch (pos_auxi_der(0)) {
568  case 0 :
569  // Dr nothing
570  break ;
571  case 1 :
572  // Dtheta
573  // Different cases of the source index
574  switch (pos_auxi_der(ind_sum+1)) {
575  case 0 :
576  //Dtheta S_r
577  pos_so.set(ind_sum) = 1 ;
578  auxi_der.set(pos_auxi_der).set_domain(dom) -=
579  ((*so.der_t)(pos_so)(dom) -(*so.val_t)(pos_so)(dom)*(*rad_term_eq->der_t)()(dom)/(*rad_term_eq->val_t)()(dom)).div_r() ;
580  break ;
581  case 1 :
582  // Dtheta S_theta
583  pos_so.set(ind_sum) = 0 ;
584  auxi_der.set(pos_auxi_der).set_domain(dom) +=
585  ((*so.der_t)(pos_so)(dom) -(*so.val_t)(pos_so)(dom)*(*rad_term_eq->der_t)()(dom)/(*rad_term_eq->val_t)()(dom)).div_r() ;
586  break ;
587  case 2 :
588  //Dtheta S_phi
589  break ;
590  default :
591  cerr << "Bad indice in Domain_shell_outer_adapted::connection_spher" << endl ;
592  abort() ;
593  }
594  break ;
595  case 2 :
596  // Dphi
597  // Different cases of the source index
598  switch (pos_auxi_der(ind_sum+1)) {
599  case 0 :
600  //Dphi S_r
601  pos_so.set(ind_sum) = 2 ;
602  auxi_der.set(pos_auxi_der).set_domain(dom) -=
603  ((*so.der_t)(pos_so)(dom) -(*so.val_t)(pos_so)(dom)*(*rad_term_eq->der_t)()(dom)/(*rad_term_eq->val_t)()(dom)).div_r() ;
604  break ;
605  case 1 :
606  // Dphi S_theta
607  pos_so.set(ind_sum) = 2 ;
608  auxi_der.set(pos_auxi_der).set_domain(dom) -=
609  ((*so.der_t)(pos_so)(dom) -(*so.val_t)(pos_so)(dom)*(*rad_term_eq->der_t)()(dom)/(*rad_term_eq->val_t)()(dom)).div_r().mult_cos_theta().div_sin_theta() ;
610  break ;
611  case 2 :
612  //Dphi S_phi
613  pos_so.set(ind_sum) = 0 ;
614  auxi_der.set(pos_auxi_der).set_domain(dom) +=
615  ((*so.der_t)(pos_so)(dom) -(*so.val_t)(pos_so)(dom)*(*rad_term_eq->der_t)()(dom)/(*rad_term_eq->val_t)()(dom)).div_r() ;
616  pos_so.set(ind_sum) = 1 ;
617  auxi_der.set(pos_auxi_der).set_domain(dom) +=
618  ((*so.der_t)(pos_so)(dom) -(*so.val_t)(pos_so)(dom)*(*rad_term_eq->der_t)()(dom)/(*rad_term_eq->val_t)()(dom)).div_r().mult_cos_theta().div_sin_theta() ;
619  break ;
620  default :
621  cerr << "Bad indice in Domain_shell_outer_adapted::connection_spher" << endl ;
622  abort() ;
623  }
624  break ;
625  default :
626  cerr << "Bad indice in Domain_shell_outer_adapted::connection_spher" << endl ;
627  abort() ;
628  }
629  }
630  while (pos_auxi_der.inc()) ;
631  }
632 
633  return Term_eq (dom, auxi_val, auxi_der) ;
634  }
635 }
636 
638  int dom = so.get_dom() ;
639  assert (dom == num_dom) ;
640 
641  int valence = so.val_t->get_valence() ;
642 
643  Array<int> type_ind (valence+1) ;
644  type_ind.set(0) = COV ;
645  for (int i=1 ; i<valence+1 ; i++)
646  type_ind.set(i) = so.val_t->get_index_type(i-1) ;
647 
648  Base_tensor basis (so.get_val_t().get_space()) ;
649  basis.set_basis(dom) = CARTESIAN_BASIS ;
650 
651  Term_eq comp_r (derive_r(so)) ;
652  Term_eq comp_t (derive_t(so) / (*rad_term_eq)) ;
654 
655  Tensor val_res (so.get_val_t().get_space(), valence+1, type_ind, basis) ;
656  {
657  Index pos_so (*so.val_t) ;
658  Index pos_res (val_res) ;
659  do {
660  for (int i=1 ; i<valence+1 ; i++)
661  pos_res.set(i) = pos_so(i-1) ;
662 
663  Val_domain auxi (mult_sin_theta((*comp_r.val_t)(pos_so)(num_dom)) + mult_cos_theta((*comp_t.val_t)(pos_so)(num_dom))) ;
664  // X part
665  pos_res.set(0) = 0 ;
666  val_res.set(pos_res).set_domain(num_dom) = mult_cos_phi(auxi) - mult_sin_phi((*comp_p.val_t)(pos_so)(num_dom)) ;
667  // Y part
668  pos_res.set(0) = 1 ;
669  val_res.set(pos_res).set_domain(num_dom) = mult_sin_phi(auxi) + mult_cos_phi((*comp_p.val_t)(pos_so)(num_dom)) ;
670  // Z part
671  pos_res.set(0) = 2 ;
672  val_res.set(pos_res).set_domain(num_dom) =
673  mult_cos_theta((*comp_r.val_t)(pos_so)(num_dom)) - mult_sin_theta((*comp_t.val_t)(pos_so)(num_dom)) ;
674  }
675  while (pos_so.inc()) ;
676 }
677 
678  bool doder = ((so.der_t==0x0) || (rad_term_eq->der_t==0x0)) ? false : true ;
679  if (doder) {
680  Tensor der_res (so.get_val_t().get_space(), valence+1, type_ind, basis) ;
681  Index pos_so (*so.der_t) ;
682  Index pos_res (val_res) ;
683  do {
684  for (int i=1 ; i<valence+1 ; i++)
685  pos_res.set(i) = pos_so(i-1) ;
686 
687  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))) ;
688  // X part
689  pos_res.set(0) = 0 ;
690  der_res.set(pos_res).set_domain(num_dom) = mult_cos_phi(auxi) - mult_sin_phi((*comp_p.der_t)(pos_so)(num_dom)) ;
691  // Y part
692  pos_res.set(0) = 1 ;
693  der_res.set(pos_res).set_domain(num_dom) = mult_sin_phi(auxi) + mult_cos_phi((*comp_p.der_t)(pos_so)(num_dom)) ;
694  // Z part
695  pos_res.set(0) = 2 ;
696  der_res.set(pos_res).set_domain(num_dom) =
697  mult_cos_theta((*comp_r.der_t)(pos_so)(num_dom)) - mult_sin_theta((*comp_t.der_t)(pos_so)(num_dom)) ;
698  }
699  while (pos_so.inc()) ;
700 
701  return Term_eq (num_dom, val_res, der_res) ;
702  }
703  else
704  return Term_eq (num_dom, val_res) ;
705 }
706 
707 const Term_eq* Domain_shell_outer_adapted::give_normal (int bound, int tipe) const {
708  assert (bound==OUTER_BC) ;
709  switch (tipe) {
710  case CARTESIAN_BASIS :
711  if (normal_cart==0x0)
712  do_normal_cart() ;
713  return normal_cart ;
714 
715  case SPHERICAL_BASIS :
716  if (normal_spher==0x0)
717  do_normal_spher() ;
718  return normal_spher ;
719 
720  default:
721  cerr << "Unknown type of tensorial basis in Domain_shell_inner_adapted::give_normal" << endl ;
722  abort() ;
723  }
724 }
725 
726 double integral_1d (int, const Array<double>&) ;
728 
729  int dom = target.get_dom() ;
730  // Check it is a tensor
731 #ifndef REMOVE_ALL_CHECKS
732  if (target.type_data != TERM_T) {
733  cerr << "Ope_int_volume only defined with respect for a tensor" << endl ;
734  abort() ;
735  }
736 
737  if (target.val_t->get_n_comp() != 1) {
738  cerr << "Ope_int_volume only defined with respect to a scalar" << endl ;
739  abort() ;
740  }
741 #endif
742 
744 
745  // The value
746  Array<int> ind (target.val_t->indices(0)) ;
747  Val_domain value ((*integrant.val_t)(ind)(dom)) ;
748  double resval = 0 ;
749  if (value.check_if_zero())
750  resval= 0. ;
751  else {
752  int baset = (*value.get_base().bases_1d[1])(0) ;
753  assert(baset==SIN_ODD) ;
754  Index pos (nbr_coefs) ;
755  for (int j=0 ; j<nbr_coefs(1) ; j++) {
756  pos.set(1) = j ;
757  int baser = (*value.get_base().bases_1d[0]) (j, 0) ;
758  assert (baser==CHEB) ;
759 
760  Array<double> cf (nbr_coefs(0)) ;
761  for (int i=0 ; i<nbr_coefs(0) ; i++) {
762  pos.set(0) = i ;
763  cf.set(i) = value.get_coef(pos) ;
764  }
765  resval += 2./(2.*j+1) * integral_1d(CHEB, cf) ;
766  }
767  resval *= 2*M_PI ;
768  }
769 
770  // The der
771  if (integrant.der_t!=0x0) {
772  Val_domain derval ((*integrant.der_t)(ind)(dom)) ;
773  double resder = 0 ;
774  if (derval.check_if_zero())
775  resder= 0. ;
776  else {
777  int baset = (*derval.get_base().bases_1d[1]) (0) ;
778  assert(baset==SIN_ODD) ;
779  Index pos (nbr_coefs) ;
780  for (int j=0 ; j<nbr_coefs(1) ; j++) {
781  pos.set(1) = j ;
782  int baser = (*derval.get_base().bases_1d[0]) (j, 0) ;
783  assert (baser==CHEB) ;
784 
785  Array<double> cf (nbr_coefs(0)) ;
786  for (int i=0 ; i<nbr_coefs(0) ; i++) {
787  pos.set(0) = i ;
788  cf.set(i) = derval.get_coef(pos) ;
789  }
790  resder += 2./(2.*j+1) * integral_1d(CHEB, cf) ;
791  }
792  resder *= 2*M_PI ;
793  }
794  return Term_eq (dom, resval, resder) ;
795  }
796  else {
797  return Term_eq (dom, resval) ;
798  }
799 }
800 
801 
803  cerr << "Domain_shell_outer_adapted::integ_term_eq not implemented yet" << endl ;
804  abort() ;
805 }
806 }
807 
reference set(const Index &pos)
Read/write of an element.
Definition: array.hpp:186
Bases_container bases_1d
Arrays containing the various basis of decomposition.
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 lap_term_eq(const Term_eq &, int) const
Returns the flat 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.
virtual Val_domain mult_sin_phi(const Val_domain &) const
Multiplication by .
virtual Term_eq integ_volume_term_eq(const Term_eq &) const
Volume integral of a Term_eq.
virtual Term_eq dr_term_eq(const Term_eq &) const
Radial derivative of a Term_eq.
Term_eq * outer_radius_term_eq
Pointer on the outer boundary , as a Term_eq.
Definition: adapted.hpp:375
Term_eq * dp_rad_term_eq
Pointer on the Term_eq containing the .
Definition: adapted.hpp:400
virtual Val_domain mult_cos_phi(const Val_domain &) const
Multiplication by .
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...
Term_eq flat_grad_spher(const Term_eq &) const
Computes the flat gradient of a field, in orthonormal spherical coordinates.
void do_normal_cart() const
Computes the normal wrt the inner boundary, in Cartesian coordinates.
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...
Term_eq derive_t(const Term_eq &so) const
Computes .
Term_eq derive_r(const Term_eq &so) const
Computes .
Term_eq * normal_spher
Pointer on the Term_eq containing the normal vector to the outer boundary, in orthonormal spherical c...
Definition: adapted.hpp:380
Term_eq * normal_cart
Pointer on the Term_eq containing the normal vector to the outer boundary, in Cartesian coordinates.
Definition: adapted.hpp:384
virtual Term_eq connection_spher(const Term_eq &) const
Computes the part of the gradient involving the connections, in spherical orthonormal coordinates.
Term_eq * der_rad_term_eq
Pointer on the Term_eq containing the .
Definition: adapted.hpp:392
virtual Val_domain mult_sin_theta(const Val_domain &) const
Multiplication by .
Term_eq * rad_term_eq
Pointer on the Term_eq containing the radius.
Definition: adapted.hpp:388
virtual Term_eq integ_term_eq(const Term_eq &, int) const
Surface integral of a Term_eq.
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.
void do_normal_spher() const
Computes the normal wrt the inner boundary, in orthonormal spherical coordinates.
Term_eq * dt_rad_term_eq
Pointer on the Term_eq containing the .
Definition: adapted.hpp:396
virtual const Term_eq * give_normal(int, int) const
Returns the vector normal to a surface.
virtual Term_eq mult_r_term_eq(const Term_eq &) const
Multiplication by of a Term_eq.
virtual Val_domain div_r(const Val_domain &) const
Division by .
virtual Val_domain mult_cos_theta(const Val_domain &) const
Multiplication by .
const Space & sp
The corresponding Space ; required for updating fields whene the mapping changes.
Definition: adapted.hpp:373
Term_eq derive_p(const Term_eq &so) const
Computes .
virtual Val_domain mult_sin_theta(const Val_domain &) const
Multiplication by .
Definition: domain.cpp:1236
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
Dim_array nbr_coefs
Number of coefficients.
Definition: space.hpp:66
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
const int type_data
Flag describing the type of data :
Definition: term_eq.hpp:75
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
Val_domain div_sin_theta() const
Division by .
Val_domain mult_cos_theta() const
Multiplication by .
Base_spectral & set_base()
Sets the basis of decomposition.
Definition: val_domain.hpp:126
Array< double > get_coef() const
Definition: val_domain.hpp:136
void allocate_conf()
Allocates the values in the configuration space and destroys the values in the coefficients space.
Definition: val_domain.cpp:209
const Base_spectral & get_base() const
Returns the basis of decomposition.
Definition: val_domain.hpp:122
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.