Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT > Class Template Reference

Mesh for 2D Bretherton problem – based on single layer mesh. Templated by spine-ified Navier-Stokes element type (e.g. SpineElement<QCrouzeixRaviartElement<2> > and the corresponding interface element (e.g. SpineLineFluidInterfaceElement<SpineElement<QCrouzeixRaviartElement<2> > > More...

#include <bretherton_spine_mesh.template.h>

Inheritance diagram for oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >:
oomph::SingleLayerSpineMesh< ELEMENT > oomph::RectangularQuadMesh< ELEMENT >

Public Member Functions

 BrethertonSpineMesh (const unsigned &nx1, const unsigned &nx2, const unsigned &nx3, const unsigned &nh, const unsigned &nhalf, const double &h, GeomObject *lower_wall_pt, GeomObject *upper_wall_pt, const double &zeta_start, const double &zeta_transition_start, const double &zeta_transition_end, const double &zeta_end, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 Constructor. Arguments: More...
 
FiniteElement *& interface_element_pt (const unsigned long &i)
 Access functions for pointers to interface elements. More...
 
unsigned long ninterface_element () const
 Number of elements on interface. More...
 
FiniteElement *& bulk_element_pt (const unsigned long &i)
 Access functions for pointers to elements in bulk. More...
 
unsigned long nbulk () const
 Number of elements in bulk. More...
 
unsigned nfree_surface_spines ()
 Number of free-surface spines (i.e. excluding the dummy spines in the channel region) More...
 
double find_distance_to_free_surface (GeomObject *const &fs_geom_object_pt, Vector< double > &initial_zeta, const Vector< double > &spine_base, const Vector< double > &spine_end)
 Recalculate the spine lengths after repositioning. More...
 
void reposition_spines (const double &zeta_lo_transition_start, const double &zeta_lo_transition_end, const double &zeta_up_transition_start, const double &zeta_up_transition_end)
 Reposition the spines in response to changes in geometry. More...
 
void pin_all_spines ()
 Pin all spines so the mesh can be used for computation without free surfaces. More...
 
void spine_node_update (SpineNode *spine_node_pt)
 General node update function implements pure virtual function defined in SpineMesh base class and performs specific update actions, depending on the node update fct id stored for each node. More...
 
ELEMENT * control_element_pt ()
 Pointer to control element (just under the symmetry line near the bubble tip, so the bubble tip is located at s=[1.0,1.0] in this element. More...
 
double spine_centre_fraction () const
 Read the value of the spine centre's vertical fraction. More...
 
void set_spine_centre_fraction_pt (double *const &fraction_pt)
 Set the pointer to the spine centre's vertial fraction. More...
 
- Public Member Functions inherited from oomph::SingleLayerSpineMesh< ELEMENT >
 SingleLayerSpineMesh (const unsigned &nx, const unsigned &ny, const double &lx, const double &h, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 Constructor: Pass number of elements in x-direction, number of elements in y-direction, axial length, height of layer, and pointer to timestepper (defaults to Steady timestepper) More...
 
 SingleLayerSpineMesh (const unsigned &nx, const unsigned &ny, const double &lx, const double &h, const bool &periodic_in_x, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 Constructor: Pass number of elements in x-direction, number of elements in y-direction, axial length, height of layer, a boolean flag to make the mesh periodic in the x-direction, and a pointer to timestepper (defaults to Steady timestepper) More...
 
- Public Member Functions inherited from oomph::RectangularQuadMesh< ELEMENT >
 RectangularQuadMesh (const unsigned &nx, const unsigned &ny, const double &lx, const double &ly, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 Simple constructor: nx: number of elements in x direction; ny: number of elements in y direction; lx, length of domain in x direction (0,lx); ly, length of domain in y direction (0,ly) Also pass pointer to timestepper (defaults to Steady) More...
 
 RectangularQuadMesh (const unsigned &nx, const unsigned &ny, const double &xmin, const double &xmax, const double &ymin, const double &ymax, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 Constructor that allows the specification of minimum and maximum values of x and y coordinates. More...
 
 RectangularQuadMesh (const unsigned &nx, const unsigned &ny, const double &lx, const double &ly, const bool &periodic_in_x, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 Simple constructor: nx: number of elements in x direction; ny: number of elements in y direction; lx, length of domain in x direction (0,lx); ly, length of domain in y direction (0,ly) Boolean flag specifies if the mesh is periodic in the x-direction. Also pass pointer to timestepper (defaults to Steady) More...
 
 RectangularQuadMesh (const unsigned &nx, const unsigned &ny, const double &xmin, const double &xmax, const double &ymin, const double &ymax, const bool &periodic_in_x, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 Constructor that allows the specification of minimum and maximum values of x and y coordinates. Boolean flag specifies if the mesh is periodic in the x-direction. More...
 
const unsigned & nx () const
 Return number of elements in x direction. More...
 
const unsigned & ny () const
 Return number of elements in y direction. More...
 
const double x_min () const
 Return the minimum value of x coordinate. More...
 
const double x_max () const
 Return the maximum value of x coordinate. More...
 
const double y_min () const
 Return the minimum value of y coordinate. More...
 
const double y_max () const
 Return the maximum value of y coordinate. More...
 
virtual void element_reorder ()
 Reorder the elements: By default they are ordered in "horizontal" layers (increasing in x, then in y). This function changes this to an ordering in the vertical direction (y first, then x). This is more efficient if a frontal solver is used and the mesh has more elements in the x than the y direction. Can be overloaded in specific derived meshes. More...
 
virtual double x_spacing_function (unsigned xelement, unsigned xnode, unsigned yelement, unsigned ynode)
 Return the value of the x-coordinate at the node given by the local node number (xnode, ynode) in the element (xelement,yelement). The description is in a "psudeo" two-dimensional coordinate system, so the range of xelement is [0,Nx-1], yelement is [0,Ny-1], and that of xnode and ynode is [0,Np-1]. The default is to return nodes that are equally spaced in the x coodinate. More...
 
virtual double y_spacing_function (unsigned xelement, unsigned xnode, unsigned yelement, unsigned ynode)
 Return the value of the y-coordinate at the node given by the local node number (xnode, ynode) in the element (xelement,yelement). The description is in a "psudeo" two-dimensional coordinate system, so the range of xelement is [0,Nx-1], yelement is [0,Ny-1], and that of xnode and ynode is [0,Np-1]. The default it to return nodes that are equally spaced in the y coordinate. More...
 

Protected Member Functions

void spine_node_update_film_lower (SpineNode *spine_node_pt)
 Update function for the deposited film region in the lower part of the domain: Vertical spines. More...
 
void spine_node_update_horizontal_transition_lower (SpineNode *spine_node_pt)
 Update function for the horizontal transitition region in the lower part of the domain: Spine points from wall to origin. More...
 
void spine_node_update_vertical_transition_lower (SpineNode *spine_node_pt)
 Update function for the vertical transitition region in the lower part of the domain: Spine points to origin. More...
 
void spine_node_update_vertical_transition_upper (SpineNode *spine_node_pt)
 Update function for the vertical transitition region in the upper part of the domain: Spine points to origin. More...
 
void spine_node_update_horizontal_transition_upper (SpineNode *spine_node_pt)
 Update function for the horizontal transitition region in the upper part of the domain: Spine points towards origin. More...
 
void spine_node_update_film_upper (SpineNode *spine_node_pt)
 Update function for the deposited film region in the upper part of the domain: Vertical spines. More...
 
void spine_node_update_channel (SpineNode *spine_node_pt)
 Update function for the nodes in the channel region ahead of the finger tip: Nodes are evenly distributed along vertical lines between the top and bottom walls. More...
 
void initial_element_reorder ()
 Initial reordering elements of the elements – before the channel mesh is added. Vertical stacks of elements, each topped by their interface element – this is (currently) identical to the version in the SingleLayerSpineMesh but it's important that the element reordering is maintained in exactly this form for the rest of the mesh generation process to work properly. Therefore we keep a copy of the function in here. More...
 
- Protected Member Functions inherited from oomph::SingleLayerSpineMesh< ELEMENT >
virtual void build_single_layer_mesh (TimeStepper *time_stepper_pt)
 Helper function to actually build the single-layer spine mesh (called from various constructors) More...
 
- Protected Member Functions inherited from oomph::RectangularQuadMesh< ELEMENT >
void build_mesh (TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 Generic mesh construction function: contains all the hard work. More...
 
 RectangularQuadMesh (const unsigned &nx, const unsigned &ny, const double &xmin, const double &xmax, const double &ymin, const double &ymax, const bool &periodic_in_x, const bool &build, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 Constructor that allows the specification of minimum and maximum values of x and y coordinates and does not build the mesh This is intend to be used in derived classes that overload the spacing functions. THis is scheduled to be changed, however. The reason why this MUST be done is because the virtual spacing functions cannot be called in the base constructur, because they will not have been overloaded yet!! More...
 

Protected Attributes

unsigned Nx1
 Number of elements along wall in deposited film region. More...
 
unsigned Nx2
 Number of elements along wall in horizontal transition region. More...
 
unsigned Nx3
 Number of elements along wall in channel region. More...
 
unsigned Nhalf
 Number of elements in vertical transition region (there are twice as many elements across the channel) More...
 
unsigned Nh
 Number of elements across the deposited film. More...
 
double H
 Thickness of deposited film. More...
 
GeomObject * Upper_wall_pt
 Pointer to geometric object that represents the upper wall. More...
 
GeomObject * Lower_wall_pt
 Pointer to geometric object that represents the lower wall. More...
 
double Zeta_start
 Start coordinate on wall. More...
 
double Zeta_end
 Wall coordinate of end of liquid filled region (inflow) More...
 
double Zeta_transition_start
 Wall coordinate of start of the transition region. More...
 
double Zeta_transition_end
 Wall coordinate of end of transition region. More...
 
double * Spine_centre_fraction_pt
 Pointer to vertical fraction of the spine centre. More...
 
double Default_spine_centre_fraction
 Default spine fraction. More...
 
ELEMENT * Control_element_pt
 Pointer to control element (just under the symmetry line near the bubble tip; the bubble tip is located at s=[1.0,1.0] in this element. More...
 
Vector< FiniteElement * > Bulk_element_pt
 Vector of pointers to element in the fluid layer. More...
 
Vector< FiniteElement * > Interface_element_pt
 Vector of pointers to interface elements. More...
 
- Protected Attributes inherited from oomph::RectangularQuadMesh< ELEMENT >
unsigned Nx
 Nx: number of elements in x-direction. More...
 
unsigned Ny
 Ny: number of elements in y-direction. More...
 
unsigned Np
 Np: number of (linear) points in the element. More...
 
double Xmin
 Minimum value of x coordinate. More...
 
double Xmax
 Maximum value of x coordinate. More...
 
double Ymin
 Minimum value of y coordinate. More...
 
double Ymax
 Maximum value of y coordinate. More...
 
bool Xperiodic
 Boolean variable used to determine whether the mesh is periodic in the x-direction. More...
 

Detailed Description

template<class ELEMENT, class INTERFACE_ELEMENT>
class oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >

Mesh for 2D Bretherton problem – based on single layer mesh. Templated by spine-ified Navier-Stokes element type (e.g. SpineElement<QCrouzeixRaviartElement<2> > and the corresponding interface element (e.g. SpineLineFluidInterfaceElement<SpineElement<QCrouzeixRaviartElement<2> > >

Definition at line 45 of file bretherton_spine_mesh.template.h.

Constructor & Destructor Documentation

◆ BrethertonSpineMesh()

template<class ELEMENT , class INTERFACE_ELEMENT >
oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::BrethertonSpineMesh ( const unsigned &  nx1,
const unsigned &  nx2,
const unsigned &  nx3,
const unsigned &  nh,
const unsigned &  nhalf,
const double &  h,
GeomObject *  lower_wall_pt,
GeomObject *  upper_wall_pt,
const double &  zeta_start,
const double &  zeta_transition_start,
const double &  zeta_transition_end,
const double &  zeta_end,
TimeStepper *  time_stepper_pt = &Mesh::Default_TimeStepper 
)

Constructor. Arguments:

  • nx1: Number of elements along wall in deposited film region
  • nx2: Number of elements along wall in horizontal transition region
  • nx3: Number of elements along wall in channel region
  • nhalf: Number of elements in vertical transition region (there are twice as many elements across the channel)
  • nh: Number of elements across the deposited film
  • h: Thickness of deposited film
  • zeta0: Start coordinate on wall
  • zeta1: Wall coordinate of start of transition region
  • zeta2: Wall coordinate of end of liquid filled region (inflow)
  • lower_wall_pt: Pointer to geometric object that represents the lower wall
  • upper_wall_pt: Pointer to geometric object that represents the upper wall
  • time_stepper_pt: Pointer to timestepper; defaults to Steady
  • nx1: Number of elements along wall in deposited film region
  • nx2: Number of elements along wall in horizontal transition region
  • nx3: Number of elements along wall in channel region
  • nhalf: Number of elements in vertical transition region (there are twice as many elements across the channel)
  • nh: Number of elements across the deposited film
  • h: Thickness of deposited film
  • zeta0: Start coordinate on wall
  • zeta1: Wall coordinate of start of transition region
  • zeta2: Wall coordinate of end of liquid filled region (inflow)
  • lower_wall_pt: Pointer to geometric object that represents the lower wall
  • upper_wall_pt: Pointer to geometric object that represents the upper wall
  • time_stepper_pt: Pointer to timestepper; defaults to Static

Definition at line 58 of file bretherton_spine_mesh.template.cc.

References oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Bulk_element_pt, oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Control_element_pt, oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::H, oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::initial_element_reorder(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Interface_element_pt, oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Lower_wall_pt, oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Nx3, oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_centre_fraction(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Upper_wall_pt, oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Zeta_end, oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Zeta_start, oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Zeta_transition_end, and oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Zeta_transition_start.

Member Function Documentation

◆ bulk_element_pt()

template<class ELEMENT , class INTERFACE_ELEMENT >
FiniteElement*& oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::bulk_element_pt ( const unsigned long &  i)
inline

Access functions for pointers to elements in bulk.

Definition at line 93 of file bretherton_spine_mesh.template.h.

References oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Bulk_element_pt.

◆ control_element_pt()

template<class ELEMENT , class INTERFACE_ELEMENT >
ELEMENT* oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::control_element_pt ( )
inline

Pointer to control element (just under the symmetry line near the bubble tip, so the bubble tip is located at s=[1.0,1.0] in this element.

Definition at line 187 of file bretherton_spine_mesh.template.h.

References oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Control_element_pt.

◆ find_distance_to_free_surface()

template<class ELEMENT , class INTERFACE_ELEMENT >
double oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::find_distance_to_free_surface ( GeomObject *const &  fs_geom_object_pt,
Vector< double > &  initial_zeta,
const Vector< double > &  spine_base,
const Vector< double > &  spine_end 
)

Recalculate the spine lengths after repositioning.

Calculate the distance from the spine base to the free surface, i.e. the spine height.

Definition at line 1065 of file bretherton_spine_mesh.template.cc.

◆ initial_element_reorder()

template<class ELEMENT , class INTERFACE_ELEMENT >
void oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::initial_element_reorder
protected

Initial reordering elements of the elements – before the channel mesh is added. Vertical stacks of elements, each topped by their interface element – this is (currently) identical to the version in the SingleLayerSpineMesh but it's important that the element reordering is maintained in exactly this form for the rest of the mesh generation process to work properly. Therefore we keep a copy of the function in here.

Reorder elements: Vertical stacks of elements, each topped by their interface element – this is (currently) identical to the version in the SingleLayerSpineMesh but it's important that element reordering is maintained in exactly this form so to be on the safe side, we move the function in here.

Definition at line 1028 of file bretherton_spine_mesh.template.cc.

Referenced by oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::BrethertonSpineMesh().

◆ interface_element_pt()

template<class ELEMENT , class INTERFACE_ELEMENT >
FiniteElement*& oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::interface_element_pt ( const unsigned long &  i)
inline

Access functions for pointers to interface elements.

Definition at line 81 of file bretherton_spine_mesh.template.h.

References oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Interface_element_pt.

◆ nbulk()

template<class ELEMENT , class INTERFACE_ELEMENT >
unsigned long oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::nbulk ( ) const
inline

Number of elements in bulk.

Definition at line 99 of file bretherton_spine_mesh.template.h.

References oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Bulk_element_pt.

◆ nfree_surface_spines()

template<class ELEMENT , class INTERFACE_ELEMENT >
unsigned oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::nfree_surface_spines ( )
inline

◆ ninterface_element()

template<class ELEMENT , class INTERFACE_ELEMENT >
unsigned long oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::ninterface_element ( ) const
inline

Number of elements on interface.

Definition at line 87 of file bretherton_spine_mesh.template.h.

References oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Interface_element_pt.

◆ pin_all_spines()

template<class ELEMENT , class INTERFACE_ELEMENT >
void oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::pin_all_spines ( )
inline

Pin all spines so the mesh can be used for computation without free surfaces.

Definition at line 128 of file bretherton_spine_mesh.template.h.

◆ reposition_spines()

template<class ELEMENT , class INTERFACE_ELEMENT >
void oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::reposition_spines ( const double &  zeta_lo_transition_start,
const double &  zeta_lo_transition_end,
const double &  zeta_up_transition_start,
const double &  zeta_up_transition_end 
)

Reposition the spines in response to changes in geometry.

Reposition the spines that emenate from the lower wall.

Definition at line 1158 of file bretherton_spine_mesh.template.cc.

◆ set_spine_centre_fraction_pt()

template<class ELEMENT , class INTERFACE_ELEMENT >
void oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::set_spine_centre_fraction_pt ( double *const &  fraction_pt)
inline

Set the pointer to the spine centre's vertial fraction.

Definition at line 200 of file bretherton_spine_mesh.template.h.

References oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Spine_centre_fraction_pt.

◆ spine_centre_fraction()

template<class ELEMENT , class INTERFACE_ELEMENT >
double oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_centre_fraction ( ) const
inline

◆ spine_node_update()

template<class ELEMENT , class INTERFACE_ELEMENT >
void oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update ( SpineNode *  spine_node_pt)
inlinevirtual

◆ spine_node_update_channel()

template<class ELEMENT , class INTERFACE_ELEMENT >
void oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_channel ( SpineNode *  spine_node_pt)
inlineprotected

Update function for the nodes in the channel region ahead of the finger tip: Nodes are evenly distributed along vertical lines between the top and bottom walls.

Definition at line 464 of file bretherton_spine_mesh.template.h.

Referenced by oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update().

◆ spine_node_update_film_lower()

template<class ELEMENT , class INTERFACE_ELEMENT >
void oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_film_lower ( SpineNode *  spine_node_pt)
inlineprotected

Update function for the deposited film region in the lower part of the domain: Vertical spines.

Definition at line 208 of file bretherton_spine_mesh.template.h.

Referenced by oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update().

◆ spine_node_update_film_upper()

template<class ELEMENT , class INTERFACE_ELEMENT >
void oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_film_upper ( SpineNode *  spine_node_pt)
inlineprotected

Update function for the deposited film region in the upper part of the domain: Vertical spines.

Definition at line 440 of file bretherton_spine_mesh.template.h.

Referenced by oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update().

◆ spine_node_update_horizontal_transition_lower()

template<class ELEMENT , class INTERFACE_ELEMENT >
void oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_horizontal_transition_lower ( SpineNode *  spine_node_pt)
inlineprotected

Update function for the horizontal transitition region in the lower part of the domain: Spine points from wall to origin.

Definition at line 231 of file bretherton_spine_mesh.template.h.

References oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_centre_fraction().

Referenced by oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update().

◆ spine_node_update_horizontal_transition_upper()

template<class ELEMENT , class INTERFACE_ELEMENT >
void oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_horizontal_transition_upper ( SpineNode *  spine_node_pt)
inlineprotected

Update function for the horizontal transitition region in the upper part of the domain: Spine points towards origin.

Definition at line 393 of file bretherton_spine_mesh.template.h.

References oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_centre_fraction().

Referenced by oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update().

◆ spine_node_update_vertical_transition_lower()

template<class ELEMENT , class INTERFACE_ELEMENT >
void oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_vertical_transition_lower ( SpineNode *  spine_node_pt)
inlineprotected

Update function for the vertical transitition region in the lower part of the domain: Spine points to origin.

Definition at line 277 of file bretherton_spine_mesh.template.h.

References oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_centre_fraction().

Referenced by oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update().

◆ spine_node_update_vertical_transition_upper()

template<class ELEMENT , class INTERFACE_ELEMENT >
void oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update_vertical_transition_upper ( SpineNode *  spine_node_pt)
inlineprotected

Update function for the vertical transitition region in the upper part of the domain: Spine points to origin.

Definition at line 335 of file bretherton_spine_mesh.template.h.

References oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_centre_fraction().

Referenced by oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::spine_node_update().

Member Data Documentation

◆ Bulk_element_pt

template<class ELEMENT , class INTERFACE_ELEMENT >
Vector<FiniteElement*> oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Bulk_element_pt
protected

◆ Control_element_pt

template<class ELEMENT , class INTERFACE_ELEMENT >
ELEMENT* oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Control_element_pt
protected

Pointer to control element (just under the symmetry line near the bubble tip; the bubble tip is located at s=[1.0,1.0] in this element.

Definition at line 540 of file bretherton_spine_mesh.template.h.

Referenced by oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::BrethertonSpineMesh(), and oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::control_element_pt().

◆ Default_spine_centre_fraction

template<class ELEMENT , class INTERFACE_ELEMENT >
double oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Default_spine_centre_fraction
protected

Default spine fraction.

Definition at line 535 of file bretherton_spine_mesh.template.h.

◆ H

template<class ELEMENT , class INTERFACE_ELEMENT >
double oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::H
protected

Thickness of deposited film.

Definition at line 511 of file bretherton_spine_mesh.template.h.

Referenced by oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::BrethertonSpineMesh().

◆ Interface_element_pt

template<class ELEMENT , class INTERFACE_ELEMENT >
Vector<FiniteElement*> oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Interface_element_pt
protected

◆ Lower_wall_pt

template<class ELEMENT , class INTERFACE_ELEMENT >
GeomObject* oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Lower_wall_pt
protected

Pointer to geometric object that represents the lower wall.

Definition at line 517 of file bretherton_spine_mesh.template.h.

Referenced by oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::BrethertonSpineMesh().

◆ Nh

template<class ELEMENT , class INTERFACE_ELEMENT >
unsigned oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Nh
protected

Number of elements across the deposited film.

Definition at line 508 of file bretherton_spine_mesh.template.h.

◆ Nhalf

template<class ELEMENT , class INTERFACE_ELEMENT >
unsigned oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Nhalf
protected

Number of elements in vertical transition region (there are twice as many elements across the channel)

Definition at line 505 of file bretherton_spine_mesh.template.h.

Referenced by oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::nfree_surface_spines().

◆ Nx1

template<class ELEMENT , class INTERFACE_ELEMENT >
unsigned oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Nx1
protected

Number of elements along wall in deposited film region.

Definition at line 495 of file bretherton_spine_mesh.template.h.

Referenced by oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::nfree_surface_spines().

◆ Nx2

template<class ELEMENT , class INTERFACE_ELEMENT >
unsigned oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Nx2
protected

Number of elements along wall in horizontal transition region.

Definition at line 498 of file bretherton_spine_mesh.template.h.

Referenced by oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::nfree_surface_spines().

◆ Nx3

template<class ELEMENT , class INTERFACE_ELEMENT >
unsigned oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Nx3
protected

Number of elements along wall in channel region.

Definition at line 501 of file bretherton_spine_mesh.template.h.

Referenced by oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::BrethertonSpineMesh().

◆ Spine_centre_fraction_pt

template<class ELEMENT , class INTERFACE_ELEMENT >
double* oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Spine_centre_fraction_pt
protected

◆ Upper_wall_pt

template<class ELEMENT , class INTERFACE_ELEMENT >
GeomObject* oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Upper_wall_pt
protected

Pointer to geometric object that represents the upper wall.

Definition at line 514 of file bretherton_spine_mesh.template.h.

Referenced by oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::BrethertonSpineMesh().

◆ Zeta_end

template<class ELEMENT , class INTERFACE_ELEMENT >
double oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Zeta_end
protected

Wall coordinate of end of liquid filled region (inflow)

Definition at line 523 of file bretherton_spine_mesh.template.h.

Referenced by oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::BrethertonSpineMesh().

◆ Zeta_start

template<class ELEMENT , class INTERFACE_ELEMENT >
double oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Zeta_start
protected

◆ Zeta_transition_end

template<class ELEMENT , class INTERFACE_ELEMENT >
double oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Zeta_transition_end
protected

Wall coordinate of end of transition region.

Definition at line 529 of file bretherton_spine_mesh.template.h.

Referenced by oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::BrethertonSpineMesh().

◆ Zeta_transition_start

template<class ELEMENT , class INTERFACE_ELEMENT >
double oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::Zeta_transition_start
protected

Wall coordinate of start of the transition region.

Definition at line 526 of file bretherton_spine_mesh.template.h.

Referenced by oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::BrethertonSpineMesh().


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