KADATH
spheric_time.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 __SPHERIC_TIME_HPP_
21 #define __SPHERIC_TIME_HPP_
22 
23 #include "space.hpp"
24 
25 namespace Kadath {
26 
46 
47  private:
48  double alpha ;
49  double tmin;
50  double tmax ;
51 
52  public:
62  Domain_spheric_time_nucleus (int num, int ttype, double tmmin, double tmmax, double radius, const Dim_array& nbr) ;
69  Domain_spheric_time_nucleus (int num, FILE* ff) ;
70 
71  virtual ~Domain_spheric_time_nucleus() ;
72  virtual void save (FILE*) const ;
73 
74  private:
75  virtual void do_absol () const ;
76 
77  private:
78  virtual void set_cheb_base(Base_spectral&) const ;
79  virtual void set_legendre_base(Base_spectral&) const ;
80  virtual void do_coloc () ;
81  virtual void do_radius () const ;
82 
83  public:
84  virtual bool is_in(const Point&xx, double prec=1e-13) const ;
85  virtual const Point absol_to_num(const Point&) const;
86  virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const ;
87  virtual Base_spectral mult (const Base_spectral&, const Base_spectral&) const ;
91  double get_alpha () const {return alpha ;} ;
92 
93  public:
94  virtual Val_domain div_x (const Val_domain&) const ;
95  virtual Val_domain mult_r (const Val_domain&) const ;
96  virtual Val_domain div_r (const Val_domain&) const ;
97  virtual Val_domain laplacian (const Val_domain&, int) const ;
98  virtual Val_domain der_r (const Val_domain&) const ;
99  virtual Val_domain srdr (const Val_domain&) const ;
100  virtual Val_domain ddt (const Val_domain&) const ;
101  virtual Val_domain dt (const Val_domain&) const ;
102 
103  virtual double val_boundary (int, const Val_domain&, const Index&) const ;
104  virtual void find_other_dom (int, int, int&, int&) const ;
105  virtual Val_domain der_normal (const Val_domain&, int) const ;
106 
107  virtual int nbr_unknowns (const Tensor&, int) const ;
114  int nbr_unknowns_val_domain (const Val_domain& so) const ;
115  virtual Array<int> nbr_conditions_array (const Tensor&, int, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
123  int nbr_conditions_val_domain_array (const Val_domain& so, const Array<int>& order) const ;
124  virtual Array<int> nbr_conditions_boundary_array (const Tensor&, int, int, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
134  int nbr_conditions_val_domain_boundary_array (const Val_domain& eq, int bound, const Array<int>& order) const ;
135  virtual void export_tau_array (const Tensor&, int, const Array<int>&, Array<double>&, int&, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
145  void export_tau_val_domain_array (const Val_domain& eq, const Array<int>& order, Array<double>& res, int& pos_res, int ncond) const ;
146  virtual void export_tau_boundary_array (const Tensor&, int, int, const Array<int>&, Array<double>&, int&, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
157  void export_tau_val_domain_boundary_array (const Val_domain& eq, int bound, const Array<int>& order, Array<double>& res, int& pos_res, int ncond) const ;
158  virtual void affecte_tau (Tensor&, int, const Array<double>&, int&) const ;
166  void affecte_tau_val_domain (Val_domain& so, const Array<double>& cf, int& pos_cf) const ;
167  virtual void affecte_tau_one_coef (Tensor&, int, int, int&) const ;
175  void affecte_tau_one_coef_val_domain (Val_domain& so, int cc, int& pos_cf) const ;
176 
177 public:
178  virtual ostream& print (ostream& o) const ;
179 } ;
180 
200 
201  private:
202  double alpha ;
203  double beta ;
204  double tmin ;
205  double tmax ;
206 
207  public:
218  Domain_spheric_time_shell (int num, int ttype, double tmmin, double tmmax , double r_int, double r_ext, const Dim_array& nbr) ;
225  Domain_spheric_time_shell (int num, FILE* ff) ;
226 
227  virtual ~Domain_spheric_time_shell() ;
228  virtual void save (FILE*) const ;
229 
230  private:
231  virtual void do_absol () const ;
232 
233  private :
234  virtual void set_cheb_base(Base_spectral&) const ;
235  virtual void set_legendre_base(Base_spectral&) const ;
236  virtual void do_coloc () ;
237  virtual void do_radius () const ;
238 
239  public:
240  virtual bool is_in(const Point& xx, double prec=1e-13) const ;
241  virtual const Point absol_to_num(const Point&) const;
242  virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const ;
243 
244  virtual Base_spectral mult (const Base_spectral&, const Base_spectral&) const ;
248  double get_alpha () const {return alpha ;} ;
249 
250  public:
251  virtual Val_domain mult_r (const Val_domain&) const ;
252  virtual Val_domain div_r (const Val_domain&) const ;
253  virtual Val_domain laplacian (const Val_domain&, int) const ;
254  virtual Val_domain der_r (const Val_domain&) const ;
255  virtual Val_domain ddt (const Val_domain&) const ;
256  virtual Val_domain dt (const Val_domain&) const ;
257 
258  virtual void find_other_dom (int, int, int&, int&) const ;
259  virtual Val_domain der_normal (const Val_domain&, int) const ;
260  virtual double val_boundary (int, const Val_domain&, const Index&) const ;
261 
262  virtual int nbr_unknowns (const Tensor&, int) const ;
269  int nbr_unknowns_val_domain (const Val_domain& so) const ;
270  virtual Array<int> nbr_conditions_array (const Tensor&, int, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
278  int nbr_conditions_val_domain_array (const Val_domain& so, const Array<int>& order) const ;
279  virtual Array<int> nbr_conditions_boundary_array (const Tensor&, int, int, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
289  int nbr_conditions_val_domain_boundary_array (const Val_domain& eq, int bound, const Array<int>& order) const ;
290  virtual void export_tau_array (const Tensor&, int, const Array<int>&, Array<double>&, int&, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
300  void export_tau_val_domain_array (const Val_domain& eq, const Array<int>& order, Array<double>& res, int& pos_res, int ncond) const ;
301  virtual void export_tau_boundary_array (const Tensor&, int, int, const Array<int>&, Array<double>&, int&, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
312  void export_tau_val_domain_boundary_array (const Val_domain& eq, int bound, const Array<int>& order, Array<double>& res, int& pos_res, int ncond) const ;
313  virtual void affecte_tau (Tensor&, int, const Array<double>&, int&) const ;
321  void affecte_tau_val_domain (Val_domain& so, const Array<double>& cf, int& pos_cf) const ;
322  virtual void affecte_tau_one_coef (Tensor&, int, int, int&) const ;
330  void affecte_tau_one_coef_val_domain (Val_domain& so, int cc, int& pos_cf) const ;
331 
332 public:
333  virtual ostream& print (ostream& o) const ;
334 } ;
335 
355 
356  private:
357  double alpha ;
358  double tmin ;
359  double tmax ;
360 
361  public:
371  Domain_spheric_time_compact (int num, int ttype, double tmmin, double tmmax , double r_int, const Dim_array& nbr) ;
378  Domain_spheric_time_compact (int num, FILE* ff) ;
379 
380  virtual ~Domain_spheric_time_compact() ;
381  virtual void save (FILE*) const ;
382 
383  private:
384  virtual void do_absol () const ;
385 
386  private :
387  virtual void set_cheb_base(Base_spectral&) const ;
388  virtual void set_legendre_base(Base_spectral&) const ;
389  virtual void do_coloc () ;
390  virtual void do_radius () const ;
391 
392  public:
393  virtual bool is_in(const Point& xx, double prec=1e-13) const ;
394  virtual const Point absol_to_num(const Point&) const;
395  virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const ;
396 
397  virtual Base_spectral mult (const Base_spectral&, const Base_spectral&) const ;
401  double get_alpha () const {return alpha ;} ;
402 
403  public:
404  virtual Val_domain mult_xm1 (const Val_domain&) const ;
405  virtual Val_domain div_xm1 (const Val_domain&) const ;
406  virtual Val_domain mult_r (const Val_domain&) const ;
407  virtual Val_domain div_r (const Val_domain&) const ;
408  virtual Val_domain der_r (const Val_domain&) const ;
409  virtual Val_domain ddt (const Val_domain&) const ;
410  virtual Val_domain dt (const Val_domain&) const ;
411  virtual Val_domain laplacian (const Val_domain&, int) const ;
412 
413  virtual void find_other_dom (int, int, int&, int&) const ;
414  virtual Val_domain der_normal (const Val_domain&, int) const ;
415  virtual double val_boundary (int, const Val_domain&, const Index&) const ;
416 
417  virtual int nbr_unknowns (const Tensor&, int) const ;
424  int nbr_unknowns_val_domain (const Val_domain& so) const ;
425  virtual Array<int> nbr_conditions_array (const Tensor&, int, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
433  int nbr_conditions_val_domain_array (const Val_domain& so, const Array<int>& order) const ;
434  virtual Array<int> nbr_conditions_boundary_array (const Tensor&, int, int, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
444  int nbr_conditions_val_domain_boundary_array (const Val_domain& eq, int bound, const Array<int>& order) const ;
445  virtual void export_tau_array (const Tensor&, int, const Array<int>&, Array<double>&, int&, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
455  void export_tau_val_domain_array (const Val_domain& eq, const Array<int>& order, Array<double>& res, int& pos_res, int ncond) const ;
456  virtual void export_tau_boundary_array (const Tensor&, int, int, const Array<int>&, Array<double>&, int&, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
467  void export_tau_val_domain_boundary_array (const Val_domain& eq, int bound, const Array<int>& order, Array<double>& res, int& pos_res, int ncond) const ;
468  virtual void affecte_tau (Tensor&, int, const Array<double>&, int&) const ;
476  void affecte_tau_val_domain (Val_domain& so, const Array<double>& cf, int& pos_cf) const ;
477  virtual void affecte_tau_one_coef (Tensor&, int, int, int&) const ;
485  void affecte_tau_one_coef_val_domain (Val_domain& so, int cc, int& pos_cf) const ;
486 
487 public:
488  virtual ostream& print (ostream& o) const ;
489 } ;
490 
495 class Space_spheric_time : public Space {
496 
497  protected:
498  double tmin ;
499  double tmax ;
500  bool withcompact ;
501 
502  public:
512  Space_spheric_time (int ttype, const Dim_array& nbr, const Array<double>& bounds, double tmmin, double tmmax, bool wc = false) ;
513  Space_spheric_time (FILE*) ;
514  virtual ~Space_spheric_time() ;
515  virtual void save(FILE*) const ;
516 
520  double get_tmax() const {return tmax ;} ;
521 
525  double get_tmin() const {return tmin ;} ;
526 } ;
527 }
528 #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 compactified spherical domain and a symetry with respect to the plane .
virtual void affecte_tau(Tensor &, int, const Array< double > &, int &) const
Affects some coefficients to a Tensor.
virtual Val_domain der_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
virtual const Point absol_to_num(const Point &) const
Computes the numerical coordinates from the physical ones.
virtual Val_domain div_r(const Val_domain &) const
Division by .
double alpha
Relates the numerical to the physical radii.
virtual void do_absol() const
Computes the absolute coordinates.
virtual void set_cheb_base(Base_spectral &) const
Gives the standard base for Chebyshev polynomials.
Domain_spheric_time_compact(int num, int ttype, double tmmin, double tmmax, double r_int, const Dim_array &nbr)
Standard constructor :
virtual Val_domain div_xm1(const Val_domain &) const
Division by .
virtual void affecte_tau_one_coef(Tensor &, int, int, int &) const
Sets at most one coefficient of a Tensor to 1.
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
virtual int nbr_unknowns(const Tensor &, int) const
Computes the number of true unknowns of a Tensor, in a given domain.
int nbr_conditions_val_domain_boundary_array(const Val_domain &eq, int bound, const Array< int > &order) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual void export_tau_boundary_array(const Tensor &, int, int, const Array< int > &, Array< double > &, int &, const Array< int > &, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to one tensorial one on a given boundary It makes us...
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 do_coloc()
Computes the colocation points.
virtual void save(FILE *) const
Saving function.
void affecte_tau_val_domain(Val_domain &so, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
virtual Array< int > nbr_conditions_array(const Tensor &, int, const Array< 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 export_tau_array(const Tensor &, int, const Array< int > &, Array< double > &, int &, const Array< int > &, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to one tensorial one in the bulk.
void export_tau_val_domain_boundary_array(const Val_domain &eq, int bound, const Array< int > &order, 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 find_other_dom(int, int, int &, int &) const
Gives the informations corresponding the a touching neighboring domain.
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
void export_tau_val_domain_array(const Val_domain &eq, const Array< int > &order, Array< double > &res, int &pos_res, int ncond) const
Exports a residual equation in the bulk.
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 laplacian(const Val_domain &, int) const
Computes the ordinary flat Laplacian for a scalar field with an harmonic index m.
virtual void set_legendre_base(Base_spectral &) const
Gives the standard base for Legendre polynomials.
virtual Array< int > nbr_conditions_boundary_array(const Tensor &, int, int, const Array< 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 double val_boundary(int, const Val_domain &, const Index &) const
Computes the value of a field at a boundary.
virtual Val_domain mult_xm1(const Val_domain &) const
Multiplication by .
int nbr_unknowns_val_domain(const Val_domain &so) const
Computes the number of true unknowns of a Val_domain.
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 Val_domain ddt(const Val_domain &) const
Compute the second derivative with respect to of a scalar field.
virtual Val_domain der_r(const Val_domain &) const
Compute the radial derivative of a scalar field.
virtual Val_domain mult_r(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.
int nbr_conditions_val_domain_array(const Val_domain &so, const Array< int > &order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
Class for a 2-dimensional spherical domain containing the origin and a symetry with respect to the pl...
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 affecte_tau_one_coef(Tensor &, int, int, int &) const
Sets at most one coefficient of a Tensor to 1.
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
virtual const Point absol_to_num(const Point &) const
Computes the numerical coordinates from the physical ones.
Domain_spheric_time_nucleus(int num, int ttype, double tmmin, double tmmax, double radius, const Dim_array &nbr)
Standard constructor :
virtual Val_domain dt(const Val_domain &) const
Compute the derivative with respect to of a scalar field.
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
virtual Val_domain ddt(const Val_domain &) const
Compute the second derivative with respect to of a scalar field.
virtual void do_absol() const
Computes the absolute coordinates.
virtual void set_legendre_base(Base_spectral &) const
Gives the standard base for Legendre polynomials.
void export_tau_val_domain_boundary_array(const Val_domain &eq, int bound, const Array< int > &order, 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 int nbr_unknowns(const Tensor &, int) const
Computes the number of true unknowns of a Tensor, in a given domain.
virtual Val_domain div_x(const Val_domain &) const
Division by .
virtual void set_cheb_base(Base_spectral &) const
Gives the standard base for Chebyshev polynomials.
virtual Base_spectral mult(const Base_spectral &, const Base_spectral &) const
Method for the multiplication of two Base_spectral.
virtual void export_tau_array(const Tensor &, int, const Array< int > &, Array< double > &, int &, const Array< int > &, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to one tensorial one in the bulk.
virtual void affecte_tau(Tensor &, int, const Array< double > &, int &) const
Affects some coefficients to a Tensor.
virtual Array< int > nbr_conditions_array(const Tensor &, int, const Array< 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 find_other_dom(int, int, int &, int &) const
Gives the informations corresponding the a touching neighboring domain.
virtual Val_domain srdr(const Val_domain &) const
Compute the of a scalar field .
virtual Array< int > nbr_conditions_boundary_array(const Tensor &, int, int, const Array< 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 export_tau_boundary_array(const Tensor &, int, int, const Array< int > &, Array< double > &, int &, const Array< int > &, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to one tensorial one on a given boundary It makes us...
double alpha
Relates the numerical to the physical radii.
int nbr_conditions_val_domain_array(const Val_domain &so, const Array< int > &order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
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.
int nbr_unknowns_val_domain(const Val_domain &so) const
Computes the number of true unknowns of a Val_domain.
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...
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 save(FILE *) const
Saving function.
virtual Val_domain der_r(const Val_domain &) const
Compute the radial derivative of a scalar field.
virtual Val_domain mult_r(const Val_domain &) const
Multiplication by .
virtual Val_domain div_r(const Val_domain &) const
Division by .
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_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
int nbr_conditions_val_domain_boundary_array(const Val_domain &eq, int bound, const Array< int > &order) const
Computes number of discretized equations associated with a given equation on a boundary.
void export_tau_val_domain_array(const Val_domain &eq, const Array< int > &order, Array< double > &res, int &pos_res, int ncond) const
Exports a residual equation in the bulk.
virtual void do_radius() const
Computes the generalized radius.
Class for a 2-dimensional spherical domain bounded between two finite radii and a symetry with respec...
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_spheric_time_shell(int num, int ttype, double tmmin, double tmmax, double r_int, double r_ext, const Dim_array &nbr)
Standard constructor :
virtual void set_legendre_base(Base_spectral &) const
Gives the standard base for Legendre polynomials.
int nbr_conditions_val_domain_array(const Val_domain &so, const Array< int > &order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
double beta
Relates the numerical to the physical radii.
virtual int nbr_unknowns(const Tensor &, int) const
Computes the number of true unknowns of a Tensor, in a given domain.
int nbr_unknowns_val_domain(const Val_domain &so) const
Computes the number of true unknowns of a Val_domain.
virtual Val_domain ddt(const Val_domain &) const
Compute the second derivative with respect to of a scalar field.
virtual void do_radius() const
Computes the generalized radius.
virtual void export_tau_array(const Tensor &, int, const Array< int > &, Array< double > &, int &, const Array< int > &, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to one tensorial one in the bulk.
virtual Val_domain der_r(const Val_domain &) const
Compute the radial derivative of a scalar field.
void export_tau_val_domain_array(const Val_domain &eq, const Array< 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.
int nbr_conditions_val_domain_boundary_array(const Val_domain &eq, int bound, const Array< int > &order) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual Val_domain div_r(const Val_domain &) const
Division by .
virtual void do_coloc()
Computes the colocation points.
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 find_other_dom(int, int, int &, int &) const
Gives the informations corresponding the a touching neighboring domain.
virtual Val_domain der_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
virtual void export_tau_boundary_array(const Tensor &, int, int, const Array< int > &, Array< double > &, int &, const Array< int > &, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to one tensorial one on a given boundary It makes us...
virtual Val_domain mult_r(const Val_domain &) const
Multiplication by .
virtual Val_domain laplacian(const Val_domain &, int) const
Computes the ordinary flat Laplacian for a scalar field with an harmonic index m.
virtual const Point absol_to_num(const Point &) const
Computes the numerical coordinates from the physical ones.
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 Val_domain dt(const Val_domain &) const
Compute the derivative with respect to of a scalar field.
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 val_boundary(int, const Val_domain &, const Index &) const
Computes the value of a field at a boundary.
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
virtual void affecte_tau_one_coef(Tensor &, int, int, int &) const
Sets at most one coefficient of a Tensor to 1.
void export_tau_val_domain_boundary_array(const Val_domain &eq, int bound, const Array< int > &order, 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_absol() const
Computes the absolute coordinates.
virtual Base_spectral mult(const Base_spectral &, const Base_spectral &) const
Method for the multiplication of two Base_spectral.
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
virtual Array< int > nbr_conditions_boundary_array(const Tensor &, int, int, const Array< 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(Base_spectral &) const
Gives the standard base for Chebyshev polynomials.
virtual Array< int > nbr_conditions_array(const Tensor &, int, const Array< 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.
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 2-dimensioanl domains : spherically symmetric and time depe...
double get_tmax() const
Returns the final time.
bool withcompact
Indicator of the presence of a compactified Domain.
Space_spheric_time(int ttype, const Dim_array &nbr, const Array< double > &bounds, double tmmin, double tmmax, bool wc=false)
Standard constructor.
double get_tmin() const
Returns the initial time.
double tmax
Final time .
double tmin
Initial time .
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
Class for storing the basis of decompositions of a field and its values on both the configuration and...
Definition: val_domain.hpp:69