Functions | Variables
oomph::FSI_functions Namespace Reference

Namespace for "global" FSI functions. More...

Functions

void apply_no_slip_on_moving_wall (Node *node_pt)
 Apply no-slip condition for N.St. on a moving wall node u = St dR/dt, where the Strouhal number St = a/(UT) is defined by FSI_functions::Strouhal_for_no_slip and is initialised to 1.0. Note: This requires the x,y,[z] velocity components to be stored in nodal values 0,1,[2]. This is the default for all currently existing Navier-Stokes elements. If you use any others, use this function at your own risk. More...
 
template<class FLUID_ELEMENT , unsigned DIM_FLUID>
void setup_fluid_load_info_for_solid_elements (Problem *problem_pt, Vector< unsigned > &boundary_in_fluid_mesh, Mesh *const &fluid_mesh_pt, Vector< Mesh * > &solid_mesh_pt, const unsigned &face=0)
 Set up the information that the FSIWallElements in the specified solid mesh require to obtain the fluid loading from the adjacent fluid elements in the specified fluid mesh. The parameter b specifies the boundary in the fluid mesh that is adjacent to the solid mesh. The template parameters specify the type of the fluid element and their spatial dimension. The optional final argument, face, identifies the face of the FSIWallElements that is exposed to the fluid. face defaults to 0, indicating that the front is loaded along the specified fluid mesh boundary. Set it to 1 to set up the FSI lookup schemes for fluid loading along the "back" of the FSIWallElements. This routine uses the procedures in the Multi_domain_functions namespace to set up the interaction by locating the adjacent (source) elements for each integration point of each solid element. More...
 
template<class FLUID_ELEMENT , unsigned DIM_FLUID>
void setup_fluid_load_info_for_solid_elements (Problem *problem_pt, const unsigned &boundary_in_fluid_mesh, Mesh *const &fluid_mesh_pt, Mesh *const &solid_mesh_pt, const unsigned &face=0)
 Set up the information that the FSIWallElements in the specified solid mesh require to obtain the fluid loading from the adjacent fluid elements in the specified fluid mesh. The parameter b specifies the boundary in the fluid mesh that is adjacent to the solid mesh. The template parameters specify the type of the fluid element and their spatial dimension. The optional final argument, face, identifies the face of the FSIWallElements that is exposed to the fluid. face defaults to 0, indicating that the front is loaded along the specified fluid mesh boundary. Set it to 1 to set up the FSI lookup schemes for fluid loading along the "back" of the FSIWallElements. This routine uses the procedures in the Multi_domain_functions namespace to set up the interaction by locating the adjacent (source) elements for each integration point of each solid element. More...
 
template<class SOLID_ELEMENT , unsigned DIM_SOLID>
void setup_solid_elements_for_displacement_bc (Problem *problem_pt, const Vector< unsigned > &b_solid_fsi, Mesh *const &solid_mesh_pt, Vector< Mesh * > &lagrange_multiplier_mesh_pt)
 Setup multi-domain interaction required for imposition of solid displacements onto the pseudo-solid fluid mesh by Lagrange multipliers: This function locates the bulk solid elements next to boundary b_solid_fsi (the FSI boundary) in the solid mesh pointed to by Solid_mesh_pt. The deformation of these elements drives the deformation of the pseudo-solid fluid mesh via the Lagrange multiplier elements stored in lagrange_multiplier_mesh_pt. The template parameters specify the type of the bulk solid elements and their spatial dimension. More...
 
template<class SOLID_ELEMENT , unsigned DIM_SOLID>
void setup_solid_elements_for_displacement_bc (Problem *problem_pt, const unsigned &b_solid_fsi, Mesh *const &solid_mesh_pt, Mesh *const &lagrange_multiplier_mesh_pt)
 Setup multi-domain interaction required for imposition of solid displacements onto the pseudo-solid fluid mesh by Lagrange multipliers: This function locates the bulk solid elements next to boundary b_solid_fsi (the FSI boundary) in the solid mesh pointed to by Solid_mesh_pt. The deformation of these elements drives the deformation of the pseudo-solid fluid mesh via the Lagrange multiplier elements stored in l lagrange_multiplier_mesh_pt. The template parameters specify the type of the bulk solid elements and their spatial dimension. More...
 
template<class NODE >
void doc_fsi (Mesh *fluid_mesh_pt, SolidMesh *wall_mesh_pt, DocInfo &doc_info)
 Doc FSI: More...
 

Variables

