Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes | Static Private Attributes | Friends | List of all members
oomph::FluidInterfaceElement Class Referenceabstract

Base class establishing common interfaces and functions for all Navier-Stokes-like fluid interface elements. Namely, elements that represent either a free surface or an interface between two fluids that have distinct momentum-like equation for each velocity component. More...

#include <interface_elements.h>

Inheritance diagram for oomph::FluidInterfaceElement:
oomph::ElasticUpdateFluidInterfaceElement< FluidInterfaceElement, SurfaceDerivatives, ELEMENT > oomph::ElasticUpdateFluidInterfaceElement< FluidInterfaceElement, AxisymmetricDerivatives, ELEMENT > oomph::ElasticUpdateFluidInterfaceElement< FluidInterfaceElement, LineDerivatives, ELEMENT > oomph::SpineUpdateFluidInterfaceElement< FluidInterfaceElement, SurfaceDerivatives, ELEMENT > oomph::SpineUpdateFluidInterfaceElement< FluidInterfaceElement, AxisymmetricDerivatives, ELEMENT > oomph::SpineUpdateFluidInterfaceElement< FluidInterfaceElement, LineDerivatives, ELEMENT > oomph::SurfactantTransportInterfaceElement oomph::ElasticSurfaceFluidInterfaceElement< ELEMENT > oomph::ElasticAxisymmetricFluidInterfaceElement< ELEMENT > oomph::ElasticLineFluidInterfaceElement< ELEMENT > oomph::SpineSurfaceFluidInterfaceElement< ELEMENT > oomph::SpineAxisymmetricFluidInterfaceElement< ELEMENT > oomph::SpineLineFluidInterfaceElement< ELEMENT > oomph::ElasticUpdateFluidInterfaceElement< SurfactantTransportInterfaceElement, AxisymmetricDerivatives, ELEMENT > oomph::SpineUpdateFluidInterfaceElement< SurfactantTransportInterfaceElement, LineDerivatives, ELEMENT > oomph::SpineUpdateFluidInterfaceElement< SurfactantTransportInterfaceElement, AxisymmetricDerivatives, ELEMENT > oomph::SpineUpdateFluidInterfaceElement< SurfactantTransportInterfaceElement, SurfaceDerivatives, ELEMENT >

Public Member Functions

 FluidInterfaceElement ()
 Constructor, set the default values of the booleans and pointers (null) More...
 
virtual double sigma (const Vector< double > &s_local)
 Virtual function that specifies the non-dimensional surface tension as a function of local position within the element. The default behaviour is a constant surface tension of value 1.0 This function can be overloaded in more specialised elements to incorporate variations in surface tension. More...
 
void fill_in_contribution_to_residuals (Vector< double > &residuals)
 Calculate the residuals by calling the generic residual contribution. More...
 
const double & ca () const
 The value of the Capillary number. More...
 
double *& ca_pt ()
 Pointer to the Capillary number. More...
 
const double & st () const
 The value of the Strouhal number. More...
 
double *& st_pt ()
 The pointer to the Strouhal number. More...
 
double u (const unsigned &j, const unsigned &i)
 Return the i-th velocity component at local node j. More...
 
double interpolated_u (const Vector< double > &s, const unsigned &i)
 Calculate the i-th velocity component at the local coordinate s. More...
 
double pext () const
 Return the value of the external pressure. More...
 
void set_external_pressure_data (Data *external_pressure_data_pt)
 Set the Data that contains the single pressure value that specifies the "external pressure" for the interface/free-surface. Setting this only makes sense if the interface is, in fact, a free surface (well, an interface to another inviscid fluid if you want to be picky). More...
 
void set_external_pressure_data (Data *external_pressure_data_pt, const unsigned &index_of_external_pressure_value)
 Set the Data that contains the pressure value that specifies the "external pressure" for the interface/free-surface. Setting this only makes sense if the interface is, in fact, a free surface (well, an interface to another inviscid fluid if you want to be picky). Second argument specifies the index of the pressure value within the Data object. More...
 
virtual FluidInterfaceBoundingElementmake_bounding_element (const int &face_index)
 Create a bounding element e.g. to apply a contact angle boundary condition. More...
 
