Hermes2D
3.0
|
Caches precalculated shape function values. More...
#include <precalc.h>
Public Member Functions | |
SpaceType | get_space_type () const |
Returns type of space. | |
PrecalcShapeset (Shapeset *shapeset) | |
Constructs a standard (master) precalculated shapeset class. More... | |
virtual | ~PrecalcShapeset () |
Destructor. | |
virtual void | set_active_shape (int index) |
Public Member Functions inherited from Hermes::Hermes2D::Function< double > | |
Function () | |
Default constructor. | |
virtual | ~Function () |
unsigned char | get_num_components () const |
Returns the number of components of the function being represented by the class. | |
virtual const double * | get_fn_values (int component=0) const |
Returns function values. More... | |
virtual const double * | get_dx_values (int component=0) const |
Returns the x partial derivative. More... | |
virtual const double * | get_dy_values (int component=0) const |
Returns the y partial derivative. More... | |
double * | deep_copy_array (int component=0, int item=0) const |
Returns function values. More... | |
Quad2D * | get_quad_2d () const |
Returns the current quadrature points. | |
void | set_quad_order (unsigned short order, unsigned short mask=H2D_FN_DEFAULT) |
virtual const double * | get_values (int component, int item) const |
virtual int | get_fn_order () const |
Returns the polynomial degree of the function being represented by the class. | |
virtual void | push_transform (int son) |
virtual void | pop_transform () |
virtual void | set_active_element (Element *e) |
Sets the active element. | |
virtual void | set_transform (uint64_t idx) |
Public Member Functions inherited from Hermes::Hermes2D::Transformable | |
Element * | get_active_element () const |
uint64_t | get_transform () const |
Protected Member Functions | |
virtual void | set_quad_2d (Quad2D *quad_2d) |
Selects the quadrature points in which the function will be evaluated. More... | |
virtual void | free () |
Frees all precalculated tables. | |
int | get_active_shape () const |
Returns the index of the active shape (can be negative if the shape is constrained). | |
Shapeset * | get_shapeset () const |
Returns a pointer to the shapeset which is being precalculated. | |
virtual unsigned short | get_edge_fn_order (int edge) |
Returns the polynomial order of the active shape function on given edge. | |
virtual void | precalculate (unsigned short order, unsigned short mask) |
precalculates the current function at the current integration points. | |
void | update_max_index () |
Protected Member Functions inherited from Hermes::Hermes2D::Function< double > | |
virtual void | reset_transform () |
Empties the stack, loads identity transform. | |
virtual void | force_transform (uint64_t sub_idx, Trf *ctm) |
For internal use only. | |
virtual int | get_edge_fn_order (unsigned char edge) const |
Returns the polynomial degree of the function at given edge. To be overridden in derived classes. More... | |
void | invalidate_values () |
Protected Member Functions inherited from Hermes::Hermes2D::Transformable | |
double | get_transform_jacobian () const |
Trf * | get_ctm () const |
unsigned int | get_depth () const |
Protected Attributes | |
Shapeset * | shapeset |
int | index |
unsigned short | max_index [H2D_NUM_MODES] |
double2 | ref_points [H2D_MAX_INTEGRATION_POINTS_COUNT] |
Transformed points to the reference domain, used by precalculate. | |
Protected Attributes inherited from Hermes::Hermes2D::Function< double > | |
double | values [H2D_MAX_SOLUTION_COMPONENTS][H2D_NUM_FUNCTION_VALUES][H2D_MAX_INTEGRATION_POINTS_COUNT] |
The data. | |
bool | values_valid |
Flag that the data are not 'dirty'. | |
int | order |
Current function polynomial order. | |
unsigned char | num_components |
Number of vector components. | |
Quad2D * | quads [H2D_MAX_QUADRATURES] |
List of available quadratures. | |
int | cur_quad |
Active quadrature (index into 'quads') | |
Protected Attributes inherited from Hermes::Hermes2D::Transformable | |
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. | |
Additional Inherited Members | |
Static Public Attributes inherited from Hermes::Hermes2D::Transformable | |
static const unsigned int | H2D_MAX_TRN_LEVEL = 15 |
If this changes, NeighborSearch::H2D_MAX_NEIGHBORS must change too. | |
Static Protected Member Functions inherited from Hermes::Hermes2D::Transformable | |
static void | push_transforms (std::set< Transformable * > &transformables, int son) |
static void | pop_transforms (std::set< Transformable * > &transformables) |
Static Protected Attributes inherited from Hermes::Hermes2D::Function< double > | |
static int | idx2mask [H2D_NUM_FUNCTION_VALUES][2] |
Index to mask table. | |
Static Protected Attributes inherited from Hermes::Hermes2D::Transformable | |
static const uint64_t | H2D_MAX_IDX = (1ULL << 3 * H2D_MAX_TRN_LEVEL) - 1 |
Caches precalculated shape function values.
PrecalcShapeset is a cache of precalculated shape function values.
Hermes::Hermes2D::PrecalcShapeset::PrecalcShapeset | ( | Shapeset * | shapeset | ) |
Constructs a standard (master) precalculated shapeset class.
shapeset[in] | Pointer to the shapeset to be precalculated. |
Definition at line 29 of file precalc.cpp.
|
virtual |
Activates a shape function given by its index. The values of the shape function can then be obtained by setting the required integration rule order by calling set_quad_order() and after that calling get_values(), get_dx_values(), etc.
index[in] | Shape index. |
Definition at line 50 of file precalc.cpp.
Referenced by Hermes::Hermes2D::CurvMapStatic::calculate_bubble_projection_matrix().
|
protectedvirtual |
Selects the quadrature points in which the function will be evaluated.
It is possible to switch back and forth between different quadrature points: no precalculated values are freed. The standard quadrature is always selected by default already.
quad_2d[in] | The quadrature points. |
Reimplemented from Hermes::Hermes2D::Function< double >.
Definition at line 45 of file precalc.cpp.
Referenced by PrecalcShapeset(), and Hermes::Hermes2D::CurvMap::update_refmap_coeffs().