double Strouhal_for_no_slip = 1.0
 Strouhal number St = a/(UT) for application of no slip condition. Initialised to 1.0. More...
 

Detailed Description

Namespace for "global" FSI functions.

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

Function Documentation

◆ apply_no_slip_on_moving_wall()

void oomph::FSI_functions::apply_no_slip_on_moving_wall ( Node node_pt)

Apply no-slip condition for N.St. on a moving wall node u = St dR/dt, where the Strouhal number St = a/(UT) is defined by FSI_functions::Strouhal_for_no_slip and is initialised to 1.0. Note: This requires the x,y,[z] velocity components to be stored in nodal values 0,1,[2]. This is the default for all currently existing Navier-Stokes elements. If you use any others, use this function at your own risk.

Apply no-slip condition for N.St. on a moving wall node, u = St dR/dt, where the Strouhal number St = a/(UT) is defined by FSI_functions::Strouhal_for_no_slip and is initialised to 1.0. Note: This requires the x,y,[z] velocity components to be stored in nodal values 0,1,[2]. This is the default for all currently existing Navier-Stokes elements. If you use any others, use this function at your own risk.

Definition at line 48 of file fsi.cc.

References oomph::Node::dposition_dt(), i, oomph::Node::ndim(), oomph::Data::set_value(), and Strouhal_for_no_slip.

◆ doc_fsi()

template<class NODE >
void oomph::FSI_functions::doc_fsi ( Mesh fluid_mesh_pt,
SolidMesh wall_mesh_pt,
DocInfo doc_info 
)

Doc FSI:

  1. Which Data values affect the traction onto the FSIWallElements and what type of Data is it stored in ([fluid-]nodes, internal Data [in fluid elements] or SolidNodes?
  2. Which SolidNodes affect the node update of the fluid nodes?

Output is in tecplot readable form: Use fs1.mcr and fsi2.mcr (or straightforward modifications thereof), stored in doc/interaction/fsi_collapsible_channel/nondist_figures to process.

Pass pointer to fluid and solid meshes and pointer to the DocInfo object that specifies the directory and the overall step number. Template parameter specifies the type of node that is used to implement the node update strategy.

Definition at line 724 of file fsi.h.

References oomph::FiniteElement::dim(), oomph::DocInfo::directory(), e, oomph::Mesh::element_pt(), oomph::ElementWithExternalElement::external_element_local_coord(), oomph::ElementWithExternalElement::external_element_pt(), oomph::Mesh::external_halo_element_pt(), oomph::Mesh::external_halo_proc(), oomph::ElementWithExternalElement::external_interaction_field_data_pt(), oomph::ElementWithExternalElement::external_interaction_geometric_data_pt(), oomph::Mesh::finite_element_pt(), i, oomph::FiniteElement::integral_pt(), oomph::GeneralisedElement::internal_data_pt(), oomph::FiniteElement::interpolated_x(), oomph::GeneralisedElement::is_halo(), oomph::Integral::knot(), oomph::Node::ndim(), oomph::Mesh::nelement(), oomph::Mesh::nexternal_halo_element(), oomph::GeneralisedElement::ninternal_data(), oomph::Mesh::nnode(), oomph::FiniteElement::nodal_dimension(), oomph::Mesh::node_pt(), oomph::SolidMesh::node_pt(), oomph::DocInfo::number(), oomph::Data::nvalue(), oomph::Integral::nweight(), oomph::FSIWallElement::only_front_is_loaded_by_fluid(), oomph::Mesh::output(), s, oomph::FiniteElement::s_max(), oomph::FiniteElement::s_min(), oomph::SolidNode::variable_position_pt(), and oomph::Node::x().

◆ setup_fluid_load_info_for_solid_elements() [1/2]

template<class FLUID_ELEMENT , unsigned DIM_FLUID>
void oomph::FSI_functions::setup_fluid_load_info_for_solid_elements ( Problem problem_pt,
const unsigned &  boundary_in_fluid_mesh,
Mesh *const &  fluid_mesh_pt,
Mesh *const &  solid_mesh_pt,
const unsigned &  face = 0 
)

Set up the information that the FSIWallElements in the specified solid mesh require to obtain the fluid loading from the adjacent fluid elements in the specified fluid mesh. The parameter b specifies the boundary in the fluid mesh that is adjacent to the solid mesh. The template parameters specify the type of the fluid element and their spatial dimension. The optional final argument, face, identifies the face of the FSIWallElements that is exposed to the fluid. face defaults to 0, indicating that the front is loaded along the specified fluid mesh boundary. Set it to 1 to set up the FSI lookup schemes for fluid loading along the "back" of the FSIWallElements. This routine uses the procedures in the Multi_domain_functions namespace to set up the interaction by locating the adjacent (source) elements for each integration point of each solid element.

Definition at line 632 of file fsi.h.

References oomph::Multi_domain_functions::setup_bulk_elements_adjacent_to_face_mesh().

◆ setup_fluid_load_info_for_solid_elements() [2/2]

template<class FLUID_ELEMENT , unsigned DIM_FLUID>
void oomph::FSI_functions::setup_fluid_load_info_for_solid_elements ( Problem problem_pt,
Vector< unsigned > &  boundary_in_fluid_mesh,
Mesh *const &  fluid_mesh_pt,
Vector< Mesh * > &  solid_mesh_pt,
const unsigned &  face = 0 
)

Set up the information that the FSIWallElements in the specified solid mesh require to obtain the fluid loading from the adjacent fluid elements in the specified fluid mesh. The parameter b specifies the boundary in the fluid mesh that is adjacent to the solid mesh. The template parameters specify the type of the fluid element and their spatial dimension. The optional final argument, face, identifies the face of the FSIWallElements that is exposed to the fluid. face defaults to 0, indicating that the front is loaded along the specified fluid mesh boundary. Set it to 1 to set up the FSI lookup schemes for fluid loading along the "back" of the FSIWallElements. This routine uses the procedures in the Multi_domain_functions namespace to set up the interaction by locating the adjacent (source) elements for each integration point of each solid element.

This is the vector based version it works simultaneously on fluid fsi boundaries identified in the vector boundary_in_fluid_mesh and the corresponding solid meshes in solid_mesh_pt.

Definition at line 601 of file fsi.h.

References oomph::Multi_domain_functions::setup_bulk_elements_adjacent_to_face_mesh().

◆ setup_solid_elements_for_displacement_bc() [1/2]

template<class SOLID_ELEMENT , unsigned DIM_SOLID>
void oomph::FSI_functions::setup_solid_elements_for_displacement_bc ( Problem problem_pt,
const unsigned &  b_solid_fsi,
Mesh *const &  solid_mesh_pt,
Mesh *const &  lagrange_multiplier_mesh_pt 
)

Setup multi-domain interaction required for imposition of solid displacements onto the pseudo-solid fluid mesh by Lagrange multipliers: This function locates the bulk solid elements next to boundary b_solid_fsi (the FSI boundary) in the solid mesh pointed to by Solid_mesh_pt. The deformation of these elements drives the deformation of the pseudo-solid fluid mesh via the Lagrange multiplier elements stored in l lagrange_multiplier_mesh_pt. The template parameters specify the type of the bulk solid elements and their spatial dimension.

Definition at line 692 of file fsi.h.

References oomph::Multi_domain_functions::setup_bulk_elements_adjacent_to_face_mesh().

◆ setup_solid_elements_for_displacement_bc() [2/2]

template<class SOLID_ELEMENT , unsigned DIM_SOLID>
void oomph::FSI_functions::setup_solid_elements_for_displacement_bc ( Problem problem_pt,
const Vector< unsigned > &  b_solid_fsi,
Mesh *const &  solid_mesh_pt,
Vector< Mesh * > &  lagrange_multiplier_mesh_pt 
)

Setup multi-domain interaction required for imposition of solid displacements onto the pseudo-solid fluid mesh by Lagrange multipliers: This function locates the bulk solid elements next to boundary b_solid_fsi (the FSI boundary) in the solid mesh pointed to by Solid_mesh_pt. The deformation of these elements drives the deformation of the pseudo-solid fluid mesh via the Lagrange multiplier elements stored in lagrange_multiplier_mesh_pt. The template parameters specify the type of the bulk solid elements and their spatial dimension.

This is the vector based version it works simultaneously on solid fsi boundaries identified in the vector b_solid_fsi and the corresponding Lagrange multiplier meshes in lagrange_multiplier_mesh_pt.

Definition at line 665 of file fsi.h.

References oomph::Multi_domain_functions::setup_bulk_elements_adjacent_to_face_mesh().

Variable Documentation

◆ Strouhal_for_no_slip

double oomph::FSI_functions::Strouhal_for_no_slip = 1.0

Strouhal number St = a/(UT) for application of no slip condition. Initialised to 1.0.

Definition at line 39 of file fsi.cc.

Referenced by apply_no_slip_on_moving_wall().