![]() |
|
#include <pseudosolid_node_update_elements.h>
Public Member Functions | |
PseudoSolidNodeUpdateElement () | |
void | describe_local_dofs (std::ostream &out, const std::string ¤t_string) const |
void | compute_norm (double &el_norm) |
unsigned | required_nvalue (const unsigned &n) const |
The required number of values is the sum of the two. More... | |
int | solid_p_nodal_index () const |
void | fill_in_contribution_to_residuals (Vector< double > &residuals) |
void | fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian) |
void | fill_in_contribution_to_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix) |
void | evaluate_shape_derivs_by_direct_fd () |
Evaluate shape derivatives by direct finite differencing. More... | |
void | evaluate_shape_derivs_by_chain_rule () |
Evaluate shape derivatives by chain rule. More... | |
void | fill_in_shape_derivatives (DenseMatrix< double > &jacobian) |
void | fill_in_shape_derivatives_by_fd (DenseMatrix< double > &jacobian) |
void | identify_geometric_data (std::set< Data * > &geometric_data_pt) |
void | output (std::ostream &outfile) |
Overload the output function: Call that of the basic element. More... | |
void | output (std::ostream &outfile, const unsigned &n_p) |
void | output (FILE *file_pt) |
Overload the output function: Call that of the basic element. More... | |
void | output (FILE *file_pt, const unsigned &n_p) |
Output function is just the same as the basic equations. More... | |
unsigned | num_Z2_flux_terms () |
void | compute_exact_Z2_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_flux_pt, double &error, double &norm) |
void | get_Z2_flux (const Vector< double > &s, Vector< double > &flux) |
unsigned | nvertex_node () const |
Number of vertex nodes in the element. More... | |
Node * | vertex_node_pt (const unsigned &j) const |
Pointer to the j-th vertex node in the element. More... | |
unsigned | nrecovery_order () |
unsigned | ndof_types () const |
unsigned | nbasic_dof_types () const |
unsigned | nsolid_dof_types () const |
void | get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const |
Private Attributes | |
bool | Shape_derivs_by_direct_fd |
Boolean flag to indicate shape derivative method. More... | |
A templated class that permits combination two different element types, for the solution of problems in deforming domains. The first template paremter BASIC is the standard element and the second SOLID solves the equations that are used to control the mesh deformation.
|
inline |
Constructor, call the BASIC and SOLID elements' constructors and set the "density" parameter for solid element to zero
References oomph::PseudoSolidHelper::Zero.
|
inline |
Plot the error when compared against a given exact flux. Also calculates the norm of the error and that of the exact flux. Use version in BASIC element
References calibrate::error.
|
inline |
|
inline |
Function to describe the local dofs of the element. The ostream specifies the output stream to which the description is written; the string stores the currently assembled output that is ultimately written to the output stream by Data::describe_dofs(...); it is typically built up incrementally as we descend through the call hierarchy of this function when called from Problem::describe_dofs(...)
References out().
|
inline |
Evaluate shape derivatives by chain rule.
References oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::Shape_derivs_by_direct_fd.
|
inline |
Evaluate shape derivatives by direct finite differencing.
References oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::Shape_derivs_by_direct_fd.
|
inline |
Final override for jacobian function: Contributions are included from both the underlying element types
References oomph::fill_in_contribution_to_jacobian(), and oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_shape_derivatives().
|
inline |
Final override for mass matrix function: contributions are included from both the underlying element types
References oomph::fill_in_contribution_to_jacobian_and_mass_matrix(), and oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_shape_derivatives().
|
inline |
Final override for the residuals function. Contributions are added from both underlying element types
References oomph::fill_in_contribution_to_residuals().
|
inline |
Fill in the shape derivatives of the BASIC equations w.r.t. the solid position dofs
References oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_shape_derivatives_by_fd(), i, k, n, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::Shape_derivs_by_direct_fd.
Referenced by oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_contribution_to_jacobian(), and oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_contribution_to_jacobian_and_mass_matrix().
|
inline |
Fill in the derivatives of the BASIC equations w.r.t. the solid position dofs
References oomph::fill_in_contribution_to_residuals(), i, k, m, and n.
Referenced by oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_shape_derivatives().
|
inline |
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contains the global equation number of the unknown, while the second one contains the number of the "DOF type" that this unknown is associated with. This method combines the get_dof_numbers_for_unknowns(...) method for the BASIC and SOLID elements. The basic elements retain their DOF type numbering and the SOLID elements DOF type numbers are incremented by nbasic_dof_types().
References oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::nbasic_dof_types().
|
inline |
'Flux' vector for Z2 error estimation: Error estimation is based on error in BASIC element
References ProblemParameters::flux(), get_Z2_flux(), and s.
|
inline |
Specify Data that affects the geometry of the element by adding the position Data to the set that's passed in. (This functionality is required in FSI problems; set is used to avoid double counting).
References j.
|
inline |
return the number of DOF types associated with the BASIC elements in this combined element
Referenced by oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::get_dof_numbers_for_unknowns().
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
Overload the output function: Call that of the basic element.
References output().
|
inline |
Output function is just the same as the basic equations.
References output().
|
inline |
Overload the output function: Call that of the basic element.
References output().
|
inline |
Output function: Plot at n_p plot points using the basic element's output function
References output().
|
inline |
|
inline |
We assume that the solid stuff is stored at the end of the nodes, i.e. its index is the number of continuously interplated values in the BASIC equations.
References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.
|
inline |
|
private |
Boolean flag to indicate shape derivative method.
Referenced by oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::evaluate_shape_derivs_by_chain_rule(), oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::evaluate_shape_derivs_by_direct_fd(), and oomph::PseudoSolidNodeUpdateElement< BASIC, SOLID >::fill_in_shape_derivatives().