KADATH
domain_polar_inner_adapted.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 
29 namespace Kadath {
30 void coef_1d (int, Array<double>&) ;
31 void coef_i_1d (int, Array<double>&) ;
32 int der_1d (int, Array<double>&) ;
33 
34 // Standard constructor
35 Domain_polar_shell_inner_adapted::Domain_polar_shell_inner_adapted (const Space& sss, int num, int ttype, double rin, double rout, const Point& cr, const Dim_array& nbr) :
36  Domain(num, ttype, nbr), sp(sss), outer_radius(rout), center(cr) {
37  assert (nbr.get_ndim()==2) ;
38  assert (cr.get_ndim()==2) ;
39 
40 
41  inner_radius_term_eq = nullptr ;
42  rad_term_eq = nullptr ;
43  der_rad_term_eq = nullptr ;
44  dt_rad_term_eq = nullptr ;
45  normal_spher = nullptr ;
46  normal_cart = nullptr ;
47 
48  do_coloc() ;
49  inner_radius = new Val_domain (this) ;
50  *inner_radius = rin ;
52 }
53 
54 Domain_polar_shell_inner_adapted::Domain_polar_shell_inner_adapted (const Space& sss, int num, int ttype, const Val_domain& rin, double rout, const Point& cr, const Dim_array& nbr) :
55  Domain(num, ttype, nbr), sp(sss), outer_radius(rout), center(cr) {
56  assert (nbr.get_ndim()==2) ;
57  assert (cr.get_ndim()==2) ;
58 
59 
60  inner_radius_term_eq = nullptr ;
61  rad_term_eq = nullptr ;
62  der_rad_term_eq = nullptr ;
63  dt_rad_term_eq = nullptr ;
64  normal_spher = nullptr ;
65  normal_cart = nullptr ;
66 
67  do_coloc() ;
68  inner_radius = new Val_domain(rin) ;
69 }
70 
71 // Constructor by copy
73  outer_radius (so.outer_radius), center(so.center) {
74 
76  if (so.inner_radius_term_eq != nullptr)
78  if (so.rad_term_eq !=nullptr)
79  rad_term_eq = new Term_eq (*so.rad_term_eq) ;
80  if (so.der_rad_term_eq !=nullptr)
82  if (so.dt_rad_term_eq !=nullptr)
84  if (so.normal_spher !=nullptr)
85  normal_spher = new Term_eq (*so.normal_spher) ;
86  if (so.normal_cart !=nullptr)
87  normal_cart = new Term_eq (*so.normal_cart) ;
88 }
89 
90 Domain_polar_shell_inner_adapted::Domain_polar_shell_inner_adapted (const Space& sss, int num, FILE* fd) : Domain(num, fd), sp(sss), center(fd) {
91  fread_be (&outer_radius, sizeof(double), 1, fd) ;
92  inner_radius = new Val_domain(this, fd) ;
93 
94  inner_radius_term_eq = nullptr ;
95  rad_term_eq = nullptr ;
96  der_rad_term_eq = nullptr ;
97  dt_rad_term_eq = nullptr ;
98  normal_spher = nullptr ;
99  normal_cart = nullptr ;
100 
101  do_coloc() ;
102 }
103 
105  for (int i=0 ; i<2 ; i++)
106  assert (coloc[i] != nullptr) ;
107  assert (radius == nullptr) ;
108  radius = new Val_domain(this) ;
109  radius->allocate_conf() ;
110  Index index (nbr_points) ;
111  do
112  radius->set(index) = (outer_radius - (*inner_radius)(index))/2.* ((*coloc[0])(index(0))) + (outer_radius + (*inner_radius)(index))/2. ;
113  while (index.inc()) ;
114  radius->std_r_base() ;
115 }
116 
117 // Destructor
118 Domain_polar_shell_inner_adapted::~Domain_polar_shell_inner_adapted() {
119 
120  delete inner_radius ;
121  if (inner_radius_term_eq != nullptr)
122  delete inner_radius_term_eq ;
123  if (rad_term_eq != nullptr)
124  delete rad_term_eq ;
125  if (der_rad_term_eq != nullptr)
126  delete der_rad_term_eq ;
127  if (normal_spher != nullptr)
128  delete normal_spher ;
129  if (normal_cart != nullptr)
130  delete normal_cart ;
131  if (dt_rad_term_eq != nullptr)
132  delete dt_rad_term_eq ;
133 }
134 
136  safe_delete(inner_radius_term_eq);
137  safe_delete(rad_term_eq);
138  safe_delete(der_rad_term_eq);
139  safe_delete(normal_spher);
140  safe_delete(normal_cart);
141  safe_delete(dt_rad_term_eq);
142 }
143 
145  return nbr_coefs(1) ;
146 }
147 
149 
150  if (inner_radius_term_eq != nullptr)
151  delete inner_radius_term_eq ;
152  Scalar val (sp) ;
154 
155  inner_radius_term_eq = new Term_eq (num_dom, val) ;
156  update() ;
157 }
158 
159 void Domain_polar_shell_inner_adapted::affecte_coef(int& conte, int cc, bool& found) const {
160 
161  Val_domain auxi (this) ;
162  auxi.std_base() ;
163  auxi.set_in_coef() ;
164  auxi.allocate_coef() ;
165  *auxi.cf = 0 ;
166 
167  found = false ;
168 
169 
170  for (int j=0 ; j<nbr_coefs(1) ; j++) {
171  if (conte==cc) {
172  Index pos_cf (nbr_coefs) ;
173  pos_cf.set(0) = 0 ;
174  pos_cf.set(1) = j ;
175  auxi.cf->set(pos_cf) = 1 ;
176  found = true ;
177  }
178  conte ++ ;
179  }
180 
181  if (found) {
182  Scalar auxi_scal (sp) ;
183  auxi_scal.set_domain(num_dom) = auxi ;
184  inner_radius_term_eq->set_der_t(auxi_scal) ;
185  }
186  else {
188  }
189  update() ;
190 }
191 
192 void Domain_polar_shell_inner_adapted::xx_to_vars_from_adapted(Val_domain& new_inner_radius, const Array<double>& xx, int& pos) const {
193 
194  new_inner_radius.allocate_coef() ;
195  *new_inner_radius.cf = 0 ;
196 
197  Index pos_cf (nbr_coefs) ;
198  pos_cf.set(0) = 0 ;
199 
200  for (int j=0 ; j<nbr_coefs(1) ; j++) {
201  pos_cf.set(1)= j ;
202  new_inner_radius.cf->set(pos_cf) -= xx(pos) ;
203  pos ++ ;
204  }
205 
206  new_inner_radius.set_base() = inner_radius->get_base() ;
207 }
208 
210 
211  *inner_radius += cor ;
212  for (int l=0 ; l<ndim ; l++) {
213  if (absol[l] !=nullptr) delete absol[l] ;
214  if (cart[l] !=nullptr) delete cart[l] ;
215  if (cart_surr[l] !=nullptr) delete cart_surr[l] ;
216  absol[l] = nullptr ;
217  cart_surr[l] = nullptr ;
218  cart[l]= nullptr ;
219  }
220  if (radius !=nullptr)
221  delete radius ;
222  radius = nullptr ;
223  update() ;
224 }
225 
227 
228  *inner_radius = cor ;
229  for (int l=0 ; l<ndim ; l++) {
230  if (absol[l] !=nullptr) delete absol[l] ;
231  if (cart[l] !=nullptr) delete cart[l] ;
232  if (cart_surr[l] !=nullptr) delete cart_surr[l] ;
233  absol[l] = nullptr ;
234  cart_surr[l] = nullptr ;
235  cart[l]= nullptr ;
236  }
237  if (radius !=nullptr)
238  delete radius ;
239  radius = nullptr ;
240  vars_to_terms() ;
241 }
242 
243 void Domain_polar_shell_inner_adapted::update_variable (const Val_domain& cor_inner_radius, const Scalar& old, Scalar& res) const {
244 
245  Val_domain dr (old(num_dom).der_r()) ;
246  if (dr.check_if_zero())
247  res.set_domain(num_dom) = 0 ;
248  else {
249 
250  Index pos(nbr_points) ;
252  do {
253  res.set_domain(num_dom).set(pos) = dr(pos) * (1-(*coloc[0])(pos(0))) / 2. ;
254  }
255  while (pos.inc()) ;
256 
257  res.set_domain(num_dom) = cor_inner_radius * res(num_dom) + old(num_dom) ;
258  res.set_domain(num_dom).set_base() = old(num_dom).get_base() ;
259  }
260 }
261 
262 
264 
265  Val_domain auxi (this) ;
266  auxi.std_base() ;
267  auxi.set_in_coef() ;
268  auxi.allocate_coef() ;
269  *auxi.cf = 0 ;
270 
271  Index pos_cf (nbr_coefs) ;
272  pos_cf.set(0) = 0 ;
273 
274  for (int j=0 ; j<nbr_coefs(1) ; j++) {
275  pos_cf.set(1)= j ;
276  auxi.cf->set(pos_cf) = xx(pos) ;
277  pos ++ ;
278  }
279 
280  Scalar auxi_scal (sp) ;
281  auxi_scal.set_domain(num_dom) = auxi ;
282  inner_radius_term_eq->set_der_t(auxi_scal) ;
283  update() ;
284 }
285 
286 
288 
289  if (rad_term_eq != nullptr)
290  delete rad_term_eq ;
291  if (der_rad_term_eq != nullptr)
292  delete der_rad_term_eq ;
293  if (dt_rad_term_eq != nullptr)
294  delete dt_rad_term_eq ;
295 
296  // Computation of rad_term_eq
297  Scalar val_res (sp) ;
298  val_res.set_domain(num_dom).allocate_conf() ;
299  Index index (nbr_points) ;
300  do {
301  val_res.set_domain(num_dom).set(index) = (outer_radius - ((*inner_radius_term_eq->val_t)()(num_dom))(index))/2. * ((*coloc[0])(index(0))) +
302  (outer_radius + ((*inner_radius_term_eq->val_t)()(num_dom))(index))/2. ;
303  }
304  while (index.inc()) ;
305  val_res.set_domain(num_dom).std_r_base() ;
306 
307  bool doder = (inner_radius_term_eq->der_t ==nullptr) ? false : true ;
308  if (doder) {
309  Scalar der_res (sp) ;
310  if ((*inner_radius_term_eq->der_t)()(num_dom).check_if_zero()) {
311  der_res.set_domain(num_dom).set_zero() ;
312  }
313  else {
314  der_res.set_domain(num_dom).allocate_conf() ;
315  index.set_start() ;
316  do {
317  der_res.set_domain(num_dom).set(index) = -((*inner_radius_term_eq->der_t)()(num_dom))(index)/2. * ((*coloc[0])(index(0))) +
318  ((*inner_radius_term_eq->der_t)()(num_dom))(index)/2. ;
319  }
320  while (index.inc()) ;
321  der_res.set_domain(num_dom).std_r_base() ;
322  }
323  rad_term_eq = new Term_eq (num_dom, val_res, der_res) ;
324  }
325  else
326  rad_term_eq = new Term_eq (num_dom, val_res) ;
327 
328  // Computation of der_rad_term_eq which is dr / dxi
329  val_res.set_domain(num_dom) = (outer_radius - (*inner_radius_term_eq->val_t)()(num_dom)) / 2. ;
330  if (doder) {
331  Scalar der_res (sp) ;
332  der_res.set_domain(num_dom) = - (*inner_radius_term_eq->der_t)()(num_dom) / 2. ;
333  der_rad_term_eq = new Term_eq (num_dom, val_res, der_res) ;
334  }
335  else
336  der_rad_term_eq = new Term_eq (num_dom, val_res) ;
337 
338  // Computation of dt_rad_term_eq which is dr / d theta
339  val_res.set_domain(num_dom) = (*rad_term_eq->val_t)()(num_dom).der_var(2) ;
340  if (doder) {
341  Scalar der_res (sp) ;
342  der_res.set_domain(num_dom) = (*rad_term_eq->der_t)()(num_dom).der_var(2) ;
343  dt_rad_term_eq = new Term_eq (num_dom, val_res, der_res) ;
344  }
345  else
346  dt_rad_term_eq = new Term_eq (num_dom, val_res) ;
347 
348  do_normal_cart() ;
349 }
350 
352  nbr_points.save(fd) ;
353  nbr_coefs.save(fd) ;
354  fwrite_be (&ndim, sizeof(int), 1, fd) ;
355  fwrite_be (&type_base, sizeof(int), 1, fd) ;
356  center.save(fd) ;
357  fwrite_be (&outer_radius, sizeof(double), 1, fd) ;
358  inner_radius->save(fd) ;
359 }
360 
361 ostream& Domain_polar_shell_inner_adapted::print (ostream& o) const {
362  o << "Adapted polar shell on the inside boundary" << endl ;
363  o << "Center = " << center << endl ;
364  o << "Nbr pts = " << nbr_points << endl ;
365  o << "Outer radius " << outer_radius << endl ;
366  o << "Inner radius " << endl ;
367  o << *inner_radius << endl ;
368  o << endl ;
369  return o ;
370 }
371 
372 
373 
375  cerr << "Domain_polar_shell_inner_adapted::der_normal not implemeted" << endl ;
376  abort() ;
377 }
378 
379 // Computes the cartesian coordinates
381  for (int i=0 ; i<2 ; i++)
382  assert (coloc[i] != nullptr) ;
383  for (int i=0 ; i<2 ; i++)
384  assert (absol[i] == nullptr) ;
385  for (int i=0 ; i<2 ; i++) {
386  absol[i] = new Val_domain(this) ;
387  absol[i]->allocate_conf() ;
388  }
389  Index index (nbr_points) ;
390  do {
391 
392  double rr = (outer_radius - (*inner_radius)(index))/2. * ((*coloc[0])(index(0))) +
393  (outer_radius + (*inner_radius)(index))/2. ;
394  absol[0]->set(index) = rr *
395  sin((*coloc[1])(index(1))) + center(1);
396  absol[1]->set(index) = rr * cos((*coloc[1])(index(1))) + center(2) ;
397  }
398  while (index.inc()) ;
399 
400 }
401 
402 // Computes the cartesian coordinates
404  for (int i=0 ; i<2 ; i++)
405  assert (coloc[i] != nullptr) ;
406  for (int i=0 ; i<2 ; i++)
407  assert (cart[i] == nullptr) ;
408  for (int i=0 ; i<2 ; i++) {
409  cart[i] = new Val_domain(this) ;
410  cart[i]->allocate_conf() ;
411  }
412  Index index (nbr_points) ;
413  do {
414  double rr = (outer_radius - (*inner_radius)(index))/2. * ((*coloc[0])(index(0))) +
415  (outer_radius + (*inner_radius)(index))/2. ;
416  cart[0]->set(index) = rr *
417  sin((*coloc[1])(index(1))) + center(1);
418  cart[1]->set(index) = rr * cos((*coloc[1])(index(1))) + center(2) ;
419  }
420  while (index.inc()) ;
421 
422 }
423 
424 
425 // Computes the cartesian coordinates over the radius
427  for (int i=0 ; i<2 ; i++)
428  assert (coloc[i] != nullptr) ;
429  for (int i=0 ; i<2 ; i++)
430  assert (cart_surr[i] == nullptr) ;
431  for (int i=0 ; i<2 ; i++) {
432  cart_surr[i] = new Val_domain(this) ;
433  cart_surr[i]->allocate_conf() ;
434  }
435  Index index (nbr_points) ;
436  do {
437  cart_surr[0]->set(index) = sin((*coloc[1])(index(1))) ;
438  cart_surr[1]->set(index) = cos((*coloc[1])(index(1))) ;
439  }
440  while (index.inc()) ;
441 
442 }
443 
444 
445 // Is a point inside this domain ?
446 bool Domain_polar_shell_inner_adapted::is_in (const Point& xx, double prec) const {
447 
448  assert (xx.get_ndim()==2) ;
449  Point num (absol_to_num(xx)) ;
450  bool res = ((num(1)>-1-prec) && (num(1)<1+prec)) ? true : false ;
451  return res ;
452 }
453 
454 
455 // Convert absolute coordinates to numerical ones
457  Point num(2) ;
458 
459  double rho_loc = fabs(abs(1) - center(1)) ;
460  double z_loc = abs(2) - center(2) ;
461  double air = sqrt(rho_loc*rho_loc+z_loc*z_loc) ;
462 
463  if (rho_loc==0) {
464  // On the axis ?
465  num.set(2) = (z_loc>=0) ? 0 : M_PI ;
466  }
467  else {
468  num.set(2) = atan(rho_loc/z_loc) ;
469  }
470 
471  if (num(2) <0)
472  num.set(2) = M_PI + num(2) ;
473 
474  // Get the boundary for those angles
475  num.set(1) = 1 ;
476  inner_radius->coef() ;
477  double inner = inner_radius->get_base().summation(num, inner_radius->get_coef()) ;
478  num.set(1) = (2./(outer_radius-inner)) * (air - (outer_radius + inner)/2.) ;
479 
480  return num ;
481 }
482 
483 
484 double coloc_leg(int, int) ;
486 
487  switch (type_base) {
488  case CHEB_TYPE:
490  del_deriv() ;
491  for (int i=0 ; i<ndim ; i++)
492  coloc[i] = new Array<double> (nbr_points(i)) ;
493  for (int i=0 ; i<nbr_points(0) ; i++)
494  coloc[0]->set(i) = -cos(M_PI*i/(nbr_points(0)-1)) ;
495  for (int j=0 ; j<nbr_points(1) ; j++)
496  coloc[1]->set(j) = M_PI/2.*j/(nbr_points(1)-1) ;
497  break ;
498  case LEG_TYPE:
500  del_deriv() ;
501  for (int i=0 ; i<ndim ; i++)
502  coloc[i] = new Array<double> (nbr_points(i)) ;
503  for (int i=0 ; i<nbr_points(0) ; i++)
504  coloc[0]->set(i) = coloc_leg(i, nbr_points(0)) ;
505  for (int j=0 ; j<nbr_points(1) ; j++)
506  coloc[1]->set(j) = M_PI/2.*j/(nbr_points(1)-1) ;
507  break ;
508  default :
509  cerr << "Unknown type of basis in Domain_polar_shell_inner_adapted::do_coloc" << endl ;
510  abort() ;
511  }
512 }
513 
515  set_cheb_base_with_m(base, 0) ;
516 }
517 
518 // standard base for a anti-symetric function in z, using Chebyshev
520  set_anti_cheb_base_with_m(base, 0) ;
521 }
522 
523 // standard base for a symetric function in z, using Legendre
525  set_legendre_base_with_m(base, 0) ;
526 }
527 
528 // standard base for a anti-symetric function in z, using Legendre
531 }
532 
533 
534 // standard base for a symetric function in z, using Chebyshev
536 
537  assert (type_base == CHEB_TYPE) ;
538 
539  base.allocate(nbr_coefs) ;
540 
541  base.def =true ;
542 
543  base.bases_1d[1]->set(0) = (m%2==0) ? COS_EVEN : SIN_ODD ;
544  for (int j=0 ; j<nbr_coefs(1) ; j++)
545  base.bases_1d[0]->set(j) = CHEB ;
546 }
547 
548 // standard base for a anti-symetric function in z, using Chebyshev
550 
551  assert (type_base == CHEB_TYPE) ;
552 
553  base.allocate(nbr_coefs) ;
554 
555  base.def =true ;
556 
557  base.bases_1d[1]->set(0) = (m%2==0) ? COS_ODD : SIN_EVEN ;
558  for (int j=0 ; j<nbr_coefs(1) ; j++)
559  base.bases_1d[0]->set(j) = CHEB ;
560 }
561 
562 // standard base for a symetric function in z, using Legendre
564 
565  assert (type_base == LEG_TYPE) ;
566 
567  base.allocate(nbr_coefs) ;
568 
569  base.def =true ;
570 
571  base.bases_1d[1]->set(0) = (m%2==0) ? COS_EVEN : SIN_ODD ;
572  for (int j=0 ; j<nbr_coefs(1) ; j++)
573  base.bases_1d[0]->set(j) = LEG ;
574 }
575 
576 // standard base for a anti-symetric function in z, using Legendre
578 
579  assert (type_base == LEG_TYPE) ;
580 
581  base.allocate(nbr_coefs) ;
582 
583  base.def =true ;
584 
585  base.bases_1d[1]->set(0) = (m%2==0) ? COS_ODD : SIN_EVEN ;
586  for (int j=0 ; j<nbr_coefs(1) ; j++)
587  base.bases_1d[0]->set(j) = LEG ;
588 }
589 
591  set_cheb_base_with_m(base, 0) ;
592 }
593 
595  set_legendre_base_with_m(base, 0) ;
596 }
597 
598 // Computes the derivatives with respect to XYZ as a function of the numerical ones.
599 void Domain_polar_shell_inner_adapted::do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const {
600  Val_domain rr (get_radius()) ;
601  Val_domain dr (*der_var[0] / (*der_rad_term_eq->val_t)()(num_dom)) ;
602  Val_domain dtsr ((*der_var[1] - rr.der_var(2) * dr) / rr);
603  dtsr.set_base() = der_var[1]->get_base() ;
604 
605  // d/dx :
606  Val_domain sintdr (dr.mult_sin_theta()) ;
607  Val_domain costdtsr (dtsr.mult_cos_theta()) ;
608 
609  der_abs[0] = new Val_domain (sintdr+costdtsr) ;
610 
611  // d/dz :
612  der_abs[1] = new Val_domain (dr.mult_cos_theta() - dtsr.mult_sin_theta()) ;
613 }
614 
615 // Rules for the multiplication of two basis.
617 
618  assert (a.ndim==2) ;
619  assert (b.ndim==2) ;
620 
621  Base_spectral res(2) ;
622  bool res_def = true ;
623 
624  if (!a.def)
625  res_def=false ;
626  if (!b.def)
627  res_def=false ;
628 
629  if (res_def) {
630 
631  // Bases in theta :
632  res.bases_1d[1] = new Array<int> (a.bases_1d[1]->get_dimensions()) ;
633  switch ((*a.bases_1d[1])(0)) {
634  case COS_EVEN:
635  switch ((*b.bases_1d[1])(0)) {
636  case COS_EVEN:
637  res.bases_1d[1]->set(0) = COS_EVEN ;
638  break ;
639  case COS_ODD:
640  res.bases_1d[1]->set(0) = COS_ODD ;
641  break ;
642  case SIN_EVEN:
643  res.bases_1d[1]->set(0) = SIN_EVEN ;
644  break ;
645  case SIN_ODD:
646  res.bases_1d[1]->set(0) = SIN_ODD ;
647  break ;
648  default:
649  res_def = false ;
650  break ;
651  }
652  break ;
653  case COS_ODD:
654  switch ((*b.bases_1d[1])(0)) {
655  case COS_EVEN:
656  res.bases_1d[1]->set(0) = COS_ODD ;
657  break ;
658  case COS_ODD:
659  res.bases_1d[1]->set(0) = COS_EVEN ;
660  break ;
661  case SIN_EVEN:
662  res.bases_1d[1]->set(0) = SIN_ODD ;
663  break ;
664  case SIN_ODD:
665  res.bases_1d[1]->set(0) = SIN_EVEN ;
666  break ;
667  default:
668  res_def = false ;
669  break ;
670  }
671  break ;
672  case SIN_EVEN:
673  switch ((*b.bases_1d[1])(0)) {
674  case COS_EVEN:
675  res.bases_1d[1]->set(0) = SIN_EVEN ;
676  break ;
677  case COS_ODD:
678  res.bases_1d[1]->set(0) = SIN_ODD ;
679  break ;
680  case SIN_EVEN:
681  res.bases_1d[1]->set(0) = COS_EVEN ;
682  break ;
683  case SIN_ODD:
684  res.bases_1d[1]->set(0) = COS_ODD ;
685  break ;
686  default:
687  res_def = false ;
688  break ;
689  }
690  break ;
691  case SIN_ODD:
692  switch ((*b.bases_1d[1])(0)) {
693  case COS_EVEN:
694  res.bases_1d[1]->set(0) = SIN_ODD ;
695  break ;
696  case COS_ODD:
697  res.bases_1d[1]->set(0) = SIN_EVEN ;
698  break ;
699  case SIN_EVEN:
700  res.bases_1d[1]->set(0) = COS_ODD ;
701  break ;
702  case SIN_ODD:
703  res.bases_1d[1]->set(0) = COS_EVEN ;
704  break ;
705  default:
706  res_def = false ;
707  break ;
708  }
709  break ;
710  default:
711  res_def = false ;
712  break ;
713  }
714 
715 
716 
717  // Base in r :
718  Index index_0 (a.bases_1d[0]->get_dimensions()) ;
719  res.bases_1d[0] = new Array<int> (a.bases_1d[0]->get_dimensions()) ;
720  do {
721  switch ((*a.bases_1d[0])(index_0)) {
722  case CHEB:
723  switch ((*b.bases_1d[0])(index_0)) {
724  case CHEB:
725  res.bases_1d[0]->set(index_0) = CHEB ;
726  break ;
727  default:
728  res_def = false ;
729  break ;
730  }
731  break ;
732  case LEG:
733  switch ((*b.bases_1d[0])(index_0)) {
734  case LEG:
735  res.bases_1d[0]->set(index_0) = LEG ;
736  break ;
737  default:
738  res_def = false ;
739  break ;
740  }
741  break ;
742  default:
743  res_def = false ;
744  break ;
745  }
746  }
747  while (index_0.inc()) ;
748  }
749 
750  if (!res_def)
751  for (int dim=0 ; dim<a.ndim ; dim++)
752  if (res.bases_1d[dim]!= nullptr) {
753  delete res.bases_1d[dim] ;
754  res.bases_1d[dim] = nullptr ;
755  }
756  res.def = res_def ;
757  return res ;
758 }
759 
760 
761 
762 void Domain_polar_shell_inner_adapted::update_constante (const Val_domain& cor_inner_radius, const Scalar& old, Scalar& res) const {
763 
764  Point MM(2) ;
765  Index pos(nbr_points) ;
767  do {
768 
769  double rr = (outer_radius - cor_inner_radius(pos) - (*inner_radius)(pos)) * (*coloc[0])(pos(0)) / 2.
770  + (outer_radius + cor_inner_radius(pos) + (*inner_radius)(pos)) / 2.;
771  double theta = (*coloc[1])(pos(1)) ;
772 
773  MM.set(1) = rr*sin(theta) + center(1) ;
774  MM.set(2) = rr*cos(theta) + center(2) ;
775 
776  res.set_domain(num_dom).set(pos) = old.val_point(MM, +1) ;
777 
778  }
779  while (pos.inc()) ;
780 
781  res.set_domain(num_dom).set_base() = old(num_dom).get_base() ;
782 }
783 
784 
786  int res = -1 ;
787  if (strcmp(p,"R ")==0)
788  res = 0 ;
789  if (strcmp(p,"T ")==0)
790  res = 1 ;
791  return res ;
792 }
793 }
794 
reference set(const Index &pos)
Read/write of an element.
Definition: array.hpp:186
Class for storing the basis of decompositions of a field.
Bases_container bases_1d
Arrays containing the various basis of decomposition.
double summation(const Point &num, const Array< double > &tab) const
Computes the spectral summation.
void allocate(const Dim_array &nbr_coefs)
Allocates the various arrays, for a given number of coefficients.
bool def
true if the Base_spectral is defined and false otherwise.
int ndim
Number of dimensions.
Class for storing the dimensions of an array.
Definition: dim_array.hpp:34
int get_ndim() const
Returns the number of dimensions.
Definition: dim_array.hpp:63
void save(FILE *) const
Save function.
Definition: dim_array.cpp:32
Class for a spherical-like domain, having a symmetry with respect to the plane .
virtual void set_anti_cheb_base(Base_spectral &) const
Gives the base using Chebyshev polynomials, for functions antisymetric with respect to .
virtual Val_domain der_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
void set_mapping(const Val_domain &so) const
Affects the inner radius.
virtual void do_cart_surr() const
Computes the Cartesian coordinates over the radius.
virtual const Point absol_to_num(const Point &) const
Computes the numerical coordinates from the physical ones.
Term_eq * der_rad_term_eq
Pointer on the Term_eq containing the .
virtual void xx_to_ders_from_adapted(const Array< double > &, int &) const
Affects the derivative part of variable a Domain from a set of values.
virtual void update_constante(const Val_domain &, const Scalar &, Scalar &) const
Update the value of a scalar, after the shape of the Domain has been changed by the system.
virtual void set_cheb_base_with_m(Base_spectral &, int m) const
Gives the standard base using Chebyshev polynomials.
virtual void update_variable(const Val_domain &, const Scalar &, Scalar &) const
Update the value of a scalar, after the shape of the Domain has been changed by the system.
Point center
Absolute coordinates of the center.
virtual void update_mapping(const Val_domain &)
Updates the variables parts of the Domain.
void del_deriv() override
Destroys the derivated members (like coloc, cart and radius), when changing the type of colocation po...
void do_normal_cart() const
Computes the normal wrt the inner boundary, in Cartesian coordinates.
Val_domain * inner_radius
Pointer on the inner boundary , as a Val_domain.
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.
virtual void set_legendre_base(Base_spectral &) const
Gives the standard base for Legendre polynomials.
Term_eq * inner_radius_term_eq
Pointer on the inner boundary , as a Term_eq.
virtual void vars_to_terms() const
The Term_eq describing the variable shape of the Domain are updated.
virtual void set_legendre_r_base(Base_spectral &) const
Gives the base using odd Legendre polynomials$ for the radius.
virtual void set_legendre_base_with_m(Base_spectral &, int m) const
Gives the stnadard base using Legendre polynomials.
virtual Val_domain der_r(const Val_domain &) const
Compute the radial derivative of a scalar field.
void update() const
Updates all the quantities that depend on the inner radius (like the normal vectors).
const Space & sp
The corresponding Space ; required for updating fields whene the mapping changes.
Term_eq * normal_cart
Pointer on the Term_eq containing the normal vector to the inner boundary, in Cartesian coordinates.
virtual void save(FILE *) const
Saving function.
virtual void xx_to_vars_from_adapted(Val_domain &, const Array< double > &, int &) const
Computes the new boundary of a Domain from a set of values.
virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const
Computes the derivative with respect to the absolute Cartesian coordinates from the derivative with r...
virtual void set_anti_cheb_base_with_m(Base_spectral &, int m) const
Gives the base using Chebyshev polynomials, for functions antisymetric with respect to .
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
virtual void do_radius() const
Computes the generalized radius.
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
virtual void set_cheb_r_base(Base_spectral &) const
Gives the base using odd Chebyshev polynomials$ for the radius.
virtual void do_absol() const
Computes the absolute coordinates.
virtual void set_anti_legendre_base_with_m(Base_spectral &, int m) const
Gives the base using Legendre polynomials, for functions antisymetric with respect to .
virtual int give_place_var(char *) const
Translates a name of a coordinate into its corresponding numerical name.
virtual void affecte_coef(int &, int, bool &) const
The variation of the functions describing the shape of the Domain are affected from the unknowns of t...
virtual void set_anti_legendre_base(Base_spectral &) const
Gives the base using Legendre polynomials, for functions antisymetric with respect to .
virtual int nbr_unknowns_from_adapted() const
Gives the number of unknowns coming from the variable shape of the domain.
Domain_polar_shell_inner_adapted(const Space &sp, int num, int ttype, double rin, double rout, const Point &cr, const Dim_array &nbr)
Constructor :
virtual void do_coloc()
Computes the colocation points.
Term_eq * normal_spher
Pointer on the Term_eq containing the normal vector to the inner boundary, in orthonormal spherical c...
virtual Base_spectral mult(const Base_spectral &, const Base_spectral &) const
Method for the multiplication of two Base_spectral.
virtual void set_cheb_base(Base_spectral &) const
Gives the standard base for Chebyshev polynomials.
virtual void do_cart() const
Computes the Cartesian coordinates.
Abstract class that implements the fonctionnalities common to all the type of domains.
Definition: space.hpp:60
Val_domain * radius
The generalized radius.
Definition: space.hpp:78
Memory_mapped_array< Val_domain * > cart
Cartesian coordinates.
Definition: space.hpp:77
Memory_mapped_array< Val_domain * > absol
Asbolute coordinates (if defined ; usually Cartesian-like)
Definition: space.hpp:76
int ndim
Number of dimensions.
Definition: space.hpp:64
Memory_mapped_array< Val_domain * > cart_surr
Cartesian coordinates divided by the radius.
Definition: space.hpp:79
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
Val_domain const & get_radius() const
Returns the generalized radius.
Definition: space.hpp:1465
Dim_array nbr_points
Number of colocation points.
Definition: space.hpp:65
int type_base
Type of colocation point :
Definition: space.hpp:73
Memory_mapped_array< Array< double > * > coloc
Colocation points in each dimension (stored in ndim 1d- arrays)
Definition: space.hpp:75
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
void set_start()
Sets the position to zero in all dimensions.
Definition: index.hpp:88
bool inc(int increm, int var=0)
Increments the position of the Index.
Definition: index.hpp:99
The class Point is used to store the coordinates of a point.
Definition: point.hpp:30
void save(FILE *) const
Saving function.
Definition: point.cpp:33
const int & get_ndim() const
Returns the number of dimensions.
Definition: point.hpp:51
double & set(int i)
Read/write of a coordinate.
Definition: point.hpp:47
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
double val_point(const Point &xxx, int sens=-1) const
Computes the value of the field at a given point, by doing the spectral summation.
Definition: scalar.cpp:62
The Space class is an ensemble of domains describing the whole space of the computation.
Definition: space.hpp:1362
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
void set_der_zero()
Sets the variation of the approriate type to zero.
Definition: term_eq.cpp:187
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
void save(FILE *) const
Saving on a file.
Definition: val_domain.cpp:110
void set_in_coef()
Destroys the values in the configuration space.
Definition: val_domain.cpp:203
void set_zero()
Sets the Val_domain to zero (logical state to zero and arrays destroyed).
Definition: val_domain.cpp:223
void allocate_coef()
Allocates the values in the coefficient space and destroys the values in the configuration space.
Definition: val_domain.cpp:216
bool check_if_zero() const
Check whether the logical state is zero or not.
Definition: val_domain.hpp:142
Val_domain mult_sin_theta() const
Multiplication by .
Array< double > * cf
Pointer on the Array of the values in the coefficients space.
Definition: val_domain.hpp:77
void std_r_base()
Sets the basis for the radius.
Definition: val_domain.cpp:262
double & set(const Index &pos)
Read/write the value of the field in the configuration space.
Definition: val_domain.cpp:171
void std_base()
Sets the standard basis of decomposition.
Definition: val_domain.cpp:246
Val_domain mult_cos_theta() const
Multiplication by .
void coef() const
Computes the coefficients.
Definition: val_domain.cpp:622
Val_domain der_var(int i) const
Computes the derivative with respect to a numerical coordinate.
Definition: val_domain.cpp:670
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