//////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////// More...
#include <specific_node_update_interface_elements.h>
Public Member Functions | |
ElasticUpdateFluidInterfaceElement (FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0) | |
Constructor, pass a pointer to the bulk element and the face index of the bulk element to which the element is to be attached to. The optional identifier can be used to distinguish the additional nodal value (Lagr mult) created by this element from those created by other FaceElements. More... | |
double | zeta_nodal (const unsigned &n, const unsigned &k, const unsigned &i) const |
The "global" intrinsic coordinate of the element when viewed as part of a geometric object should be given by the FaceElement representation, by default. More... | |
double & | lagrange (const unsigned &n) |
Return the lagrange multiplier at local node n. More... | |
void | fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian) |
Fill in contribution to residuals and Jacobian. More... | |
void | output (std::ostream &outfile) |
Overload the output function. More... | |
void | output (std::ostream &outfile, const unsigned &n_plot) |
Output the element. More... | |
void | output (FILE *file_pt) |
Overload the C-style output function. More... | |
void | output (FILE *file_pt, const unsigned &n_plot) |
C-style Output function. More... | |
void | add_additional_residual_contributions_interface (Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag, const Shape &psif, const DShape &dpsifds, const DShape &dpsifdS, const DShape &dpsifdS_div, const Vector< double > &s, const Vector< double > &interpolated_x, const Vector< double > &interpolated_n, const double &W, const double &J) |
Helper function to calculate the additional contributions to be added at each integration point. This deals with Lagrange multiplier contribution, as well as any additional contributions by the other equations. More... | |
virtual FluidInterfaceBoundingElement * | make_bounding_element (const int &face_index) |
Create an "bounding" element (here actually a 2D line element of type ElasticLineFluidInterfaceBoundingElement<ELEMENT> that allows the application of a contact angle boundary condition on the the specified face. More... | |
Protected Member Functions | |
double | compute_surface_derivatives (const Shape &psi, const DShape &dpsids, const DenseMatrix< double > &interpolated_t, const Vector< double > &interpolated_x, DShape &surface_gradient, DShape &surface_divergence) |
Fill in the specific surface derivative calculations by calling the appropriate function from the derivative class. More... | |
Private Member Functions | |
unsigned | lagrange_index (const unsigned &n) |
Return the index at which the lagrange multiplier is stored at the n-th node. More... | |
int | kinematic_local_eqn (const unsigned &n) |
Equation number of the kinematic BC associated with node j. (This is the equation for the Lagrange multiplier) More... | |
void | hijack_kinematic_conditions (const Vector< unsigned > &bulk_node_number) |
Hijacking the kinematic condition corresponds to hijacking the variables associated with the Lagrange multipliers that are assigned on construction of this element. More... | |
Private Attributes | |
Vector< unsigned > | Lagrange_index |
Storage for the location of the Lagrange multiplier (If other additional values have been added we need to add the Lagrange multiplier at the end) More... | |
//////////////////////////////////////////////////////////////////// ////////////////////////////////////////////////////////////////////
Generic Elastic node update interface template class that can be combined with a given surface equations class and surface derivative class to provide a concrete implementation of any surface element that uses elastic node updates.
Definition at line 672 of file specific_node_update_interface_elements.h.
|
inline |
Constructor, pass a pointer to the bulk element and the face index of the bulk element to which the element is to be attached to. The optional identifier can be used to distinguish the additional nodal value (Lagr mult) created by this element from those created by other FaceElements.
Definition at line 739 of file specific_node_update_interface_elements.h.
References oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::Lagrange_index, oomph::FluidInterfaceAdditionalValues< ELEMENT >::nadditional_values(), and oomph::FluidInterfaceAdditionalValues< ELEMENT >::setup_equation_indices().
|
inline |
Helper function to calculate the additional contributions to be added at each integration point. This deals with Lagrange multiplier contribution, as well as any additional contributions by the other equations.
Definition at line 877 of file specific_node_update_interface_elements.h.
References oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::kinematic_local_eqn(), and oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::lagrange().
|
inlineprotected |
Fill in the specific surface derivative calculations by calling the appropriate function from the derivative class.
Definition at line 716 of file specific_node_update_interface_elements.h.
|
inline |
Fill in contribution to residuals and Jacobian.
Definition at line 837 of file specific_node_update_interface_elements.h.
|
inlineprivate |
Hijacking the kinematic condition corresponds to hijacking the variables associated with the Lagrange multipliers that are assigned on construction of this element.
Definition at line 701 of file specific_node_update_interface_elements.h.
References oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::lagrange_index().
|
inlineprivate |
Equation number of the kinematic BC associated with node j. (This is the equation for the Lagrange multiplier)
Definition at line 692 of file specific_node_update_interface_elements.h.
References oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::lagrange_index().
|
inline |
Return the lagrange multiplier at local node n.
Definition at line 830 of file specific_node_update_interface_elements.h.
References oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::lagrange_index().
|
inlineprivate |
Return the index at which the lagrange multiplier is stored at the n-th node.
Definition at line 685 of file specific_node_update_interface_elements.h.
Referenced by oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::hijack_kinematic_conditions(), oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::kinematic_local_eqn(), and oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::lagrange().
|
inlinevirtual |
Create an "bounding" element (here actually a 2D line element of type ElasticLineFluidInterfaceBoundingElement<ELEMENT> that allows the application of a contact angle boundary condition on the the specified face.
Definition at line 960 of file specific_node_update_interface_elements.h.
|
inline |
Overload the C-style output function.
Definition at line 861 of file specific_node_update_interface_elements.h.
|
inline |
C-style Output function.
Definition at line 867 of file specific_node_update_interface_elements.h.
|
inline |
Overload the output function.
Definition at line 849 of file specific_node_update_interface_elements.h.
|
inline |
Output the element.
Definition at line 855 of file specific_node_update_interface_elements.h.
|
inline |
The "global" intrinsic coordinate of the element when viewed as part of a geometric object should be given by the FaceElement representation, by default.
Definition at line 822 of file specific_node_update_interface_elements.h.
|
private |
Storage for the location of the Lagrange multiplier (If other additional values have been added we need to add the Lagrange multiplier at the end)
Definition at line 681 of file specific_node_update_interface_elements.h.
Referenced by oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::ElasticUpdateFluidInterfaceElement().