20 #ifndef __SPHERIC_HPP_
21 #define __SPHERIC_HPP_
32 class Domain_shell_log ;
33 class Domain_shell_surr ;
91 virtual void save (FILE*)
const ;
173 virtual bool is_in(
const Point&xx,
double prec=1e-13)
const ;
524 virtual ostream&
print (ostream& o)
const ;
582 virtual void save (FILE*)
const ;
587 virtual void do_cart ()
const ;
665 virtual bool is_in(
const Point& xx,
double prec=1e-13)
const ;
944 int type_exception,
const Val_domain& exception)
const ;
959 int type_exception,
const Val_domain& exception)
const ;
972 virtual ostream&
print (ostream& o)
const ;
1032 virtual void save(FILE*)
const ;
1037 virtual void do_cart ()
const ;
1120 virtual bool is_in(
const Point& xx,
double prec=1e-13)
const ;
1343 virtual ostream&
print (ostream& o)
const ;
1384 virtual void save(FILE*)
const ;
1530 virtual void save (FILE*)
const ;
1535 virtual void do_cart ()
const ;
1548 virtual bool is_in(
const Point& xx,
double prec=1e-13)
const ;
1589 virtual ostream&
print (ostream& o)
const ;
1620 virtual void save (FILE*)
const ;
1625 virtual void do_cart ()
const ;
1638 virtual bool is_in(
const Point& xx,
double prec=1e-13)
const ;
1679 virtual ostream&
print (ostream& o)
const ;
Class for storing the basis of decompositions of a field.
Class for storing the dimensions of an array.
Class for a spherical compactified domain and a symmetry with respect to the plane .
virtual Tensor change_basis_spher_to_cart(int dd, const Tensor &) const
Changes the tensorial basis from spherical to Cartesian in a given domain.
int nbr_unknowns_val_domain(const Val_domain &so, int mlim) const
Computes the number of true unknowns of a Val_domain.
virtual double val_boundary(int, const Val_domain &, const Index &) const
Computes the value of a field at a boundary.
virtual void set_legendre_r_base(Base_spectral &) const
Gives the base using odd Legendre polynomials$ for the radius.
void affecte_tau_one_coef_val_domain(Val_domain &so, int mlim, int cc, int &pos_cf) const
Sets at most one coefficient of a Val_domain to 1.
void affecte_tau_val_domain_mquant(Val_domain &so, int mquant, 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 Term_eq derive_flat_cart(int, char, const Term_eq &, const Metric *) const
Computes the flat derivative of a Term_eq, in Cartesian coordinates.
virtual Val_domain div_xm1(const Val_domain &) const
Division by .
virtual void set_cheb_base(Base_spectral &so) const
Sets the base to the standard one for Chebyshev polynomials.
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 Val_domain mult_sin_phi(const Val_domain &) const
Multiplication by .
int nbr_unknowns_val_domain_mquant(const Val_domain &so, int mquant) const
Computes the number of true unknowns of a Val_domain.
virtual void set_val_inf(Val_domain &, double) const
Sets the value at infinity of a Val_domain : not implemented for this type of Domain.
virtual Val_domain der_r_rtwo(const Val_domain &) const
Compute the radial derivative multiplied by of a scalar field.
virtual const Point absol_to_num_bound(const Point &, int) const
Computes the numerical coordinates from the physical ones for a point lying on a boundary.
virtual void save(FILE *) const
Saving function.
virtual void do_cart() const
Computes the Cartesian coordinates.
virtual void set_cheb_base_r_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the radial component of a vector in the MTZ setting.
virtual void set_legendre_base_t_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector in the MTZ context.
virtual Val_domain mult_cos_theta(const Val_domain &) const
Multiplication by .
double alpha
Relates the numerical to the physical radii.
virtual Val_domain div_cos_theta(const Val_domain &) const
Division by .
virtual void find_other_dom(int, int, int &, int &) const
Gives the informations corresponding the a touching neighboring domain.
virtual Val_domain div_sin_theta(const Val_domain &) const
Division by .
virtual Val_domain mult_r(const Val_domain &) const
Multiplication by .
Point center
Absolute coordinates of the center.
int nbr_conditions_val_domain(const Val_domain &eq, int mlim, int order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual Val_domain dt(const Val_domain &) const
Compute the derivative with respect to of a scalar field.
Domain_compact(int num, int ttype, double r_int, const Point &cr, const Dim_array &nbr)
Standard constructor :
virtual void do_cart_surr() const
Computes the Cartesian coordinates over the radius.
virtual const Point absol_to_num(const Point &xxx) const
Computes the numerical coordinates from the physical ones.
virtual Val_domain mult_sin_theta(const Val_domain &) const
Multiplication by .
virtual Val_domain der_partial_var(const Val_domain &, int) const
Partial derivative with respect to a coordinate.
virtual void set_cheb_base_p_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector.
double get_alpha() const
Returns the center ;.
int nbr_conditions_val_domain_boundary(const Val_domain &eq, int mlim) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual void set_legendre_base_r_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the radial component of a vector.
int nbr_conditions_val_domain_boundary_mquant(const Val_domain &eq, int mquant) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual void set_legendre_base_r_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the radial component of a vector in the MTZ context.
void export_tau_val_domain(const Val_domain &eq, int mlim, 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 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 set_legendre_base_p_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector in the MTZ context.
virtual Val_domain mult_cos_phi(const Val_domain &) const
Multiplication by .
void export_tau_val_domain_mquant(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 void affecte_tau(Tensor &, int, const Array< double > &, int &) const
Affects some coefficients to a Tensor.
void affecte_tau_val_domain(Val_domain &so, int mlim, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
virtual Tensor change_basis_cart_to_spher(int dd, const Tensor &) const
Changes the tensorial basis from Cartsian to spherical in a given domain.
void export_tau_val_domain_boundary_mquant(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 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_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector in the MTZ setting.
int nbr_conditions_val_domain_mquant(const Val_domain &eq, int mquant, int order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual void set_cheb_base_t_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector in the MTZ setting.
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 &so) const
Sets the base to the standard one for Legendre polynomials for functions antisymetric in The bases a...
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.
virtual int give_place_var(char *) const
Translates a name of a coordinate into its corresponding numerical name.
virtual void set_anti_cheb_base(Base_spectral &so) const
Sets the base to the standard one for Chebyshev polynomials for functions antisymetric in The bases ...
virtual int nbr_unknowns(const Tensor &, int) const
Computes the number of true unknowns of a Tensor, in a given domain.
virtual double integ_volume(const Val_domain &so) const
Volume integral.
virtual void set_legendre_base(Base_spectral &so) const
Sets the base to the standard one for Legendre polynomials.
virtual Val_domain mult_xm1(const Val_domain &) const
Multiplication by .
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside 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.
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 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 Point get_center() const
Returns the center.
virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const
Computes the derivative with respect to the absolute Cartesian coordinates from the derivative with r...
virtual void set_legendre_base_t_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
void affecte_tau_one_coef_val_domain_mquant(Val_domain &so, int mquant, int cc, int &pos_cf) const
Sets at most one coefficient of a Val_domain to 1.
virtual void do_which_points_boundary(int, const Base_spectral &, Index **, int) const
Lists all the indices corresponding to true collocation points on 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 set_cheb_r_base(Base_spectral &) const
Gives the base using odd Chebyshev polynomials$ for the radius.
virtual double integ(const Val_domain &, int) const
Surface integral on a given boundary.
virtual void do_radius() const
Computes the generalized radius.
virtual Val_domain div_r(const Val_domain &) const
Division by .
void export_tau_val_domain_boundary(const Val_domain &eq, int mlim, 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 der_r(const Val_domain &) const
Compute the radial derivative of a scalar field.
virtual int nbr_points_boundary(int, const Base_spectral &) const
Computes the number of relevant collocation points on a boundary.
virtual void set_cheb_base_with_m(Base_spectral &so, int m) const
Gives the standard base using Chebyshev polynomials.
virtual void set_legendre_base_p_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
Class for a spherical domain containing the origin and a symmetry with respect to the plane .
void export_tau_val_domain_boundary_vr(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_radius() const
Computes the generalized radius.
void export_tau_val_domain(const Val_domain &eq, int mlim, int llim, int order, Array< double > &res, int &pos_res, int ncond) const
Exports a residual equation in the bulk.
void export_tau_val_domain_boundary(const Val_domain &eq, int mlim, 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 ...
Point center
Absolute coordinates of the center.
virtual Val_domain der_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
virtual void affecte_tau(Tensor &, int, const Array< double > &, int &) const
Affects some coefficients to a Tensor.
virtual void do_cart() const
Computes the Cartesian coordinates.
virtual void set_cheb_r_base(Base_spectral &) const
Gives the base using odd Chebyshev polynomials$ for the radius.
virtual Val_domain div_x(const Val_domain &) const
Division by .
virtual Val_domain der_r(const Val_domain &) const
Compute the radial derivative of a scalar field.
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
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_vr(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.
int nbr_conditions_val_domain_boundary_vt(const Val_domain &eq) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual void set_legendre_base_p_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
int nbr_unknowns_val_domain_vp(const Val_domain &so) const
Computes the number of true unknowns of a Val_domain.
void export_tau_val_domain_vt(const Val_domain &eq, int order, Array< double > &res, int &pos_res, int ncond) const
Exports a residual equation in the bulk.
void affecte_tau_val_domain_vp(Val_domain &so, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
virtual Val_domain div_r(const Val_domain &) const
Division by .
int nbr_unknowns_val_domain_vt(const Val_domain &so) const
Computes the number of true unknowns of a Val_domain.
void affecte_tau_val_domain_vt(Val_domain &so, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
virtual double integ(const Val_domain &so, int bound) const
Surface integral on a given boundary.
virtual Tensor change_basis_cart_to_spher(int dd, const Tensor &) const
Changes the tensorial basis from Cartsian to spherical in a given domain.
virtual Tensor change_basis_spher_to_cart(int dd, const Tensor &) const
Changes the tensorial basis from spherical to Cartesian in a given domain.
virtual void set_legendre_base(Base_spectral &so) const
Sets the base to the standard one for Legendre polynomials.
virtual Val_domain srdr(const Val_domain &) const
Compute the of a scalar field .
virtual const Point absol_to_num_bound(const Point &, int) const
Computes the numerical coordinates from the physical ones for a point lying on a boundary.
virtual void do_cart_surr() const
Computes the Cartesian coordinates over the radius.
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 &so) const
Sets the base to the standard one for Chebyshev polynomials for functions antisymetric in The bases ...
void affecte_tau_one_coef_val_domain(Val_domain &so, int mlim, int llim, 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, int mlim) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual void set_cheb_base(Base_spectral &so) const
Sets the base to the standard one for Chebyshev polynomials.
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_vp(Val_domain &so, int cc, int &pos_cf) const
Sets at most one coefficient of a Val_domain to 1.
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_vt(const Val_domain &so, int order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
void affecte_tau_val_domain_vr(Val_domain &so, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
virtual void set_legendre_base_t_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector in the MTZ context.
int nbr_conditions_val_domain_boundary_vr(const Val_domain &eq) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual void set_cheb_base_with_m(Base_spectral &so, int m) const
Gives the standard base using Chebyshev polynomials.
virtual void set_cheb_base_r_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the radial component of a vector.
virtual Val_domain div_sin_theta(const Val_domain &) const
Division by .
virtual Val_domain ddp(const Val_domain &) const
Compute the second derivative with respect to of a scalar field.
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_vp(const Val_domain &eq, int order, Array< double > &res, int &pos_res, int ncond) const
Exports a residual equation in the bulk.
virtual void set_legendre_base_r_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the radial component of a vector.
int nbr_conditions_val_domain(const Val_domain &so, int mlim, int llim, int order) 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 Val_domain mult_sin_theta(const Val_domain &) const
Multiplication by .
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.
Domain_nucleus(int num, int ttype, double radius, const Point &cr, const Dim_array &nbr)
Standard constructor :
virtual void save(FILE *) const
Saving function.
virtual void set_cheb_base_p_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector.
void export_tau_val_domain_boundary_vp(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 Point get_center() const
Returns the center.
virtual Val_domain mult_cos_theta(const Val_domain &) const
Multiplication by .
int nbr_conditions_val_domain_boundary_vp(const Val_domain &eq) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual double integ_volume(const Val_domain &so) const
Volume integral.
void export_tau_val_domain_boundary_vt(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_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the radial component of a vector in the MTZ context.
int nbr_conditions_val_domain_vr(const Val_domain &so, int order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual Term_eq derive_flat_cart(int, char, const Term_eq &, const Metric *) const
Computes the flat derivative of a Term_eq, in Cartesian coordinates.
virtual Val_domain div_cos_theta(const Val_domain &) const
Division by .
virtual const Point absol_to_num(const Point &xxx) const
Computes the numerical coordinates from the physical ones.
virtual Val_domain mult_sin_phi(const Val_domain &) const
Multiplication by .
virtual void set_legendre_r_base(Base_spectral &) const
Gives the base using odd Legendre polynomials$ for the radius.
virtual void do_absol() const
Computes the absolute coordinates.
virtual Val_domain mult_r(const Val_domain &) const
Multiplication by .
int nbr_unknowns_val_domain_vr(const Val_domain &so) const
Computes the number of true unknowns of a Val_domain.
double alpha
Relates the numerical to the physical radii.
virtual void do_coloc()
Computes the colocation points.
void affecte_tau_one_coef_val_domain_vr(Val_domain &so, int cc, int &pos_cf) const
Sets at most one coefficient of a Val_domain to 1.
virtual int nbr_points_boundary(int, const Base_spectral &) const
Computes the number of relevant collocation points on a boundary.
int nbr_unknowns_val_domain(const Val_domain &so, int mlim, int llim) const
Computes the number of true unknowns of a Val_domain.
virtual double get_rmax() const
Returns the maximum radius.
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 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_t_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector in the MTZ setting.
virtual Val_domain laplacian2(const Val_domain &, int) const
Computes the ordinary flat 2dè- Laplacian for a scalar field with an harmonic index m.
void affecte_tau_one_coef_val_domain_vt(Val_domain &so, int cc, int &pos_cf) const
Sets at most one coefficient of a Val_domain to 1.
virtual Val_domain mult_cos_phi(const Val_domain &) const
Multiplication by .
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_legendre_base_p_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector in the MTZ context.
void affecte_tau_val_domain(Val_domain &so, int mlim, int llim, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
int nbr_conditions_val_domain_vp(const Val_domain &so, int order) 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 Val_domain div_1mx2(const Val_domain &) const
Division by .
virtual void do_which_points_boundary(int, const Base_spectral &, Index **, int) const
Lists all the indices corresponding to true collocation points on a boundary.
virtual void set_cheb_base_p_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector in the MTZ setting.
virtual 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 set_cheb_base_r_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the radial component of a vector in the MTZ setting.
virtual void set_anti_legendre_base(Base_spectral &so) const
Sets the base to the standard one for Legendre polynomials for functions antisymetric in The bases a...
Class for a spherical shell and a symmetry with respect to the plane .
virtual Val_domain div_r(const Val_domain &) const
Division by .
virtual const Point absol_to_num_bound(const Point &, int) const
Computes the numerical coordinates from the physical ones for a point lying on a boundary.
virtual void 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 Val_domain der_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
virtual void do_radius() const
Computes the generalized radius.
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 &xxx) const
Computes the numerical coordinates from the physical ones.
virtual double get_rmin() const
Returns the minimum radius.
virtual double get_rmax() const
Returns the maximum radius.
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
virtual void save(FILE *) const
Saving function.
virtual void do_cart() const
Computes the Cartesian coordinates.
Class for a spherical shell and a symmetry with respect to the plane .
virtual void do_cart() const
Computes the Cartesian coordinates.
virtual const Point absol_to_num(const Point &xxx) const
Computes the numerical coordinates from the physical ones.
virtual const Point absol_to_num_bound(const Point &, int) const
Computes the numerical coordinates from the physical ones for a point lying on a boundary.
virtual Val_domain der_r(const Val_domain &) const
Compute the radial derivative of a scalar field.
virtual void do_absol() const
Computes the absolute coordinates.
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
virtual Val_domain mult_r(const Val_domain &) const
Multiplication by .
virtual double get_rmin() const
Returns the minimum radius.
virtual Val_domain div_r(const Val_domain &) const
Division by .
virtual void save(FILE *) const
Saving function.
virtual double get_rmax() const
Returns the maximum radius.
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_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
virtual void do_radius() const
Computes the generalized radius.
Class for a spherical shell and a symmetry with respect to the plane .
virtual Term_eq der_harmonics_sym(const Term_eq &, const Term_eq &, int, Term_eq(*f)(const Space &, int, int, const Term_eq &, const Param &), const Param &, const Array< double > &) const
Fit, spherical harmonic by spherical harmonic, for the radial derivative of a symmetric function.
virtual Val_domain div_xp1(const Val_domain &) const
Division by .
virtual Val_domain div_r(const Val_domain &) const
Division by .
void export_tau_val_domain_boundary_mquant(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 Term_eq derive_flat_cart(int, char, const Term_eq &, const Metric *) const
Computes the flat derivative of a Term_eq, in Cartesian coordinates.
virtual void affecte_tau(Tensor &, int, const Array< double > &, int &) const
Affects some coefficients to a Tensor.
virtual Val_domain mult_cos_theta(const Val_domain &) const
Multiplication by .
int nbr_unknowns_val_domain_mquant(const Val_domain &so, int mquant) const
Computes the number of true unknowns of a Val_domain.
virtual void set_legendre_base_t_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
void affecte_tau_one_coef_val_domain(Val_domain &so, int mlim, int cc, int &pos_cf) const
Sets at most one coefficient of a Val_domain to 1.
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_with_m(Base_spectral &so, int m) const
Gives the standard base using Chebyshev polynomials.
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
virtual Val_domain div_1mrsL(const Val_domain &) const
Division by .
virtual double integ_volume(const Val_domain &so) const
Volume integral.
virtual void affecte_tau_one_coef(Tensor &, int, int, int &) const
Sets at most one coefficient of a Tensor to 1.
int nbr_conditions_val_domain_boundary_mquant(const Val_domain &so, int mquant) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual double integ(const Val_domain &so, int bound) const
Surface integral on a given boundary.
double beta
Relates the numerical to the physical radii.
void export_tau_val_domain_boundary_exception(const Val_domain &eq, int mlim, int bound, Array< double > &res, int &pos_res, int ncond, const Param ¶m, int type_exception, const Val_domain &exception) const
Exports all the residual equations corresponding to one tensorial one on a given boundary,...
virtual void set_cheb_base_r_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the radial component of a vector in the MTZ setting.
virtual void set_cheb_base_p_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector in the MTZ setting.
virtual 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.
void export_tau_val_domain(const Val_domain &eq, int mlim, int order, Array< double > &res, int &pos_res, int ncond) const
Exports a residual equation in the bulk.
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 Val_domain ddp(const Val_domain &) const
Compute the second derivative with respect to of a scalar field.
virtual void set_cheb_base_t_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector in the MTZ setting.
virtual Val_domain der_partial_var(const Val_domain &, int) const
Partial derivative with respect to a coordinate.
virtual Val_domain der_r(const Val_domain &) const
Compute the radial derivative of a scalar field.
double alpha
Relates the numerical to the physical radii.
virtual double multipoles_asym(int, int, int, const Val_domain &, const Array< double > &) const
Extraction of a given multipole, at some boundary, for a anti-symmetric scalar function.
virtual void set_cheb_r_base(Base_spectral &) const
Gives the base using odd Chebyshev polynomials$ for the radius.
virtual Val_domain div_cos_theta(const Val_domain &) const
Division by .
virtual void do_coloc()
Computes the colocation points.
virtual void save(FILE *) const
Saving function.
virtual int give_place_var(char *) const
Translates a name of a coordinate into its corresponding numerical name.
virtual Term_eq der_multipoles_sym(int, int, int, const Term_eq &, const Array< double > &) const
Extraction of a given multipole, at some boundary, for the radial derivative a symmetric scalar funct...
virtual Term_eq harmonics_sym(const Term_eq &, const Term_eq &, int, Term_eq(*f)(const Space &, int, int, const Term_eq &, const Param &), const Param &, const Array< double > &) const
Fit, spherical harmonic by spherical harmonic, for a symmetric function.
Domain_shell(int num, int ttype, double r_int, double r_ext, const Point &cr, const Dim_array &nbr)
Standard constructor :
virtual void set_legendre_base_r_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the radial component of a vector.
virtual Val_domain div_xm1(const Val_domain &) const
Division by .
virtual void do_which_points_boundary(int, const Base_spectral &, Index **, int) const
Lists all the indices corresponding to true collocation points on a boundary.
void export_tau_val_domain_mquant(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 void set_cheb_base_p_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector.
virtual void set_legendre_base_t_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector in the MTZ context.
Val_domain fitschwarz(const Val_domain &, int) const
Fit some field with a decay (Val_domain version).
virtual Val_domain mult_cos_phi(const Val_domain &) const
Multiplication by .
virtual Term_eq der_harmonics_asym(const Term_eq &, const Term_eq &, int, Term_eq(*f)(const Space &, int, int, const Term_eq &, const Param &), const Param &, const Array< double > &) const
Fit, spherical harmonic by spherical harmonic, for the radial derivative of an anti-symmetric functio...
virtual Term_eq der_multipoles_asym(int, int, int, const Term_eq &, const Array< double > &) const
Extraction of a given multipole, at some boundary, for the radial derivative of an anti-symmetric sca...
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 set_legendre_base(Base_spectral &so) const
Sets the base to the standard one for Legendre polynomials.
virtual Point get_center() const
Returns the center.
virtual void set_cheb_base(Base_spectral &so) const
Sets the base to the standard one for Chebyshev polynomials.
virtual Term_eq der_radial_part_asym(const Space &, int, int, const Term_eq &, Term_eq(*f)(const Space &, int, int, const Term_eq &, const Param &), const Param &) const
Gives some radial fit for a given multipole, intended for the radial derivative of an anti-symmetric ...
int nbr_conditions_val_domain_mquant(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 set_anti_legendre_base(Base_spectral &so) const
Sets the base to the standard one for Legendre polynomials for functions antisymetric in The bases a...
virtual double multipoles_sym(int, int, int, const Val_domain &, const Array< double > &) const
Extraction of a given multipole, at some boundary, for a symmetric scalar function.
Term_eq bc_waves(const Term_eq &gamma, const Term_eq &omega) const
Gives an matching of the spatial metric, based on homogeneous solutions of outgoing waves.
virtual void set_cheb_base_r_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the radial component of a vector.
int nbr_conditions_val_domain_boundary(const Val_domain &so, int mlim) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual double get_rmin() const
Returns the minimum radius.
virtual Term_eq der_radial_part_sym(const Space &, int, int, const Term_eq &, Term_eq(*f)(const Space &, int, int, const Term_eq &, const Param ¶m), const Param &) const
Gives some radial fit for a given multipole, intended for the radial derivative of a symmetric scalar...
void affecte_tau_val_domain(Val_domain &so, int mlim, 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.
virtual Term_eq radial_part_asym(const Space &, int, int, const Term_eq &, Term_eq(*f)(const Space &, int, int, const Term_eq &, const Param &), const Param &) const
Gives some radial fit for a given multipole, intended for anti-symmetric scalar function.
virtual void set_cheb_base_t_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector.
virtual const Point absol_to_num(const Point &xxx) const
Computes the numerical coordinates from the physical ones.
int nbr_unknowns_val_domain(const Val_domain &so, int mlim) const
Computes the number of true unknowns of a Val_domain.
virtual Val_domain mult_sin_phi(const Val_domain &) const
Multiplication by .
virtual Term_eq harmonics_asym(const Term_eq &, const Term_eq &, int, Term_eq(*f)(const Space &, int, int, const Term_eq &, const Param &), const Param &, const Array< double > &) const
Fit, spherical harmonic by spherical harmonic, for an anti-symmetric function.
virtual void set_anti_cheb_base(Base_spectral &so) const
Sets the base to the standard one for Chebyshev polynomials for functions antisymetric in The bases ...
void export_tau_val_domain_boundary_exception_mquant(const Val_domain &eq, int mquant, int bound, Array< double > &res, int &pos_res, int ncond, const Param ¶m, int type_exception, const Val_domain &exception) const
Exports all the residual equations corresponding to one tensorial one on a given boundary,...
virtual void do_absol() const
Computes the absolute coordinates.
virtual double val_boundary(int, const Val_domain &, const Index &) const
Computes the value of a field at a boundary.
virtual Val_domain laplacian2(const Val_domain &, int) const
Computes the ordinary flat 2dè- Laplacian for a scalar field with an harmonic index m.
virtual Tensor change_basis_cart_to_spher(int dd, const Tensor &) const
Changes the tensorial basis from Cartsian to spherical in a given domain.
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
int nbr_conditions_val_domain(const Val_domain &so, int mlim, int order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual Val_domain mult_r(const Val_domain &) const
Multiplication by .
Point center
Absolute coordinates of the center.
virtual int nbr_points_boundary(int, const Base_spectral &) const
Computes the number of relevant collocation points on a boundary.
virtual void set_legendre_base_p_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const
Computes the derivative with respect to the absolute Cartesian coordinates from the derivative with r...
void export_tau_val_domain_boundary(const Val_domain &eq, int mlim, 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_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the radial component of a vector in the MTZ context.
virtual void set_legendre_r_base(Base_spectral &) const
Gives the base using odd Legendre polynomials$ for the radius.
virtual Val_domain mult_1mrsL(const Val_domain &) const
Multiplication by .
virtual Term_eq radial_part_sym(const Space &, int, int, const Term_eq &, Term_eq(*f)(const Space &, int, int, const Term_eq &, const Param &), const Param &) const
Gives some radial fit for a given multipole, intended for symmetric scalar function.
void affecte_tau_val_domain_mquant(Val_domain &so, int mquant, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
virtual double get_rmax() const
Returns the maximum radius.
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 const Point absol_to_num_bound(const Point &, int) const
Computes the numerical coordinates from the physical ones for a point lying on a boundary.
virtual void set_legendre_base_p_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector in the MTZ context.
virtual void 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 do_cart() const
Computes the Cartesian coordinates.
virtual Val_domain div_1mx2(const Val_domain &) const
Division by .
void affecte_tau_one_coef_val_domain_mquant(Val_domain &so, int mquant, int cc, int &pos_cf) const
Sets at most one coefficient of a Val_domain to 1.
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 Tensor change_basis_spher_to_cart(int dd, const Tensor &) const
Changes the tensorial basis from spherical to Cartesian in a given 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 Term_eq derive_flat_mtz(int, char, const Term_eq &, const Metric *) const
Computes the flat derivative of a Term_eq, in spherical coordinates where the constant radii sections...
virtual Val_domain mult_sin_theta(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 do_radius() const
Computes the generalized radius.
virtual Val_domain mult_xm1(const Val_domain &) const
Multiplication by .
virtual void do_cart_surr() const
Computes the Cartesian coordinates over the radius.
Abstract class that implements the fonctionnalities common to all the type of domains.
Val_domain * radius
The generalized radius.
Class that gives the position inside a multi-dimensional Array.
Purely abstract class for metric handling.
The class Point is used to store the coordinates of a point.
The class Scalar does not really implements scalars in the mathematical sense but rather tensorial co...
The Space_spheric class fills the space with one nucleus, several shells and a compactified domain,...
virtual Array< int > get_indices_matching_non_std(int, int) const
Gives the number of the other domains, touching a given boundary.
void add_eq(System_of_eqs &syst, const char *eq, int nused=-1, Array< int > **pused=0x0)
Adds a bulk equation in all the domains of a given system (equation is assumed to be second order)
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.
double int_inf(const Scalar &) const
Computes the surface integral at infinity.
void add_eq_ori(System_of_eqs &syst, const char *eq)
Adds an equation being the value of some field at the origin.
Space_spheric(int ttype, const Point &cr, const Dim_array &nbr, const Array< double > &bounds, bool withzec=true)
Standard constructor.
void add_eq_mode_mid(System_of_eqs &syst, const char *f, int itarget, int jtarget, int ktarget)
Adds an equation saying that one coefficient of a field is zero (in the first shell)
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_full(System_of_eqs &syst, const char *eq, int nused=-1, Array< int > **pused=0x0)
Adds a bulk equation in all the domains of a given system (equation is assumed to be 0th order)
void add_outer_bc(System_of_eqs &syst, const char *eq, int nused=-1, Array< int > **pused=0x0)
Sets a boundary condition at the outer radius of the compactified domain.
virtual void save(FILE *) const
Saving function.
void add_matching(System_of_eqs &syst, const char *eq, int nused=-1, Array< int > **pused=0x0)
Adds a matching condition, at all the interface present in a given system.
virtual ~Space_spheric()
Destructor
void add_inner_bc(System_of_eqs &syst, const char *eq, int nused=-1, Array< int > **pused=0x0)
Sets a boundary condition at the inner radius of the innermost shell.
void add_eq_one_side(System_of_eqs &syst, const char *eq, int nused=-1, Array< int > **pused=0x0)
Adds a bulk equation in all the domains of a given system (equation is assumed to be first order)
The Space class is an ensemble of domains describing the whole space of the computation.
Class used to describe and solve a system of equations.
This class is intended to describe the manage objects appearing in the equations.
Class for storing the basis of decompositions of a field and its values on both the configuration and...