KADATH
domain_shell.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 
22 #include "utilities.hpp"
23 #include "spheric.hpp"
24 #include "point.hpp"
25 #include "val_domain.hpp"
26 namespace Kadath {
27 // Standard constructor
28 Domain_shell::Domain_shell (int num, int ttype, double rint, double rext, const Point& cr, const Dim_array& nbr) : Domain(num, ttype, nbr),
29  alpha((rext-rint)/2.), beta((rext+rint)/2.), center(cr) {
30 
31  assert (nbr.get_ndim()==3) ;
32  assert (cr.get_ndim()==3) ; // Affectation de type_point :
33  do_coloc() ;
34 
35 }
36 
37 
38 // Constructor by copy
39 Domain_shell::Domain_shell (const Domain_shell& so) : Domain(so), alpha(so.alpha),
40  beta(so.beta), center(so.center){}
41 
42 
43 
44 Domain_shell::Domain_shell (int num, FILE* fd) : Domain(num, fd), center(fd) {
45  fread_be (&alpha, sizeof(double), 1, fd) ;
46  fread_be (&beta, sizeof(double), 1, fd) ;
47  do_coloc() ;
48 }
49 
50 // Destructor
51 Domain_shell::~Domain_shell() {}
52 
53 void Domain_shell::save (FILE* fd) const {
54  nbr_points.save(fd) ;
55  nbr_coefs.save(fd) ;
56  fwrite_be (&ndim, sizeof(int), 1, fd) ;
57  fwrite_be (&type_base, sizeof(int), 1, fd) ;
58  center.save(fd) ;
59  fwrite_be (&alpha, sizeof(double), 1, fd) ;
60  fwrite_be (&beta, sizeof(double), 1, fd) ;
61 }
62 
63 ostream& Domain_shell::print (ostream& o) const {
64  o << "Shell" << endl ;
65  o << beta-alpha << " < r < " << beta+alpha << endl ;
66  o << "Center = " << center << endl ;
67  o << "Nbr pts = " << nbr_points << endl ;
68  o << endl ;
69  return o ;
70 }
71 
72 
73 
74 Val_domain Domain_shell::der_normal (const Val_domain& so, int bound) const {
75 
76  Val_domain res (so.der_var(1)) ;
77  switch (bound) {
78  case OUTER_BC :
79  res /= alpha ;
80  break ;
81  case INNER_BC :
82  res /= alpha ;
83  break ;
84  default:
85  cerr << "Unknown boundary case in Domain_shell::der_normal" << endl ;
86  abort() ;
87  }
88 return res ;
89 }
90 
91 // Computes the Cartesian coordinates
92 void Domain_shell::do_absol () const {
93  for (int i=0 ; i<3 ; i++)
94  assert (coloc[i] != 0x0) ;
95  for (int i=0 ; i<3 ; i++)
96  assert (absol[i] == 0x0) ;
97  for (int i=0 ; i<3 ; i++) {
98  absol[i] = new Val_domain(this) ;
99  absol[i]->allocate_conf() ;
100  }
101  Index index (nbr_points) ;
102 
103  do {
104  absol[0]->set(index) = (alpha* ((*coloc[0])(index(0))) +beta)*
105  sin((*coloc[1])(index(1)))*cos((*coloc[2])(index(2))) + center(1) ;
106  absol[1]->set(index) = (alpha* ((*coloc[0])(index(0))) +beta) *
107  sin((*coloc[1])(index(1)))*sin((*coloc[2])(index(2))) + center(2) ;
108  absol[2]->set(index) = (alpha* ((*coloc[0])(index(0))) + beta) * cos((*coloc[1])(index(1))) + center(3) ;
109  }
110  while (index.inc()) ;
111 }
112 
113 // Computes the radius
114 void Domain_shell::do_radius () const {
115  for (int i=0 ; i<3 ; i++)
116  assert (coloc[i] != 0x0) ;
117  assert (radius == 0x0) ;
118  radius = new Val_domain(this) ;
119  radius->allocate_conf() ;
120  Index index (nbr_points) ;
121  do
122  radius->set(index) = alpha* ((*coloc[0])(index(0))) + beta ;
123  while (index.inc()) ;
124 }
125 
126 // Computes the Cartesian coordinates
127 void Domain_shell::do_cart () const {
128  for (int i=0 ; i<3 ; i++)
129  assert (coloc[i] != 0x0) ;
130  for (int i=0 ; i<3 ; i++)
131  assert (cart[i] == 0x0) ;
132  for (int i=0 ; i<3 ; i++) {
133  cart[i] = new Val_domain(this) ;
134  cart[i]->allocate_conf() ;
135  }
136  Index index (nbr_points) ;
137 
138  do {
139  cart[0]->set(index) = (alpha* ((*coloc[0])(index(0))) +beta)*
140  sin((*coloc[1])(index(1)))*cos((*coloc[2])(index(2))) + center(1) ;
141  cart[1]->set(index) = (alpha* ((*coloc[0])(index(0))) +beta) *
142  sin((*coloc[1])(index(1)))*sin((*coloc[2])(index(2))) + center(2) ;
143  cart[2]->set(index) = (alpha* ((*coloc[0])(index(0))) + beta) * cos((*coloc[1])(index(1))) + center(3) ;
144  }
145  while (index.inc()) ;
146 }
147 
148 // Computes the Cartesian coordinates
150  for (int i=0 ; i<3 ; i++)
151  assert (coloc[i] != 0x0) ;
152  for (int i=0 ; i<3 ; i++)
153  assert (cart_surr[i] == 0x0) ;
154  for (int i=0 ; i<3 ; i++) {
155  cart_surr[i] = new Val_domain(this) ;
156  cart_surr[i]->allocate_conf() ;
157  }
158  Index index (nbr_points) ;
159 
160  do {
161  cart_surr[0]->set(index) = sin((*coloc[1])(index(1)))*cos((*coloc[2])(index(2))) ;
162  cart_surr[1]->set(index) = sin((*coloc[1])(index(1)))*sin((*coloc[2])(index(2))) ;
163  cart_surr[2]->set(index) = cos((*coloc[1])(index(1))) ;
164  }
165  while (index.inc()) ;
166 }
167 
168 // Check if a point is inside this domain
169 bool Domain_shell::is_in (const Point& xx, double prec) const {
170 
171  assert (xx.get_ndim()==3) ;
172 
173  double x_loc = xx(1) - center(1) ;
174  double y_loc = xx(2) - center(2) ;
175  double z_loc = xx(3) - center(3) ;
176  double air_loc = sqrt (x_loc*x_loc + y_loc*y_loc + z_loc*z_loc) ;
177 
178  bool res = ((air_loc / (alpha+beta) - 1 <= prec) && (air_loc/ (beta-alpha) - 1 >= -prec)) ? true : false ;
179  return res ;
180 }
181 
182 //Computes the numerical coordinates, as a function of the absolute ones.
183 const Point Domain_shell::absol_to_num(const Point& abs) const {
184 
185  assert (is_in(abs)) ;
186  Point num(3) ;
187 
188  double x_loc = abs(1) - center(1) ;
189  double y_loc = abs(2) - center(2) ;
190  double z_loc = abs(3) - center(3) ;
191  double air = sqrt(x_loc*x_loc+y_loc*y_loc+z_loc*z_loc) ;
192  num.set(1) = (air-beta)/alpha ;
193  double rho = sqrt(x_loc*x_loc+y_loc*y_loc) ;
194 
195  if (rho==0) {
196  // On the axis ?
197  num.set(2) = (z_loc>=0) ? 0 : M_PI ;
198  num.set(3) = 0 ;
199  }
200  else {
201  num.set(2) = atan(rho/z_loc) ;
202  num.set(3) = atan2 (y_loc, x_loc) ;
203  }
204 
205  if (num(2) <0)
206  num.set(2) = M_PI + num(2) ;
207 
208  return num ;
209 }
210 
211 // Convert absolute coordinates to numerical ones
212 const Point Domain_shell::absol_to_num_bound(const Point& abs, int bound) const {
213 
214  assert ((bound==OUTER_BC) || (bound==INNER_BC)) ;
215  assert (is_in(abs, 1e-3)) ;
216  Point num(3) ;
217 
218  double x_loc = abs(1) - center(1) ;
219  double y_loc = abs(2) - center(2) ;
220  double z_loc = abs(3) - center(3) ;
221 
222  switch (bound) {
223  case INNER_BC:
224  num.set(1) = -1 ;
225  break ;
226  case OUTER_BC:
227  num.set(1) = 1 ;
228  break ;
229  default:
230  cerr << "unknown boundary in Domain_shell::absol_to_num" << endl ;
231  abort() ;
232  }
233 
234  double rho = sqrt(x_loc*x_loc+y_loc*y_loc) ;
235 
236  if (rho==0) {
237  // Sur l'axe
238  num.set(2) = (z_loc>=0) ? 0 : M_PI ;
239  num.set(3) = 0 ;
240  }
241  else {
242  num.set(2) = atan(rho/z_loc) ;
243  num.set(3) = atan2 (y_loc, x_loc) ;
244  }
245 
246  if (num(2) <0)
247  num.set(2) = M_PI + num(2) ;
248 
249  return num ;
250 }
251 double coloc_leg (int, int) ;
253 
254  switch (type_base) {
255  case CHEB_TYPE:
257  nbr_coefs.set(2) += 2 ;
258  del_deriv() ;
259  for (int i=0 ; i<ndim ; i++)
260  coloc[i] = new Array<double> (nbr_points(i)) ;
261  for (int i=0 ; i<nbr_points(0) ; i++)
262  coloc[0]->set(i) = -cos(M_PI*i/(nbr_points(0)-1)) ;
263  for (int j=0 ; j<nbr_points(1) ; j++)
264  coloc[1]->set(j) = M_PI/2.*j/(nbr_points(1)-1) ;
265  for (int k=0 ; k<nbr_points(2) ; k++)
266  coloc[2]->set(k) = M_PI*2.*k/nbr_points(2) ;
267  break ;
268  case LEG_TYPE:
270  nbr_coefs.set(2) += 2 ;
271  del_deriv() ;
272  for (int i=0 ; i<ndim ; i++)
273  coloc[i] = new Array<double> (nbr_points(i)) ;
274  for (int i=0 ; i<nbr_points(0) ; i++)
275  coloc[0]->set(i) = coloc_leg(i, nbr_points(0)) ;
276  for (int j=0 ; j<nbr_points(1) ; j++)
277  coloc[1]->set(j) = M_PI/2.*j/(nbr_points(1)-1) ;
278  for (int k=0 ; k<nbr_points(2) ; k++)
279  coloc[2]->set(k) = M_PI*2.*k/nbr_points(2) ;
280  break ;
281  default :
282  cerr << "Unknown type of basis in Domain_shell::do_coloc" << endl ;
283  abort() ;
284  }
285 }
286 
287 // standard base for a symetric function in z, using Chebyshev
289 
290  int m ;
291 
292  assert (type_base == CHEB_TYPE) ;
293 
294  base.allocate(nbr_coefs) ;
295 
296  base.def =true ;
297  base.bases_1d[2]->set(0) = COSSIN ;
298 
299  Index index (base.bases_1d[0]->get_dimensions()) ;
300 
301  for (int k=0 ; k<nbr_coefs(2) ; k++) {
302  m = (k%2==0) ? k/2 : (k-1)/2 ;
303  base.bases_1d[1]->set(k) = (m%2==0) ? COS_EVEN : SIN_ODD ;
304  for (int j=0 ; j<nbr_coefs(1) ; j++) {
305  index.set(0) = j ; index.set(1) = k ;
306  base.bases_1d[0]->set(index) = CHEB ;
307  }
308  }
309 }
310 
311 void Domain_shell::set_cheb_base_with_m(Base_spectral& base, int mquant) const {
312 
313  int m ;
314 
315  assert (type_base == CHEB_TYPE) ;
316 
317  base.allocate(nbr_coefs) ;
318 
319  if (mquant%2==0) {
320  base.def =true ;
321  base.bases_1d[2]->set(0) = COSSIN ;
322 
323  Index index (base.bases_1d[0]->get_dimensions()) ;
324 
325  for (int k=0 ; k<nbr_coefs(2) ; k++) {
326  m = (k%2==0) ? k/2 : (k-1)/2 ;
327  base.bases_1d[1]->set(k) = (m%2==0) ? COS_EVEN : SIN_ODD ;
328  for (int j=0 ; j<nbr_coefs(1) ; j++) {
329  index.set(0) = j ; index.set(1) = k ;
330  base.bases_1d[0]->set(index) = CHEB ;
331  }
332  }
333  }
334  else {
335  base.def =true ;
336  base.bases_1d[2]->set(0) = COSSIN ;
337 
338  Index index (base.bases_1d[0]->get_dimensions()) ;
339 
340  for (int k=0 ; k<nbr_coefs(2) ; k++) {
341  m = (k%2==0) ? k/2 : (k-1)/2 ;
342  base.bases_1d[1]->set(k) = (m%2==0) ? SIN_ODD : COS_EVEN ;
343  for (int j=0 ; j<nbr_coefs(1) ; j++) {
344  index.set(0) = j ; index.set(1) = k ;
345  base.bases_1d[0]->set(index) = CHEB ;
346  }
347  }
348  }
349 }
350 
352  set_cheb_base(base) ;
353 }
354 
356  set_legendre_base(base) ;
357 }
358 
359 // standard base for a anti-symetric function in z, using Chebyshev
361 
362  int m ;
363 
364  assert (type_base == CHEB_TYPE) ;
365 
366  base.allocate(nbr_coefs) ;
367 
368  base.def =true ;
369  base.bases_1d[2]->set(0) = COSSIN ;
370 
371  Index index (base.bases_1d[0]->get_dimensions()) ;
372 
373  for (int k=0 ; k<nbr_coefs(2) ; k++) {
374  m = (k%2==0) ? k/2 : (k-1)/2 ;
375  base.bases_1d[1]->set(k) = (m%2==0) ? COS_ODD : SIN_EVEN ;
376  for (int j=0 ; j<nbr_coefs(1) ; j++) {
377  index.set(0) = j ; index.set(1) = k ;
378  base.bases_1d[0]->set(index) = CHEB ;
379  }
380  }
381 }
382 
384 
385  int m ;
386 
387  assert (type_base == CHEB_TYPE) ;
388 
389  base.allocate(nbr_coefs) ;
390 
391  base.def =true ;
392  base.bases_1d[2]->set(0) = COSSIN ;
393 
394  Index index (base.bases_1d[0]->get_dimensions()) ;
395 
396  for (int k=0 ; k<nbr_coefs(2) ; k++) {
397  m = (k%2==0) ? k/2 : (k-1)/2 ;
398  base.bases_1d[1]->set(k) = (m%2==0) ? COS_EVEN : SIN_ODD ;
399  for (int j=0 ; j<nbr_coefs(1) ; j++) {
400  index.set(0) = j ; index.set(1) = k ;
401  base.bases_1d[0]->set(index) = CHEB ;
402  }
403  }
404 }
405 
407 
408  int m ;
409 
410  assert (type_base == CHEB_TYPE) ;
411 
412  base.allocate(nbr_coefs) ;
413 
414  base.def =true ;
415  base.bases_1d[2]->set(0) = COSSIN ;
416 
417  Index index (base.bases_1d[0]->get_dimensions()) ;
418 
419  for (int k=0 ; k<nbr_coefs(2) ; k++) {
420  m = (k%2==0) ? k/2 : (k-1)/2 ;
421  base.bases_1d[1]->set(k) = (m%2==0) ? SIN_EVEN : COS_ODD ;
422  for (int j=0 ; j<nbr_coefs(1) ; j++) {
423  index.set(0) = j ; index.set(1) = k ;
424  base.bases_1d[0]->set(index) = CHEB ;
425  }
426  }
427 }
428 
430 
431  int m ;
432 
433  assert (type_base == CHEB_TYPE) ;
434 
435  base.allocate(nbr_coefs) ;
436 
437  base.def =true ;
438  base.bases_1d[2]->set(0) = COSSIN ;
439 
440  Index index (base.bases_1d[0]->get_dimensions()) ;
441 
442  for (int k=0 ; k<nbr_coefs(2) ; k++) {
443  m = (k%2==0) ? k/2 : (k-1)/2 ;
444  base.bases_1d[1]->set(k) = (m%2==0) ? SIN_ODD : COS_EVEN ;
445  for (int j=0 ; j<nbr_coefs(1) ; j++) {
446  index.set(0) = j ; index.set(1) = k ;
447  base.bases_1d[0]->set(index) = CHEB ;
448  }
449  }
450 }
452 
453  int m ;
454 
455  assert (type_base == CHEB_TYPE) ;
456 
457  base.allocate(nbr_coefs) ;
458 
459  base.def =true ;
460  base.bases_1d[2]->set(0) = COSSIN ;
461 
462  Index index (base.bases_1d[0]->get_dimensions()) ;
463 
464  for (int k=0 ; k<nbr_coefs(2) ; k++) {
465  m = (k%2==0) ? k/2 : (k-1)/2 ;
466  base.bases_1d[1]->set(k) = (m%2==0) ? COS_EVEN : SIN_ODD ;
467  for (int j=0 ; j<nbr_coefs(1) ; j++) {
468  index.set(0) = j ; index.set(1) = k ;
469  base.bases_1d[0]->set(index) = CHEB ;
470  }
471  }
472 }
473 
475 
476  int m ;
477 
478  assert (type_base == CHEB_TYPE) ;
479 
480  base.allocate(nbr_coefs) ;
481 
482  base.def =true ;
483  base.bases_1d[2]->set(0) = COSSIN ;
484 
485  Index index (base.bases_1d[0]->get_dimensions()) ;
486 
487  for (int k=0 ; k<nbr_coefs(2) ; k++) {
488  m = (k%2==0) ? k/2 : (k-1)/2 ;
489  base.bases_1d[1]->set(k) = (m%2==0) ? SIN_ODD : COS_EVEN ;
490  for (int j=0 ; j<nbr_coefs(1) ; j++) {
491  index.set(0) = j ; index.set(1) = k ;
492  base.bases_1d[0]->set(index) = CHEB ;
493  }
494  }
495 }
496 
498 
499  int m ;
500 
501  assert (type_base == CHEB_TYPE) ;
502 
503  base.allocate(nbr_coefs) ;
504 
505  base.def =true ;
506  base.bases_1d[2]->set(0) = COSSIN ;
507 
508  Index index (base.bases_1d[0]->get_dimensions()) ;
509 
510  for (int k=0 ; k<nbr_coefs(2) ; k++) {
511  m = (k%2==0) ? k/2 : (k-1)/2 ;
512  base.bases_1d[1]->set(k) = (m%2==0) ? SIN_EVEN: COS_ODD ;
513  for (int j=0 ; j<nbr_coefs(1) ; j++) {
514  index.set(0) = j ; index.set(1) = k ;
515  base.bases_1d[0]->set(index) = CHEB ;
516  }
517  }
518 }
519 
520 // standard base for a symetric function in z, using Legendre
522 
523  int m ;
524 
525  assert (type_base == LEG_TYPE) ;
526  base.allocate(nbr_coefs) ;
527 
528  base.def = true ;
529  base.bases_1d[2]->set(0) = COSSIN ;
530 
531  Index index (base.bases_1d[0]->get_dimensions()) ;
532 
533  for (int k=0 ; k<nbr_coefs(2) ; k++) {
534  m = (k%2==0) ? k/2 : (k-1)/2 ;
535  base.bases_1d[1]->set(k) = (m%2==0) ? COS_EVEN : SIN_ODD ;
536  for (int j=0 ; j<nbr_coefs(1) ; j++) {
537  index.set(0) = j ; index.set(1) = k ;
538  base.bases_1d[0]->set(index) = LEG ;
539  }
540  }
541 }
542 
543 // standard base for a anti-symetric function in z, using Legendre
545 
546  int m ;
547  assert (type_base == LEG_TYPE) ;
548 
549  base.allocate(nbr_coefs) ;
550 
551  base.def = true ;
552  base.bases_1d[2]->set(0) = COSSIN ;
553 
554  Index index (base.bases_1d[0]->get_dimensions()) ;
555 
556  for (int k=0 ; k<nbr_coefs(2) ; k++) {
557  m = (k%2==0) ? k/2 : (k-1)/2 ;
558  base.bases_1d[1]->set(k) = (m%2==0) ? COS_ODD : SIN_EVEN ;
559  for (int j=0 ; j<nbr_coefs(1) ; j++) {
560  index.set(0) = j ; index.set(1) = k ;
561  base.bases_1d[0]->set(index) = LEG ;
562  }
563  }
564 }
565 
567 
568  int m ;
569 
570  assert (type_base == LEG_TYPE) ;
571 
572  base.allocate(nbr_coefs) ;
573 
574  base.def =true ;
575  base.bases_1d[2]->set(0) = COSSIN ;
576 
577  Index index (base.bases_1d[0]->get_dimensions()) ;
578 
579  for (int k=0 ; k<nbr_coefs(2) ; k++) {
580  m = (k%2==0) ? k/2 : (k-1)/2 ;
581  base.bases_1d[1]->set(k) = (m%2==0) ? COS_EVEN : SIN_ODD ;
582  for (int j=0 ; j<nbr_coefs(1) ; j++) {
583  index.set(0) = j ; index.set(1) = k ;
584  base.bases_1d[0]->set(index) = LEG ;
585  }
586  }
587 }
588 
590 
591  int m ;
592 
593  assert (type_base == LEG_TYPE) ;
594 
595  base.allocate(nbr_coefs) ;
596 
597  base.def =true ;
598  base.bases_1d[2]->set(0) = COSSIN ;
599 
600  Index index (base.bases_1d[0]->get_dimensions()) ;
601 
602  for (int k=0 ; k<nbr_coefs(2) ; k++) {
603  m = (k%2==0) ? k/2 : (k-1)/2 ;
604  base.bases_1d[1]->set(k) = (m%2==0) ? SIN_EVEN : COS_ODD ;
605  for (int j=0 ; j<nbr_coefs(1) ; j++) {
606  index.set(0) = j ; index.set(1) = k ;
607  base.bases_1d[0]->set(index) = LEG ;
608  }
609  }
610 }
611 
613 
614  int m ;
615 
616  assert (type_base == LEG_TYPE) ;
617 
618  base.allocate(nbr_coefs) ;
619 
620  base.def =true ;
621  base.bases_1d[2]->set(0) = COSSIN ;
622 
623  Index index (base.bases_1d[0]->get_dimensions()) ;
624 
625  for (int k=0 ; k<nbr_coefs(2) ; k++) {
626  m = (k%2==0) ? k/2 : (k-1)/2 ;
627  base.bases_1d[1]->set(k) = (m%2==0) ? SIN_ODD : COS_EVEN ;
628  for (int j=0 ; j<nbr_coefs(1) ; j++) {
629  index.set(0) = j ; index.set(1) = k ;
630  base.bases_1d[0]->set(index) = LEG ;
631  }
632  }
633 }
634 
636 
637  int m ;
638 
639  assert (type_base == LEG_TYPE) ;
640 
641  base.allocate(nbr_coefs) ;
642 
643  base.def =true ;
644  base.bases_1d[2]->set(0) = COSSIN ;
645 
646  Index index (base.bases_1d[0]->get_dimensions()) ;
647 
648  for (int k=0 ; k<nbr_coefs(2) ; k++) {
649  m = (k%2==0) ? k/2 : (k-1)/2 ;
650  base.bases_1d[1]->set(k) = (m%2==0) ? COS_EVEN : SIN_ODD ;
651  for (int j=0 ; j<nbr_coefs(1) ; j++) {
652  index.set(0) = j ; index.set(1) = k ;
653  base.bases_1d[0]->set(index) = LEG ;
654  }
655  }
656 }
657 
659 
660  int m ;
661 
662  assert (type_base == LEG_TYPE) ;
663 
664  base.allocate(nbr_coefs) ;
665 
666  base.def =true ;
667  base.bases_1d[2]->set(0) = COSSIN ;
668 
669  Index index (base.bases_1d[0]->get_dimensions()) ;
670 
671  for (int k=0 ; k<nbr_coefs(2) ; k++) {
672  m = (k%2==0) ? k/2 : (k-1)/2 ;
673  base.bases_1d[1]->set(k) = (m%2==0) ? SIN_ODD : COS_EVEN ;
674  for (int j=0 ; j<nbr_coefs(1) ; j++) {
675  index.set(0) = j ; index.set(1) = k ;
676  base.bases_1d[0]->set(index) = LEG ;
677  }
678  }
679 }
680 
682 
683  int m ;
684 
685  assert (type_base == LEG_TYPE) ;
686 
687  base.allocate(nbr_coefs) ;
688 
689  base.def =true ;
690  base.bases_1d[2]->set(0) = COSSIN ;
691 
692  Index index (base.bases_1d[0]->get_dimensions()) ;
693 
694  for (int k=0 ; k<nbr_coefs(2) ; k++) {
695  m = (k%2==0) ? k/2 : (k-1)/2 ;
696  base.bases_1d[1]->set(k) = (m%2==0) ? SIN_EVEN : COS_ODD ;
697  for (int j=0 ; j<nbr_coefs(1) ; j++) {
698  index.set(0) = j ; index.set(1) = k ;
699  base.bases_1d[0]->set(index) = LEG ;
700  }
701  }
702 }
703 
704 // Computes the derivatives with respect to XYZ as a function of the numerical ones.
705 void Domain_shell::do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const {
706 
707  // d/dx :
708  Val_domain sintdr (der_var[0]->mult_sin_theta()/alpha) ;
709  Val_domain dtsr (*der_var[1]/get_radius()) ;
710  dtsr.set_base() = der_var[1]->get_base() ;
711  Val_domain dpsr (*der_var[2]/get_radius()) ;
712  dpsr.set_base() = der_var[2]->get_base() ;
713  Val_domain costdtsr (dtsr.mult_cos_theta()) ;
714 
715  Val_domain dpsrssint (dpsr.div_sin_theta()) ;
716  der_abs[0] = new Val_domain ((sintdr+costdtsr).mult_cos_phi() - dpsrssint.mult_sin_phi()) ;
717 
718  // d/dy :
719  der_abs[1] = new Val_domain ((sintdr+costdtsr).mult_sin_phi() + dpsrssint.mult_cos_phi()) ;
720  // d/dz :
721  der_abs[2] = new Val_domain (der_var[0]->mult_cos_theta()/alpha - dtsr.mult_sin_theta()) ;
722 }
723 
724 // Rules for the multiplication of two basis.
726 
727  assert (a.ndim==3) ;
728  assert (b.ndim==3) ;
729 
730  Base_spectral res(3) ;
731  bool res_def = true ;
732 
733  if (!a.def)
734  res_def=false ;
735  if (!b.def)
736  res_def=false ;
737 
738  if (res_def) {
739 
740  // Base in phi :
741  res.bases_1d[2] = new Array<int> (a.bases_1d[2]->get_dimensions()) ;
742  switch ((*a.bases_1d[2])(0)) {
743  case COSSIN:
744  switch ((*b.bases_1d[2])(0)) {
745  case COSSIN:
746  res.bases_1d[2]->set(0) = COSSIN ;
747  break ;
748  default:
749  res_def = false ;
750  break ;
751  }
752  break ;
753  default:
754  res_def = false ;
755  break ;
756  }
757 
758  // Bases in theta :
759  // On check l'alternance :
760  Index index_1 (a.bases_1d[1]->get_dimensions()) ;
761  res.bases_1d[1] = new Array<int> (a.bases_1d[1]->get_dimensions()) ;
762  do {
763  switch ((*a.bases_1d[1])(index_1)) {
764  case COS_EVEN:
765  switch ((*b.bases_1d[1])(index_1)) {
766  case COS_EVEN:
767  res.bases_1d[1]->set(index_1) = (index_1(0)%4<2) ? COS_EVEN : SIN_ODD ;
768  break ;
769  case COS_ODD:
770  res.bases_1d[1]->set(index_1) = (index_1(0)%4<2) ? COS_ODD : SIN_EVEN ;
771  break ;
772  case SIN_EVEN:
773  res.bases_1d[1]->set(index_1) = (index_1(0)%4<2) ? SIN_EVEN : COS_ODD ;
774  break ;
775  case SIN_ODD:
776  res.bases_1d[1]->set(index_1) = (index_1(0)%4<2) ? SIN_ODD : COS_EVEN ;
777  break ;
778  default:
779  res_def = false ;
780  break ;
781  }
782  break ;
783  case COS_ODD:
784  switch ((*b.bases_1d[1])(index_1)) {
785  case COS_EVEN:
786  res.bases_1d[1]->set(index_1) = (index_1(0)%4<2) ? COS_ODD : SIN_EVEN ;
787  break ;
788  case COS_ODD:
789  res.bases_1d[1]->set(index_1) = (index_1(0)%4<2) ? COS_EVEN : SIN_ODD ;
790  break ;
791  case SIN_EVEN:
792  res.bases_1d[1]->set(index_1) = (index_1(0)%4<2) ? SIN_ODD : COS_EVEN ;
793  break ;
794  case SIN_ODD:
795  res.bases_1d[1]->set(index_1) = (index_1(0)%4<2) ? SIN_EVEN : COS_ODD ;
796  break ;
797  default:
798  res_def = false ;
799  break ;
800  }
801  break ;
802  case SIN_EVEN:
803  switch ((*b.bases_1d[1])(index_1)) {
804  case COS_EVEN:
805  res.bases_1d[1]->set(index_1) = (index_1(0)%4<2) ? SIN_EVEN : COS_ODD ;
806  break ;
807  case COS_ODD:
808  res.bases_1d[1]->set(index_1) = (index_1(0)%4<2) ? SIN_ODD : COS_EVEN ;
809  break ;
810  case SIN_EVEN:
811  res.bases_1d[1]->set(index_1) = (index_1(0)%4<2) ? COS_EVEN : SIN_ODD ;
812  break ;
813  case SIN_ODD:
814  res.bases_1d[1]->set(index_1) = (index_1(0)%4<2) ? COS_ODD : SIN_EVEN ;
815  break ;
816  default:
817  res_def = false ;
818  break ;
819  }
820  break ;
821  case SIN_ODD:
822  switch ((*b.bases_1d[1])(index_1)) {
823  case COS_EVEN:
824  res.bases_1d[1]->set(index_1) = (index_1(0)%4<2) ? SIN_ODD : COS_EVEN ;
825  break ;
826  case COS_ODD:
827  res.bases_1d[1]->set(index_1) = (index_1(0)%4<2) ? SIN_EVEN : COS_ODD ;
828  break ;
829  case SIN_EVEN:
830  res.bases_1d[1]->set(index_1) = (index_1(0)%4<2) ? COS_ODD : SIN_EVEN ;
831  break ;
832  case SIN_ODD:
833  res.bases_1d[1]->set(index_1) = (index_1(0)%4<2) ? COS_EVEN : SIN_ODD ;
834  break ;
835  default:
836  res_def = false ;
837  break ;
838  }
839  break ;
840  default:
841  res_def = false ;
842  break ;
843  }
844  }
845  while (index_1.inc()) ;
846 
847 
848  // Base in r :
849  Index index_0 (a.bases_1d[0]->get_dimensions()) ;
850  res.bases_1d[0] = new Array<int> (a.bases_1d[0]->get_dimensions()) ;
851  do {
852  switch ((*a.bases_1d[0])(index_0)) {
853  case CHEB:
854  switch ((*b.bases_1d[0])(index_0)) {
855  case CHEB:
856  res.bases_1d[0]->set(index_0) = CHEB ;
857  break ;
858  default:
859  res_def = false ;
860  break ;
861  }
862  break ;
863  case LEG:
864  switch ((*b.bases_1d[0])(index_0)) {
865  case LEG:
866  res.bases_1d[0]->set(index_0) = LEG ;
867  break ;
868  default:
869  res_def = false ;
870  break ;
871  }
872  break ;
873  default:
874  res_def = false ;
875  break ;
876  }
877  }
878  while (index_0.inc()) ;
879  }
880 
881  if (!res_def)
882  for (int dim=0 ; dim<a.ndim ; dim++)
883  if (res.bases_1d[dim]!= 0x0) {
884  delete res.bases_1d[dim] ;
885  res.bases_1d[dim] = 0x0 ;
886  }
887  res.def = res_def ;
888  return res ;
889 }
890 
891 int Domain_shell::give_place_var (char* p) const {
892  int res = -1 ;
893  if (strcmp(p,"R ")==0)
894  res = 0 ;
895  if (strcmp(p,"T ")==0)
896  res = 1 ;
897  if (strcmp(p,"P ")==0)
898  res = 2 ;
899  return res ;
900 }
901 
902 double Domain_shell::integ(const Val_domain& so, int bound) const
903 {
904  Val_domain rrso (mult_r(mult_r(mult_sin_theta(so)))) ;
905 
906  double res = 0 ;
907  if (!so.check_if_zero())
908  {
909 
910  int baset = (*rrso.get_base().bases_1d[1]) (0) ;
911  Index pcf (nbr_coefs) ;
912  switch (baset) {
913  case COS_ODD :
914  break ;
915  case SIN_EVEN :
916  break ;
917  case COS_EVEN : {
918  res += M_PI*val_boundary(bound, rrso, pcf) ;
919  break ;
920  }
921  case SIN_ODD : {
922  for (int j=0 ; j<nbr_coefs(1) ; j++) {
923  pcf.set(1) = j ;
924  res += 2./(2*double(j)+1) * val_boundary(bound, rrso, pcf) ;
925  }
926  break ;
927  }
928 
929  default :
930  cerr << "Case not yet implemented in Domain_shell::integ" << endl ;
931  abort() ;
932  }
933  res *= 2*M_PI ;
934  }
935  return res ;
936 }
937 }
Class for storing the basis of decompositions of a field.
Bases_container bases_1d
Arrays containing the various basis of decomposition.
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
int & set(int i)
Read/write of the size of a given dimension.
Definition: dim_array.hpp:54
void save(FILE *) const
Save function.
Definition: dim_array.cpp:32
Class for a spherical shell and a symmetry with respect to the plane .
Definition: spheric.hpp:555
virtual Val_domain mult_cos_theta(const Val_domain &) const
Multiplication by .
virtual void set_legendre_base_t_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
virtual void set_cheb_base_with_m(Base_spectral &so, int m) const
Gives the standard base using Chebyshev polynomials.
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
virtual double integ(const Val_domain &so, int bound) const
Surface integral on a given boundary.
double beta
Relates the numerical to the physical radii.
Definition: spheric.hpp:559
virtual void set_cheb_base_r_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the radial component of a vector in the MTZ setting.
virtual void set_cheb_base_p_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector in the MTZ setting.
virtual Base_spectral mult(const Base_spectral &, const Base_spectral &) const
Method for the multiplication of two Base_spectral.
virtual void set_cheb_base_t_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector in the MTZ setting.
double alpha
Relates the numerical to the physical radii.
Definition: spheric.hpp:558
virtual void set_cheb_r_base(Base_spectral &) const
Gives the base using odd Chebyshev polynomials$ for the radius.
virtual void do_coloc()
Computes the colocation points.
virtual void save(FILE *) const
Saving function.
virtual int give_place_var(char *) const
Translates a name of a coordinate into its corresponding numerical name.
Domain_shell(int num, int ttype, double r_int, double r_ext, const Point &cr, const Dim_array &nbr)
Standard constructor :
virtual void set_legendre_base_r_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the radial component of a vector.
virtual void set_cheb_base_p_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector.
virtual void set_legendre_base_t_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector in the MTZ context.
virtual Val_domain mult_cos_phi(const Val_domain &) const
Multiplication by .
virtual Val_domain der_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
virtual void set_legendre_base(Base_spectral &so) const
Sets the base to the standard one for Legendre polynomials.
virtual void set_cheb_base(Base_spectral &so) const
Sets the base to the standard one for Chebyshev polynomials.
virtual void set_anti_legendre_base(Base_spectral &so) const
Sets the base to the standard one for Legendre polynomials for functions antisymetric in The bases a...
virtual void set_cheb_base_r_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the radial component of a vector.
virtual void set_cheb_base_t_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector.
virtual const Point absol_to_num(const Point &xxx) const
Computes the numerical coordinates from the physical ones.
virtual Val_domain mult_sin_phi(const Val_domain &) const
Multiplication by .
virtual void set_anti_cheb_base(Base_spectral &so) const
Sets the base to the standard one for Chebyshev polynomials for functions antisymetric in The bases ...
virtual void do_absol() const
Computes the absolute coordinates.
virtual double val_boundary(int, const Val_domain &, const Index &) const
Computes the value of a field at a boundary.
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
virtual Val_domain mult_r(const Val_domain &) const
Multiplication by .
Point center
Absolute coordinates of the center.
Definition: spheric.hpp:560
virtual void set_legendre_base_p_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
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_legendre_base_r_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the radial component of a vector in the MTZ context.
virtual void set_legendre_r_base(Base_spectral &) const
Gives the base using odd Legendre polynomials$ for the radius.
virtual const Point absol_to_num_bound(const Point &, int) const
Computes the numerical coordinates from the physical ones for a point lying on a boundary.
virtual void set_legendre_base_p_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector in the MTZ context.
virtual void do_cart() const
Computes the Cartesian coordinates.
virtual Val_domain mult_sin_theta(const Val_domain &) const
Multiplication by .
virtual void do_radius() const
Computes the generalized radius.
virtual void do_cart_surr() const
Computes the Cartesian coordinates over the radius.
Abstract class that implements the fonctionnalities common to all the type of domains.
Definition: space.hpp:60
virtual void del_deriv()
Destroys the derivated members (like coloc, cart and radius), when changing the type of colocation po...
Definition: domain.cpp:77
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
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
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
Class for storing the basis of decompositions of a field and its values on both the configuration and...
Definition: val_domain.hpp:69
Val_domain mult_sin_phi() const
Multiplication by .
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 .
Val_domain mult_cos_phi() const
Multiplication by .
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 .
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
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