virtual void hijack_kinematic_conditions (const Vector< unsigned > &bulk_node_number)=0
 Hijack the kinematic condition at the node numbers passed in the vector. The node numbers correspond to the local numbers of nodes in the associated bulk element. This is required so that contact-angle conditions can be applied by the FluidInterfaceBoundingElements. More...
 
void output (std::ostream &outfile)
 Overload the output function. More...
 
void output (std::ostream &outfile, const unsigned &n_plot)
 Output function. 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...
 

Protected Member Functions

virtual int kinematic_local_eqn (const unsigned &n)=0
 Access function that returns the local equation number for the (scalar) kinematic equation associated with the j-th local node. This must be overloaded by specific interface elements and depends on the method for handing the free-surface deformation. More...
 
int pext_local_eqn ()
 Access function for the local equation number that corresponds to the external pressure. More...
 
virtual void fill_in_generic_residual_contribution_interface (Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
 Helper function to calculate the residuals and (if flag==1) the Jacobian of the equations. This is implemented generically using the surface divergence information that is overloaded in each element i.e. axisymmetric, two- or three-dimensional. More...
 
virtual double compute_surface_derivatives (const Shape &psi, const DShape &dpsids, const DenseMatrix< double > &interpolated_t, const Vector< double > &interpolated_x, DShape &dpsidS, DShape &dpsidS_div)=0
 Compute the surface gradient and surface divergence operators given the shape functions, derivatives, tangent vectors and position. All derivatives and tangent vectors should be formed with respect to the local coordinates. More...
 
virtual 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 the resisuals and Jacobian that arise from specific node update strategies. This is called within the integration loop over the element (for efficiency) and therefore requires a fairly large number of input parameters: More...
 

Protected Attributes

Vector< unsigned > U_index_interface
 Nodal index at which the i-th velocity component is stored. More...
 
int External_data_number_of_external_pressure
 The Data that contains the external pressure is stored as external Data for the element. Which external Data item is it? (int so it can be initialised to -1, indicating that external pressure hasn't been set). More...
 
Data * Pext_data_pt
 Pointer to the Data item that stores the external pressure. More...
 
unsigned Index_of_external_pressure_value
 Which of the values in Pext_data_pt stores the external pressure. More...
 

Private Attributes

double * Ca_pt
 Pointer to the Capillary number. More...
 
double * St_pt
 Pointer to the Strouhal number. More...
 

Static Private Attributes

static double Default_Physical_Constant_Value = 1.0
 Default value for physical constants. More...
 

Friends

class FluidInterfaceBoundingElement
 

Detailed Description

Base class establishing common interfaces and functions for all Navier-Stokes-like fluid interface elements. Namely, elements that represent either a free surface or an interface between two fluids that have distinct momentum-like equation for each velocity component.

Definition at line 320 of file interface_elements.h.

Constructor & Destructor Documentation

◆ FluidInterfaceElement()

oomph::FluidInterfaceElement::FluidInterfaceElement ( )
inline

Constructor, set the default values of the booleans and pointers (null)

Definition at line 444 of file interface_elements.h.

References Ca_pt, Default_Physical_Constant_Value, and St_pt.

Member Function Documentation

◆ add_additional_residual_contributions_interface()

virtual void oomph::FluidInterfaceElement::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 
)
inlineprotectedvirtual

Helper function to calculate the additional contributions to the resisuals and Jacobian that arise from specific node update strategies. This is called within the integration loop over the element (for efficiency) and therefore requires a fairly large number of input parameters:

  • the velocity shape functions and their derivatives w.r.t. the local coordinates
  • the surface gradient and divergence of the velocity shape functions
  • The local and Eulerian coordinates,
  • the outer unit normal,
  • the integration weight from the integration scheme
  • the Jacobian of the mapping between the local and global coordinates along the element. (Note that in the axisymmmetric case this includes the r term)!

Reimplemented in oomph::SurfactantTransportInterfaceElement, oomph::ElasticUpdateFluidInterfaceElement< FluidInterfaceElement, SurfaceDerivatives, ELEMENT >, oomph::ElasticUpdateFluidInterfaceElement< SurfactantTransportInterfaceElement, AxisymmetricDerivatives, ELEMENT >, oomph::ElasticUpdateFluidInterfaceElement< FluidInterfaceElement, AxisymmetricDerivatives, ELEMENT >, oomph::ElasticUpdateFluidInterfaceElement< FluidInterfaceElement, LineDerivatives, ELEMENT >, oomph::SpineUpdateFluidInterfaceElement< SurfactantTransportInterfaceElement, LineDerivatives, ELEMENT >, oomph::SpineUpdateFluidInterfaceElement< FluidInterfaceElement, SurfaceDerivatives, ELEMENT >, oomph::SpineUpdateFluidInterfaceElement< SurfactantTransportInterfaceElement, AxisymmetricDerivatives, ELEMENT >, oomph::SpineUpdateFluidInterfaceElement< SurfactantTransportInterfaceElement, SurfaceDerivatives, ELEMENT >, oomph::SpineUpdateFluidInterfaceElement< FluidInterfaceElement, AxisymmetricDerivatives, ELEMENT >, and oomph::SpineUpdateFluidInterfaceElement< FluidInterfaceElement, LineDerivatives, ELEMENT >.

Definition at line 426 of file interface_elements.h.

Referenced by fill_in_generic_residual_contribution_interface().

◆ ca()

const double& oomph::FluidInterfaceElement::ca ( ) const
inline

◆ ca_pt()

double*& oomph::FluidInterfaceElement::ca_pt ( )
inline

Pointer to the Capillary number.

Definition at line 492 of file interface_elements.h.

References Ca_pt.

◆ compute_surface_derivatives()

virtual double oomph::FluidInterfaceElement::compute_surface_derivatives ( const Shape &  psi,
const DShape &  dpsids,
const DenseMatrix< double > &  interpolated_t,
const Vector< double > &  interpolated_x,
DShape &  dpsidS,
DShape &  dpsidS_div 
)
protectedpure virtual

Compute the surface gradient and surface divergence operators given the shape functions, derivatives, tangent vectors and position. All derivatives and tangent vectors should be formed with respect to the local coordinates.

Return the jacobian of the surface, as well as the dpsidS, and dpsidS_div objects.

This is the only function that needs to be overloaded to specify different geometries.

In order to compute the surface gradient of a scalar function one needs only compute the sum over the nodes of dpsidS(l,i) * nodal_value(l,scalar_index) To compute the surface divergence of a vector quantity one computes a sum over nodes and coordinate directions dpsidS_div(l,i) * nodal_value(l,vector_index[i]) In Cartesian cordinates the two surface derivatives are the same, but in Axisymmetric coordinates they are not!

Implemented in oomph::ElasticUpdateFluidInterfaceElement< FluidInterfaceElement, SurfaceDerivatives, ELEMENT >, oomph::ElasticUpdateFluidInterfaceElement< SurfactantTransportInterfaceElement, AxisymmetricDerivatives, ELEMENT >, oomph::ElasticUpdateFluidInterfaceElement< FluidInterfaceElement, AxisymmetricDerivatives, ELEMENT >, oomph::ElasticUpdateFluidInterfaceElement< FluidInterfaceElement, LineDerivatives, ELEMENT >, oomph::SpineUpdateFluidInterfaceElement< SurfactantTransportInterfaceElement, LineDerivatives, ELEMENT >, oomph::SpineUpdateFluidInterfaceElement< FluidInterfaceElement, SurfaceDerivatives, ELEMENT >, oomph::SpineUpdateFluidInterfaceElement< SurfactantTransportInterfaceElement, AxisymmetricDerivatives, ELEMENT >, oomph::SpineUpdateFluidInterfaceElement< SurfactantTransportInterfaceElement, SurfaceDerivatives, ELEMENT >, oomph::SpineUpdateFluidInterfaceElement< FluidInterfaceElement, AxisymmetricDerivatives, ELEMENT >, and oomph::SpineUpdateFluidInterfaceElement< FluidInterfaceElement, LineDerivatives, ELEMENT >.

Referenced by fill_in_generic_residual_contribution_interface(), and oomph::SurfactantTransportInterfaceElement::integrate_c().

◆ fill_in_contribution_to_residuals()

void oomph::FluidInterfaceElement::fill_in_contribution_to_residuals ( Vector< double > &  residuals)
inline

Calculate the residuals by calling the generic residual contribution.

Definition at line 464 of file interface_elements.h.

References fill_in_generic_residual_contribution_interface().

◆ fill_in_generic_residual_contribution_interface()

void oomph::FluidInterfaceElement::fill_in_generic_residual_contribution_interface ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian,
unsigned  flag 
)
protectedvirtual

