KADATH
oned.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 __ONED_HPP_
21 #define __ONED_HPP_
22 
23 #include "space.hpp"
24 #include "val_domain.hpp"
25 
26 namespace Kadath {
40 class Domain_oned_ori : public Domain {
41 
42  private:
43  double alpha ;
44 
45  public:
53  Domain_oned_ori (int num, int ttype, double radius, const Dim_array& nbr) ;
54  Domain_oned_ori (const Domain_oned_ori& so) ;
60  Domain_oned_ori (int num, FILE* ff) ;
61 
62  virtual ~Domain_oned_ori() ;
63  virtual void save (FILE*) const ;
64 
65  private:
66  virtual void do_absol () const ;
67 
68  private:
69  virtual void set_cheb_base(Base_spectral&) const ;
70  virtual void set_legendre_base(Base_spectral&) const ;
71  virtual void set_cheb_base_odd(Base_spectral&) const ;
72  virtual void set_legendre_base_odd(Base_spectral&) const ;
73  virtual void do_coloc () ;
74  virtual int give_place_var (char*) const ;
75  virtual void do_radius () const ;
76 
77  public:
78  virtual bool is_in(const Point&xx, double prec=1e-13) const ;
79  virtual const Point absol_to_num(const Point&) const;
80  virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const ;
81  virtual Base_spectral mult (const Base_spectral&, const Base_spectral&) const ;
82 
83  virtual double val_boundary (int, const Val_domain&, const Index&) const ;
84  virtual void find_other_dom (int, int, int&, int&) const ;
85  virtual Val_domain der_normal (const Val_domain&, int) const ;
86  virtual Val_domain der_partial_var (const Val_domain&, int) const ;
87  virtual Val_domain div_r (const Val_domain&) const ;
88  virtual int nbr_unknowns (const Tensor&, int) const ;
95  int nbr_unknowns_val_domain (const Val_domain& so) const ;
96  virtual Array<int> nbr_conditions (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
97 
105  int nbr_conditions_val_domain (const Val_domain& so, int order) const ;
106  virtual Array<int> nbr_conditions_boundary (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
107 
115  int nbr_conditions_val_domain_boundary (const Val_domain& eq) const ;
116  virtual void export_tau (const Tensor&, int, int, Array<double>&, int&, const Array<int>&,int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
117 
127  void export_tau_val_domain (const Val_domain& eq, int order, Array<double>& res, int& pos_res, int ncond) const ;
128  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 ;
129 
139  void export_tau_val_domain_boundary (const Val_domain& eq, int bound, Array<double>& res, int& pos_res, int ncond) const ;
140  virtual void affecte_tau (Tensor&, int, const Array<double>&, int&) const ;
141 
149  void affecte_tau_val_domain (Val_domain& so, const Array<double>& cf, int& pos_cf) const ;
150  virtual void affecte_tau_one_coef (Tensor&, int, int, int&) const ;
158  void affecte_tau_one_coef_val_domain (Val_domain& so, int cc, int& pos_cf) const ;
159 
160  virtual double integrale(const Val_domain&) const ;
161 
162 public:
163  virtual ostream& print (ostream& o) const ;
164 } ;
165 
179 class Domain_oned_qcq : public Domain {
180 
181  private:
182  double alpha ;
183  double beta ;
184 
185  public:
194  Domain_oned_qcq (int num, int ttype, double x_int, double x_ext, const Dim_array& nbr) ;
195  Domain_oned_qcq (const Domain_oned_qcq& so) ;
201  Domain_oned_qcq (int num, FILE* ff) ;
202 
203  virtual ~Domain_oned_qcq() ;
204  virtual void save (FILE*) const ;
205 
206  private:
207  virtual void do_absol () const ;
208 
209  private :
210  virtual void set_cheb_base(Base_spectral&) const ;
211  virtual void set_legendre_base(Base_spectral&) const ;
212  virtual void set_cheb_base_odd(Base_spectral&) const ;
213  virtual void set_legendre_base_odd(Base_spectral&) const ;
214 
215  virtual void do_coloc () ;
216  virtual int give_place_var (char*) const ;
217  virtual void do_radius () const ;
218 
219  public:
220  virtual bool is_in(const Point& xx, double prec=1e-13) const ;
221  virtual const Point absol_to_num(const Point&) const;
222  virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const ;
223 
224  virtual Base_spectral mult (const Base_spectral&, const Base_spectral&) const ;
225 
226  public:
227  virtual void find_other_dom (int, int, int&, int&) const ;
228  virtual Val_domain der_normal (const Val_domain&, int) const ;
229  virtual double val_boundary (int, const Val_domain&, const Index&) const ;
230  virtual Val_domain der_partial_var (const Val_domain&, int) const ;
231 
232  virtual int nbr_unknowns (const Tensor&, int) const ;
239  int nbr_unknowns_val_domain (const Val_domain& so) const ;
240  virtual Array<int> nbr_conditions (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
248  int nbr_conditions_val_domain (const Val_domain& so, int order) const ;
249  virtual Array<int> nbr_conditions_boundary (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
257  int nbr_conditions_val_domain_boundary (const Val_domain& eq) const ;
258  virtual void export_tau (const Tensor&, int, int, Array<double>&, int&, const Array<int>&,int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
268  void export_tau_val_domain (const Val_domain& eq, int order, Array<double>& res, int& pos_res, int ncond) const ;
269  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 ;
279  void export_tau_val_domain_boundary (const Val_domain& eq, int bound, Array<double>& res, int& pos_res, int ncond) const ;
280  virtual void affecte_tau (Tensor&, int, const Array<double>&, int&) const ;
288  void affecte_tau_val_domain (Val_domain& so, const Array<double>& cf, int& pos_cf) const ;
289  virtual void affecte_tau_one_coef (Tensor&, int, int, int&) const ;
297  void affecte_tau_one_coef_val_domain (Val_domain& so, int cc, int& pos_cf) const ;
298 
299  virtual Val_domain div_xp1 (const Val_domain&) const ;
300  virtual double integrale(const Val_domain&) const ;
301 
302 public:
303  virtual ostream& print (ostream& o) const ;
304 } ;
305 
319 class Domain_oned_inf : public Domain {
320 
321  private:
322  double alpha ;
323 
324  public:
332  Domain_oned_inf (int num, int ttype, double radius, const Dim_array& nbr) ;
333  Domain_oned_inf (const Domain_oned_inf& so) ;
339  Domain_oned_inf (int num, FILE* ff) ;
340 
341  virtual ~Domain_oned_inf() ;
342  virtual void save(FILE*) const ;
343 
344  private:
345  virtual void set_cheb_base(Base_spectral&) const ;
346  virtual void set_legendre_base(Base_spectral&) const ;
347  virtual void set_cheb_base_odd(Base_spectral&) const ;
348  virtual void set_legendre_base_odd(Base_spectral&) const ;
349 
350  virtual void do_coloc() ;
351 
352  virtual void do_absol () const ;
353  virtual int give_place_var (char*) const ;
354  virtual void do_radius () const ;
355 
356  public:
360  double get_alpha() const {return alpha ;} ;
361 
362  virtual bool is_in(const Point& xx, double prec=1e-13) const ;
363  virtual const Point absol_to_num(const Point&) const;
364  virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const ;
365  virtual Base_spectral mult (const Base_spectral&, const Base_spectral&) const ;
366 
367  public:
368 
369  virtual void set_val_inf (Val_domain& so, double xx) const ;
370 
371  virtual void find_other_dom (int, int, int&, int&) const ;
372  virtual Val_domain der_normal (const Val_domain&, int) const ;
373  virtual double val_boundary (int, const Val_domain&, const Index&) const ;
374  virtual Val_domain der_partial_var (const Val_domain&, int) const ;
375  virtual Val_domain mult_xm1 (const Val_domain&) const ;
376  virtual Val_domain div_xm1 (const Val_domain&) const ;
377  virtual Val_domain mult_x (const Val_domain&) const ;
378  virtual Val_domain mult_r (const Val_domain& so) const {return mult_x(so) ;} ;
379 
380  virtual int nbr_unknowns (const Tensor&, int) const ;
387  int nbr_unknowns_val_domain (const Val_domain& so) const ;
388  virtual Array<int> nbr_conditions (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
396  int nbr_conditions_val_domain (const Val_domain& so, int order) const ;
397  virtual Array<int> nbr_conditions_boundary (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
405  int nbr_conditions_val_domain_boundary (const Val_domain& eq) const ;
406  virtual void export_tau (const Tensor&, int, int, Array<double>&, int&, const Array<int>&,int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
416  void export_tau_val_domain (const Val_domain& eq, int order, Array<double>& res, int& pos_res, int ncond) const ;
417  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 ;
427  void export_tau_val_domain_boundary (const Val_domain& eq, int bound, Array<double>& res, int& pos_res, int ncond) const ;
428  virtual void affecte_tau (Tensor&, int, const Array<double>&, int&) const ;
436  void affecte_tau_val_domain (Val_domain& so, const Array<double>& cf, int& pos_cf) const ;
437  virtual void affecte_tau_one_coef (Tensor&, int, int, int&) const ;
445  void affecte_tau_one_coef_val_domain (Val_domain& so, int cc, int& pos_cf) const ;
446 
447  virtual double integrale(const Val_domain&) const ;
448 
449 public:
450  virtual ostream& print (ostream& o) const ;
451 } ;
452 
453 
458 class Space_oned : public Space {
459  public:
466  Space_oned (int ttype, const Dim_array& nbr, const Array<double>& bounds) ;
467  Space_oned (FILE*) ;
468  virtual ~Space_oned() ;
469  virtual void save(FILE*) const ;
470 
476  void add_eq_ori (System_of_eqs& syst, const char* eq) ;
477 } ;
478 }
479 #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 1-dimensional compactified spherical domain.
Definition: oned.hpp:319
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 Val_domain mult_xm1(const Val_domain &) const
Multiplication by .
virtual int nbr_unknowns(const Tensor &, int) const
Computes the number of true unknowns of a Tensor, in a given domain.
virtual void set_legendre_base_odd(Base_spectral &) const
Gives the base using odd Legendre polynomials$.
virtual int give_place_var(char *) const
Translates a name of a coordinate into its corresponding numerical name.
virtual void save(FILE *) const
Saving function.
Domain_oned_inf(int num, int ttype, double radius, const Dim_array &nbr)
Standard constructor :
virtual const Point absol_to_num(const Point &) const
Computes the numerical coordinates from the physical ones.
double get_alpha() const
Returns the of the mapping.
Definition: oned.hpp:360
virtual void set_cheb_base_odd(Base_spectral &) const
Gives the base using odd Chebyshev 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 Val_domain div_xm1(const Val_domain &) const
Division by .
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
double alpha
Relates the numerical radius to the physical one.
Definition: oned.hpp:322
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 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_partial_var(const Val_domain &, int) const
Partial derivative with respect to a coordinate.
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 do_coloc()
Computes the colocation points.
virtual Val_domain der_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
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 mult_r(const Val_domain &so) const
Multiplication by .
Definition: oned.hpp:378
virtual void find_other_dom(int, int, int &, int &) const
Gives the informations corresponding the a touching neighboring domain.
virtual void do_absol() const
Computes the absolute coordinates.
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.
virtual Val_domain mult_x(const Val_domain &) const
Multiplication by .
virtual void affecte_tau(Tensor &, int, const Array< double > &, int &) const
Affects some coefficients to a Tensor.
virtual void do_radius() const
Computes the generalized radius.
virtual void set_legendre_base(Base_spectral &) const
Gives the standard base for Legendre polynomials.
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.
virtual Base_spectral mult(const Base_spectral &, const Base_spectral &) const
Method for the multiplication of two Base_spectral.
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.
int nbr_conditions_val_domain_boundary(const Val_domain &eq) const
Computes number of discretized equations associated with a given equation on a boundary.
int nbr_unknowns_val_domain(const Val_domain &so) const
Computes the number of true unknowns of a Val_domain.
virtual void affecte_tau_one_coef(Tensor &, int, int, int &) const
Sets at most one coefficient of a Tensor to 1.
virtual void set_cheb_base(Base_spectral &) const
Gives the standard base for Chebyshev polynomials.
virtual double integrale(const Val_domain &) const
Volume integral.
void affecte_tau_val_domain(Val_domain &so, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
virtual double val_boundary(int, const Val_domain &, const Index &) const
Computes the value of a field at a boundary.
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.
Class for a 1-dimensional spherical domain containing the origin.
Definition: oned.hpp:40
virtual void do_absol() const
Computes the absolute coordinates.
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 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 const Point absol_to_num(const Point &) const
Computes the numerical coordinates from the physical ones.
virtual Base_spectral mult(const Base_spectral &, const Base_spectral &) const
Method for the multiplication of two Base_spectral.
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.
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, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
virtual void save(FILE *) const
Saving function.
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 bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
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 void do_radius() const
Computes the generalized radius.
virtual void do_coloc()
Computes the colocation points.
virtual void set_legendre_base(Base_spectral &) const
Gives the standard base for Legendre polynomials.
virtual void set_cheb_base_odd(Base_spectral &) const
Gives the base using odd Chebyshev polynomials$.
double alpha
Relates the numerical radius to the physical one.
Definition: oned.hpp:43
int nbr_unknowns_val_domain(const Val_domain &so) const
Computes the number of true unknowns of a Val_domain.
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
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 Val_domain der_partial_var(const Val_domain &, int) const
Partial derivative with respect to a coordinate.
virtual void affecte_tau(Tensor &, int, const Array< double > &, int &) const
Affects some coefficients to a Tensor.
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 double integrale(const Val_domain &) const
Volume integral.
virtual Val_domain div_r(const Val_domain &) const
Division by .
virtual double val_boundary(int, const Val_domain &, const Index &) const
Computes the value of a field at a boundary.
virtual void affecte_tau_one_coef(Tensor &, int, int, int &) const
Sets at most one coefficient of a Tensor to 1.
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 set_cheb_base(Base_spectral &) const
Gives the standard base for Chebyshev polynomials.
Domain_oned_ori(int num, int ttype, double radius, const Dim_array &nbr)
Standard constructor :
virtual void set_legendre_base_odd(Base_spectral &) const
Gives the base using odd Legendre polynomials$.
virtual void find_other_dom(int, int, int &, int &) const
Gives the informations corresponding the a touching neighboring domain.
virtual ~Domain_oned_ori()
Destructor.
virtual Val_domain der_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
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 ...
Class for a 1-dimensional spherical domain bounded between two raii.
Definition: oned.hpp:179
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 der_partial_var(const Val_domain &, int) const
Partial derivative with respect to a coordinate.
virtual void do_coloc()
Computes the colocation points.
virtual double val_boundary(int, const Val_domain &, const Index &) const
Computes the value of a field at a boundary.
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.
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(Tensor &, int, const Array< double > &, int &) const
Affects some coefficients to a Tensor.
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 double integrale(const Val_domain &) const
Volume integral.
virtual void set_cheb_base_odd(Base_spectral &) const
Gives the base using odd Chebyshev polynomials$.
virtual void find_other_dom(int, int, int &, int &) const
Gives the informations corresponding the a touching neighboring domain.
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
virtual int give_place_var(char *) const
Translates a name of a coordinate into its corresponding numerical name.
virtual void set_cheb_base(Base_spectral &) const
Gives the standard base for Chebyshev polynomials.
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 const Point absol_to_num(const Point &) const
Computes the numerical coordinates from the physical ones.
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.
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 save(FILE *) const
Saving function.
virtual int nbr_unknowns(const Tensor &, int) const
Computes the number of true unknowns of a Tensor, in a given domain.
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.
Domain_oned_qcq(int num, int ttype, double x_int, double x_ext, const Dim_array &nbr)
Standard constructor :
virtual void set_legendre_base_odd(Base_spectral &) const
Gives the base using odd Legendre polynomials$.
int nbr_unknowns_val_domain(const Val_domain &so) const
Computes the number of true unknowns of a Val_domain.
void affecte_tau_val_domain(Val_domain &so, 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) const
Computes number of discretized equations associated with a given equation on a boundary.
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 Val_domain div_xp1(const Val_domain &) const
Division by .
double beta
Relates the numerical radius the physical one.
Definition: oned.hpp:183
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 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 ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
virtual Base_spectral mult(const Base_spectral &, const Base_spectral &) const
Method for the multiplication of two Base_spectral.
double alpha
Relates the numerical radius to the physical one.
Definition: oned.hpp:182
virtual void do_absol() const
Computes the absolute coordinates.
Abstract class that implements the fonctionnalities common to all the type of domains.
Definition: space.hpp:60
Val_domain * radius
The generalized radius.
Definition: space.hpp:78
Class that gives the position inside a multi-dimensional Array.
Definition: index.hpp:38
The class Point is used to store the coordinates of a point.
Definition: point.hpp:30
The Space_oned class fills the space with 1-dimensional domains (spherical ones), intended for spheri...
Definition: oned.hpp:458
void add_eq_ori(System_of_eqs &syst, const char *eq)
Adds an equation being the value of some field at the origin.
Definition: oned_add_eq.cpp:25
virtual void save(FILE *) const
Saving function.
Definition: space_oned.cpp:60
Space_oned(int ttype, const Dim_array &nbr, const Array< double > &bounds)
Standard constructor.
Definition: space_oned.cpp:26
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