Hermes2D
2.0
|
Represents the solution of a PDE.
More...
#include <solution.h>
Public Member Functions | |
Solution (const Mesh *mesh) | |
Solution (Space< Scalar > *s, Vector< Scalar > *coeff_vec) | |
Solution (Space< Scalar > *s, Scalar *coeff_vec) | |
std::string | getClassName () const |
State querying helpers. | |
void | assign (Solution *sln) |
Solution & | operator= (Solution &sln) |
virtual void | copy (const Solution< Scalar > *sln) |
void | set_dirichlet_lift (const Space< Scalar > *space, PrecalcShapeset *pss=NULL) |
Sets solution equal to Dirichlet lift only, solution vector = 0. | |
virtual void | save (const char *filename) const |
void | load (const char *filename, Space< Scalar > *space) |
Scalar | get_ref_value (Element *e, double xi1, double xi2, int component=0, int item=0) |
Scalar | get_ref_value_transformed (Element *e, double xi1, double xi2, int a, int b) |
virtual Func< Scalar > * | get_pt_value (double x, double y) |
void | multiply (Scalar coef) |
Multiplies the function represented by this class by the given coefficient. | |
SolutionType | get_type () const |
Returns solution type. | |
SpaceType | get_space_type () const |
virtual void | set_active_element (Element *e) |
Internal. | |
virtual MeshFunction< Scalar > * | clone () const |
template<> | |
Solution (const Mesh *mesh) | |
template<> | |
Solution (const Mesh *mesh) | |
template<> | |
Solution (Space< double > *s, Vector< double > *coeff_vec) | |
template<> | |
Solution (Space< std::complex< double > > *s, Vector< std::complex< double > > *coeff_vec) | |
template<> | |
Solution (Space< double > *s, double *coeff_vec) | |
template<> | |
Solution (Space< std::complex< double > > *s, std::complex< double > *coeff_vec) | |
template<> | |
void | save (const char *filename) const |
template<> | |
void | save (const char *filename) const |
template<> | |
void | load (const char *filename, Space< double > *space) |
template<> | |
void | load (const char *filename, Space< std::complex< double > > *space) |
![]() | |
MeshFunction (const Mesh *mesh) | |
virtual bool | isOkay () const |
State querying helpers. | |
virtual void | reinit () |
virtual void | set_quad_2d (Quad2D *quad_2d) |
Selects the quadrature points in which the function will be evaluated. More... | |
const Mesh * | get_mesh () const |
RefMap * | get_refmap (bool update=true) |
void | set_refmap (RefMap *refmap_to_set) |
virtual void | handle_overflow_idx () |
virtual void | push_transform (int son) |
See Transformable::push_transform. | |
virtual void | pop_transform () |
![]() | |
Function () | |
Default constructor. | |
virtual | ~Function () |
int | get_num_components () const |
Returns the number of components of the function being represented by the class. | |
Scalar * | get_fn_values (int component=0) |
Returns function values. More... | |
Scalar * | get_dx_values (int component=0) |
Returns the x partial derivative. More... | |
Scalar * | get_dy_values (int component=0) |
Returns the y partial derivative. More... | |
void | get_dx_dy_values (Scalar *&dx, Scalar *&dy, int component=0) |
Returns both x and y partial derivatives. This function provides the both often-used dx and dy values in one call. More... | |
Scalar * | get_dxx_values (int component=0) |
Returns the second x partial derivative. More... | |
Scalar * | get_dyy_values (int component=0) |
Returns the second y partial derivative. More... | |
Scalar * | get_dxy_values (int component=0) |
Returns the second mixed derivative. More... | |
Quad2D * | get_quad_2d () const |
Returns the current quadrature points. | |
void | set_quad_order (unsigned int order, int mask=H2D_FN_DEFAULT) |
Scalar * | get_values (int a, int b) |
int | get_fn_order () const |
Returns the polynomial degree of the function being represented by the class. | |
![]() | |
Element * | get_active_element () const |
void | set_transform (uint64_t idx) |
uint64_t | get_transform () const |
![]() | |
void | check () const |
Method to handle the state. | |
![]() | |
XMLParsing () | |
Constructor. | |
void | set_validation (bool to_set) |
Set to validate / not to validate. | |
Static Public Member Functions | |
static void | vector_to_solutions (const Scalar *solution_vector, Hermes::vector< const Space< Scalar > * > spaces, Hermes::vector< Solution< Scalar > * > solutions, Hermes::vector< bool > add_dir_lift=Hermes::vector< bool >(), Hermes::vector< int > start_indices=Hermes::vector< int >()) |
Passes solution components calculated from solution vector as Solutions. | |
static void | vector_to_solution (const Scalar *solution_vector, const Space< Scalar > *space, Solution< Scalar > *solution, bool add_dir_lift=true, int start_index=0) |
static void | vector_to_solutions (const Vector< Scalar > *vec, Hermes::vector< const Space< Scalar > * > spaces, Hermes::vector< Solution< Scalar > * > solutions, Hermes::vector< bool > add_dir_lift=Hermes::vector< bool >(), Hermes::vector< int > start_indices=Hermes::vector< int >()) |
static void | vector_to_solutions_common_dir_lift (const Vector< Scalar > *vec, Hermes::vector< const Space< Scalar > * > spaces, Hermes::vector< Solution< Scalar > * > solutions, bool add_dir_lift=false) |
static void | vector_to_solutions_common_dir_lift (const Scalar *solution_vector, Hermes::vector< const Space< Scalar > * > spaces, Hermes::vector< Solution< Scalar > * > solutions, bool add_dir_lift=false) |
static void | vector_to_solution (const Vector< Scalar > *vec, const Space< Scalar > *space, Solution< Scalar > *solution, bool add_dir_lift=true, int start_index=0) |
static void | vector_to_solutions (const Scalar *solution_vector, Hermes::vector< const Space< Scalar > * > spaces, Hermes::vector< Solution< Scalar > * > solutions, Hermes::vector< PrecalcShapeset * > pss, Hermes::vector< bool > add_dir_lift=Hermes::vector< bool >(), Hermes::vector< int > start_indices=Hermes::vector< int >()) |
static void | vector_to_solution (const Scalar *solution_vector, const Space< Scalar > *space, Solution< Scalar > *solution, PrecalcShapeset *pss, bool add_dir_lift=true, int start_index=0) |
static void | set_static_verbose_output (bool verbose) |
Protected Member Functions | |
virtual int | get_edge_fn_order (int edge) |
void | enable_transform (bool enable=true) |
virtual void | init () |
virtual void | free () |
Frees all precalculated tables. | |
virtual void | set_coeff_vector (const Space< Scalar > *space, const Vector< Scalar > *vec, bool add_dir_lift, int start_index) |
Converts a coefficient vector into a Solution. | |
virtual void | set_coeff_vector (const Space< Scalar > *space, PrecalcShapeset *pss, const Scalar *coeffs, bool add_dir_lift, int start_index) |
virtual void | set_coeff_vector (const Space< Scalar > *space, const Scalar *coeffs, bool add_dir_lift, int start_index) |
void | transform_values (int order, struct Function< Scalar >::Node *node, int newmask, int oldmask, int np) |
virtual void | precalculate (int order, int mask) |
precalculates the current function at the current integration points. | |
double ** | calc_mono_matrix (int o, int *&perm) |
void | init_dxdy_buffer () |
void | free_tables () |
template<> | |
void | init () |
template<> | |
void | init () |
template<> | |
void | free () |
Frees all precalculated tables. | |
template<> | |
void | free () |
Frees all precalculated tables. | |
![]() | |
void | force_transform (MeshFunction< Scalar > *mf) |
void | update_refmap () |
void | force_transform (uint64_t sub_idx, Trf *ctm) |
![]() | |
virtual int | get_edge_fn_order (int edge) const |
Returns the polynomial degree of the function at given edge. To be overridden in derived classes. More... | |
void | update_nodes_ptr () |
void | force_transform (uint64_t sub_idx, Trf *ctm) |
For internal use only. | |
Node * | new_node (int mask, int num_points) |
allocates a new Node structure | |
void | replace_cur_node (Node *node) |
![]() | |
void | reset_transform () |
Empties the stack, loads identity transform. | |
double | get_transform_jacobian () const |
Trf * | get_ctm () const |
unsigned int | get_depth () const |
Protected Attributes | |
SolutionType | sln_type |
SpaceType | space_type |
bool | transform |
std::map< uint64_t, LightArray < struct Function< Scalar > ::Node * > * > * | tables [H2D_MAX_QUADRATURES][H2D_SOLUTION_ELEMENT_CACHE_SIZE] |
Element * | elems [H2D_MAX_QUADRATURES][H2D_SOLUTION_ELEMENT_CACHE_SIZE] |
int | cur_elem |
int | oldest [H2D_SOLUTION_ELEMENT_CACHE_SIZE] |
Scalar * | mono_coeffs |
monomial coefficient array | |
int * | elem_coeffs [H2D_MAX_SOLUTION_COMPONENTS] |
int * | elem_orders |
Stored element orders in the mathematical sense. The polynomial degree of the highest basis function + increments due to the element shape, etc. . | |
int | num_coeffs |
int | num_elems |
int | num_dofs |
Scalar * | dxdy_coeffs [H2D_MAX_SOLUTION_COMPONENTS][6] |
Scalar * | dxdy_buffer |
Element * | e_last |
last visited element when getting solution values at specific points | |
![]() | |
ElementMode2D | mode |
const Mesh * | mesh |
RefMap * | refmap |
![]() | |
int | order |
current function polynomial order | |
int | num_components |
number of vector components | |
std::map< uint64_t, LightArray < Node * > * > * | sub_tables |
Table of Node tables, for each possible transformation there can be a different Node table. | |
LightArray< Node * > * | nodes |
Table of nodes. | |
Node * | cur_node |
Current Node. | |
LightArray< Node * > * | overflow_nodes |
Nodes for the overflow sub-element transformation. | |
Quad2D * | quads [8] |
list of available quadratures | |
int | cur_quad |
active quadrature (index into 'quads') | |
int | total_mem |
total memory in bytes used by the tables | |
int | max_mem |
peak memory usage | |
![]() | |
Element * | element |
The active element. | |
Trf * | ctm |
Current sub-element transformation matrix. | |
uint64_t | sub_idx |
Sub-element transformation index. | |
Trf | stack [21] |
Transformation matrix stack. | |
unsigned int | top |
Stack top. | |
![]() | |
bool | validate |
Internal. | |
Static Protected Attributes | |
static bool | static_verbose_output = false |
Friends | |
class | RefMap |
template<typename T > | |
class | KellyTypeAdapt |
template<typename T > | |
class | CalculationContinuity |
template<typename T > | |
class | OGProjection |
template<typename T > | |
class | OGProjectionNOX |
template<typename T > | |
class | Adapt |
template<typename T > | |
class | Func |
template<typename T > | |
class | DiscontinuousFunc |
template<typename T > | |
class | DiscreteProblem |
template<typename T > | |
class | DiscreteProblemLinear |
template<typename T > | |
class | NeighborSearch |
template<typename T > | |
class | RefinementSelectors::ProjBasedSelector |
template<typename T > | |
class | RefinementSelectors::H1ProjBasedSelector |
template<typename T > | |
class | RefinementSelectors::L2ProjBasedSelector |
template<typename T > | |
class | RefinementSelectors::HcurlProjBasedSelector |
template<typename T > | |
HERMES_API Func< T > * | init_fn (Solution< T > *fu, const int order) |
Additional Inherited Members | |
![]() | |
static void | check_params (int component, Node *cur_node, int num_components) |
static void | check_table (int component, Node *cur_node, int n, const char *msg) |
Represents the solution of a PDE.
The Solution class represents the solution of a PDE. Given a space and a solution vector, it calculates the appropriate linear combination of basis functions at the specified element and integration points.
The higher-order solution on elements is best calculated not as a linear combination of shape functions (the usual approach), but as a linear combination of monomials. This has the advantage that no shape function table calculation and look-ups are necessary (except for the conversion of the solution coefficients), which means that visualization and multi-mesh calculations are much faster (all the push_transforms and table searches take the most time when evaluating the solution).
The linear combination of monomials can be calculated using the Horner's scheme, which requires the same number of multiplications as the calculation of the linear combination of shape functions. However, sub-element transforms are trivial and cheap. Moreover, after the solution on all elements is expressed as a combination of monomials, the Space can be forgotten. This is comfortable for the user, since the Solution class acts as a self-contained unit, internally containing just a copy of the mesh and a table of monomial coefficients. It is also very straight-forward to obtain all derivatives of a solution defined in this way. Finally, it is possible to store the solution on the disk easily (no need to store the Space, which is difficult).
The following is an example of the set of monomials for a cubic quad and a cubic triangle. (Note that these are actually the definitions of the polynomial spaces on these elements.)
[ x^3*y^3 x^2*y^3 x*y^3 y^3 ] [ y^3 ] [ x^3*y^2 x^2*y^2 x*y^2 y^2 ] [ x*y^2 y^2 ] [ x^3*y x^2*y x*y y ] [ x^2*y x*y y ] [ x^3 x^2 x 1 ] [ x^3 x^2 x 1 ]
The number of monomials is (p + 1)^2 for quads and (p + 1)*(p + 2)/2 for triangles, where 'p' is the polynomial degree.
|
protected |
Enables or disables transformation of the solution derivatives (H1 case) or values (vector (Hcurl) case). This means H2D_FN_DX_0 and H2D_FN_DY_0 or H2D_FN_VAL_0 and H2D_FN_VAL_1 will or will not be returned premultiplied by the reference mapping matrix. The default is enabled (true).
Definition at line 801 of file solution.cpp.
Referenced by Hermes::Hermes2D::RefinementSelectors::ProjBasedSelector< Scalar >::calc_projection_errors().
|
virtual |
Returns solution value or derivatives at the physical domain point (x, y). 'item' controls the returned value: H2D_FN_VAL_0, H2D_FN_VAL_1, H2D_FN_DX_0, H2D_FN_DX_1, H2D_FN_DY_0, .... NOTE: This function should be used for postprocessing only, it is not effective enough for calculations. Since it searches for an element sequentinally, it is extremelly slow. Prefer Solution::get_ref_value if possible.
Implements Hermes::Hermes2D::MeshFunction< Scalar >.
Definition at line 1627 of file solution.cpp.
Scalar Hermes::Hermes2D::Solution< Scalar >::get_ref_value | ( | Element * | e, |
double | xi1, | ||
double | xi2, | ||
int | component = 0 , |
||
int | item = 0 |
||
) |
Returns solution value or derivatives at element e, in its reference domain point (xi1, xi2). 'item' controls the returned value: 0 = value, 1 = dx, 2 = dy, 3 = dxx, 4 = dyy, 5 = dxy. NOTE: This function should be used for postprocessing only, it is not effective enough for calculations.
Definition at line 1544 of file solution.cpp.
Scalar Hermes::Hermes2D::Solution< Scalar >::get_ref_value_transformed | ( | Element * | e, |
double | xi1, | ||
double | xi2, | ||
int | a, | ||
int | b | ||
) |
Returns solution value or derivatives (correctly transformed) at element e, in its reference domain point (xi1, xi2). 'item' controls the returned value: 0 = value, 1 = dx, 2 = dy, 3 = dxx, 4 = dyy, 5 = dxy. NOTE: This function should be used for postprocessing only, it is not effective enough for calculations.
Definition at line 1566 of file solution.cpp.
void Hermes::Hermes2D::Solution< Scalar >::load | ( | const char * | filename, |
Space< Scalar > * | space | ||
) |
Loads the solution from a file previously created by Solution::save(). This completely restores the solution in the memory.
Referenced by Hermes::Hermes2D::CalculationContinuity< Scalar >::Record::load_solution().
|
virtual |
Saves the complete solution (i.e., including the internal copy of the mesh and element orders) to an XML file.
Reimplemented in Hermes::Hermes2D::ZeroSolutionVector< Scalar >, Hermes::Hermes2D::ZeroSolutionVector< Scalar >, Hermes::Hermes2D::ConstantSolutionVector< Scalar >, Hermes::Hermes2D::ConstantSolutionVector< Scalar >, Hermes::Hermes2D::ZeroSolution< Scalar >, Hermes::Hermes2D::ZeroSolutionVector< Scalar >, Hermes::Hermes2D::ConstantSolutionVector< Scalar >, Hermes::Hermes2D::ZeroSolution< Scalar >, Hermes::Hermes2D::ZeroSolution< Scalar >, Hermes::Hermes2D::ConstantSolution< Scalar >, Hermes::Hermes2D::ConstantSolution< Scalar >, and Hermes::Hermes2D::ConstantSolution< Scalar >.
Referenced by Hermes::Hermes2D::ConstantSolution< Scalar >::save(), Hermes::Hermes2D::ZeroSolution< Scalar >::save(), Hermes::Hermes2D::ConstantSolutionVector< Scalar >::save(), Hermes::Hermes2D::ZeroSolutionVector< Scalar >::save(), and Hermes::Hermes2D::CalculationContinuity< Scalar >::Record::save_solution().
|
protected |
array of pointers into mono_coeffs
Definition at line 213 of file solution.h.
|
protected |
Precalculated tables for last four used elements. There is a 2-layer structure of the precalculated tables. The first (the lowest) one is the layer where mapping of integral orders to Function::Node takes place. See function.h for details. The second one is the layer with mapping of sub-element transformation to a table from the lowest layer. The highest layer (in contrast to the PrecalcShapeset class) is represented here only by this array.
Definition at line 207 of file solution.h.