Helper function to calculate the residuals and (if flag==1) the Jacobian of the equations. This is implemented generically using the surface divergence information that is overloaded in each element i.e. axisymmetric, two- or three-dimensional.

Calculate the contribution to the residuals from the interface implemented generically with geometric information to be added from the specific elements.

Definition at line 471 of file interface_elements.cc.

References add_additional_residual_contributions_interface(), ca(), compute_surface_derivatives(), interpolated_u(), kinematic_local_eqn(), pext(), Pext_data_pt, pext_local_eqn(), sigma(), st(), u(), and U_index_interface.

Referenced by fill_in_contribution_to_residuals().

◆ hijack_kinematic_conditions()

virtual void oomph::FluidInterfaceElement::hijack_kinematic_conditions ( const Vector< unsigned > &  bulk_node_number)
pure virtual

◆ interpolated_u()

double oomph::FluidInterfaceElement::interpolated_u ( const Vector< double > &  s,
const unsigned &  i 
)

Calculate the i-th velocity component at the local coordinate s.

Calculate the i-th velocity component at local coordinate s.

Definition at line 442 of file interface_elements.cc.

References u().

Referenced by oomph::SurfactantTransportInterfaceElement::add_additional_residual_contributions_interface(), fill_in_generic_residual_contribution_interface(), output(), and oomph::SurfactantTransportInterfaceElement::output().

