KADATH
polar.hpp
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 #ifndef __POLAR_HPP_
21 #define __POLAR_HPP_
22 
23 #include "space.hpp"
24 namespace Kadath {
25 
26 
45 class Domain_polar_nucleus : public Domain {
46 
47  private:
48  double alpha ;
50 
51  public:
60  Domain_polar_nucleus (int num, int ttype, double radius, const Point& cr, const Dim_array& nbr) ;
67  Domain_polar_nucleus (int num, FILE* ff) ;
68 
69  virtual ~Domain_polar_nucleus() ;
70  virtual void save (FILE*) const ;
71 
72  private:
73  virtual void do_absol () const ;
74  virtual void do_radius () const ;
75  virtual void do_cart () const ;
76  virtual void do_cart_surr () const ;
77 
78  private:
79  virtual void set_cheb_base(Base_spectral&) const ;
80  virtual void set_legendre_base(Base_spectral&) const ;
81  virtual void set_anti_cheb_base(Base_spectral&) const ;
82  virtual void set_anti_legendre_base(Base_spectral&) const ;
83  virtual void set_cheb_base_with_m(Base_spectral&, int m) const ;
84  virtual void set_legendre_base_with_m(Base_spectral&, int m) const ;
85  virtual void set_anti_cheb_base_with_m(Base_spectral&, int m) const ;
86  virtual void set_anti_legendre_base_with_m(Base_spectral&, int m) const ;
87  virtual void do_coloc () ;
88  virtual int give_place_var (char*) const ;
89 
90  public:
91  virtual Point get_center () const {return center ;} ;
92  virtual bool is_in(const Point&xx, double prec=1e-13) const ;
93  virtual const Point absol_to_num(const Point&) const;
94  virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const ;
95  virtual Base_spectral mult (const Base_spectral&, const Base_spectral&) const ;
96 
97  public:
98 
99  virtual Val_domain mult_cos_theta (const Val_domain&) const ;
100  virtual Val_domain mult_sin_theta (const Val_domain&) const ;
101  virtual Val_domain div_sin_theta (const Val_domain&) const ;
102  virtual Val_domain div_x (const Val_domain&) const ;
103  virtual Val_domain mult_r (const Val_domain&) const ;
104  virtual Val_domain div_r (const Val_domain&) const ;
105  virtual Val_domain laplacian (const Val_domain&, int) const ;
106  virtual Val_domain laplacian2 (const Val_domain&, int) const ;
107  virtual Val_domain der_r (const Val_domain&) const ;
108  virtual Val_domain dt (const Val_domain&) const ;
109  virtual Val_domain srdr (const Val_domain&) const ;
110  virtual double integrale (const Val_domain&) const ;
111  virtual double integ_volume (const Val_domain&) const ;
112 
113 
114  virtual double val_boundary (int, const Val_domain&, const Index&) const ;
115  virtual void find_other_dom (int, int, int&, int&) const ;
116  virtual Val_domain der_normal (const Val_domain&, int) const ;
117 
118  virtual int nbr_unknowns (const Tensor&, int) const ;
127  int nbr_unknowns_val_domain (const Val_domain& so, int mquant, int llim) const ;
128  virtual Array<int> nbr_conditions (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
138  int nbr_conditions_val_domain (const Val_domain& so, int mquant, int llim, int order) const ;
139  virtual Array<int> nbr_conditions_boundary (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
148  int nbr_conditions_val_domain_boundary (const Val_domain& eq, int mquant) const ;
149  virtual void export_tau (const Tensor&, int, int, Array<double>&, int&, const Array<int>&,int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
161  void export_tau_val_domain (const Val_domain& eq, int mquant, int llim, int order, Array<double>& res, int& pos_res, int ncond) const ;
162  virtual void export_tau_boundary (const Tensor&, int, int, Array<double>&, int&, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
173  void export_tau_val_domain_boundary (const Val_domain& eq, int mquant, int bound, Array<double>& res, int& pos_res, int ncond) const ;
174  virtual void affecte_tau (Tensor&, int, const Array<double>&, int&) const ;
184  void affecte_tau_val_domain (Val_domain& so, int mquant, int llim, const Array<double>& cf, int& pos_cf) const ;
185  virtual void affecte_tau_one_coef (Tensor&, int, int, int&) const ;
195  void affecte_tau_one_coef_val_domain (Val_domain& so, int mquant, int llim, int cc, int& pos_cf) const ;
196 
197 public:
198  virtual ostream& print (ostream& o) const ;
199 } ;
200 
201 
221 class Domain_polar_shell : public Domain {
222 
223  private:
224  double alpha ;
225  double beta ;
227 
228  public:
238  Domain_polar_shell (int num, int ttype, double r_int, double r_ext, const Point& cr, const Dim_array& nbr) ;
245  Domain_polar_shell (int num, FILE* ff) ;
246 
247  virtual ~Domain_polar_shell() ;
248  virtual void save (FILE*) const ;
249 
250  private:
251  virtual void do_absol () const ;
252  virtual void do_radius () const ;
253  virtual void do_cart () const ;
254  virtual void do_cart_surr () const ;
255 
256  private :
257  virtual void set_cheb_base(Base_spectral&) const ;
258  virtual void set_legendre_base(Base_spectral&) const ;
259  virtual void set_anti_cheb_base(Base_spectral&) const ;
260  virtual void set_anti_legendre_base(Base_spectral&) const ;
261  virtual void set_cheb_base_with_m(Base_spectral&, int m) const ;
262  virtual void set_legendre_base_with_m(Base_spectral&, int m) const ;
263  virtual void set_anti_cheb_base_with_m(Base_spectral&, int m) const ;
264  virtual void set_anti_legendre_base_with_m(Base_spectral&, int m) const ;
265 
266  virtual void do_coloc () ;
267  virtual int give_place_var (char*) const ;
268  public:
269  virtual Point get_center () const {return center ;} ;
270  virtual bool is_in(const Point& xx, double prec=1e-13) const ;
271  virtual const Point absol_to_num(const Point&) const;
272  virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const ;
273 
274  virtual Base_spectral mult (const Base_spectral&, const Base_spectral&) const ;
275 
276  public:
277 
278  virtual Val_domain mult_cos_theta (const Val_domain&) const ;
279  virtual Val_domain mult_sin_theta (const Val_domain&) const ;
280  virtual Val_domain div_sin_theta (const Val_domain&) const ;
281  virtual Val_domain mult_r (const Val_domain&) const ;
282  virtual Val_domain div_r (const Val_domain&) const ;
283  virtual Val_domain laplacian (const Val_domain&, int) const ;
284  virtual Val_domain laplacian2 (const Val_domain&, int) const ;
285  virtual Val_domain der_r (const Val_domain&) const ;
286  virtual Val_domain dt (const Val_domain&) const ;
287  virtual Val_domain div_xp1 (const Val_domain&) const ;
288  virtual double integrale (const Val_domain&) const ;
289  virtual double integ_volume (const Val_domain&) const ;
290 
291  virtual void find_other_dom (int, int, int&, int&) const ;
292  virtual Val_domain der_normal (const Val_domain&, int) const ;
293  virtual double val_boundary (int, const Val_domain&, const Index&) const ;
294 
295  virtual int nbr_unknowns (const Tensor&, int) const ;
303  int nbr_unknowns_val_domain (const Val_domain& so, int mquant) const ;
304  virtual Array<int> nbr_conditions (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
313  int nbr_conditions_val_domain (const Val_domain& so, int mquant, int order) const ;
314  virtual Array<int> nbr_conditions_boundary (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
323  int nbr_conditions_val_domain_boundary (const Val_domain& eq, int mquant) const ;
324  virtual void export_tau (const Tensor&, int, int, Array<double>&, int&, const Array<int>&,int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
335  void export_tau_val_domain (const Val_domain& eq, int mquant, int order, Array<double>& res, int& pos_res, int ncond) const ;
336  virtual void export_tau_boundary (const Tensor&, int, int, Array<double>&, int&, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
347  void export_tau_val_domain_boundary (const Val_domain& eq, int mquant, int bound, Array<double>& res, int& pos_res, int ncond) const ;
348  virtual void affecte_tau (Tensor&, int, const Array<double>&, int&) const ;
357  void affecte_tau_val_domain (Val_domain& so, int mquant, const Array<double>& cf, int& pos_cf) const ;
358  virtual void affecte_tau_one_coef (Tensor&, int, int, int&) const ;
367  void affecte_tau_one_coef_val_domain (Val_domain& so, int mquant, int cc, int& pos_cf) const ;
381  void export_tau_val_domain_boundary_exception_mquant (const Val_domain& so, int mquant, int bound, Array<double>& res, int& pos_res, int ncond, const Param& param,
382  int type_exception, const Val_domain& exception) const ;
383  virtual void export_tau_boundary_exception (const Tensor&, int, int, Array<double>&, int&, const Array<int>&, const Param&, int,
384  const Tensor&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
385 public:
386  virtual ostream& print (ostream& o) const ;
387 } ;
388 
408 class Domain_polar_compact : public Domain {
409 
410  private:
411  double alpha ;
413 
414  public:
423  Domain_polar_compact (int num, int ttype, double r_int, const Point& cr, const Dim_array& nbr) ;
430  Domain_polar_compact (int num, FILE* ff) ;
431 
432  virtual ~Domain_polar_compact() ;
433  virtual void save(FILE*) const ;
434 
435  private:
436  virtual void do_radius () const ;
437 
438  private:
439  virtual void set_cheb_base(Base_spectral&) const ;
440  virtual void set_legendre_base(Base_spectral&) const ;
441  virtual void set_anti_cheb_base(Base_spectral&) const ;
442  virtual void set_anti_legendre_base(Base_spectral&) const ;
443  virtual void set_cheb_base_with_m(Base_spectral&, int m) const ;
444  virtual void set_legendre_base_with_m(Base_spectral&, int m) const ;
445  virtual void set_anti_cheb_base_with_m(Base_spectral&, int m) const ;
446  virtual void set_anti_legendre_base_with_m(Base_spectral&, int m) const ;
447  virtual void do_coloc() ;
448 
449  virtual void do_absol () const ;
450  virtual void do_cart () const ;
451  virtual void do_cart_surr () const ;
452  virtual int give_place_var (char*) const ;
453  public:
454  virtual Point get_center () const {return center ;} ;
458  double get_alpha() const {return alpha ;} ;
459 
460 
461  virtual bool is_in(const Point& xx, double prec=1e-13) const ;
462  virtual const Point absol_to_num(const Point&) const;
463  virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const ;
464 
465  virtual Base_spectral mult (const Base_spectral&, const Base_spectral&) const ;
466 
467  public:
468 
469  virtual Val_domain mult_cos_theta (const Val_domain&) const ;
470  virtual Val_domain mult_sin_theta (const Val_domain&) const ;
471  virtual Val_domain div_sin_theta (const Val_domain&) const ;
472  virtual Val_domain div_xm1 (const Val_domain&) const ;
473 
474  virtual Val_domain mult_xm1 (const Val_domain&) const ;
475  virtual Val_domain mult_r (const Val_domain&) const ;
476  virtual Val_domain div_r (const Val_domain&) const ;
477  virtual Val_domain laplacian (const Val_domain&, int) const ;
478  virtual Val_domain laplacian2 (const Val_domain&, int) const ;
479  virtual Val_domain der_r (const Val_domain&) const ;
480  virtual Val_domain der_r_rtwo (const Val_domain&) const ;
481  virtual Val_domain dt (const Val_domain&) const ;
482  virtual Val_domain div_xp1 (const Val_domain&) const ;
483  virtual double integrale (const Val_domain&) const ;
484  virtual double integ_volume (const Val_domain&) const ;
485  virtual void set_val_inf (Val_domain& so, double xx) const ;
486 
487  virtual void find_other_dom (int, int, int&, int&) const ;
488  virtual Val_domain der_normal (const Val_domain&, int) const ;
489  virtual double integ (const Val_domain&, int) const ;
490  virtual double val_boundary (int, const Val_domain&, const Index&) const ;
491 
492  virtual int nbr_unknowns (const Tensor&, int) const ;
500  int nbr_unknowns_val_domain (const Val_domain& so, int mquant) const ;
501  virtual Array<int> nbr_conditions (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
510  int nbr_conditions_val_domain (const Val_domain& eq, int mquant, int order) const ;
511  virtual Array<int> nbr_conditions_boundary (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
520  int nbr_conditions_val_domain_boundary (const Val_domain& eq, int mquant) const ;
521  virtual void export_tau (const Tensor&, int, int, Array<double>&, int&, const Array<int>&,int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
532  void export_tau_val_domain (const Val_domain& eq, int mquant, int order, Array<double>& res, int& pos_res, int ncond) const ;
533  virtual void export_tau_boundary (const Tensor&, int, int, Array<double>&, int&, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
544  void export_tau_val_domain_boundary (const Val_domain& eq, int mquant, int bound, Array<double>& res, int& pos_res, int ncond) const ;
545  virtual void affecte_tau (Tensor&, int, const Array<double>&, int&) const ;
554  void affecte_tau_val_domain (Val_domain& so, int mquant, const Array<double>& cf, int& pos_cf) const ;
555  virtual void affecte_tau_one_coef (Tensor&, int, int, int&) const ;
564  void affecte_tau_one_coef_val_domain (Val_domain& so, int mquant, int cc, int& pos_cf) const ;
565 
566 public:
567  virtual ostream& print (ostream& o) const ;
568 } ;
569 
570 
575 class Space_polar : public Space {
576  public:
584  Space_polar (int ttype, const Point& cr, const Dim_array& nbr, const Array<double>& bounds) ;
585  Space_polar (FILE*) ;
586  virtual ~Space_polar() ;
587  virtual void save(FILE*) const ;
588 
597  void add_inner_bc (System_of_eqs& syst, const char* eq, int nused=-1, Array<int>** pused=0x0) const ;
598 
606  void add_outer_bc (System_of_eqs& syst, const char* eq, int nused=-1, Array<int>** pused=0x0) const ;
607 
615  void add_eq (System_of_eqs& syst, const char* eq, int nused=-1, Array<int>** pused=0x0) const ;
616 
624  void add_eq_full (System_of_eqs& syst, const char* eq, int nused=-1, Array<int>** pused=0x0) const ;
625 
633  void add_matching (System_of_eqs& syst, const char* eq, int nused=-1, Array<int>** pused=0x0) const ;
634 
644  void add_eq (System_of_eqs& syst, const char* eq, const char* rac, const char* rac_der, int nused=-1, Array<int>** pused=0x0) const ;
645 
651  void add_eq_int_volume (System_of_eqs& syst, const char* eq) ;
652 
658  void add_eq_int_inf (System_of_eqs& syst, const char* eq) ;
659 
665  void add_eq_int_inner (System_of_eqs& syst, const char* eq) ;
666 
675  void add_eq_mode (System_of_eqs& syst, const char* f, int domtarget, int itarget, int jtarget) ;
676 
677 
683  void add_eq_ori (System_of_eqs& syst, const char* eq) ;
684 
691  void add_eq_point (System_of_eqs& syst, const Point& pp, const char* eq) ;
692 
698  double int_inf (const Scalar& so) const ;
699 } ;
700 }
701 #endif
Class for storing the basis of decompositions of a field.
Class for storing the dimensions of an array.
Definition: dim_array.hpp:34
Class for a 2-dimensional spherical shell and a symmetry with respect to the plane .
Definition: polar.hpp:408
Point center
Absolute coordinates of the center.
Definition: polar.hpp:412
virtual void export_tau_boundary(const Tensor &, int, int, Array< double > &, int &, const Array< int > &, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to a tensorial one on a given boundary It makes use ...
virtual void save(FILE *) const
Saving function.
virtual void affecte_tau(Tensor &, int, const Array< double > &, int &) const
Affects some coefficients to a Tensor.
virtual Base_spectral mult(const Base_spectral &, const Base_spectral &) const
Method for the multiplication of two Base_spectral.
virtual Val_domain mult_sin_theta(const Val_domain &) const
Multiplication by .
virtual void set_legendre_base(Base_spectral &) const
Gives the standard base for Legendre polynomials.
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
virtual void set_legendre_base_with_m(Base_spectral &, int m) const
Gives the stnadard base using Legendre polynomials.
virtual void set_val_inf(Val_domain &so, double xx) const
Sets the value at infinity of a Val_domain : not implemented for this type of Domain.
int nbr_unknowns_val_domain(const Val_domain &so, int mquant) const
Computes the number of true unknowns of a Val_domain.
void export_tau_val_domain(const Val_domain &eq, int mquant, int order, Array< double > &res, int &pos_res, int ncond) const
Exports a residual equation in the bulk.
virtual double integ_volume(const Val_domain &) const
Volume integral.
virtual Val_domain div_r(const Val_domain &) const
Division by .
virtual Val_domain der_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
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 int nbr_unknowns(const Tensor &, int) const
Computes the number of true unknowns of a Tensor, in a given domain.
virtual Val_domain laplacian(const Val_domain &, int) const
Computes the ordinary flat Laplacian for a scalar field with an harmonic index m.
virtual Val_domain der_r_rtwo(const Val_domain &) const
Compute the radial derivative multiplied by of a scalar field.
virtual void affecte_tau_one_coef(Tensor &, int, int, int &) const
Sets at most one coefficient of a Tensor to 1.
double alpha
Relates the numerical to the physical radii.
Definition: polar.hpp:411
virtual double integrale(const Val_domain &) const
Volume integral.
virtual void do_absol() const
Computes the absolute coordinates.
virtual Array< int > nbr_conditions_boundary(const Tensor &, int, int, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Computes number of discretized equations associated with a given tensorial equation on a boundary.
virtual void set_anti_legendre_base(Base_spectral &) const
Gives the base using Legendre polynomials, for functions antisymetric with respect to .
void affecte_tau_val_domain(Val_domain &so, int mquant, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
int nbr_conditions_val_domain_boundary(const Val_domain &eq, int mquant) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual void do_cart() const
Computes the Cartesian coordinates.
virtual Point get_center() const
Returns the center.
Definition: polar.hpp:454
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
virtual void do_coloc()
Computes the colocation points.
Domain_polar_compact(int num, int ttype, double r_int, const Point &cr, const Dim_array &nbr)
Standard constructor :
virtual Val_domain mult_cos_theta(const Val_domain &) const
Multiplication by .
virtual Val_domain div_xp1(const Val_domain &) const
Division by .
virtual const Point absol_to_num(const Point &) const
Computes the numerical coordinates from the physical ones.
virtual double val_boundary(int, const Val_domain &, const Index &) const
Computes the value of a field at a boundary.
virtual void do_cart_surr() const
Computes the Cartesian coordinates over the radius.
virtual Val_domain div_sin_theta(const Val_domain &) const
Division by .
virtual double integ(const Val_domain &, int) const
Surface integral on a given boundary.
void affecte_tau_one_coef_val_domain(Val_domain &so, int mquant, int cc, int &pos_cf) const
Sets at most one coefficient of a Val_domain to 1.
virtual void find_other_dom(int, int, int &, int &) const
Gives the informations corresponding the a touching neighboring domain.
virtual int give_place_var(char *) const
Translates a name of a coordinate into its corresponding numerical name.
virtual Array< int > nbr_conditions(const Tensor &, int, int, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
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_legendre_base_with_m(Base_spectral &, int m) const
Gives the base using Legendre polynomials, for functions antisymetric with respect to .
virtual Val_domain div_xm1(const Val_domain &) const
Division by .
virtual void set_cheb_base_with_m(Base_spectral &, int m) const
Gives the standard base using Chebyshev polynomials.
virtual Val_domain laplacian2(const Val_domain &, int) const
Computes the ordinary flat 2dè- Laplacian for a scalar field with an harmonic index m.
double get_alpha() const
Returns the of the mapping.
Definition: polar.hpp:458
virtual void export_tau(const Tensor &, int, int, Array< double > &, int &, const Array< int > &, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to a tensorial one in the bulk.
void export_tau_val_domain_boundary(const Val_domain &eq, int mquant, int bound, Array< double > &res, int &pos_res, int ncond) const
Exports all the residual equations corresponding to a tensorial one on a given boundary It makes use ...
virtual Val_domain dt(const Val_domain &) const
Compute the derivative with respect to of a scalar field.
virtual void do_radius() const
Computes the generalized radius.
virtual Val_domain mult_r(const Val_domain &) const
Multiplication by .
virtual Val_domain mult_xm1(const Val_domain &) const
Multiplication by .
virtual void set_cheb_base(Base_spectral &) const
Gives the standard base for Chebyshev polynomials.
int nbr_conditions_val_domain(const Val_domain &eq, int mquant, int order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual Val_domain der_r(const Val_domain &) const
Compute the radial derivative of a scalar field.
virtual void set_anti_cheb_base(Base_spectral &) const
Gives the base using Chebyshev polynomials, for functions antisymetric with respect to .
Class for a 2-dimensional spherical domain containing the origin and a symetry with respect to the pl...
Definition: polar.hpp:45
virtual Val_domain der_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
virtual void export_tau(const Tensor &, int, int, Array< double > &, int &, const Array< int > &, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to a tensorial one in the bulk.
Domain_polar_nucleus(int num, int ttype, double radius, const Point &cr, const Dim_array &nbr)
Standard constructor :
int nbr_unknowns_val_domain(const Val_domain &so, int mquant, int llim) const
Computes the number of true unknowns of a Val_domain.
virtual Val_domain div_sin_theta(const Val_domain &) const
Division by .
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 void set_legendre_base_with_m(Base_spectral &, int m) const
Gives the stnadard base using Legendre polynomials.
virtual void set_cheb_base(Base_spectral &) const
Gives the standard base for Chebyshev polynomials.
virtual void export_tau_boundary(const Tensor &, int, int, Array< double > &, int &, const Array< int > &, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to a tensorial one on a given boundary It makes use ...
virtual void affecte_tau(Tensor &, int, const Array< double > &, int &) const
Affects some coefficients to a Tensor.
virtual Val_domain der_r(const Val_domain &) const
Compute the radial derivative of a scalar field.
virtual void do_cart_surr() const
Computes the Cartesian coordinates over the radius.
virtual Val_domain div_r(const Val_domain &) const
Division by .
int nbr_conditions_val_domain(const Val_domain &so, int mquant, int llim, int order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual void do_absol() const
Computes the absolute coordinates.
virtual Val_domain laplacian2(const Val_domain &, int) const
Computes the ordinary flat 2dè- Laplacian for a scalar field with an harmonic index m.
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 Val_domain laplacian(const Val_domain &, int) const
Computes the ordinary flat Laplacian for a scalar field with an harmonic index m.
void affecte_tau_one_coef_val_domain(Val_domain &so, int mquant, int llim, int cc, int &pos_cf) const
Sets at most one coefficient of a Val_domain to 1.
virtual double val_boundary(int, const Val_domain &, const Index &) const
Computes the value of a field at a boundary.
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 Base_spectral mult(const Base_spectral &, const Base_spectral &) const
Method for the multiplication of two Base_spectral.
virtual int give_place_var(char *) const
Translates a name of a coordinate into its corresponding numerical name.
virtual Val_domain div_x(const Val_domain &) const
Division by .
virtual double integrale(const Val_domain &) const
Volume integral.
Point center
Absolute coordinates of the center.
Definition: polar.hpp:49
virtual void do_coloc()
Computes the colocation points.
virtual Val_domain mult_sin_theta(const Val_domain &) const
Multiplication by .
void export_tau_val_domain_boundary(const Val_domain &eq, int mquant, int bound, Array< double > &res, int &pos_res, int ncond) const
Exports all the residual equations corresponding to a tensorial one on a given boundary It makes use ...
int nbr_conditions_val_domain_boundary(const Val_domain &eq, int mquant) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual void set_cheb_base_with_m(Base_spectral &, int m) const
Gives the standard base using Chebyshev polynomials.
void export_tau_val_domain(const Val_domain &eq, int mquant, int llim, int order, Array< double > &res, int &pos_res, int ncond) const
Exports a residual equation in the bulk.
virtual double integ_volume(const Val_domain &) const
Volume integral.
virtual Val_domain mult_r(const Val_domain &) const
Multiplication by .
virtual Point get_center() const
Returns the center.
Definition: polar.hpp:91
virtual void set_anti_legendre_base(Base_spectral &) const
Gives the base using Legendre polynomials, for functions antisymetric with respect to .
virtual void set_legendre_base(Base_spectral &) const
Gives the standard base for Legendre polynomials.
virtual Array< int > nbr_conditions(const Tensor &, int, int, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual void do_cart() const
Computes the Cartesian coordinates.
virtual void set_anti_cheb_base(Base_spectral &) const
Gives the base using Chebyshev polynomials, for functions antisymetric with respect to .
virtual Array< int > nbr_conditions_boundary(const Tensor &, int, int, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Computes number of discretized equations associated with a given tensorial equation on a boundary.
virtual const Point absol_to_num(const Point &) const
Computes the numerical coordinates from the physical ones.
virtual void affecte_tau_one_coef(Tensor &, int, int, int &) const
Sets at most one coefficient of a Tensor to 1.
virtual Val_domain srdr(const Val_domain &) const
Compute the of a scalar field .
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
virtual void find_other_dom(int, int, int &, int &) const
Gives the informations corresponding the a touching neighboring domain.
virtual void save(FILE *) const
Saving function.
virtual Val_domain dt(const Val_domain &) const
Compute the derivative with respect to of a scalar field.
virtual void do_radius() const
Computes the generalized radius.
virtual int nbr_unknowns(const Tensor &, int) const
Computes the number of true unknowns of a Tensor, in a given domain.
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
virtual Val_domain mult_cos_theta(const Val_domain &) const
Multiplication by .
double alpha
Relates the numerical to the physical radii.
Definition: polar.hpp:48
void affecte_tau_val_domain(Val_domain &so, int mquant, int llim, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
Class for a 2-dimensional spherical shell and a symmetry with respect to the plane .
Definition: polar.hpp:221
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 void find_other_dom(int, int, int &, int &) const
Gives the informations corresponding the a touching neighboring domain.
virtual Array< int > nbr_conditions(const Tensor &, int, int, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual void affecte_tau_one_coef(Tensor &, int, int, int &) const
Sets at most one coefficient of a Tensor to 1.
virtual void set_legendre_base_with_m(Base_spectral &, int m) const
Gives the stnadard base using Legendre polynomials.
virtual Val_domain mult_r(const Val_domain &) const
Multiplication by .
virtual Val_domain mult_sin_theta(const Val_domain &) const
Multiplication by .
virtual void do_cart_surr() const
Computes the Cartesian coordinates over the radius.
double beta
Relates the numerical to the physical radii.
Definition: polar.hpp:225
virtual void set_anti_cheb_base(Base_spectral &) const
Gives the base using Chebyshev polynomials, for functions antisymetric with respect to .
virtual double val_boundary(int, const Val_domain &, const Index &) const
Computes the value of a field at a boundary.
virtual Val_domain laplacian2(const Val_domain &, int) const
Computes the ordinary flat 2dè- Laplacian for a scalar field with an harmonic index m.
virtual void export_tau(const Tensor &, int, int, Array< double > &, int &, const Array< int > &, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to a tensorial one in the bulk.
virtual Val_domain mult_cos_theta(const Val_domain &) const
Multiplication by .
Domain_polar_shell(int num, int ttype, double r_int, double r_ext, const Point &cr, const Dim_array &nbr)
Standard constructor :
virtual void affecte_tau(Tensor &, int, const Array< double > &, int &) const
Affects some coefficients to a Tensor.
virtual Base_spectral mult(const Base_spectral &, const Base_spectral &) const
Method for the multiplication of two Base_spectral.
virtual double integ_volume(const Val_domain &) const
Volume integral.
void export_tau_val_domain(const Val_domain &eq, int mquant, int order, Array< double > &res, int &pos_res, int ncond) const
Exports a residual equation in the bulk.
double alpha
Relates the numerical to the physical radii.
Definition: polar.hpp:224
virtual Array< int > nbr_conditions_boundary(const Tensor &, int, int, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Computes number of discretized equations associated with a given tensorial equation on a boundary.
int nbr_conditions_val_domain_boundary(const Val_domain &eq, int mquant) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual Point get_center() const
Returns the center.
Definition: polar.hpp:269
int nbr_conditions_val_domain(const Val_domain &so, int mquant, int order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual void save(FILE *) const
Saving function.
void export_tau_val_domain_boundary_exception_mquant(const Val_domain &so, int mquant, int bound, Array< double > &res, int &pos_res, int ncond, const Param &param, int type_exception, const Val_domain &exception) const
Exports all the residual equations corresponding to a tensorial one on a given boundary but for one c...
virtual Val_domain div_sin_theta(const Val_domain &) const
Division by .
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 Val_domain dt(const Val_domain &) const
Compute the derivative with respect to of a scalar field.
virtual const Point absol_to_num(const Point &) const
Computes the numerical coordinates from the physical ones.
virtual Val_domain der_r(const Val_domain &) const
Compute the radial derivative of a scalar field.
Point center
Absolute coordinates of the center.
Definition: polar.hpp:226
virtual void do_radius() const
Computes the generalized radius.
virtual int nbr_unknowns(const Tensor &, int) const
Computes the number of true unknowns of a Tensor, in a given domain.
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
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 Val_domain div_r(const Val_domain &) const
Division by .
void affecte_tau_one_coef_val_domain(Val_domain &so, int mquant, int cc, int &pos_cf) const
Sets at most one coefficient of a Val_domain to 1.
virtual void set_cheb_base(Base_spectral &) const
Gives the standard base for Chebyshev polynomials.
virtual Val_domain div_xp1(const Val_domain &) const
Division by .
void export_tau_val_domain_boundary(const Val_domain &eq, int mquant, int bound, Array< double > &res, int &pos_res, int ncond) const
Exports all the residual equations corresponding to a tensorial one on a given boundary It makes use ...
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
virtual int give_place_var(char *) const
Translates a name of a coordinate into its corresponding numerical name.
virtual Val_domain laplacian(const Val_domain &, int) const
Computes the ordinary flat Laplacian for a scalar field with an harmonic index m.
virtual void set_cheb_base_with_m(Base_spectral &, int m) const
Gives the standard base using Chebyshev polynomials.
virtual void export_tau_boundary(const Tensor &, int, int, Array< double > &, int &, const Array< int > &, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to a tensorial one on a given boundary It makes use ...
int nbr_unknowns_val_domain(const Val_domain &so, int mquant) const
Computes the number of true unknowns of a Val_domain.
virtual void set_legendre_base(Base_spectral &) const
Gives the standard base for Legendre polynomials.
virtual Val_domain der_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
virtual void do_cart() const
Computes the Cartesian coordinates.
virtual void export_tau_boundary_exception(const Tensor &, int, int, Array< double > &, int &, const Array< int > &, const Param &, int, const Tensor &, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to one tensorial one on a given boundary,...
virtual void do_coloc()
Computes the colocation points.
virtual void set_anti_legendre_base(Base_spectral &) const
Gives the base using Legendre polynomials, for functions antisymetric with respect to .
virtual void do_absol() const
Computes the absolute coordinates.
virtual double integrale(const Val_domain &) const
Volume integral.
void affecte_tau_val_domain(Val_domain &so, int mquant, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
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
Class that gives the position inside a multi-dimensional Array.
Definition: index.hpp:38
Parameter storage.
Definition: param.hpp:30
The class Point is used to store the coordinates of a point.
Definition: point.hpp:30
The class Scalar does not really implements scalars in the mathematical sense but rather tensorial co...
Definition: scalar.hpp:67
The Space_polar class fills the space with one polar nucleus, several polar shells and a compactified...
Definition: polar.hpp:575
void add_outer_bc(System_of_eqs &syst, const char *eq, int nused=-1, Array< int > **pused=0x0) const
Sets a boundary condition at the outer radius of the compactified polar domain.
double int_inf(const Scalar &so) const
Computes the surface integral at infinity.
Definition: space_polar.cpp:71
void add_inner_bc(System_of_eqs &syst, const char *eq, int nused=-1, Array< int > **pused=0x0) const
Sets a boundary condition at the inner radius of the innermost polar shell.
void add_eq_point(System_of_eqs &syst, const Point &pp, const char *eq)
Adds an equation being the value of some field at a given point.
void add_eq(System_of_eqs &syst, const char *eq, int nused=-1, Array< int > **pused=0x0) const
Adds a bulk equation in all the domains of a given system (equation is assumed to be second order)
void add_eq_full(System_of_eqs &syst, const char *eq, int nused=-1, Array< int > **pused=0x0) const
Adds a bulk equation in all the domains of a given system (equation is assumed to be 0th order)
void add_eq_mode(System_of_eqs &syst, const char *f, int domtarget, int itarget, int jtarget)
Adds an equation saying that one coefficient of a field is zero in a given domain.
void add_eq_int_volume(System_of_eqs &syst, const char *eq)
Adds an equation being a volume integral in the whole computational space.
void add_eq_int_inf(System_of_eqs &syst, const char *eq)
Adds an equation being a surface integral at infinity.
void add_matching(System_of_eqs &syst, const char *eq, int nused=-1, Array< int > **pused=0x0) const
Adds a matching condition, at all the interface present in a given system.
virtual void save(FILE *) const
Saving function.
Definition: space_polar.cpp:62
Space_polar(int ttype, const Point &cr, const Dim_array &nbr, const Array< double > &bounds)
Standard constructor.
Definition: space_polar.cpp:27
void add_eq_int_inner(System_of_eqs &syst, const char *eq)
Adds an equation being a surface integral in the innermost shell.
void add_eq_ori(System_of_eqs &syst, const char *eq)
Adds an equation being the value of some field at the origin.
The Space class is an ensemble of domains describing the whole space of the computation.
Definition: space.hpp:1362
Class used to describe and solve a system of equations.
Tensor handling.
Definition: tensor.hpp:149
Class for storing the basis of decompositions of a field and its values on both the configuration and...
Definition: val_domain.hpp:69