29#ifndef OOMPH_SOLID_TRACTION_ELEMENTS_HEADER
30#define OOMPH_SOLID_TRACTION_ELEMENTS_HEADER
34#include <oomph-lib-config.h>
48 namespace SolidTractionElementHelper
75 template<
class ELEMENT>
135 if (element_pt->
dim() == 3)
145 "This face element will not work correctly if nodes are "
146 "hanging.\nUse the refineable version instead. ",
233 for (
unsigned i = 0;
i <
n_dim;
i++)
239 for (
unsigned i = 0;
i <
n_dim;
i++)
245 for (
unsigned i = 0;
i <
n_dim;
i++)
285 template<
class ELEMENT>
298 this->interpolated_xi(
s, xi);
315 template<
class ELEMENT>
363 for (
unsigned i = 0;
i <
n_dim;
i++)
369 interpolated_xi[
i] +=
370 this->lagrangian_position_gen(
l, bulk_position_type(
k),
i) *
375 for (
unsigned j = 0;
j <
n_dim - 1;
j++)
387 for (
unsigned i = 0;
i <
n_dim - 1;
i++)
389 for (
unsigned j = 0;
j <
n_dim - 1;
j++)
394 for (
unsigned k = 0;
k <
n_dim;
k++)
413 Adet = A(0, 0) * A(1, 1) - A(0, 1) * A(1, 0);
417 "Wrong dimension in SolidTractionElement",
418 "SolidTractionElement::fill_in_contribution_to_residuals()",
440 for (
unsigned i = 0;
i <
n_dim;
i++)
442 local_eqn = this->position_local_eqn(
l, bulk_position_type(
k),
i);
470 template<
class ELEMENT>
537 template<
class ELEMENT>
550 if (n_position_type != 1)
553 "RefineableSolidTractionElement only works for n_position_type=1",
598 for (
unsigned i = 0;
i <
n_dim;
i++)
604 interpolated_xi[
i] +=
605 this->lagrangian_position_gen(
l, this->bulk_position_type(
k),
i) *
610 for (
unsigned j = 0;
j <
n_dim - 1;
j++)
622 for (
unsigned i = 0;
i <
n_dim - 1;
i++)
624 for (
unsigned j = 0;
j <
n_dim - 1;
j++)
629 for (
unsigned k = 0;
k <
n_dim;
k++)
648 Adet = A(0, 0) * A(1, 1) - A(0, 1) * A(1, 0);
652 "Wrong dimension in RefineableSolidTractionElement",
715 for (
unsigned i = 0;
i <
n_dim;
i++)
718 this->position_local_eqn(
l,
k,
i);
730 for (
unsigned i = 0;
i <
n_dim;
i++)
777 template<
class ELEMENT,
unsigned DIM>
867 for (
unsigned i = 0;
i <
DIM;
i++)
881 outfile <<
"ZONE" << std::endl;
925 for (
unsigned i = 0;
i <
DIM;
i++)
929 for (
unsigned i = 0;
i <
DIM;
i++)
933 for (
unsigned i = 0;
i <
DIM;
i++)
949 throw OomphLibError(
"It doesn't make sense to specify an external "
950 "traction in an FSI context",
984 template<
class ELEMENT,
unsigned DIM>
1041 template<
class ELEMENT,
unsigned DIM>
1103 template<
class ELEMENT>
1117 const unsigned&
id = 0,
1142 if (element_pt->
dim() == 3)
1152 "This face element will not work correctly if nodes are "
1153 "hanging\nUse the refineable version instead. ",
1168 "ImposeDisplacementByLagrangeMultiplierElement cannot (currently) "
1169 "be used with elements that have generalised positional dofs",
1177 unsigned dim = element_pt->
dim();
1236 std::ostringstream error_message;
1237 error_message <<
"About to wipe external data for "
1238 <<
"ImposeDisplacementByLagrangeMultiplierElement.\n"
1239 <<
"I noted that N_assigned_geom_data = "
1243 <<
"so we're going to wipe some external data that\n"
1244 <<
"is not geometric Data of the GeomObject that\n"
1245 <<
"specifies the desired boundary shape.\n"
1253 for (
unsigned i = 0;
i < n_geom_data;
i++)
1284 std::ostringstream error_message;
1285 error_message <<
"About to wipe external data for "
1286 <<
"ImposeDisplacementByLagrangeMultiplierElement.\n"
1287 <<
"I noted that N_assigned_geom_data = "
1291 <<
"so we're going to wipe some external data that\n"
1292 <<
"is not geometric Data of the GeomObject that\n"
1293 <<
"specifies the desired boundary shape.\n"
1392 for (
unsigned i = 0;
i < n_dof;
i++)
1395 for (
unsigned j = 0;
j <
n_dof;
j++)
1397 jacobian(
i,
j) *= -1.0;
1413 if (n_position_type != 1)
1416 "ImposeDisplacementByLagrangeMultiplierElement cannot (currently) be "
1417 "used with elements that have generalised positional dofs",
1460 bnod_pt->index_of_first_value_assigned_by_face_element(
Id);
1463 for (
unsigned i = 0;
i <
dim_el + 1;
i++)
1486 for (
unsigned i = 0;
i <
dim_el + 1;
i++)
1490 for (
unsigned i = 0;
i <
dim_el + 1;
i++)
1494 for (
unsigned i = 0;
i <
dim_el + 1;
i++)
1523 if (n_position_type != 1)
1526 "ImposeDisplacementByLagrangeMultiplierElement cannot (currently) be "
1527 "used with elements that have generalised positional dofs",
1577 bnod_pt->index_of_first_value_assigned_by_face_element(
Id);
1580 for (
unsigned i = 0;
i <
dim_el + 1;
i++)
1615 for (
unsigned k = 0;
k <
dim_el + 1;
k++)
1632 adet = a(0, 0) * a(1, 1) - a(0, 1) * a(1, 0);
1638 "fill_in_generic_contribution_to_residuals_displ_lagr_multiplier",
1639 "ImposeDisplacementByLagrangeMultiplierElement::fill_in_generic_"
1640 "contribution_to_residuals_displ_lagr_multiplier()",
1662 for (
unsigned i = 0;
i <
dim_el + 1;
i++)
1679 const unsigned&
flag)
1685 if (n_position_type != 1)
1688 "ImposeDisplacementByLagrangeMultiplierElement cannot (currently) be "
1689 "used with elements that have generalised positional dofs",
1735 bnod_pt->index_of_first_value_assigned_by_face_element(
Id);
1738 for (
unsigned i = 0;
i <
dim_el + 1;
i++)
1773 for (
unsigned k = 0;
k <
dim_el + 1;
k++)
1790 adet = a(0, 0) * a(1, 1) - a(0, 1) * a(1, 0);
1796 "fill_in_generic_contribution_to_residuals_displ_lagr_multiplier",
1797 "ImposeDisplacementByLagrangeMultiplierElement::fill_in_generic_"
1798 "contribution_to_residuals_displ_lagr_multiplier()",
1820 for (
unsigned i = 0;
i <
dim_el + 1;
i++)
1834 bnod_pt->index_of_first_value_assigned_by_face_element(
Id) +
i);
1880 bnode_pt->index_of_first_value_assigned_by_face_element(
1903 return this->
dim() + 1;
1924 for (
unsigned i = 0;
i < dim_el + 1;
i++)
1935 j,
bnod_pt->index_of_first_value_assigned_by_face_element(
Id) +
i);
2008 template<
class ELEMENT>
2023 const unsigned&
id = 0)
2066 const unsigned&
flag)
2072 if (n_position_type != 1)
2075 "RefineableImposeDisplacementByLagrangeMultiplierElement \n cannot "
2076 "(currently) be used with elements that have generalised\n "
2077 "positional dofs. Upgrade should be straightforward though the code "
2078 "is in a \n bit of state, with generalised degrees of freedom "
2079 "sometimes half taken into \n account, sometimes completely "
2132 bnod_pt->index_of_first_value_assigned_by_face_element(this->
Id);
2135 for (
unsigned i = 0;
i < dim_el + 1;
i++)
2170 for (
unsigned k = 0;
k <
dim_el + 1;
k++)
2187 adet = a(0, 0) * a(1, 1) - a(0, 1) * a(1, 0);
2193 "refineable_fill_in_generic_contribution_to_residuals_displ_lagr_"
2195 "RefineableImposeDisplacementByLagrangeMultiplierElement::"
2196 "refineablefill_in_generic_contribution_to_residuals_displ_lagr_"
2258 for (
unsigned i = 0;
i <
dim_el + 1;
i++)
2272 bnod_pt->index_of_first_value_assigned_by_face_element(
2289 bnod_pt->index_of_first_value_assigned_by_face_element(
2466 ->index_of_first_value_assigned_by_face_element(
2484 ->index_of_first_value_assigned_by_face_element(
2533 template<
class ELEMENT>
2563 const unsigned&
id = 0,
2582 if (element_pt->
dim() == 3)
2592 "This face element will not work correctly if nodes are "
2593 "hanging\nUse the refineable version instead. ",
2607 throw OomphLibError(
"FSIImposeDisplacementByLagrangeMultiplierElement"
2608 " cannot (currently) be used with elements that "
2609 "have generalised positional dofs",
2617 unsigned dim = element_pt->
dim();
2666 for (
unsigned i = 0;
i < n_dof;
i++)
2669 for (
unsigned j = 0;
j <
n_dof;
j++)
2671 jacobian(
i,
j) *= -1.0;
2687 if (n_position_type != 1)
2690 "FSIImposeDisplacementByLagrangeMultiplierElement cannot (currently) "
2691 "be used with elements that have generalised positional dofs",
2734 bnod_pt->index_of_first_value_assigned_by_face_element(
Id);
2737 for (
unsigned i = 0;
i <
dim_el + 1;
i++)
2756 for (
unsigned i = 0;
i <
dim_el + 1;
i++)
2760 for (
unsigned i = 0;
i <
dim_el + 1;
i++)
2783 const unsigned&
flag)
2789 if (n_position_type != 1)
2792 "FSIImposeDisplacementByLagrangeMultiplierElement cannot (currently) "
2793 "be used with elements that have generalised positional dofs",
2839 bnod_pt->index_of_first_value_assigned_by_face_element(
Id);
2842 for (
unsigned i = 0;
i <
dim_el + 1;
i++)
2863 for (
unsigned k = 0;
k <
dim_el + 1;
k++)
2880 adet = a(0, 0) * a(1, 1) - a(0, 1) * a(1, 0);
2886 "fill_in_generic_contribution_to_residuals_displ_lagr_multiplier",
2887 "FSIImposeDisplacementByLagrangeMultiplierElement::fill_in_"
2888 "generic_contribution_to_residuals_fsi_displ_lagr_multiplier()",
2910 for (
unsigned i = 0;
i <
dim_el + 1;
i++)
2924 bnod_pt->index_of_first_value_assigned_by_face_element(
Id) +
i);
2970 bnode_pt->index_of_first_value_assigned_by_face_element(
2992 return this->
dim() + 1;
3013 for (
unsigned i = 0;
i < dim_el + 1;
i++)
3024 j,
bnod_pt->index_of_first_value_assigned_by_face_element(
Id) +
i);
3068 template<
class ELEMENT>
3094 const unsigned&
id = 0)
3136 const unsigned&
flag)
3142 if (n_position_type != 1)
3145 "RefineableImposeDisplacementByLagrangeMultiplierElement \n cannot "
3146 "(currently) be used with elements that have generalised\n "
3147 "positional dofs. Upgrade should be straightforward though the code "
3148 "is in a \n bit of state, with generalised degrees of freedom "
3149 "sometimes half taken into \n account, sometimes completely "
3202 bnod_pt->index_of_first_value_assigned_by_face_element(this->
Id);
3205 for (
unsigned i = 0;
i < dim_el + 1;
i++)
3227 for (
unsigned k = 0;
k <
dim_el + 1;
k++)
3244 adet = a(0, 0) * a(1, 1) - a(0, 1) * a(1, 0);
3250 "refineable_fill_in_generic_contribution_to_residuals_displ_lagr_"
3252 "RefineableImposeDisplacementByLagrangeMultiplierElement::"
3253 "refineablefill_in_generic_contribution_to_residuals_displ_lagr_"
3316 for (
unsigned i = 0;
i <
dim_el + 1;
i++)
3330 bnod_pt->index_of_first_value_assigned_by_face_element(
3347 bnod_pt->index_of_first_value_assigned_by_face_element(
3524 ->index_of_first_value_assigned_by_face_element(
3542 ->index_of_first_value_assigned_by_face_element(
A class that contains the information required by Nodes that are located on Mesh boundaries....
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
unsigned nvalue() const
Return number of values stored in data object (incl pinned ones).
This is a base class for all elements that require external sources (e.g. FSI, multi-domain problems ...
Vector< double > & external_element_local_coord(const unsigned &interaction_index, const unsigned &ipt)
Access function to get source element's local coords for specified interaction index at specified int...
void set_ninteraction(const unsigned &n_interaction)
Set the number of interactions in the element This function is usually called in the specific element...
void describe_local_dofs(std::ostream &out, const std::string &curr_string) const
Function to describe the local dofs of the element. The ostream specifies the output stream to which ...
void fill_in_jacobian_from_external_interaction_by_fd(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Calculate the contributions to the jacobian from all external interaction degrees of freedom (geometr...
FiniteElement *& external_element_pt(const unsigned &interaction_index, const unsigned &ipt)
Access function to source element for specified interaction index at specified integration point.
A class for elements that allow the imposition of a displacement constraint for bulk solid elements v...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
FSIImposeDisplacementByLagrangeMultiplierElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0, const bool &called_from_refineable_constructor=false)
Constructor takes a "bulk" element and the index that identifies which face the FaceElement is suppos...
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned > > &dof_lookup_list) const
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contain...
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Fill in contribution to Mass matrix and Jacobian. There is no contributiont to mass matrix so simply ...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the residuals.
void fill_in_generic_contribution_to_residuals_fsi_displ_lagr_multiplier(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Helper function to compute the residuals and, if flag==1, the Jacobian.
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into: Just the soli...
void describe_local_dofs(std::ostream &out, const std::string ¤t_string) const
Function to describe the local dofs of the elements. The ostream specifies the output stream to which...
void output(std::ostream &outfile)
Output function.
SolidTractionElement "upgraded" to a FSIWallElement (and thus, by inheritance, a GeomObject),...
void set_normal_pointing_into_fluid()
Set the normal computed by underlying FaceElement to point into the fluid.
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into: Just the soli...
void dposition_dlagrangian_at_local_coordinate(const Vector< double > &s, DenseMatrix< double > &drdxi) const
Derivative of position vector w.r.t. the SolidFiniteElement's Lagrangian coordinates; evaluated at cu...
FSISolidTractionElement(FiniteElement *const &element_pt, const int &face_index, const bool &called_from_refineable_constructor=false)
Constructor: Create element as FSIWallElement (and thus, by inheritance, a GeomObject) with DIM-1 Lag...
~FSISolidTractionElement()
Destructor: empty.
virtual void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Final overload... Forwards to the version in the FSIWallElement.
virtual void get_traction(const unsigned &intpt, const Vector< double > &xi, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction)
Get the load vector: Pass number of the integration point, Lagr. coordinate, Eulerian coordinate (nei...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: Note we can only output the traction at Gauss points so n_plot is actually ignored.
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned > > &dof_lookup_list) const
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contain...
void set_normal_pointing_out_of_fluid()
Set the normal computed by underlying FaceElement to point out of the fluid.
bool Normal_points_into_fluid
Boolean flag to indicate whether the normal is directed into the fluid.
virtual void(*&)(const Vector< double > &xi, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction) traction_fct_pt()
Broken overloaded reference to the traction function pointer. It doesn't make sense to specify an ext...
This is a base class for all SolidFiniteElements that participate in FSI computations....
void setup_fsi_wall_element(const unsigned &nlagr_solid, const unsigned &ndim_fluid)
Setup: Assign storage – pass the Eulerian dimension of the "adjacent" fluid elements and the number o...
void fluid_load_vector(const unsigned &intpt, const Vector< double > &N, Vector< double > &load)
Get FE Jacobian by systematic finite differencing w.r.t. nodal positition Data, internal and external...
FaceElements are elements that coincide with the faces of higher-dimensional "bulk" elements....
bool Boundary_number_in_bulk_mesh_has_been_set
Has the Boundary_number_in_bulk_mesh been set? Only included if compiled with PARANOID switched on.
unsigned Boundary_number_in_bulk_mesh
The boundary number in the bulk mesh to which this element is attached.
int & face_index()
Index of the face (a number that uniquely identifies the face in the element)
const unsigned & boundary_number_in_bulk_mesh() const
Broken assignment operator.
void outer_unit_normal(const Vector< double > &s, Vector< double > &unit_normal) const
Compute outer unit normal at the specified local coordinate.
void add_additional_values(const Vector< unsigned > &nadditional_values, const unsigned &id)
Helper function adding additional values for the unknowns associated with the FaceElement....
double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s. Overloaded to get information from bulk...
FiniteElement *& bulk_element_pt()
Pointer to higher-dimensional "bulk" element.
FaceGeometry class definition: This policy class is used to allow construction of face elements that ...
A general Finite Element class.
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
void position(const Vector< double > &zeta, Vector< double > &r) const
Return the parametrised position of the FiniteElement in its incarnation as a GeomObject,...
virtual void dshape_local_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsids) const
Return the geometric shape function and its derivative w.r.t. the local coordinates at the ipt-th int...
virtual void output(std::ostream &outfile)
Output the element data — typically the values at the nodes in a format suitable for post-processing.
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction")
unsigned nnodal_position_type() const
Return the number of coordinate types that the element requires to interpolate the geometry between t...
double size() const
Calculate the size of the element (length, area, volume,...) in Eulerian computational coordinates....
virtual void shape(const Vector< double > &s, Shape &psi) const =0
Calculate the geometric shape functions at local coordinate s. This function must be overloaded for e...
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
double nodal_position_gen(const unsigned &n, const unsigned &k, const unsigned &i) const
Return the value of the k-th type of the i-th positional variable at the local node n.
unsigned nnode() const
Return the number of nodes.
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Get cector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
virtual unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction")
double nodal_position(const unsigned &n, const unsigned &i) const
Return the i-th coordinate at local node n. If the node is hanging, the appropriate interpolation is ...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
virtual void build_face_element(const int &face_index, FaceElement *face_element_pt)
Function for building a lower dimensional FaceElement on the specified face of the FiniteElement....
virtual void describe_nodal_local_dofs(std::ostream &out, const std::string ¤t_string) const
Function to describe the local dofs of the element[s]. The ostream specifies the output stream to whi...
unsigned nodal_dimension() const
Return the required Eulerian dimension of the nodes in this element.
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Add tecplot zone "footer" to output stream (when plotting nplot points in each "coordinate direction"...
unsigned ngeom_data() const
A standard FiniteElement is fixed, so there are no geometric data when viewed in its GeomObject incar...
virtual void shape_at_knot(const unsigned &ipt, Shape &psi) const
Return the geometric shape function at the ipt-th integration point.
bool has_hanging_nodes() const
Return boolean to indicate if any of the element's nodes are geometrically hanging.
void fill_in_jacobian_from_external_by_fd(Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
Calculate the contributions to the jacobian from the external degrees of freedom using finite differe...
unsigned nexternal_data() const
Return the number of external data objects.
bool is_halo() const
Is this element a halo?
unsigned ndof() const
Return the number of equations/dofs in the element.
unsigned long eqn_number(const unsigned &ieqn_local) const
Return the global equation number corresponding to the ieqn_local-th local equation number.
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
void flush_external_data()
Flush all external data.
unsigned add_external_data(Data *const &data_pt, const bool &fd=true)
Add a (pointer to an) external data object to the element and return its index (i....
A geometric object is an object that provides a parametrised description of its shape via the functio...
virtual void position(const Vector< double > &zeta, Vector< double > &r) const =0
Parametrised position on object at current time: r(zeta).
virtual unsigned ngeom_data() const
How many items of Data does the shape of the object depend on? This is implemented as a broken virtua...
virtual Data * geom_data_pt(const unsigned &j)
Return pointer to the j-th Data item that the object's shape depends on. This is implemented as a bro...
virtual void locate_zeta(const Vector< double > &zeta, GeomObject *&sub_geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
A geometric object may be composed of may sub-objects (e.g. a finite-element representation of a boun...
Class that contains data for hanging nodes.
Node *const & master_node_pt(const unsigned &i) const
Return a pointer to the i-th master node.
A class for elements that allow the imposition of a displacement constraint for "bulk" solid elements...
Vector< Vector< double > > Zeta_sub_geom_object
Storage for local coordinates in sub-GeomObjects at integration points.
unsigned N_assigned_geom_data
Bool to record the number of geom Data that has been assigned to external data (we're keeping a recor...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the residuals.
Vector< GeomObject * > Sub_geom_object_pt
Storage for sub-GeomObject at the integration points.
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into: We only label...
void output(std::ostream &outfile)
Output function.
void set_boundary_shape_geom_object_pt(GeomObject *boundary_shape_geom_object_pt, const unsigned &boundary_number_in_bulk_mesh)
Set GeomObject that specifies the prescribed boundary displacement; GeomObject is assumed to be param...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
double square_of_l2_norm_of_error()
Compute square of L2 norm of error between prescribed and actual displacement.
GeomObject * Boundary_shape_geom_object_pt
GeomObject that specifies the prescribed boundary displacement; GeomObject is assumed to be parametri...
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned > > &dof_lookup_list) const
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contain...
GeomObject * boundary_shape_geom_object_pt() const
Access to GeomObject that specifies the prescribed boundary displacement; GeomObject is assumed to be...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
bool Sparsify
Boolean flag to indicate if we're using geometric Data of sub-GeomObjects that make up the (possibly ...
ImposeDisplacementByLagrangeMultiplierElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0, const bool &called_from_refineable_constructor=false)
Constructor takes a "bulk" element and the index that identifies which face the FaceElement is suppos...
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Fill in contribution to Mass matrix and Jacobian. There is no contributiont to mass matrix so simply ...
void fill_in_generic_contribution_to_residuals_displ_lagr_multiplier(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Helper function to compute the residuals and, if flag==1, the Jacobian.
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
A base class for SolidElements that can have hanging nodes but are not refineable as such....
An OomphLibError object which should be thrown when an run-time error is encountered....
RefineableElements are FiniteElements that may be subdivided into children to provide a better local ...
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 ...
A class for elements that allow the imposition of a displacement constraint for bulk solid elements v...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
RefineableFSIImposeDisplacementByLagrangeMultiplierElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0)
Constructor takes a "bulk" element and the index that identifies which face the FaceElement is suppos...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the residuals.
void refineable_fill_in_generic_contribution_to_residuals_fsi_displ_lagr_multiplier(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Helper function to compute the residuals and, if flag==1, the Jacobian.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: Same for all nodes since it's a refineable element.
RefineableSolidTractionElement "upgraded" to a FSIWallElement (and thus, by inheritance,...
virtual void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Final overload. Get contributions from refineable solid traction element and derivatives from externa...
~RefineableFSISolidTractionElement()
Destructor: empty.
RefineableFSISolidTractionElement(FiniteElement *const &element_pt, const int &face_index)
Constructor: Create element as FSIWallElement (and thus, by inheritance, a GeomObject) with DIM-1 Lag...
A class for elements that allow the imposition of a displacement constraint for "bulk" solid elements...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the residuals.
void refineable_fill_in_generic_contribution_to_residuals_displ_lagr_multiplier(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Helper function to compute the residuals and, if flag==1, the Jacobian.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: Same for all nodes since it's a refineable element.
RefineableImposeDisplacementByLagrangeMultiplierElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0)
Constructor takes a "bulk" element and the index that identifies which face the FaceElement is suppos...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
DenseMatrix< int > & local_position_hang_eqn(Node *const &node_pt)
Access the local equation number of of hanging node variables associated with nodal positions....
A class for elements that allow the imposition of an applied traction in the principle of virtual dis...
~RefineableSolidTractionElement()
Destructor should not delete anything.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
This function returns just the residuals.
RefineableSolidTractionElement(FiniteElement *const &element_pt, const int &face_index)
Constructor, which takes a "bulk" element and the face index.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values are the same as those in the bulk element.
void refineable_fill_in_contribution_to_residuals_solid_traction(Vector< double > &residuals)
This function returns the residuals for the traction function.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
This function returns the residuals and the Jacobian.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
SolidFaceElements combine FaceElements and SolidFiniteElements and overload various functions so they...
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 ...
double interpolated_xi(const Vector< double > &s, const unsigned &i) const
Return i-th FE-interpolated Lagrangian coordinate xi[i] at local coordinate s. Overloaded from SolidF...
Data * geom_data_pt(const unsigned &j)
Return pointer to the j-th Data item that the object's shape depends on. (Redirects to the nodes' pos...
double lagrangian_position(const unsigned &n, const unsigned &i) const
Return i-th Lagrangian coordinate at local node n.
int position_local_eqn(const unsigned &n, const unsigned &k, const unsigned &j) const
Access function that returns the local equation number that corresponds to the j-th coordinate of the...
A class for elements that allow the imposition of an applied traction in the principle of virtual dis...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
void fill_in_contribution_to_residuals_solid_traction(Vector< double > &residuals)
Helper function that actually calculates the residuals.
void(*&)(const Vector< double > &xi, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction) traction_fct_pt()
Reference to the traction function pointer.
SolidTractionElement(FiniteElement *const &element_pt, const int &face_index, const bool &called_from_refineable_constructor=false)
Constructor, which takes a "bulk" element and the value of the index and its limit.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
void output(std::ostream &outfile)
Output function.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Return the residuals.
void traction(const Vector< double > &s, Vector< double > &traction)
Compute traction vector at specified local coordinate Should only be used for post-processing; ignore...
virtual void get_traction(const unsigned &intpt, const Vector< double > &xi, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction)
Get the traction vector: Pass number of integration point (dummy), Lagr. coordinate and normal vector...
void output(FILE *file_pt)
C_style output function.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
void(* Traction_fct_pt)(const Vector< double > &xi, const Vector< double > &x, const Vector< double > &n, Vector< double > &result)
Pointer to an imposed traction function. Arguments: Lagrangian coordinate; Eulerian coordinate; outer...
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
TAdvectionDiffusionReactionElement()
Constructor: Call constructors for TElement and AdvectionDiffusionReaction equations.
void Zero_traction_fct(const Vector< double > &xi, const Vector< double > &x, const Vector< double > &N, Vector< double > &load)
Default load function (zero traction)
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).