KADATH
polar_periodic.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_PERIODIC_HPP_
21 #define __POLAR_PERIODIC_HPP_
22 
23 #include "space.hpp"
24 
25 #define TO_PI 0
26 namespace Kadath {
27 
54 
55  private:
56  double alpha ;
57 
58  double ome ;
59  mutable Term_eq* ome_term_eq ;
60 
65  int type_time ;
66  double maxt ;
67 
68  public:
77  Domain_polar_periodic_nucleus (int num, int ttype, double radius, double ome, const Dim_array& nbr) ;
84  Domain_polar_periodic_nucleus (int num, FILE* ff) ;
85 
87  virtual void save (FILE*) const ;
88 
89  private:
90  virtual void do_absol () const ;
91 
92  private:
93  virtual void set_cheb_base(Base_spectral&) const ;
94  virtual void set_legendre_base(Base_spectral&) const ;
95  virtual void set_anti_cheb_base(Base_spectral&) const ;
96  virtual void set_anti_legendre_base(Base_spectral&) const ;
97  virtual void set_cheb_base_with_m(Base_spectral&, int m) const ;
98  virtual void set_legendre_base_with_m(Base_spectral&, int m) const ;
99  virtual void set_anti_cheb_base_with_m(Base_spectral&, int m) const ;
100  virtual void set_anti_legendre_base_with_m(Base_spectral&, int m) const ;
101  virtual void set_cheb_base_r_spher(Base_spectral&) const ;
102  virtual void set_cheb_base_t_spher(Base_spectral&) const ;
103  virtual void set_cheb_base_p_spher(Base_spectral&) const ;
104  virtual void set_legendre_base_r_spher(Base_spectral&) const ;
105  virtual void set_legendre_base_t_spher(Base_spectral&) const ;
106  virtual void set_legendre_base_p_spher(Base_spectral&) const ;
107 
108 
109  virtual void do_coloc () ;
110  virtual void do_radius () const ;
111 
112  public:
113  virtual bool is_in(const Point&xx, double prec=1e-13) const ;
114  virtual const Point absol_to_num(const Point&) const;
115  virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const ;
116  virtual Base_spectral mult (const Base_spectral&, const Base_spectral&) const ;
120  double get_ome() const {return ome ;} ;
121 
122  public:
123  virtual Val_domain mult_cos_theta (const Val_domain&) const ;
124  virtual Val_domain mult_sin_theta (const Val_domain&) const ;
125  virtual Val_domain div_sin_theta (const Val_domain&) const ;
126  virtual Val_domain div_x (const Val_domain&) const ;
127  virtual Val_domain mult_r (const Val_domain&) const ;
128  virtual Val_domain div_r (const Val_domain&) const ;
129  virtual Val_domain der_r (const Val_domain&) const ;
130  virtual Val_domain dt (const Val_domain&) const ;
131  virtual Val_domain srdr (const Val_domain&) const ;
132  virtual Val_domain dtime (const Val_domain&) const ;
133  virtual Val_domain mult_cos_time (const Val_domain&) const ;
134  virtual Val_domain mult_sin_time (const Val_domain&) const ;
135 
136  virtual Term_eq partial_spher (const Term_eq&) const ;
137  virtual Term_eq derive_flat_spher (int, char, const Term_eq&, const Metric*) const ;
138 
139  virtual double val_boundary (int, const Val_domain&, const Index&) const ;
140  virtual void find_other_dom (int, int, int&, int&) const ;
141  virtual Val_domain der_normal (const Val_domain&, int) const ;
142 
143 
144  // virtual int nbr_unknowns_from_adapted() const ;
145  // virtual void vars_to_terms() const ;
146  // virtual void affecte_coef (int&, int, bool&) const ;
147  // virtual void xx_to_vars_from_adapted (double, const Array<double>&, int&) const ;
148  // virtual void xx_to_ders_from_adapted (const Array<double>&, int&) const ;
149  // virtual void update_term_eq (Term_eq*) const ;
150  // virtual void update_variable (double, const Scalar&, Scalar&) const ;
151  // virtual void update_constante (double, const Scalar&, Scalar&) const ;
152  // virtual void update_mapping(double) ;
153 
157  // void update() const ;
158 
159  virtual Term_eq dtime_term_eq (const Term_eq& so) const ;
160  virtual Term_eq ddtime_term_eq (const Term_eq& so) const ;
161 
162  virtual int nbr_unknowns (const Tensor&, int) const ;
170  int nbr_unknowns_val_domain (const Val_domain& so, int llim) const ;
171 
172  virtual Array<int> nbr_conditions (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
181  int nbr_conditions_val_domain (const Val_domain& so, int llim, int order) const ;
182  virtual Array<int> nbr_conditions_boundary (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
190  int nbr_conditions_val_domain_boundary (const Val_domain& eq) const ;
191  virtual void export_tau (const Tensor&, int, int, Array<double>&, int&, const Array<int>&,int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
202  void export_tau_val_domain (const Val_domain& eq, int llim, int order, Array<double>& res, int& pos_res, int ncond) const ;
203  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 ;
213  void export_tau_val_domain_boundary (const Val_domain& eq, int bound, Array<double>& res, int& pos_res, int ncond) const ;
214  virtual void affecte_tau (Tensor&, int, const Array<double>&, int&) const ;
223  void affecte_tau_val_domain (Val_domain& so, int llim, const Array<double>& cf, int& pos_cf) const ;
224  virtual void affecte_tau_one_coef (Tensor&, int, int, int&) const ;
233  void affecte_tau_one_coef_val_domain (Val_domain& so, int llim, int cc, int& pos_cf) const ;
234 
235 
236 public:
237  virtual ostream& print (ostream& o) const ;
238 } ;
239 
266 
267  private:
268  double alpha ;
269  double beta ;
270 
271 
272  double ome ;
273  mutable Term_eq* ome_term_eq ;
274 
279  int type_time ;
280  double maxt ;
281 
282  public:
292  Domain_polar_periodic_shell (int num, int ttype, double rin, double rout, double ome, const Dim_array& nbr) ;
299  Domain_polar_periodic_shell (int num, FILE* ff) ;
300 
301  virtual ~Domain_polar_periodic_shell() ;
302  virtual void save (FILE*) const ;
303 
304  private:
305  virtual void do_absol () const ;
306 
307  private:
308  virtual void set_cheb_base(Base_spectral&) const ;
309  virtual void set_legendre_base(Base_spectral&) const ;
310  virtual void set_anti_cheb_base(Base_spectral&) const ;
311  virtual void set_anti_legendre_base(Base_spectral&) const ;
312  virtual void set_cheb_base_with_m(Base_spectral&, int m) const ;
313  virtual void set_legendre_base_with_m(Base_spectral&, int m) const ;
314  virtual void set_anti_cheb_base_with_m(Base_spectral&, int m) const ;
315  virtual void set_anti_legendre_base_with_m(Base_spectral&, int m) const ;
316  virtual void set_cheb_base_r_spher(Base_spectral&) const ;
317  virtual void set_cheb_base_t_spher(Base_spectral&) const ;
318  virtual void set_cheb_base_p_spher(Base_spectral&) const ;
319  virtual void set_legendre_base_r_spher(Base_spectral&) const ;
320  virtual void set_legendre_base_t_spher(Base_spectral&) const ;
321  virtual void set_legendre_base_p_spher(Base_spectral&) const ;
322 
323  virtual void do_coloc () ;
324  virtual void do_radius () const ;
325 
326  public:
327  virtual bool is_in(const Point&xx, double prec=1e-13) const ;
328  virtual const Point absol_to_num(const Point&) const;
329  virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const ;
330  virtual Base_spectral mult (const Base_spectral&, const Base_spectral&) const ;
331 
332  public:
333  virtual Val_domain mult_cos_theta (const Val_domain&) const ;
334  virtual Val_domain mult_sin_theta (const Val_domain&) const ;
335  virtual Val_domain div_sin_theta (const Val_domain&) const ;
336  virtual Val_domain mult_r (const Val_domain&) const ;
337  virtual Val_domain div_r (const Val_domain&) const ;
338  virtual Val_domain der_r (const Val_domain&) const ;
339  virtual Val_domain dt (const Val_domain&) const ;
340  virtual Val_domain div_1mrsL (const Val_domain&) const ;
341  virtual Val_domain dtime (const Val_domain&) const ;
342  virtual Val_domain mult_cos_time (const Val_domain&) const ;
343  virtual Val_domain mult_sin_time (const Val_domain&) const ;
344 
345  virtual Term_eq partial_spher (const Term_eq&) const ;
346  virtual Term_eq derive_flat_spher (int, char, const Term_eq&, const Metric*) const ;
347 
348  virtual double val_boundary (int, const Val_domain&, const Index&) const ;
349  virtual void find_other_dom (int, int, int&, int&) const ;
350  virtual Val_domain der_normal (const Val_domain&, int) const ;
351 
352  /* virtual int nbr_unknowns_from_adapted() const ;
353  virtual void vars_to_terms() const ;
354  virtual void affecte_coef (int&, int, bool&) const ;
355  virtual void xx_to_vars_from_adapted (double, const Array<double>&, int&) const ;
356  virtual void xx_to_ders_from_adapted (const Array<double>&, int&) const ;
357  virtual void update_term_eq (Term_eq*) const ;
358  virtual void update_variable (double, const Scalar&, Scalar&) const ;
359  virtual void update_constante (double, const Scalar&, Scalar&) const ;
360  virtual void update_mapping(double) ;
361 */
365  // void update() const ;
366 
367  virtual Term_eq dtime_term_eq (const Term_eq& so) const ;
368  virtual Term_eq ddtime_term_eq (const Term_eq& so) const ;
369 
370 
371  virtual int nbr_unknowns (const Tensor&, int) const ;
378  int nbr_unknowns_val_domain (const Val_domain& so) const ;
379 
380  virtual Array<int> nbr_conditions (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
388  int nbr_conditions_val_domain (const Val_domain& so, int order) const ;
389  virtual Array<int> nbr_conditions_boundary (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
397  int nbr_conditions_val_domain_boundary (const Val_domain& eq) const ;
398  virtual void export_tau (const Tensor&, int, int, Array<double>&, int&, const Array<int>&,int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
408  void export_tau_val_domain (const Val_domain& eq, int order, Array<double>& res, int& pos_res, int ncond) const ;
409  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 ;
419  void export_tau_val_domain_boundary (const Val_domain& eq, int bound, Array<double>& res, int& pos_res, int ncond) const ;
420  virtual void affecte_tau (Tensor&, int, const Array<double>&, int&) const ;
428  void affecte_tau_val_domain (Val_domain& so, const Array<double>& cf, int& pos_cf) const ;
429  virtual void affecte_tau_one_coef (Tensor&, int, int, int&) const ;
437  void affecte_tau_one_coef_val_domain (Val_domain& so, int cc, int& pos_cf) const ;
438 
439 
440 public:
441  virtual ostream& print (ostream& o) const ;
442 } ;
443 
444 
449 class Space_polar_periodic : public Space {
450  public:
458  Space_polar_periodic (int ttype, double omega, const Dim_array& nbr, const Array<double>& bounds) ;
459  Space_polar_periodic (FILE*) ;
460  virtual ~Space_polar_periodic() ;
461  virtual void save(FILE*) const ;
462 
463  /*virtual int nbr_unknowns_from_variable_domains() const ;
464  virtual void affecte_coef_to_variable_domains(int& , int, Array<int>&) const ;
465  virtual void xx_to_ders_variable_domains(const Array<double>&, int&) const ;
466  virtual void xx_to_vars_variable_domains(System_of_eqs*, const Array<double>&, int&) const ;
467 */
469  double get_omega() const ;
470 } ;
471 
472 }
473 #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 spherical nucleus with a symmetry in .
virtual Term_eq ddtime_term_eq(const Term_eq &so) const
Second time derivative of a Term_eq.
virtual Val_domain div_r(const Val_domain &) const
Division by .
virtual void do_radius() const
Computes the generalized radius.
virtual void set_cheb_base_with_m(Base_spectral &, int m) const
Gives the standard base using Chebyshev polynomials.
virtual void set_legendre_base_t_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
void export_tau_val_domain(const Val_domain &eq, int llim, int order, Array< double > &res, int &pos_res, int ncond) const
Exports a residual equation in the bulk.
virtual Val_domain mult_cos_theta(const Val_domain &) const
Multiplication by .
virtual Term_eq derive_flat_spher(int, char, const Term_eq &, const Metric *) const
Computes the flat derivative of a Term_eq, in spherical orthonormal coordinates.
virtual Val_domain dtime(const Val_domain &) const
Computes the time derivative of a field.
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 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...
Domain_polar_periodic_nucleus(int num, int ttype, double radius, double ome, const Dim_array &nbr)
Standard constructor :
virtual double val_boundary(int, const Val_domain &, const Index &) const
Computes the value of a field at a boundary.
virtual void affecte_tau(Tensor &, int, const Array< double > &, int &) const
Affects some coefficients to a Tensor.
virtual void set_legendre_base_with_m(Base_spectral &, int m) const
Gives the stnadard base using Legendre polynomials.
int nbr_conditions_val_domain_boundary(const Val_domain &eq) const
Computes number of discretized equations associated with a given equation on a boundary.
void export_tau_val_domain_boundary(const Val_domain &eq, 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 mult_cos_time(const Val_domain &) const
Multiplication by .
virtual Val_domain div_x(const Val_domain &) const
Division by .
virtual Val_domain dt(const Val_domain &) const
Compute the derivative with respect to of a scalar field.
virtual Term_eq dtime_term_eq(const Term_eq &so) const
Updates the quantities that depend on the frequency.
virtual Val_domain div_sin_theta(const Val_domain &) const
Division by .
virtual Base_spectral mult(const Base_spectral &, const Base_spectral &) const
Method for the multiplication of two Base_spectral.
int nbr_unknowns_val_domain(const Val_domain &so, int llim) const
Computes the number of true unknowns of a Val_domain.
virtual Val_domain der_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
int type_time
Gives the type of time periodicity.
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 void set_cheb_base_p_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector.
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
virtual void set_cheb_base(Base_spectral &) const
Gives the standard base for Chebyshev polynomials.
virtual Val_domain mult_sin_theta(const Val_domain &) const
Multiplication by .
virtual Val_domain srdr(const Val_domain &) const
Compute the of a scalar field .
virtual void set_legendre_base_r_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the radial component of a vector.
Term_eq * ome_term_eq
Pointer on the Term_eq version of the pulsation.
virtual void affecte_tau_one_coef(Tensor &, int, int, int &) const
Sets at most one coefficient of a Tensor to 1.
virtual void find_other_dom(int, int, int &, int &) const
Gives the informations corresponding the a touching neighboring domain.
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 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_coloc()
Computes the colocation points.
virtual void do_absol() const
Computes the absolute coordinates.
virtual void set_legendre_base(Base_spectral &) const
Gives the standard base for Legendre polynomials.
double get_ome() const
Returns omega.
virtual Val_domain mult_r(const Val_domain &) const
Multiplication by .
virtual Val_domain der_r(const Val_domain &) const
Compute the radial derivative of a scalar field.
virtual const Point absol_to_num(const Point &) const
Computes the numerical coordinates from the physical ones.
virtual Val_domain mult_sin_time(const Val_domain &) const
Multiplication by .
virtual void save(FILE *) const
Saving function.
virtual void set_legendre_base_p_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
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_cheb_base_r_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the radial component of a vector.
virtual void set_anti_legendre_base(Base_spectral &) const
Gives the base using Legendre polynomials, for functions antisymetric with respect to .
virtual Term_eq partial_spher(const Term_eq &) const
Computes the part of the gradient containing the partial derivative of the field, in spherical orthon...
virtual void set_cheb_base_t_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector.
virtual void set_anti_cheb_base(Base_spectral &) const
Gives the base using Chebyshev polynomials, for functions antisymetric with respect to .
double maxt
Upper bound of which is or .
void affecte_tau_one_coef_val_domain(Val_domain &so, int llim, int cc, int &pos_cf) const
Sets at most one coefficient of a Val_domain to 1.
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 .
int nbr_conditions_val_domain(const Val_domain &so, int llim, int order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual int nbr_unknowns(const Tensor &, int) const
Computes the number of true unknowns of a Tensor, in a given domain.
void affecte_tau_val_domain(Val_domain &so, int llim, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
double alpha
Relates the numerical radius to the physical one.
Class for a spherical shell with a symmetry in .
void affecte_tau_val_domain(Val_domain &so, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
virtual Val_domain dt(const Val_domain &) const
Compute the derivative with respect to of a scalar field.
Domain_polar_periodic_shell(int num, int ttype, double rin, double rout, double ome, const Dim_array &nbr)
Standard constructor :
virtual Val_domain div_1mrsL(const Val_domain &) const
Division by .
virtual void do_absol() const
Computes the absolute coordinates.
virtual Val_domain div_r(const Val_domain &) const
Division by .
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 set_legendre_base_with_m(Base_spectral &, int m) const
Gives the stnadard base using Legendre polynomials.
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
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 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 Term_eq partial_spher(const Term_eq &) const
Computes the part of the gradient containing the partial derivative of the field, in spherical orthon...
int nbr_conditions_val_domain(const Val_domain &so, int order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual double val_boundary(int, const Val_domain &, const Index &) const
Computes the value of a field at a boundary.
virtual Val_domain mult_r(const Val_domain &) const
Multiplication by .
virtual Term_eq dtime_term_eq(const Term_eq &so) const
Updates the quantities that depend on the frequency.
void export_tau_val_domain(const Val_domain &eq, int order, Array< double > &res, int &pos_res, int ncond) const
Exports a residual equation in the bulk.
int nbr_unknowns_val_domain(const Val_domain &so) const
Computes the number of true unknowns of a Val_domain.
Term_eq * ome_term_eq
Pointer on the Term_eq version of the pulsation.
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 mult_cos_time(const Val_domain &) const
Multiplication by .
virtual Val_domain mult_sin_theta(const Val_domain &) const
Multiplication by .
virtual Term_eq derive_flat_spher(int, char, const Term_eq &, const Metric *) const
Computes the flat derivative of a Term_eq, in spherical orthonormal coordinates.
virtual void save(FILE *) const
Saving function.
virtual Val_domain dtime(const Val_domain &) const
Computes the time derivative of a field.
int type_time
Gives the type of time periodicity.
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_anti_legendre_base(Base_spectral &) const
Gives the base using Legendre polynomials, for functions antisymetric with respect to .
double alpha
Relates the numerical radius to the physical one.
void affecte_tau_one_coef_val_domain(Val_domain &so, int cc, int &pos_cf) const
Sets at most one coefficient of a Val_domain to 1.
virtual void set_anti_cheb_base(Base_spectral &) const
Gives the base using Chebyshev polynomials, for functions antisymetric with respect to .
virtual void set_cheb_base(Base_spectral &) const
Gives the standard base for Chebyshev polynomials.
double beta
Relates the numerical radius to the physical one.
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_with_m(Base_spectral &, int m) const
Gives the standard base using Chebyshev polynomials.
virtual void set_cheb_base_p_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector.
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 Val_domain der_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
virtual void do_coloc()
Computes the colocation points.
virtual void set_legendre_base(Base_spectral &) const
Gives the standard base for Legendre polynomials.
virtual int nbr_unknowns(const Tensor &, int) const
Computes the number of true unknowns of a Tensor, in a given domain.
double maxt
Upper bound of which is or .
virtual void set_legendre_base_p_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
virtual Val_domain mult_cos_theta(const Val_domain &) const
Multiplication by .
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_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector.
void export_tau_val_domain_boundary(const Val_domain &eq, 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 void set_legendre_base_r_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the radial component of a vector.
virtual void find_other_dom(int, int, int &, int &) const
Gives the informations corresponding the a touching neighboring domain.
virtual void do_radius() const
Computes the generalized radius.
virtual void affecte_tau_one_coef(Tensor &, int, int, int &) const
Sets at most one coefficient of a Tensor to 1.
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 Term_eq ddtime_term_eq(const Term_eq &so) const
Second time derivative of a Term_eq.
int nbr_conditions_val_domain_boundary(const Val_domain &eq) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual Val_domain div_sin_theta(const Val_domain &) const
Division by .
virtual const Point absol_to_num(const Point &) const
Computes the numerical coordinates from the physical ones.
virtual void set_legendre_base_t_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
virtual Val_domain mult_sin_time(const Val_domain &) const
Multiplication by .
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
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.
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
Purely abstract class for metric handling.
Definition: metric.hpp:39
The class Point is used to store the coordinates of a point.
Definition: point.hpp:30
The Space_polar_periodic class fills the space with one polar nucleus and several polar shells,...
double get_omega() const
Returns omega.
Space_polar_periodic(int ttype, double omega, const Dim_array &nbr, const Array< double > &bounds)
Standard constructor.
virtual void save(FILE *) const
Saving function.
The Space class is an ensemble of domains describing the whole space of the computation.
Definition: space.hpp:1362
Tensor handling.
Definition: tensor.hpp:149
This class is intended to describe the manage objects appearing in the equations.
Definition: term_eq.hpp:62
Class for storing the basis of decompositions of a field and its values on both the configuration and...
Definition: val_domain.hpp:69