KADATH
|
The Space_spheric_adapted
class fills the space with one nucleus, one shell adapted on the outside, one shell adapted on the inside, several standard shells and a compactified domain, all centered on the same point.
More...
#include <adapted.hpp>
Public Member Functions | |
Space_spheric_adapted (int ttyp, const Point &cr, const Dim_array &nbr, const Array< double > &bounds) | |
Space_spheric_adapted (FILE *) | |
Constructor from a file. More... | |
virtual | ~Space_spheric_adapted () |
Destructor More... | |
virtual void | save (FILE *) const |
Saving function. More... | |
virtual int | nbr_unknowns_from_variable_domains () const |
Gives the number of unknowns coming from the variable shape of the domain. More... | |
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 the system. More... | |
virtual void | xx_to_ders_variable_domains (const Array< double > &, int &) const |
Update the vairable domains from a set of values. More... | |
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. More... | |
void | add_eq_ori (System_of_eqs &syst, const char *eq) |
Adds an equation being the value of some field at the origin. More... | |
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. More... | |
void | add_eq_int_inf (System_of_eqs &syst, const char *eq) |
Adds an equation being a surface integral at infinity. More... | |
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. More... | |
virtual Array< int > | get_indices_matching_non_std (int, int) const |
Gives the number of the other domains, touching a given boundary. More... | |
int | get_ndim () const |
Returns the number of dimensions. More... | |
int | get_nbr_domains () const |
Returns the number of Domains . More... | |
int | get_type_base () const |
Returns the type of basis. More... | |
const Domain * | get_domain (int i) const |
returns a pointer on the domain. More... | |
Protected Attributes | |
int | nbr_domains |
Number od Domains . More... | |
int | ndim |
Number of dimensions (should be the same for all the Domains ). More... | |
int | type_base |
Type of basis used (i.e. using either Chebyshev or Legendre polynomials). More... | |
Domain ** | domains |
Pointers on the various Domains . More... | |
The Space_spheric_adapted
class fills the space with one nucleus, one shell adapted on the outside, one shell adapted on the inside, several standard shells and a compactified domain, all centered on the same point.
Definition at line 661 of file adapted.hpp.
Kadath::Space_spheric_adapted::Space_spheric_adapted | ( | int | ttyp, |
const Point & | cr, | ||
const Dim_array & | nbr, | ||
const Array< double > & | bounds | ||
) |
Standard constructor ; all the shells are initially not deformed.
ttyp | [input] : the type of basis. |
cr | [input] : absolute coordinates of the center. |
nbr | [input] : number of points in each domain. |
bounds | [input] : radii of the various shells (and also determines the total number of domains). |
Definition at line 29 of file space_spheric_adapted.cpp.
References Kadath::Space::domains, Kadath::Array< T >::get_ndim(), Kadath::Array< T >::get_size(), Kadath::Space::nbr_domains, Kadath::Space::ndim, Kadath::Space::type_base, Kadath::Domain_shell_inner_adapted::vars_to_terms(), and Kadath::Domain_shell_outer_adapted::vars_to_terms().
Kadath::Space_spheric_adapted::Space_spheric_adapted | ( | FILE * | fd | ) |
Constructor from a file.
Definition at line 57 of file space_spheric_adapted.cpp.
References Kadath::Space::domains, Kadath::Space::nbr_domains, Kadath::Space::ndim, Kadath::Space::type_base, Kadath::Domain_shell_inner_adapted::vars_to_terms(), and Kadath::Domain_shell_outer_adapted::vars_to_terms().
|
virtual |
Destructor
Definition at line 78 of file space_spheric_adapted.cpp.
References Kadath::Domain_shell_inner_adapted::del_deriv(), Kadath::Domain_shell_outer_adapted::del_deriv(), and Kadath::Space::domains.
void Kadath::Space_spheric_adapted::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.
syst | : the System_of_eqs . |
eq | : the string describing the bulk equation. |
rac | : the string describing the first matching condition. |
rac_der | : the string describing the second matching condition. |
nused | : number of components of eq to be considered. All the components are used of it is -1. |
pused | : pointer on the indexes of the components to be considered. Not used of nused = -1 . |
Definition at line 82 of file spheric_adapted_add_eq.cpp.
References Kadath::System_of_eqs::add_eq_inside(), Kadath::System_of_eqs::add_eq_matching(), Kadath::System_of_eqs::get_dom_max(), and Kadath::System_of_eqs::get_dom_min().
void Kadath::Space_spheric_adapted::add_eq_int_inf | ( | System_of_eqs & | syst, |
const char * | eq | ||
) |
Adds an equation being a surface integral at infinity.
syst | : the System_of_eqs . |
eq | : the string describing the equation (should contain something like integ(f)=b) |
Definition at line 29 of file spheric_adapted_add_eq.cpp.
References Kadath::Space::domains, Kadath::System_of_eqs::eq_int, Kadath::System_of_eqs::give_ope(), Kadath::System_of_eqs::is_ope_bin(), Kadath::System_of_eqs::nbr_conditions, Kadath::Space::nbr_domains, and Kadath::System_of_eqs::neq_int.
void Kadath::Space_spheric_adapted::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.
syst | : the System_of_eqs . |
nz | : the integral is taken for all the Domains which number is ![]() |
eq | : the string describing the equation (should contain something like integvolume(f)=b) |
Definition at line 58 of file spheric_adapted_add_eq.cpp.
References Kadath::System_of_eqs::eq_int, Kadath::System_of_eqs::give_ope(), Kadath::System_of_eqs::is_ope_bin(), Kadath::System_of_eqs::nbr_conditions, and Kadath::System_of_eqs::neq_int.
void Kadath::Space_spheric_adapted::add_eq_ori | ( | System_of_eqs & | syst, |
const char * | eq | ||
) |
Adds an equation being the value of some field at the origin.
syst | : the System_of_eqs . |
eq | : the string describing the quantity that must be zero at the origin |
Definition at line 195 of file space_spheric_adapted.cpp.
References Kadath::System_of_eqs::add_eq_val(), and Kadath::Space::domains.
|
virtual |
The variation of the functions describing the shape of the Domain
are affected from the unknowns of the system.
conte | : current position if the unknowns. |
cc | : position of the unknown to be set. |
doms | : Array containing the index of the variable domains present in this space. |
Reimplemented from Kadath::Space.
Definition at line 98 of file space_spheric_adapted.cpp.
References Kadath::Domain::affecte_coef(), Kadath::Space::domains, and Kadath::Array< T >::set().
|
inlineinherited |
returns a pointer on the domain.
i | [input] : the index of the domain. |
Definition at line 1385 of file space.hpp.
References Kadath::Space::domains, and Kadath::Space::nbr_domains.
|
virtual |
Gives the number of the other domains, touching a given boundary.
It also gives the name of the boundary, as seen by the other domain.
dom | : the domain considered. |
bound | : the boundary. |
Reimplemented from Kadath::Space.
Definition at line 174 of file space_spheric_adapted.cpp.
References Kadath::Space::nbr_domains, and Kadath::Array< T >::set().
|
inlineinherited |
Returns the number of Domains
.
Definition at line 1375 of file space.hpp.
References Kadath::Space::nbr_domains.
|
inlineinherited |
Returns the number of dimensions.
Definition at line 1373 of file space.hpp.
References Kadath::Space::ndim.
|
inlineinherited |
Returns the type of basis.
Definition at line 1377 of file space.hpp.
References Kadath::Space::type_base.
|
virtual |
Gives the number of unknowns coming from the variable shape of the domain.
Reimplemented from Kadath::Space.
Definition at line 93 of file space_spheric_adapted.cpp.
References Kadath::Space::domains, and Kadath::Domain::nbr_unknowns_from_adapted().
|
virtual |
Saving function.
Reimplemented from Kadath::Space.
Definition at line 85 of file space_spheric_adapted.cpp.
References Kadath::Space::domains, Kadath::Space::nbr_domains, Kadath::Space::ndim, and Kadath::Space::type_base.
|
virtual |
Update the vairable domains from a set of values.
xx | : set of values used by the affectation |
conte | : current position in the values vector. |
Reimplemented from Kadath::Space.
Definition at line 117 of file space_spheric_adapted.cpp.
References Kadath::Space::domains, and Kadath::Domain::xx_to_ders_from_adapted().
|
virtual |
Update the variables of a system, from the variation of the shape of the domains.
syst | : the System_of_eqs considered. |
xx | : set of values used by the affectation |
conte | : current position in the values vector. |
Reimplemented from Kadath::Space.
Definition at line 125 of file space_spheric_adapted.cpp.
References Kadath::System_of_eqs::cst, Kadath::System_of_eqs::dom_max, Kadath::System_of_eqs::dom_min, Kadath::Space::domains, Kadath::System_of_eqs::ncst, Kadath::System_of_eqs::nvar, Kadath::Scalar::set_domain(), Kadath::Domain::update_constante(), Kadath::Domain::update_mapping(), Kadath::Domain::update_variable(), Kadath::System_of_eqs::var, and Kadath::Domain::xx_to_vars_from_adapted().
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |