Example codes
and tutorials
The (Not-so) Quick Guide
List of tutorials/demo codes
Single-Physics Problems
Poisson
Adaptivity illustrated for Poisson
Advection-Diffusion
Unsteady heat equation
Linear wave equation
The Young-Laplace equation
Navier-Stokes
Free-surface Navier-Stokes
Axisymmetric Navier-Stokes
Solid mechanics
Beam structures
Shell structures
Multi-physics Problems
Fluid-structure interaction
Boussinesq convection
Steady thermoelasticity
Methods-based example codes and tutorials
Mesh generation
Linear solvers and preconditioners
Visualisation of the results
Parallel processing
How to write a new element
How to write a new refineable element
Default nonlinear solvers -- the sequence of action functions
...
Documentation
FE theory and top-down discussion of the data structure
The (Not-so) Quick Guide
Comprehensive bottom-up discussion of the data structure
List of available structured and unstructured meshes
Linear solvers and preconditioners
Visualisation of the results
Parallel processing
Coding conventions and C++ style
Creating documentation
Optimisation - robustness vs. "raw speed"
Linear vs. nonlinear problems
Storing shape functions
Changing the default "full" integration scheme
Disabling the ALE formulation of unsteady equations
C vs. C++ output
Different sparse assembly techniques and the STL memory pool
Publications
Publications
Talks
Journal publications
Theses
Picture show
Download
Copyright
Download/installation instructions
Download page
FAQ & Contact
FAQ
Change log
Bugs and other known problems
Completeness of the library & our "To-Do List"
Contact the developers
Get involved

 


Beta release!

Please note that the library has not been "officially" released. While we continue to work on the documentation, these web pages are likely to contain broken links and documents in draft form. Please send an email to

oomph-lib AT maths DOT man DOT ac DOT uk

if you wish to be informed of the library's "official" release.

oomph::RefineableElement Class Reference

#include <refineable_elements.h>

Inheritance diagram for oomph::RefineableElement:

oomph::FiniteElement oomph::GeneralisedElement oomph::GeomObject oomph::NonRefineableElementWithHangingNodes oomph::RefineableAdvectionDiffusionEquations< DIM > oomph::RefineableAdvectionDiffusionReactionEquations< NREAGENT, DIM > oomph::RefineableAxisymmetricNavierStokesEquations oomph::RefineableGeneralisedAdvectionDiffusionEquations< DIM > oomph::RefineableLinearElasticityEquations< DIM > oomph::RefineableLinearWaveEquations< DIM > oomph::RefineableNavierStokesEquations< DIM > oomph::RefineablePoissonEquations< DIM > oomph::RefineableQElement< 2 > oomph::RefineableQElement< 3 > oomph::RefineableSolidElement oomph::RefineableUnsteadyHeatEquations< DIM > List of all members.

Detailed Description

RefineableElements are FiniteElements that may be subdivided into children to provide a better local approximation to the solution. After non-uniform refinement adjacent elements need not necessarily have nodes in common. A node that does not have a counterpart in its neighbouring element is known as a hanging node and its position and any data that it stores must be constrained to ensure inter-element continuity.

Generic data and function interfaces associated with refinement are defined in this class.

Additional data includes:

Additional functions perform the following generic tasks:

In addition, there are a number of interfaces that specify element-specific tasks. These should be overloaded in RefineableElements of particular geometric types and perform the following tasks:

In mixed element different sets of nodes are used to interpolate different unknowns. Interfaces are provided for functions that can be used to find the position of the nodes that interpolate the different unknowns. These functions are used to setup hanging node information automatically in particular elements, e.g. Taylor Hood Navier--Stokes. The default implementation assumes that the elements are isoparameteric.

Definition at line 99 of file refineable_elements.h.


Public Member Functions

 RefineableElement ()
 Constructor, calls the FiniteElement constructor and initialises the member data.
virtual ~RefineableElement ()
 RefineableElement (const RefineableElement &)
 Broken copy constructor.
void operator= (const RefineableElement &)
 Broken assignment operator.
Treetree_pt ()
 Access function: Pointer to quadtree representation of this element.
void set_tree_pt (Tree *my_tree_pt)
 Set pointer to quadtree representation of this element.
virtual unsigned required_nsons () const
 Set the number of sons that can be constructed by the element The default is none.
template<class ELEMENT>
void split (Vector< ELEMENT * > &son_pt) const
 Split the element into the number of sons to be constructed and return a vector of pointers to the sons. Elements are allocated, but they are not given any properties. The refinement level of the sons is one higher than that of the father elemern.
virtual void build (Mesh *&mesh_pt, Vector< Node * > &new_node_pt, bool &was_already_built, std::ofstream &new_nodes_file)=0
 Interface to function that builds the element: i.e. construct the nodes, assign their positions,.
void set_refinement_level (const int &refine_level)
 Set the refinement level.
unsigned refinement_level () const
 Return the Refinement level.
void select_for_refinement ()
 Select the element for refinement.
void deselect_for_refinement ()
 Deselect the element for refinement.
void select_sons_for_unrefinement ()
 Unrefinement will be performed by merging the four sons of this element.
void deselect_sons_for_unrefinement ()
 No unrefinement will be performed by merging the four sons of this element.
bool to_be_refined ()
 Has the element been selected for refinement?
bool sons_to_be_unrefined ()
 Has the element been selected for unrefinement?
virtual void rebuild_from_sons (Mesh *&mesh_pt)=0
 Rebuild the element, e.g. set internal values in line with those of the sons that have now merged.
virtual void unbuild ()
 Unbuild the element, i.e. mark the nodes that were created during its creation for possible deletion.
virtual void deactivate_element ()
 Final operations that must be performed when the element is no longer active in the mesh, but still resident in the QuadTree.
virtual bool nodes_built ()
 Return true if all the nodes have been built, false if not.
long number () const
 Element number (for debugging/plotting).
void set_number (const long &mynumber)
 Set element number (for debugging/plotting).
virtual unsigned ncont_interpolated_values () const =0
 Number of continuously interpolated values. Note: We assume that they are located at the beginning of the value_pt Vector! (Used for interpolation to son elements, for integrity check and post-processing -- we can only expect the continously interpolated values to be continous across element boundaries).
virtual void get_interpolated_values (const Vector< double > &s, Vector< double > &values)=0
 Get all continously interpolated function values in this element as a Vector. Note: Vector sets is own size to ensure that that this function can be used in black-box fashion.
virtual void get_interpolated_values (const unsigned &t, const Vector< double > &s, Vector< double > &values)=0
 Get all continously interpolated function values at previous timestep in this element as a Vector. (t=0: present; t>0: prev. timestep) Note: Vector sets is own size to ensure that that this function can be used in black-box fashion.
virtual Nodeinterpolating_node_pt (const unsigned &n, const int &value_id)
 In mixed elements, different sets of nodes are used to interpolate different unknowns. This function returns the n-th node that interpolates the value_id-th unknown. Default implementation is that all variables use the positional nodes, i.e. isoparametric elements. Note that any overloaded versions of this function MUST provide a set of nodes for the position, which always has the value_id -1.
virtual double local_one_d_fraction_of_interpolating_node (const unsigned &n1d, const unsigned &i, const int &value_id)
 Return the local one dimensional fraction of the n1d-th node in the direction of the local coordinate s[i] that is used to interpolate the value_id-th continuously interpolated unknown. Default assumes isoparametric interpolation for all unknowns.
virtual Nodeget_interpolating_node_at_local_coordinate (const Vector< double > &s, const int &value_id)
 Return a pointer to the node that interpolates the value-id-th unknown at local coordinate s. If there is not a node at that position, then return 0.
virtual unsigned ninterpolating_node (const int &value_id)
 Return the number of nodes that are used to interpolate the value_id-th unknown. Default is to assume isoparametric elements.
virtual unsigned ninterpolating_node_1d (const int &value_id)
 Return the number of nodes in a one_d direction that are used to interpolate the value_id-th unknown. Default is to assume an isoparametric mapping.
virtual void interpolating_basis (const Vector< double > &s, Shape &psi, const int &value_id) const
 Return the basis functions that are used to interpolate the value_id-th unknown. By default assume isoparameteric interpolation.
virtual void check_integrity (double &max_error)=0
 Check the integrity of the element: Continuity of positions values, etc. Essentially, check that the approximation of the functions is consistent when viewed from both sides of the element boundaries Must be overloaded for each different geometric element.
void identify_field_data_for_interactions (std::set< std::pair< Data *, unsigned > > &paired_field_data)
 The purpose of this function is to identify all possible.
void assign_nodal_local_eqn_numbers ()
 Overload the function that assigns local equation numbers for the Data stored at the nodes so that hanging data is taken into account.
virtual RefineableElementroot_element_pt ()
 Setup the local equation numbering schemes: Pointer to the root element in refinement hierarchy (must be implemented in specific elements that do refinement via tree-like refinement structure. Here we provide a default implementation that is appropriate for cases where tree-like refinement doesn't exist or if the element doesn't have root in that tree (i.e. if it's a root itself): We return "this".
virtual RefineableElementfather_element_pt () const
 Return a pointer to the father element.
void get_father_at_refinement_level (unsigned &refinement_level, RefineableElement *&father_at_reflevel_pt)
 Return a pointer to the "father" element at the specified refinement level.
virtual void further_build ()
 Further build: e.g. deal with interpolation of internal values.
virtual void setup_hanging_nodes (Vector< std::ofstream * > &output_stream)
 Mark up any hanging nodes that arise as a result of non-uniform refinement. Any hanging nodes will be documented in files addressed by the streams in the vector output_stream, if the streams are open.
virtual void further_setup_hanging_nodes ()
 Perform additional hanging node procedures for variables that are not interpolated by all nodes (e.g. lower order interpolations for the pressure in Taylor Hood).
void get_dresidual_dnodal_coordinates (RankThreeTensor< double > &dresidual_dnodal_coordinates)
 Compute derivatives of elemental residual vector with respect to nodal coordinates. Default implementation by FD can be overwritten for specific elements. dresidual_dnodal_coordinates(l,i,j) = d res(l) / dX_{ij} This version is overloaded from the version in FiniteElement and takes hanging nodes into account -- j in the above loop loops over all the nodes that actively control the shape of the element (i.e. they are non-hanging or master nodes of hanging nodes in this element).
unsigned nshape_controlling_nodes ()
 Number of shape-controlling nodes = the number of non-hanging nodes plus the number of master nodes associated with hanging nodes.
std::map< Node *, unsigned > shape_controlling_node_lookup ()
 Return lookup scheme for unique number associated with any of the nodes that actively control the shape of the element (i.e. they are either non-hanging nodes of this element or master nodes of hanging nodes.

Static Public Member Functions

static double max_integrity_tolerance ()
 Max. allowed discrepancy in element integrity check.

Protected Member Functions

void assemble_local_to_eulerian_jacobian (const DShape &dpsids, DenseMatrix< double > &jacobian) const
 Assemble the jacobian matrix for the mapping from local.
void assemble_local_to_eulerian_jacobian2 (const DShape &d2psids, DenseMatrix< double > &jacobian2) const
 Assemble the the "jacobian" matrix of second derivatives of the mapping from local to Eulerian coordinates, given the second derivatives of the shape functions w.r.t. local coordinates. Overload the standard version to use the hanging information for the Eulerian coordinates.
void assemble_eulerian_base_vectors (const DShape &dpsids, DenseMatrix< double > &interpolated_G) const
 Assemble the covariant Eulerian base vectors, assuming that the derivatives of the shape functions with respect to the local coordinates have already been constructed. Overload the standard version to account for hanging nodes.
double local_to_eulerian_mapping_diagonal (const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 Calculate the mapping from local to Eulerian coordinates given the derivatives of the shape functions w.r.t the local coordinates. assuming that the coordinates are aligned in the direction of the local coordinates, i.e. there are no cross terms and the jacobian is diagonal. This funciton returns the determinant of the jacobian, the jacobian and the inverse jacobian. Overload the standard version to take hanging info into account.
int local_hang_eqn (Node *const &node_pt, const unsigned &i)
 Access function that returns the local equation number for the hanging node variables (values stored at master nodes). The local equation number corresponds to the i-th unknown stored at the node addressed by node_pt.
void assign_hanging_local_eqn_numbers ()
 Assign the local hang eqn.
virtual void fill_in_jacobian_from_nodal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 Calculate the contributions to the jacobian from the nodal degrees of freedom using finite differences. This version is overloaded to take hanging node information into account.

Static Protected Member Functions

static void check_value_id (const int &n_continuously_interpolated_values, const int &value_id)
 Static helper function that is used to check that the value_id is in range.

Protected Attributes

TreeTree_pt
 A pointer to a general tree object.
unsigned Refine_level
 Refinement level.
bool To_be_refined
 Flag for refinement.
bool Sons_to_be_unrefined
 Flag for unrefinement.
long Number
 Global element number -- for plotting/validation purposes.

Static Protected Attributes

static double Max_integrity_tolerance = 1.0e-8
 Max. allowed discrepancy in element integrity check.

Private Attributes

std::map< Node *, int > * Local_hang_eqn
 Storage for local equation numbers of hanging node variables (values stored at master nodes). It is essential that these are indexed by a Node pointer because the Node may be internal or external to the element. local equation number = Local_hang_eqn(master_node_pt,ival).
std::map< Node *, unsigned > Shape_controlling_node_lookup
 Lookup scheme for unique number associated with any of the nodes that actively control the shape of the element (i.e. they are either non-hanging nodes of this element or master nodes of hanging nodes.

Constructor & Destructor Documentation

oomph::RefineableElement::RefineableElement (  )  [inline]

Constructor, calls the FiniteElement constructor and initialises the member data.

Definition at line 214 of file refineable_elements.h.

virtual oomph::RefineableElement::~RefineableElement (  )  [inline, virtual]

Destructor, delete the allocated storage for the hanging equations

Definition at line 221 of file refineable_elements.h.

References Local_hang_eqn.

oomph::RefineableElement::RefineableElement ( const RefineableElement  )  [inline]

Broken copy constructor.

Definition at line 224 of file refineable_elements.h.

References oomph::BrokenCopy::broken_copy().


Member Function Documentation

void oomph::RefineableElement::assemble_eulerian_base_vectors ( const DShape dpsids,
DenseMatrix< double > &  interpolated_G 
) const [protected, virtual]

Assemble the covariant Eulerian base vectors, assuming that the derivatives of the shape functions with respect to the local coordinates have already been constructed. Overload the standard version to account for hanging nodes.

Assemble the covariant Eulerian base vectors and return them in the matrix interpolated_G. The derivatives of the shape functions with respect to the local coordinate should already have been calculated before calling this function

Reimplemented from oomph::FiniteElement.

Definition at line 179 of file refineable_elements.cc.

References oomph::FiniteElement::dim(), oomph::FiniteElement::nnodal_position_type(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_dimension(), and oomph::FiniteElement::nodal_position_gen().

void oomph::RefineableElement::assemble_local_to_eulerian_jacobian ( const DShape dpsids,
DenseMatrix< double > &  jacobian 
) const [protected, virtual]

Assemble the jacobian matrix for the mapping from local.

to Eulerian coordinates, given the derivatives of the shape function w.r.t the local coordinates. Overload the standard version to use the hanging information for the Eulerian coordinates.

Reimplemented from oomph::FiniteElement.

Definition at line 75 of file refineable_elements.cc.

References oomph::FiniteElement::dim(), oomph::FiniteElement::nnodal_position_type(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_dimension(), oomph::FiniteElement::nodal_position_gen(), and OOMPH_EXCEPTION_LOCATION.

void oomph::RefineableElement::assemble_local_to_eulerian_jacobian2 ( const DShape d2psids,
DenseMatrix< double > &  jacobian2 
) const [protected, virtual]

Assemble the the "jacobian" matrix of second derivatives of the mapping from local to Eulerian coordinates, given the second derivatives of the shape functions w.r.t. local coordinates. Overload the standard version to use the hanging information for the Eulerian coordinates.

Internal function that is used to assemble the jacobian of second derivatives of the the mapping from local coordinates (s) to the eulerian coordinates (x), given the second derivatives of the shape functions.

Reimplemented from oomph::FiniteElement.

Definition at line 135 of file refineable_elements.cc.

References oomph::FiniteElement::dim(), oomph::FiniteElement::N2deriv, oomph::FiniteElement::nnodal_position_type(), oomph::FiniteElement::nnode(), and oomph::FiniteElement::nodal_position_gen().

void oomph::RefineableElement::assign_hanging_local_eqn_numbers (  )  [protected]

Assign the local hang eqn.

Definition at line 306 of file refineable_elements.cc.

References oomph::GeneralisedElement::add_global_eqn_numbers(), oomph::GeneralisedElement::eqn_number(), oomph::Node::hanging_pt(), oomph::Data::Is_pinned, oomph::GeneralisedElement::local_eqn_number(), Local_hang_eqn, ncont_interpolated_values(), oomph::GeneralisedElement::ndof(), oomph::HangInfo::nmaster(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_local_eqn(), oomph::FiniteElement::node_pt(), oomph::Data::nvalue(), OOMPH_EXCEPTION_LOCATION, and Shape_controlling_node_lookup.

Referenced by assign_nodal_local_eqn_numbers().

void oomph::RefineableElement::assign_nodal_local_eqn_numbers (  )  [inline, virtual]

Overload the function that assigns local equation numbers for the Data stored at the nodes so that hanging data is taken into account.

Reimplemented from oomph::FiniteElement.

Definition at line 433 of file refineable_elements.h.

References assign_hanging_local_eqn_numbers(), and oomph::FiniteElement::assign_nodal_local_eqn_numbers().

virtual void oomph::RefineableElement::build ( Mesh *&  mesh_pt,
Vector< Node * > &  new_node_pt,
bool &  was_already_built,
std::ofstream &  new_nodes_file 
) [pure virtual]

Interface to function that builds the element: i.e. construct the nodes, assign their positions,.

required procedures depend on the geometrical type of the element and must be implemented in specific refineable elements. Any new nodes created during the build process are returned in the vector new_node_pt.

Implemented in oomph::RefineableQElement< 3 >, oomph::RefineableSolidQElement< 3 >, oomph::NonRefineableElementWithHangingNodes, oomph::RefineableQElement< 2 >, and oomph::RefineableSolidQElement< 2 >.

virtual void oomph::RefineableElement::check_integrity ( double &  max_error  )  [pure virtual]

Check the integrity of the element: Continuity of positions values, etc. Essentially, check that the approximation of the functions is consistent when viewed from both sides of the element boundaries Must be overloaded for each different geometric element.

Implemented in oomph::RefineableQElement< 3 >, oomph::NonRefineableElementWithHangingNodes, and oomph::RefineableQElement< 2 >.

void oomph::RefineableElement::check_value_id ( const int &  n_continuously_interpolated_values,
const int &  value_id 
) [static, protected]

Static helper function that is used to check that the value_id is in range.

Helper function that is used to check that the value_id is in the range allowed by the element. The number of continuously interpolated values and the value_id must be passed as arguments.

Definition at line 47 of file refineable_elements.cc.

References OOMPH_EXCEPTION_LOCATION.

Referenced by oomph::RefineableQPVDElementWithContinuousPressure< DIM >::get_interpolating_node_at_local_coordinate(), oomph::RefineableQPVDElementWithContinuousPressure< DIM >::interpolating_basis(), oomph::RefineableQPVDElementWithContinuousPressure< DIM >::interpolating_node_pt(), oomph::RefineableQPVDElementWithContinuousPressure< DIM >::local_one_d_fraction_of_interpolating_node(), oomph::RefineableQPVDElementWithContinuousPressure< DIM >::ninterpolating_node(), and oomph::RefineableQPVDElementWithContinuousPressure< DIM >::ninterpolating_node_1d().

void oomph::RefineableElement::deactivate_element (  )  [virtual]

Final operations that must be performed when the element is no longer active in the mesh, but still resident in the QuadTree.

Deactivate the element by marking setting all local pointers to obsolete nodes to zero

Definition at line 289 of file refineable_elements.cc.

References oomph::FiniteElement::nnode(), and oomph::FiniteElement::node_pt().

Referenced by oomph::Tree::deactivate_object().

void oomph::RefineableElement::deselect_for_refinement (  )  [inline]

Deselect the element for refinement.

Definition at line 291 of file refineable_elements.h.

References To_be_refined.

void oomph::RefineableElement::deselect_sons_for_unrefinement (  )  [inline]

No unrefinement will be performed by merging the four sons of this element.

Definition at line 297 of file refineable_elements.h.

References Sons_to_be_unrefined.

Referenced by oomph::Tree::merge_sons_if_required().

virtual RefineableElement* oomph::RefineableElement::father_element_pt (  )  const [inline, virtual]

Return a pointer to the father element.

Definition at line 464 of file refineable_elements.h.

References oomph::Tree::father_pt(), oomph::Tree::object_pt(), and Tree_pt.

Referenced by oomph::RefineableUnsteadyHeatEquations< DIM >::further_build(), oomph::RefineablePVDEquationsWithPressure< DIM >::further_build(), oomph::RefineablePVDEquations< DIM >::further_build(), oomph::RefineablePoissonEquations< DIM >::further_build(), oomph::RefineableNavierStokesEquations< DIM >::further_build(), oomph::RefineableLinearWaveEquations< DIM >::further_build(), oomph::RefineableLinearElasticityEquations< DIM >::further_build(), oomph::RefineableGeneralisedAdvectionDiffusionEquations< DIM >::further_build(), oomph::RefineableSolidElement::further_build(), oomph::RefineableAxisymmetricNavierStokesEquations::further_build(), oomph::RefineableAdvectionDiffusionReactionEquations< NREAGENT, DIM >::further_build(), oomph::RefineableAdvectionDiffusionEquations< DIM >::further_build(), and unbuild().

void oomph::RefineableElement::fill_in_jacobian_from_nodal_by_fd ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian 
) [protected, virtual]

Calculate the contributions to the jacobian from the nodal degrees of freedom using finite differences. This version is overloaded to take hanging node information into account.

This function calculates the entries of Jacobian matrix, used in the Newton method, associated with the nodal degrees of freedom. This is done by finite differences to handle the general case Overload the standard case to include hanging node case

Reimplemented from oomph::FiniteElement.

Definition at line 660 of file refineable_elements.cc.

References oomph::GeneralisedElement::Default_fd_jacobian_step, oomph::GeneralisedElement::get_residuals(), oomph::Node::hanging_pt(), oomph::Node::is_hanging(), local_hang_eqn(), oomph::HangInfo::master_node_pt(), oomph::GeneralisedElement::ndof(), oomph::HangInfo::nmaster(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_local_eqn(), oomph::FiniteElement::node_pt(), oomph::FiniteElement::reset_after_nodal_fd(), oomph::FiniteElement::reset_in_nodal_fd(), oomph::FiniteElement::update_before_nodal_fd(), oomph::FiniteElement::update_in_nodal_fd(), and oomph::Data::value_pt().

virtual void oomph::RefineableElement::further_build (  )  [inline, virtual]

Further build: e.g. deal with interpolation of internal values.

Reimplemented in oomph::RefineableNavierStokesEquations< DIM >, oomph::RefineableQCrouzeixRaviartElement< DIM >, oomph::RefineableSolidElement, oomph::RefineablePoissonEquations< DIM >, oomph::RefineableAxisymmetricNavierStokesEquations, oomph::RefineableAxisymmetricQCrouzeixRaviartElement, oomph::RefineablePVDEquations< DIM >, oomph::RefineablePVDEquationsWithPressure< DIM >, oomph::RefineableQPVDElementWithPressure< DIM >, oomph::RefineableAdvectionDiffusionEquations< DIM >, oomph::RefineableGeneralisedAdvectionDiffusionEquations< DIM >, oomph::RefineableUnsteadyHeatEquations< DIM >, oomph::RefineableLinearWaveEquations< DIM >, oomph::RefineableLinearElasticityEquations< DIM >, oomph::RefineableAdvectionDiffusionReactionEquations< NREAGENT, DIM >, oomph::RefineableQCrouzeixRaviartElement< DIM >, oomph::RefineableQCrouzeixRaviartElement< DIM >, oomph::RefineableQPVDElementWithPressure< DIM >, and oomph::RefineableQPVDElementWithPressure< DIM >.

Definition at line 503 of file refineable_elements.h.

Referenced by oomph::RefineableSolidElement::further_build().

virtual void oomph::RefineableElement::further_setup_hanging_nodes (  )  [inline, virtual]

Perform additional hanging node procedures for variables that are not interpolated by all nodes (e.g. lower order interpolations for the pressure in Taylor Hood).

Reimplemented in oomph::RefineableQTaylorHoodElement< DIM >, oomph::RefineableQCrouzeixRaviartElement< DIM >, oomph::RefineableQElement< 3 >, oomph::RefineableQElement< 2 >, oomph::RefineableQPoissonElement< DIM, NNODE_1D >, oomph::RefineableQSpectralPoissonElement< DIM, NNODE_1D >, oomph::RefineableAxisymmetricQTaylorHoodElement, oomph::RefineableAxisymmetricQCrouzeixRaviartElement, oomph::RefineableQPVDElement< DIM, NNODE_1D >, oomph::RefineableQPVDElementWithPressure< DIM >, oomph::RefineableQPVDElementWithContinuousPressure< DIM >, oomph::RefineableQAdvectionDiffusionElement< DIM, NNODE_1D >, oomph::RefineableQGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D >, oomph::RefineableQSUPGAdvectionDiffusionElement< DIM, NNODE_1D >, oomph::RefineableQUnsteadyHeatElement< DIM, NNODE_1D >, oomph::RefineableQLinearWaveElement< DIM, NNODE_1D >, oomph::RefineableQLinearElasticityElement< DIM, NNODE_1D >, and oomph::RefineableQAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D >.

Definition at line 513 of file refineable_elements.h.

void oomph::RefineableElement::get_dresidual_dnodal_coordinates ( RankThreeTensor< double > &  dresidual_dnodal_coordinates  )  [virtual]

Compute derivatives of elemental residual vector with respect to nodal coordinates. Default implementation by FD can be overwritten for specific elements. dresidual_dnodal_coordinates(l,i,j) = d res(l) / dX_{ij} This version is overloaded from the version in FiniteElement and takes hanging nodes into account -- j in the above loop loops over all the nodes that actively control the shape of the element (i.e. they are non-hanging or master nodes of hanging nodes in this element).

Compute derivatives of elemental residual vector with respect to nodal coordinates. Default implementation by FD can be overwritten for specific elements. dresidual_dnodal_coordinates(l,i,j) = d res(l) / dX_{ij} This version is overloaded from the version in FiniteElement and takes hanging nodes into account -- j in the above loop loops over all the nodes that actively control the shape of the element (i.e. they are non-hanging or master nodes of hanging nodes in this element).

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::RefineableNavierStokesEquations< DIM >, and oomph::RefineablePoissonEquations< DIM >.

Definition at line 579 of file refineable_elements.cc.

References oomph::GeneralisedElement::Default_fd_jacobian_step, oomph::GeneralisedElement::get_residuals(), oomph::Node::ndim(), oomph::GeneralisedElement::ndof(), oomph::FiniteElement::nnode(), oomph::FiniteElement::node_pt(), and Shape_controlling_node_lookup.

void oomph::RefineableElement::get_father_at_refinement_level ( unsigned &  refinement_level,
RefineableElement *&  father_at_reflevel_pt 
) [inline]

Return a pointer to the "father" element at the specified refinement level.

Definition at line 479 of file refineable_elements.h.

References oomph::Tree::father_pt(), get_father_at_refinement_level(), oomph::Tree::object_pt(), refinement_level(), and Tree_pt.

Referenced by get_father_at_refinement_level().

virtual void oomph::RefineableElement::get_interpolated_values ( const unsigned &  t,
const Vector< double > &  s,
Vector< double > &  values 
) [pure virtual]

Get all continously interpolated function values at previous timestep in this element as a Vector. (t=0: present; t>0: prev. timestep) Note: Vector sets is own size to ensure that that this function can be used in black-box fashion.

Implemented in oomph::RefineableQTaylorHoodElement< DIM >, oomph::RefineableQCrouzeixRaviartElement< DIM >, oomph::NonRefineableElementWithHangingNodes, oomph::RefineablePoissonEquations< DIM >, oomph::RefineableAxisymmetricQTaylorHoodElement, oomph::RefineableAxisymmetricQCrouzeixRaviartElement, oomph::RefineablePVDEquations< DIM >, oomph::RefineablePVDEquationsWithPressure< DIM >, oomph::RefineableQPVDElementWithContinuousPressure< DIM >, oomph::RefineableAdvectionDiffusionEquations< DIM >, oomph::RefineableGeneralisedAdvectionDiffusionEquations< DIM >, oomph::RefineableUnsteadyHeatEquations< DIM >, oomph::RefineableLinearWaveEquations< DIM >, oomph::RefineableLinearElasticityEquations< DIM >, and oomph::RefineableAdvectionDiffusionReactionEquations< NREAGENT, DIM >.

virtual void oomph::RefineableElement::get_interpolated_values ( const Vector< double > &  s,
Vector< double > &  values 
) [pure virtual]

Get all continously interpolated function values in this element as a Vector. Note: Vector sets is own size to ensure that that this function can be used in black-box fashion.

Implemented in oomph::RefineableQTaylorHoodElement< DIM >, oomph::RefineableQCrouzeixRaviartElement< DIM >, oomph::NonRefineableElementWithHangingNodes, oomph::RefineablePoissonEquations< DIM >, oomph::RefineableAxisymmetricQTaylorHoodElement, oomph::RefineableAxisymmetricQCrouzeixRaviartElement, oomph::RefineablePVDEquations< DIM >, oomph::RefineablePVDEquationsWithPressure< DIM >, oomph::RefineableQPVDElementWithContinuousPressure< DIM >, oomph::RefineableAdvectionDiffusionEquations< DIM >, oomph::RefineableGeneralisedAdvectionDiffusionEquations< DIM >, oomph::RefineableUnsteadyHeatEquations< DIM >, oomph::RefineableLinearWaveEquations< DIM >, oomph::RefineableLinearElasticityEquations< DIM >, and oomph::RefineableAdvectionDiffusionReactionEquations< NREAGENT, DIM >.

virtual Node* oomph::RefineableElement::get_interpolating_node_at_local_coordinate ( const Vector< double > &  s,
const int &  value_id 
) [inline, virtual]

Return a pointer to the node that interpolates the value-id-th unknown at local coordinate s. If there is not a node at that position, then return 0.

Reimplemented in oomph::RefineableQTaylorHoodElement< DIM >, oomph::RefineableAxisymmetricQTaylorHoodElement, and oomph::RefineableQPVDElementWithContinuousPressure< DIM >.

Definition at line 388 of file refineable_elements.h.

References oomph::FiniteElement::get_node_at_local_coordinate().

Referenced by oomph::RefineableQElement< 3 >::oc_hang_helper(), and oomph::RefineableQElement< 2 >::quad_hang_helper().

void oomph::RefineableElement::identify_field_data_for_interactions ( std::set< std::pair< Data *, unsigned > > &  paired_field_data  )  [virtual]

The purpose of this function is to identify all possible.

Data that can affect the fields interpolated by the FiniteElement. This must be overloaded to include data from any hanging nodes correctly

Reimplemented from oomph::FiniteElement.

Definition at line 504 of file refineable_elements.cc.

References oomph::GeneralisedElement::ninternal_data(), and oomph::FiniteElement::nnode().

virtual void oomph::RefineableElement::interpolating_basis ( const Vector< double > &  s,
Shape psi,
const int &  value_id 
) const [inline, virtual]

Return the basis functions that are used to interpolate the value_id-th unknown. By default assume isoparameteric interpolation.

Reimplemented in oomph::RefineableQTaylorHoodElement< DIM >, oomph::RefineableAxisymmetricQTaylorHoodElement, and oomph::RefineableQPVDElementWithContinuousPressure< DIM >.

Definition at line 406 of file refineable_elements.h.

References oomph::FiniteElement::shape().

Referenced by oomph::RefineableQElement< 3 >::oc_hang_helper(), and oomph::RefineableQElement< 2 >::quad_hang_helper().

virtual Node* oomph::RefineableElement::interpolating_node_pt ( const unsigned &  n,
const int &  value_id 
) [inline, virtual]

In mixed elements, different sets of nodes are used to interpolate different unknowns. This function returns the n-th node that interpolates the value_id-th unknown. Default implementation is that all variables use the positional nodes, i.e. isoparametric elements. Note that any overloaded versions of this function MUST provide a set of nodes for the position, which always has the value_id -1.

Reimplemented in oomph::RefineableQTaylorHoodElement< DIM >, oomph::RefineableAxisymmetricQTaylorHoodElement, and oomph::RefineableQPVDElementWithContinuousPressure< DIM >.

Definition at line 371 of file refineable_elements.h.

References oomph::FiniteElement::node_pt().

Referenced by oomph::RefineableQElement< 3 >::oc_hang_helper(), and oomph::RefineableQElement< 2 >::quad_hang_helper().

int oomph::RefineableElement::local_hang_eqn ( Node *const &  node_pt,
const unsigned &  i 
) [inline, protected]

Access function that returns the local equation number for the hanging node variables (values stored at master nodes). The local equation number corresponds to the i-th unknown stored at the node addressed by node_pt.

Definition at line 181 of file refineable_elements.h.

References Local_hang_eqn, ncont_interpolated_values(), oomph::FiniteElement::node_pt(), and OOMPH_EXCEPTION_LOCATION.

Referenced by oomph::RefineableAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::RefineablePVDEquationsWithPressure< DIM >::fill_in_generic_residual_contribution_pvd_with_pressure(), and fill_in_jacobian_from_nodal_by_fd().

virtual double oomph::RefineableElement::local_one_d_fraction_of_interpolating_node ( const unsigned &  n1d,
const unsigned &  i,
const int &  value_id 
) [inline, virtual]

Return the local one dimensional fraction of the n1d-th node in the direction of the local coordinate s[i] that is used to interpolate the value_id-th continuously interpolated unknown. Default assumes isoparametric interpolation for all unknowns.

Reimplemented in oomph::RefineableQTaylorHoodElement< DIM >, oomph::RefineableAxisymmetricQTaylorHoodElement, and oomph::RefineableQPVDElementWithContinuousPressure< DIM >.

Definition at line 380 of file refineable_elements.h.

References oomph::FiniteElement::local_one_d_fraction_of_node().

double oomph::RefineableElement::local_to_eulerian_mapping_diagonal ( const DShape dpsids,
DenseMatrix< double > &  jacobian,
DenseMatrix< double > &  inverse_jacobian 
) const [protected, virtual]

Calculate the mapping from local to Eulerian coordinates given the derivatives of the shape functions w.r.t the local coordinates. assuming that the coordinates are aligned in the direction of the local coordinates, i.e. there are no cross terms and the jacobian is diagonal. This funciton returns the determinant of the jacobian, the jacobian and the inverse jacobian. Overload the standard version to take hanging info into account.

Calculate the mapping from local to eulerian coordinates assuming that the coordinates are aligned in the direction of the local coordinates, i.e. there are no cross terms and the jacobian is diagonal. The local derivatives are passed as dpsids and the jacobian and inverse jacobian are returned.

Reimplemented from oomph::FiniteElement.

Definition at line 217 of file refineable_elements.cc.

References oomph::FiniteElement::check_jacobian(), oomph::FiniteElement::dim(), oomph::FiniteElement::nnodal_position_type(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_dimension(), oomph::FiniteElement::nodal_position_gen(), and OOMPH_EXCEPTION_LOCATION.

static double oomph::RefineableElement::max_integrity_tolerance (  )  [inline, static]

Max. allowed discrepancy in element integrity check.

Definition at line 419 of file refineable_elements.h.

References Max_integrity_tolerance.

Referenced by oomph::RefineableMeshBase::adapt_mesh().

virtual unsigned oomph::RefineableElement::ncont_interpolated_values (  )  const [pure virtual]

Number of continuously interpolated values. Note: We assume that they are located at the beginning of the value_pt Vector! (Used for interpolation to son elements, for integrity check and post-processing -- we can only expect the continously interpolated values to be continous across element boundaries).

Implemented in oomph::RefineableNavierStokesFluxControlElement< ELEMENT >, oomph::RefineableQTaylorHoodElement< DIM >, oomph::RefineableQCrouzeixRaviartElement< DIM >, oomph::RefineableQPoissonElement< DIM, NNODE_1D >, oomph::RefineableQSpectralPoissonElement< DIM, NNODE_1D >, oomph::RefineableAxisymmetricQTaylorHoodElement, oomph::RefineableAxisymmetricQCrouzeixRaviartElement, oomph::RefineablePVDEquations< DIM >, oomph::RefineablePVDEquationsWithPressure< DIM >, oomph::RefineableQPVDElementWithPressure< DIM >, oomph::RefineableQPVDElementWithContinuousPressure< DIM >, oomph::RefineableQAdvectionDiffusionElement< DIM, NNODE_1D >, oomph::RefineableQGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D >, oomph::RefineableQSUPGAdvectionDiffusionElement< DIM, NNODE_1D >, oomph::RefineableQUnsteadyHeatElement< DIM, NNODE_1D >, oomph::RefineableQLinearWaveElement< DIM, NNODE_1D >, oomph::RefineableLinearElasticityEquations< DIM >, and oomph::RefineableQAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D >.

Referenced by assign_hanging_local_eqn_numbers(), oomph::RefineableQElement< 2 >::check_integrity(), and local_hang_eqn().

virtual unsigned oomph::RefineableElement::ninterpolating_node ( const int &  value_id  )  [inline, virtual]

Return the number of nodes that are used to interpolate the value_id-th unknown. Default is to assume isoparametric elements.

Reimplemented in oomph::RefineableQTaylorHoodElement< DIM >, oomph::RefineableAxisymmetricQTaylorHoodElement, and oomph::RefineableQPVDElementWithContinuousPressure< DIM >.

Definition at line 396 of file refineable_elements.h.

References oomph::FiniteElement::nnode().

Referenced by oomph::RefineableQElement< 3 >::oc_hang_helper(), and oomph::RefineableQElement< 2 >::quad_hang_helper().

virtual unsigned oomph::RefineableElement::ninterpolating_node_1d ( const int &  value_id  )  [inline, virtual]

Return the number of nodes in a one_d direction that are used to interpolate the value_id-th unknown. Default is to assume an isoparametric mapping.

Reimplemented in oomph::RefineableQTaylorHoodElement< DIM >, oomph::RefineableAxisymmetricQTaylorHoodElement, and oomph::RefineableQPVDElementWithContinuousPressure< DIM >.

Definition at line 401 of file refineable_elements.h.

References oomph::FiniteElement::nnode_1d().

virtual bool oomph::RefineableElement::nodes_built (  )  [inline, virtual]

Return true if all the nodes have been built, false if not.

Reimplemented in oomph::RefineableQSpectralElement< 3 >, and oomph::RefineableQSpectralElement< 2 >.

Definition at line 335 of file refineable_elements.h.

References oomph::FiniteElement::node_pt().

Referenced by oomph::RefineableQElement< 2 >::check_integrity(), oomph::RefineableQElement< 3 >::check_integrity(), oomph::OcTree::doc_face_neighbours(), oomph::QuadTree::doc_neighbours(), and oomph::OcTree::doc_true_edge_neighbours().

unsigned oomph::RefineableElement::nshape_controlling_nodes (  )  [inline]

Number of shape-controlling nodes = the number of non-hanging nodes plus the number of master nodes associated with hanging nodes.

Definition at line 531 of file refineable_elements.h.

References Shape_controlling_node_lookup.

Referenced by oomph::ElementWithMovingNodes::fill_in_jacobian_from_geometric_data(), oomph::ElementWithMovingNodes::get_dnodal_coordinates_dgeom_dofs(), oomph::RefineablePoissonEquations< DIM >::get_dresidual_dnodal_coordinates(), and oomph::RefineableNavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates().

long oomph::RefineableElement::number (  )  const [inline]

Element number (for debugging/plotting).

Definition at line 338 of file refineable_elements.h.

References Number.

Referenced by oomph::OcTree::doc_face_neighbours(), oomph::QuadTree::doc_neighbours(), and oomph::OcTree::doc_true_edge_neighbours().

void oomph::RefineableElement::operator= ( const RefineableElement  )  [inline]

Broken assignment operator.

Definition at line 230 of file refineable_elements.h.

References oomph::BrokenCopy::broken_assign().

virtual void oomph::RefineableElement::rebuild_from_sons ( Mesh *&  mesh_pt  )  [pure virtual]

Rebuild the element, e.g. set internal values in line with those of the sons that have now merged.

Implemented in oomph::RefineableQTaylorHoodElement< DIM >, oomph::RefineableQCrouzeixRaviartElement< DIM >, oomph::RefineableQSpectralElement< 3 >, oomph::NonRefineableElementWithHangingNodes, oomph::RefineableQSpectralElement< 2 >, oomph::RefineableQPoissonElement< DIM, NNODE_1D >, oomph::RefineableAxisymmetricQTaylorHoodElement, oomph::RefineableAxisymmetricQCrouzeixRaviartElement, oomph::RefineableQPVDElement< DIM, NNODE_1D >, oomph::RefineableQPVDElementWithPressure< DIM >, oomph::RefineableQPVDElementWithContinuousPressure< DIM >, oomph::RefineableQAdvectionDiffusionElement< DIM, NNODE_1D >, oomph::RefineableQGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D >, oomph::RefineableQSUPGAdvectionDiffusionElement< DIM, NNODE_1D >, oomph::RefineableQUnsteadyHeatElement< DIM, NNODE_1D >, oomph::RefineableQLinearWaveElement< DIM, NNODE_1D >, oomph::RefineableQLinearElasticityElement< DIM, NNODE_1D >, oomph::RefineableQAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D >, oomph::RefineableQCrouzeixRaviartElement< DIM >, oomph::RefineableQCrouzeixRaviartElement< DIM >, oomph::RefineableQPVDElementWithPressure< DIM >, and oomph::RefineableQPVDElementWithPressure< DIM >.

Referenced by oomph::Tree::merge_sons_if_required().

unsigned oomph::RefineableElement::refinement_level (  )  const [inline]

Return the Refinement level.

Definition at line 285 of file refineable_elements.h.

References Refine_level.

Referenced by get_father_at_refinement_level().

virtual unsigned oomph::RefineableElement::required_nsons (  )  const [inline, virtual]

Set the number of sons that can be constructed by the element The default is none.

Reimplemented in oomph::RefineableQElement< 3 >, and oomph::RefineableQElement< 2 >.

Definition at line 243 of file refineable_elements.h.

Referenced by split().

virtual RefineableElement* oomph::RefineableElement::root_element_pt (  )  [inline, virtual]

Setup the local equation numbering schemes: Pointer to the root element in refinement hierarchy (must be implemented in specific elements that do refinement via tree-like refinement structure. Here we provide a default implementation that is appropriate for cases where tree-like refinement doesn't exist or if the element doesn't have root in that tree (i.e. if it's a root itself): We return "this".

Definition at line 455 of file refineable_elements.h.

References oomph::Tree::object_pt(), oomph::Tree::root_pt(), and Tree_pt.

void oomph::RefineableElement::select_for_refinement (  )  [inline]

Select the element for refinement.

Definition at line 288 of file refineable_elements.h.

References To_be_refined.

void oomph::RefineableElement::select_sons_for_unrefinement (  )  [inline]

Unrefinement will be performed by merging the four sons of this element.

Definition at line 294 of file refineable_elements.h.

References Sons_to_be_unrefined.

void oomph::RefineableElement::set_number ( const long &  mynumber  )  [inline]

Set element number (for debugging/plotting).

Definition at line 341 of file refineable_elements.h.

References Number.

void oomph::RefineableElement::set_refinement_level ( const int &  refine_level  )  [inline]

Set the refinement level.

Definition at line 281 of file refineable_elements.h.

References Refine_level.

void oomph::RefineableElement::set_tree_pt ( Tree my_tree_pt  )  [inline]

Set pointer to quadtree representation of this element.

Definition at line 239 of file refineable_elements.h.

References Tree_pt.

Referenced by oomph::Tree::Tree().

virtual void oomph::RefineableElement::setup_hanging_nodes ( Vector< std::ofstream * > &  output_stream  )  [inline, virtual]

Mark up any hanging nodes that arise as a result of non-uniform refinement. Any hanging nodes will be documented in files addressed by the streams in the vector output_stream, if the streams are open.

Reimplemented in oomph::RefineableQElement< 3 >, and oomph::RefineableQElement< 2 >.

Definition at line 508 of file refineable_elements.h.

std::map<Node*,unsigned> oomph::RefineableElement::shape_controlling_node_lookup (  )  [inline]

Return lookup scheme for unique number associated with any of the nodes that actively control the shape of the element (i.e. they are either non-hanging nodes of this element or master nodes of hanging nodes.

Definition at line 540 of file refineable_elements.h.

References Shape_controlling_node_lookup.

Referenced by oomph::ElementWithMovingNodes::get_dnodal_coordinates_dgeom_dofs(), oomph::RefineablePoissonEquations< DIM >::get_dresidual_dnodal_coordinates(), and oomph::RefineableNavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates().

bool oomph::RefineableElement::sons_to_be_unrefined (  )  [inline]

Has the element been selected for unrefinement?

Definition at line 303 of file refineable_elements.h.

References Sons_to_be_unrefined.

Referenced by oomph::Tree::merge_sons_if_required().

template<class ELEMENT>
void oomph::RefineableElement::split ( Vector< ELEMENT * > &  son_pt  )  const [inline]

Split the element into the number of sons to be constructed and return a vector of pointers to the sons. Elements are allocated, but they are not given any properties. The refinement level of the sons is one higher than that of the father elemern.

Definition at line 251 of file refineable_elements.h.

References Refine_level, and required_nsons().

Referenced by oomph::Tree::split_if_required().

bool oomph::RefineableElement::to_be_refined (  )  [inline]

Has the element been selected for refinement?

Definition at line 300 of file refineable_elements.h.

References To_be_refined.

Referenced by oomph::Tree::split_if_required().

Tree* oomph::RefineableElement::tree_pt (  )  [inline]

Access function: Pointer to quadtree representation of this element.

Definition at line 236 of file refineable_elements.h.

References Tree_pt.

virtual void oomph::RefineableElement::unbuild (  )  [inline, virtual]

Unbuild the element, i.e. mark the nodes that were created during its creation for possible deletion.

Definition at line 311 of file refineable_elements.h.

References father_element_pt(), oomph::FiniteElement::get_node_number(), oomph::FiniteElement::nnode(), oomph::FiniteElement::node_pt(), and oomph::Node::set_non_obsolete().


Member Data Documentation

std::map<Node*,int>* oomph::RefineableElement::Local_hang_eqn [private]

Storage for local equation numbers of hanging node variables (values stored at master nodes). It is essential that these are indexed by a Node pointer because the Node may be internal or external to the element. local equation number = Local_hang_eqn(master_node_pt,ival).

Definition at line 168 of file refineable_elements.h.

Referenced by assign_hanging_local_eqn_numbers(), local_hang_eqn(), and ~RefineableElement().

double oomph::RefineableElement::Max_integrity_tolerance = 1.0e-8 [static, protected]

Max. allowed discrepancy in element integrity check.

Definition at line 119 of file refineable_elements.h.

Referenced by max_integrity_tolerance().

long oomph::RefineableElement::Number [protected]

Global element number -- for plotting/validation purposes.

Definition at line 116 of file refineable_elements.h.

Referenced by number(), and set_number().

unsigned oomph::RefineableElement::Refine_level [protected]

Refinement level.

Definition at line 107 of file refineable_elements.h.

Referenced by refinement_level(), set_refinement_level(), and split().

std::map<Node*,unsigned> oomph::RefineableElement::Shape_controlling_node_lookup [private]

Lookup scheme for unique number associated with any of the nodes that actively control the shape of the element (i.e. they are either non-hanging nodes of this element or master nodes of hanging nodes.

Definition at line 173 of file refineable_elements.h.

Referenced by assign_hanging_local_eqn_numbers(), get_dresidual_dnodal_coordinates(), nshape_controlling_nodes(), and shape_controlling_node_lookup().

bool oomph::RefineableElement::Sons_to_be_unrefined [protected]

Flag for unrefinement.

Definition at line 113 of file refineable_elements.h.

Referenced by deselect_sons_for_unrefinement(), select_sons_for_unrefinement(), and sons_to_be_unrefined().

bool oomph::RefineableElement::To_be_refined [protected]

Flag for refinement.

Definition at line 110 of file refineable_elements.h.

Referenced by deselect_for_refinement(), select_for_refinement(), and to_be_refined().

Tree* oomph::RefineableElement::Tree_pt [protected]

A pointer to a general tree object.

Definition at line 104 of file refineable_elements.h.

Referenced by father_element_pt(), get_father_at_refinement_level(), root_element_pt(), set_tree_pt(), and tree_pt().


The documentation for this class was generated from the following files:
Generated on Mon Aug 10 11:25:51 2009 by  doxygen 1.4.7