20 #ifndef __ADAPTED_HPP_
21 #define __ADAPTED_HPP_
24 #include "term_eq.hpp"
25 #include "spheric.hpp"
122 virtual void save (FILE*)
const ;
131 virtual void do_cart ()
const ;
181 virtual bool is_in(
const Point&xx,
double prec=1e-13)
const ;
339 virtual ostream&
print (ostream& o)
const ;
438 virtual void save (FILE*)
const ;
447 virtual void do_cart ()
const ;
494 virtual bool is_in(
const Point&xx,
double prec=1e-13)
const ;
648 virtual ostream&
print (ostream& o)
const ;
674 virtual void save(FILE*)
const ;
Class for storing the basis of decompositions of a field.
Class for storing the dimensions of an array.
Class for a spherical-like domain, having a symmetry with respect to the plane .
virtual Val_domain div_sin_theta(const Val_domain &) const
Division by .
Val_domain get_inner_radius() const
Returns the inner variable boundary.
Term_eq * der_rad_term_eq
Pointer on the Term_eq containing the .
virtual void set_legendre_r_base(Base_spectral &) const
Gives the base using odd Legendre polynomials$ for the radius.
Term_eq flat_grad_spher(const Term_eq &) const
Computes the flat gradient of a field, in orthonormal spherical coordinates.
virtual ~Domain_shell_inner_adapted()
Destructor.
virtual int nbr_unknowns(const Tensor &, int) const
Computes the number of true unknowns of a Tensor, in a given domain.
virtual void affecte_coef(int &, int, bool &) const
The variation of the functions describing the shape of the Domain are affected from the unknowns of t...
int nbr_unknowns_val_domain(const Val_domain &so, int mlim) const
Computes the number of true unknowns of a Val_domain.
virtual Base_spectral mult(const Base_spectral &, const Base_spectral &) const
Method for the multiplication of two Base_spectral.
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_legendre_base_t_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
Term_eq derive_t(const Term_eq &so) const
Computes .
virtual double val_boundary(int, const Val_domain &, const Index &) const
Computes the value of a field at a boundary.
virtual Tensor change_basis_cart_to_spher(int, const Tensor &) const
Changes the tensorial basis from Cartsian to spherical in a given domain.
void update() const
Updates all the quantities that depend on the inner radius (like the normal vectors).
Term_eq derive_p(const Term_eq &so) const
Computes .
virtual const Point absol_to_num(const Point &) const
Computes the numerical coordinates from the physical ones.
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 int give_place_var(char *) const
Translates a name of a coordinate into its corresponding numerical name.
virtual void update_mapping(const Val_domain &)
Updates the variables parts of the Domain.
virtual Term_eq connection_spher(const Term_eq &) const
Computes the part of the gradient involving the connections, in spherical orthonormal coordinates.
virtual void set_cheb_base_r_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the radial component of a vector.
virtual Term_eq integ_volume_term_eq(const Term_eq &) const
Volume integral of a Term_eq.
virtual void update_term_eq(Term_eq *) const
Update the value of a field, after the shape of the Domain has been changed by the system.
virtual Term_eq partial_spher(const Term_eq &) const
Computes the part of the gradient containing the partial derivative of the field, in spherical orthon...
virtual void set_legendre_base(Base_spectral &) const
Gives the standard base for Legendre polynomials.
virtual Val_domain div_cos_theta(const Val_domain &) const
Division by .
virtual void set_cheb_base_p_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector.
virtual Val_domain ddp(const Val_domain &) const
Compute the second derivative with respect to of a scalar field.
virtual int nbr_unknowns_from_adapted() const
Gives the number of unknowns coming from the variable shape of the domain.
Term_eq * rad_term_eq
Pointer on the Term_eq containing the radius.
virtual void do_coloc()
Computes the colocation points.
virtual void update_variable(const Val_domain &, const Scalar &, Scalar &) const
Update the value of a scalar, after the shape of the Domain has been changed by the system.
double outer_radius
The outer radius .
virtual Term_eq dr_term_eq(const Term_eq &) const
Radial derivative of a Term_eq.
virtual void xx_to_vars_from_adapted(Val_domain &, const Array< double > &, int &) const
Computes the new boundary of a Domain from a set of values.
void set_mapping(const Val_domain &so) const
Affects the inner radius.
virtual void affecte_tau(Tensor &, int, const Array< double > &, int &) const
Affects some coefficients to a Tensor.
virtual int nbr_points_boundary(int, const Base_spectral &) const
Computes the number of relevant collocation points on a boundary.
virtual Term_eq der_normal_term_eq(const Term_eq &, int) const
Returns the normal derivative of a Term_eq.
virtual Term_eq lap_term_eq(const Term_eq &, int) const
Returns the flat Laplacian of Term_eq, for a given harmonic.
virtual void set_cheb_r_base(Base_spectral &) const
Gives the base using odd Chebyshev polynomials$ for the radius.
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 set_cheb_base_t_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector.
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
virtual void set_cheb_base(Base_spectral &) const
Gives the standard base for Chebyshev polynomials.
Term_eq * normal_spher
Pointer on the Term_eq containing the normal vector to the inner boundary, in orthonormal spherical c...
void do_normal_cart() const
Computes the normal wrt the inner boundary, in Cartesian coordinates.
virtual void affecte_tau_one_coef(Tensor &, int, int, int &) const
Sets at most one coefficient of a Tensor to 1.
virtual void xx_to_ders_from_adapted(const Array< double > &, int &) const
Affects the derivative part of variable a Domain from a set of values.
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 do_which_points_boundary(int, const Base_spectral &, Index **, int) const
Lists all the indices corresponding to true collocation points on a boundary.
virtual Val_domain laplacian2(const Val_domain &, int) const
Computes the ordinary flat 2dè- Laplacian for a scalar field with an harmonic index m.
virtual void 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_r(const Val_domain &) const
Division by .
virtual Point get_center() const
Returns the center.
virtual Val_domain der_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
virtual Val_domain der_r(const Val_domain &) const
Compute the radial derivative of a scalar field.
virtual void set_anti_legendre_base(Base_spectral &) const
Gives the base using Legendre polynomials, for functions antisymetric with respect to .
Val_domain * inner_radius
Pointer on the inner boundary , as a Val_domain.
const Space & sp
The corresponding Space ; required for updating fields whene the mapping changes.
Term_eq * dt_rad_term_eq
Pointer on the Term_eq containing the .
virtual Term_eq integ_term_eq(const Term_eq &, int) const
Surface integral of a Term_eq.
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_sin_theta(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.
Term_eq * normal_cart
Pointer on the Term_eq containing the normal vector to the inner boundary, in Cartesian coordinates.
virtual Val_domain mult_cos_phi(const Val_domain &) const
Multiplication by .
void del_deriv() override
Destroys the derivated members (like coloc, cart and radius), when changing the type of colocation po...
void do_normal_spher() const
Computes the normal wrt the inner boundary, in orthonormal spherical coordinates.
Term_eq * inner_radius_term_eq
Pointer on the inner boundary , as a Term_eq.
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 Term_eq * give_normal(int, int) const
Returns the vector normal to a surface.
virtual void set_anti_cheb_base(Base_spectral &) const
Gives the base using Chebyshev polynomials, for functions antisymetric with respect to .
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 double integ_volume(const Val_domain &) const
Volume integral.
virtual void vars_to_terms() const
The Term_eq describing the variable shape of the Domain are updated.
virtual Val_domain mult_sin_phi(const Val_domain &) const
Multiplication 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 set_legendre_base_p_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
virtual Term_eq partial_cart(const Term_eq &) const
Computes the part of the gradient containing the partial derivative of the field, in Cartesian coordi...
Domain_shell_inner_adapted(const Space &sp, int num, int ttype, double rin, double rout, const Point &cr, const Dim_array &nbr)
Constructor :
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 void update_constante(const Val_domain &, const Scalar &, Scalar &) const
Update the value of a scalar, after the shape of the Domain has been changed by the system.
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
virtual Tensor change_basis_spher_to_cart(int, const Tensor &) const
Changes the tensorial basis from spherical to Cartesian in a given domain.
virtual void do_cart_surr() const
Computes the Cartesian coordinates over the radius.
virtual void do_cart() const
Computes the Cartesian coordinates.
virtual Term_eq mult_r_term_eq(const Term_eq &) const
Multiplication by of a Term_eq.
Point center
Absolute coordinates of the center.
virtual Val_domain mult_cos_theta(const Val_domain &) const
Multiplication by .
virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const
Computes the derivative with respect to the absolute Cartesian coordinates from the derivative with r...
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 void save(FILE *) const
Saving function.
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_radius() const
Computes the generalized radius.
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(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 ...
Term_eq derive_r(const Term_eq &so) const
Computes .
Term_eq * dp_rad_term_eq
Pointer on the Term_eq containing the .
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.
Class for a spherical-like domain, having a symmetry with respect to the plane .
virtual Term_eq lap_term_eq(const Term_eq &, int) const
Returns the flat Laplacian of Term_eq, for a given harmonic.
virtual Term_eq der_normal_term_eq(const Term_eq &, int) const
Returns the normal derivative of a Term_eq.
virtual void set_legendre_base_p_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
virtual Val_domain mult_sin_phi(const Val_domain &) const
Multiplication by .
void set_mapping(const Val_domain &so) const
Affects the outer radius.
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_shell_outer_adapted(const Space &sp, int num, int ttype, double rin, double rout, const Point &cr, const Dim_array &nbr)
Constructor :
virtual Val_domain div_cos_theta(const Val_domain &) const
Division by .
virtual Term_eq integ_volume_term_eq(const Term_eq &) const
Volume integral of a Term_eq.
Val_domain get_outer_radius() const
Returns the outer variable boundary.
virtual void affecte_coef(int &, int, bool &) const
The variation of the functions describing the shape of the Domain are affected from the unknowns of t...
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_absol() const
Computes the absolute coordinates.
virtual int nbr_unknowns(const Tensor &, int) const
Computes the number of true unknowns of a Tensor, in a given domain.
virtual void do_cart_surr() const
Computes the Cartesian coordinates over the radius.
virtual Term_eq dr_term_eq(const Term_eq &) const
Radial derivative of a Term_eq.
virtual void set_cheb_base(Base_spectral &) const
Gives the standard base for Chebyshev polynomials.
Term_eq * outer_radius_term_eq
Pointer on the outer boundary , as a Term_eq.
virtual void do_coloc()
Computes the colocation points.
virtual Base_spectral mult(const Base_spectral &, const Base_spectral &) const
Method for the multiplication of two Base_spectral.
virtual void affecte_tau(Tensor &, int, const Array< double > &, int &) const
Affects some coefficients to a Tensor.
virtual void set_legendre_base_r_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the radial component of a vector.
virtual double val_boundary(int, const Val_domain &, const Index &) const
Computes the value of a field at a boundary.
Term_eq * dp_rad_term_eq
Pointer on the Term_eq containing the .
void affecte_tau_one_coef_val_domain(Val_domain &so, int mlim, int cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
virtual Val_domain mult_cos_phi(const Val_domain &) const
Multiplication by .
Point center
Absolute coordinates of the center.
virtual void do_radius() const
Computes the generalized radius.
virtual Term_eq partial_spher(const Term_eq &) const
Computes the part of the gradient containing the partial derivative of the field, in spherical orthon...
virtual int give_place_var(char *) const
Translates a name of a coordinate into its corresponding numerical name.
virtual void save(FILE *) const
Saving function.
Term_eq flat_grad_spher(const Term_eq &) const
Computes the flat gradient of a field, in orthonormal spherical coordinates.
void do_normal_cart() const
Computes the normal wrt the inner boundary, in Cartesian coordinates.
virtual Term_eq partial_cart(const Term_eq &) const
Computes the part of the gradient containing the partial derivative of the field, in Cartesian coordi...
virtual void set_cheb_base_p_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector.
Term_eq derive_t(const Term_eq &so) const
Computes .
virtual void set_cheb_r_base(Base_spectral &) const
Gives the base using odd Chebyshev polynomials$ for the radius.
void update() const
Updates all the quantities that depend on the inner radius (like the normal vectors).
virtual void set_legendre_base_t_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
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.
Term_eq derive_r(const Term_eq &so) const
Computes .
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.
Term_eq * normal_spher
Pointer on the Term_eq containing the normal vector to the outer boundary, in orthonormal spherical c...
void del_deriv() override
Destroys the derivated members (like coloc, cart and radius), when changing the type of colocation po...
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 mlim, int order, Array< double > &res, int &pos_res, int ncond) const
Exports a residual equation in the bulk.
virtual Tensor change_basis_cart_to_spher(int dd, const Tensor &) const
Changes the tensorial basis from Cartsian to spherical in a given domain.
Term_eq * normal_cart
Pointer on the Term_eq containing the normal vector to the outer boundary, in Cartesian coordinates.
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 nbr_unknowns_from_adapted() const
Gives the number of unknowns coming from the variable shape of the domain.
virtual void set_cheb_base_t_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector.
virtual Term_eq connection_spher(const Term_eq &) const
Computes the part of the gradient involving the connections, in spherical orthonormal coordinates.
virtual Val_domain div_sin_theta(const Val_domain &) const
Division by .
Term_eq * der_rad_term_eq
Pointer on the Term_eq containing the .
virtual Val_domain laplacian2(const Val_domain &, int) const
Computes the ordinary flat 2dè- Laplacian for a scalar field with an harmonic index m.
virtual void export_tau_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 double integ_volume(const Val_domain &) const
Volume integral.
Val_domain * outer_radius
Pointer on the outer boundary , as a Val_domain.
virtual Val_domain mult_sin_theta(const Val_domain &) const
Multiplication by .
virtual void update_variable(const Val_domain &, const Scalar &, Scalar &) const
Update the value of a scalar, after the shape of the Domain has been changed by the system.
virtual void update_constante(const Val_domain &, const Scalar &, Scalar &) const
Update the value of a scalar, after the shape of the Domain has been changed by the system.
virtual Val_domain der_r(const Val_domain &) const
Compute the radial derivative of a scalar field.
Term_eq * rad_term_eq
Pointer on the Term_eq containing the radius.
virtual void affecte_tau_one_coef(Tensor &, int, int, int &) const
Sets at most one coefficient of a Tensor to 1.
virtual Term_eq integ_term_eq(const Term_eq &, int) const
Surface integral of a Term_eq.
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 Term_eq derive_flat_cart(int, char, const Term_eq &, const Metric *) const
Computes the flat derivative of a Term_eq, in Cartesian coordinates.
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 void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const
Computes the derivative with respect to the absolute Cartesian coordinates from the derivative with r...
virtual Val_domain laplacian(const Val_domain &, int) const
Computes the ordinary flat Laplacian for a scalar field with an harmonic index m.
virtual const Point absol_to_num(const Point &) const
Computes the numerical coordinates from the physical ones.
double inner_radius
The inner radius .
virtual void update_term_eq(Term_eq *) const
Update the value of a field, after the shape of the Domain has been changed by the system.
virtual void xx_to_ders_from_adapted(const Array< double > &, int &) const
Affects the derivative part of variable a Domain from a set of values.
virtual Tensor change_basis_spher_to_cart(int dd, const Tensor &) const
Changes the tensorial basis from spherical to Cartesian in a given domain.
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.
void do_normal_spher() const
Computes the normal wrt the inner boundary, in orthonormal spherical coordinates.
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
virtual void set_legendre_r_base(Base_spectral &) const
Gives the base using odd Legendre polynomials$ for the radius.
virtual void update_mapping(const Val_domain &)
Updates the variables parts of the Domain.
virtual void set_anti_legendre_base(Base_spectral &) const
Gives the base using Legendre polynomials, for functions antisymetric with respect to .
Term_eq * dt_rad_term_eq
Pointer on the Term_eq containing the .
virtual const Term_eq * give_normal(int, int) const
Returns the vector normal to a surface.
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 vars_to_terms() const
The Term_eq describing the variable shape of the Domain are updated.
virtual void xx_to_vars_from_adapted(Val_domain &, const Array< double > &, int &) const
Computes the new boundary of a Domain from a set of values.
virtual Val_domain ddp(const Val_domain &) const
Compute the second derivative with respect to of a scalar field.
virtual void set_anti_cheb_base(Base_spectral &) const
Gives the base using Chebyshev polynomials, for functions antisymetric with respect to .
virtual Val_domain der_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
int nbr_unknowns_val_domain(const Val_domain &so, int mlim) const
Computes the number of true unknowns of a Val_domain.
virtual Point get_center() const
Returns the center.
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 do_which_points_boundary(int, const Base_spectral &, Index **, int) const
Lists all the indices corresponding to true collocation points on a boundary.
virtual Term_eq mult_r_term_eq(const Term_eq &) const
Multiplication by of a Term_eq.
virtual Val_domain div_r(const Val_domain &) const
Division by .
virtual Val_domain mult_cos_theta(const Val_domain &) const
Multiplication by .
const Space & sp
The corresponding Space ; required for updating fields whene the mapping changes.
virtual int nbr_points_boundary(int, const Base_spectral &) const
Computes the number of relevant collocation points on a boundary.
virtual void do_cart() const
Computes the Cartesian coordinates.
Term_eq derive_p(const Term_eq &so) const
Computes .
Abstract class that implements the fonctionnalities common to all the type of domains.
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_adapted class fills the space with one shell adapted on the inside,...
Spacetime intended for binary black hole configurations in full general relativity (see constructor f...
Spacetime intended for binary black hole configurations (see constructor for details about the domain...
Spacetime intended for binary neutron stars configurations (see constructor for details about the dom...
The Space_spheric_adapted class fills the space with one nucleus, one shell adapted on the outside,...
void add_eq_ori(System_of_eqs &syst, const char *eq)
Adds an equation being the value of some field at the origin.
void add_eq(System_of_eqs &syst, const char *eq, const char *rac, const char *rac_der, int nused=-1, Array< int > **pused=0x0)
Adds a bulk equation and two matching conditions.
virtual void xx_to_ders_variable_domains(const Array< double > &, int &) const
Update the vairable domains from a set of values.
virtual ~Space_spheric_adapted()
Destructor
virtual void affecte_coef_to_variable_domains(int &, int, Array< int > &) const
The variation of the functions describing the shape of the Domain are affected from the unknowns of t...
void add_eq_int_inf(System_of_eqs &syst, const char *eq)
Adds an equation being a surface integral at infinity.
virtual Array< int > get_indices_matching_non_std(int, int) const
Gives the number of the other domains, touching a given boundary.
virtual void save(FILE *) const
Saving function.
Space_spheric_adapted(int ttyp, const Point &cr, const Dim_array &nbr, const Array< double > &bounds)
virtual int nbr_unknowns_from_variable_domains() const
Gives the number of unknowns coming from the variable shape of the domain.
virtual void xx_to_vars_variable_domains(System_of_eqs *, const Array< double > &, int &) const
Update the variables of a system, from the variation of the shape of the domains.
void add_eq_int_volume(System_of_eqs &syst, int nz, const char *eq)
Adds an equation being a volume integral in the domains below a given number.
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...