◆ kinematic_local_eqn()

virtual int oomph::FluidInterfaceElement::kinematic_local_eqn ( const unsigned &  n)
protectedpure virtual

◆ make_bounding_element()

virtual FluidInterfaceBoundingElement* oomph::FluidInterfaceElement::make_bounding_element ( const int &  face_index)
inlinevirtual

◆ output() [1/4]

void oomph::FluidInterfaceElement::output ( FILE *  file_pt)
inline

Overload the C-style output function.

Definition at line 645 of file interface_elements.h.

◆ output() [2/4]

void oomph::FluidInterfaceElement::output ( FILE *  file_pt,
const unsigned &  n_plot 
)

C-style Output function.

Overload the output function.

Definition at line 710 of file interface_elements.cc.

References interpolated_u(), and U_index_interface.

◆ output() [3/4]

void oomph::FluidInterfaceElement::output ( std::ostream &  outfile)
inline

Overload the output function.

Definition at line 636 of file interface_elements.h.

◆ output() [4/4]

void oomph::FluidInterfaceElement::output ( std::ostream &  outfile,
const unsigned &  n_plot 
)

Output function.

Overload the output functions generically.

Definition at line 672 of file interface_elements.cc.

References interpolated_u(), and U_index_interface.

◆ pext()

double oomph::FluidInterfaceElement::pext ( ) const
inline

Return the value of the external pressure.

Definition at line 519 of file interface_elements.h.

References Index_of_external_pressure_value, and Pext_data_pt.

Referenced by fill_in_generic_residual_contribution_interface().

◆ pext_local_eqn()

int oomph::FluidInterfaceElement::pext_local_eqn ( )
inlineprotected

Access function for the local equation number that corresponds to the external pressure.

Definition at line 360 of file interface_elements.h.

References External_data_number_of_external_pressure, and Index_of_external_pressure_value.

Referenced by fill_in_generic_residual_contribution_interface().

◆ set_external_pressure_data() [1/2]

void oomph::FluidInterfaceElement::set_external_pressure_data ( Data *  external_pressure_data_pt)
inline

Set the Data that contains the single pressure value that specifies the "external pressure" for the interface/free-surface. Setting this only makes sense if the interface is, in fact, a free surface (well, an interface to another inviscid fluid if you want to be picky).

Definition at line 539 of file interface_elements.h.

References External_data_number_of_external_pressure, Index_of_external_pressure_value, and Pext_data_pt.

◆ set_external_pressure_data() [2/2]

void oomph::FluidInterfaceElement::set_external_pressure_data ( Data *  external_pressure_data_pt,
const unsigned &  index_of_external_pressure_value 
)
inline

Set the Data that contains the pressure value that specifies the "external pressure" for the interface/free-surface. Setting this only makes sense if the interface is, in fact, a free surface (well, an interface to another inviscid fluid if you want to be picky). Second argument specifies the index of the pressure value within the Data object.

Definition at line 578 of file interface_elements.h.

References External_data_number_of_external_pressure, Index_of_external_pressure_value, and Pext_data_pt.

◆ sigma()

virtual double oomph::FluidInterfaceElement::sigma ( const Vector< double > &  s_local)
inlinevirtual

Virtual function that specifies the non-dimensional surface tension as a function of local position within the element. The default behaviour is a constant surface tension of value 1.0 This function can be overloaded in more specialised elements to incorporate variations in surface tension.

Reimplemented in oomph::SurfactantTransportInterfaceElement.

Definition at line 458 of file interface_elements.h.

Referenced by fill_in_generic_residual_contribution_interface().

◆ st()

const double& oomph::FluidInterfaceElement::st ( ) const
inline

The value of the Strouhal number.

Definition at line 498 of file interface_elements.h.

References St_pt.

Referenced by fill_in_generic_residual_contribution_interface().

◆ st_pt()

double*& oomph::FluidInterfaceElement::st_pt ( )
inline

The pointer to the Strouhal number.

Definition at line 504 of file interface_elements.h.

References St_pt.

◆ u()

double oomph::FluidInterfaceElement::u ( const unsigned &  j,
const unsigned &  i 
)
inline

Return the i-th velocity component at local node j.

Definition at line 510 of file interface_elements.h.

References U_index_interface.

Referenced by fill_in_generic_residual_contribution_interface(), interpolated_u(), and oomph::SurfactantTransportInterfaceElement::output().

Friends And Related Function Documentation

◆ FluidInterfaceBoundingElement

friend class FluidInterfaceBoundingElement
friend

Definition at line 323 of file interface_elements.h.

Member Data Documentation

◆ Ca_pt

double* oomph::FluidInterfaceElement::Ca_pt
private

Pointer to the Capillary number.

Definition at line 327 of file interface_elements.h.

Referenced by ca(), ca_pt(), and FluidInterfaceElement().

◆ Default_Physical_Constant_Value

double oomph::FluidInterfaceElement::Default_Physical_Constant_Value = 1.0
staticprivate

Default value for physical constants.

////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////

Default value for physical constant (static)

Definition at line 333 of file interface_elements.h.

Referenced by FluidInterfaceElement().

◆ External_data_number_of_external_pressure

int oomph::FluidInterfaceElement::External_data_number_of_external_pressure
protected

The Data that contains the external pressure is stored as external Data for the element. Which external Data item is it? (int so it can be initialised to -1, indicating that external pressure hasn't been set).

Definition at line 344 of file interface_elements.h.

Referenced by pext_local_eqn(), and set_external_pressure_data().

◆ Index_of_external_pressure_value

unsigned oomph::FluidInterfaceElement::Index_of_external_pressure_value
protected

Which of the values in Pext_data_pt stores the external pressure.

Definition at line 350 of file interface_elements.h.

Referenced by pext(), pext_local_eqn(), and set_external_pressure_data().

◆ Pext_data_pt

Data* oomph::FluidInterfaceElement::Pext_data_pt
protected

Pointer to the Data item that stores the external pressure.

Definition at line 347 of file interface_elements.h.

Referenced by fill_in_generic_residual_contribution_interface(), pext(), and set_external_pressure_data().

◆ St_pt

double* oomph::FluidInterfaceElement::St_pt
private

Pointer to the Strouhal number.

Definition at line 330 of file interface_elements.h.

Referenced by FluidInterfaceElement(), st(), and st_pt().

◆ U_index_interface

Vector<unsigned> oomph::FluidInterfaceElement::U_index_interface
protected

The documentation for this class was generated from the following files: