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

Two-layer spine mesh class derived from standard 2D mesh. The mesh contains two layers of spinified fluid elements (of type ELEMENT; e.g SpineElement<QCrouzeixRaviartElement<2>). More...

#include <two_layer_spine_mesh.template.h>

Inheritance diagram for oomph::TwoLayerSpineMesh< ELEMENT >:
oomph::RectangularQuadMesh< ELEMENT >

Public Member Functions

 TwoLayerSpineMesh (const unsigned &nx, const unsigned &ny1, const unsigned &ny2, const double &lx, const double &h1, const double &h2, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 Constructor: Pass number of elements in x-direction, number of elements in y-direction in bottom and top layer, respectively, axial length and height of top and bottom layers and pointer to timestepper (defaults to Steady timestepper) More...
 
 TwoLayerSpineMesh (const unsigned &nx, const unsigned &ny1, const unsigned &ny2, const double &lx, const double &h1, const double &h2, 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 in bottom and top layer, respectively, axial length and height of top and bottom layers, a boolean flag to make the mesh periodic in the x-direction, and pointer to timestepper (defaults to Steady timestepper) More...
 
 TwoLayerSpineMesh (const unsigned &nx, const unsigned &ny1, const unsigned &ny2, const double &lx, const double &h1, const double &h2, const bool &periodic_in_x, const bool &build_mesh, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 Constructor: Pass number of elements in x-direction, number of elements in y-direction in bottom and top layer, respectively, axial length and height of top and bottom layers, a boolean flag to make the mesh periodic in the x-direction, a boolean flag to specify whether or not to call the "build_two_layer_mesh" function, and pointer to timestepper (defaults to Steady timestepper) More...
 
FiniteElement *& upper_layer_element_pt (const unsigned long &i)
 Access functions for pointers to elements in upper layer. More...
 
FiniteElement *& lower_layer_element_pt (const unsigned long &i)
 Access functions for pointers to elements in bottom layer. More...
 
unsigned long nupper () const
 Number of elements in upper layer. More...
 
unsigned long nlower () const
 Number of elements in top layer. More...
 
FiniteElement *& interface_upper_boundary_element_pt (const unsigned long &i)
 Access functions for pointers to elements in upper layer. More...
 
FiniteElement *& interface_lower_boundary_element_pt (const unsigned long &i)
 Access functions for pointers to elements in bottom layer. More...
 
unsigned long ninterface_upper () const
 Number of elements in upper layer. More...
 
unsigned long ninterface_lower () const
 Number of elements in top layer. More...
 
int interface_upper_face_index_at_boundary (const unsigned &e)
 Index of the face of the elements next to the interface in the upper region (always -2) More...
 
int interface_lower_face_index_at_boundary (const unsigned &e)
 Index of the face of the elements next to the interface in the lower region (always 2) 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...
 
- 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...
 

Protected Member Functions

double x_spacing_function (unsigned xelement, unsigned xnode, unsigned yelement, unsigned ynode)
 The spacing function for the x co-ordinates with two regions. More...
 
double y_spacing_function (unsigned xelement, unsigned xnode, unsigned yelement, unsigned ynode)
 The spacing function for the y co-ordinates with three regions in each fluid. More...
 
void spine_node_update_lower (SpineNode *spine_node_pt)
 Update function for the lower part of the domain. More...
 
void spine_node_update_upper (SpineNode *spine_node_pt)
 Update function for the upper part of the domain. More...
 
virtual void build_two_layer_mesh (TimeStepper *time_stepper_pt)
 Helper function to actually build the two-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 Ny1
 Number of elements in lower layer. More...
 
unsigned Ny2
 Number of elements in upper layer. More...
 
double H1
 Height of the lower layer. More...
 
double H2
 Height of the upper layer. More...
 
Vector< FiniteElement * > Lower_layer_element_pt
 Vector of pointers to element in the upper layer. More...
 
Vector< FiniteElement * > Upper_layer_element_pt
 Vector of pointers to element in the lower layer. More...
 
Vector< FiniteElement * > Interface_lower_boundary_element_pt
 Vector of pointers to the elements adjacent to the interface on the lower layer. More...
 
Vector< FiniteElement * > Interface_upper_boundary_element_pt
 Vector of pointers to the element adjacent to the interface on the upper layer. 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 oomph::TwoLayerSpineMesh< ELEMENT >

Two-layer spine mesh class derived from standard 2D mesh. The mesh contains two layers of spinified fluid elements (of type ELEMENT; e.g SpineElement<QCrouzeixRaviartElement<2>).

This mesh paritions the elements into those above and below a notional interface and relabels boundaries so that the mesh is as follows

3

| | 4 | | 2

| 6 |

| | 5 | | 1

| |

0 Update information for the nodes in response to changes in spine length is given, but additional equations must be specified in order to completely specify the problem.

Definition at line 59 of file two_layer_spine_mesh.template.h.

Constructor & Destructor Documentation

◆ TwoLayerSpineMesh() [1/3]

template<class ELEMENT >
oomph::TwoLayerSpineMesh< ELEMENT >::TwoLayerSpineMesh ( const unsigned &  nx,
const unsigned &  ny1,
const unsigned &  ny2,
const double &  lx,
const double &  h1,
const double &  h2,
TimeStepper *  time_stepper_pt = &Mesh::Default_TimeStepper 
)

Constructor: Pass number of elements in x-direction, number of elements in y-direction in bottom and top layer, respectively, axial length and height of top and bottom layers and pointer to timestepper (defaults to Steady timestepper)

Constuctor for spine 2D mesh: Pass number of elements in x-direction, number of elements in y-direction in bottom and top layer, respectively, axial length and height of top and bottom layers, and pointer to timestepper (defaults to Static timestepper).

The mesh contains two layers of elements (of type ELEMENT; e.g SpineElement<QCrouzeixRaviartElement<2>) and an interfacial layer of corresponding Spine interface elements of type INTERFACE_ELEMENT, e.g. SpineLineFluidInterfaceElement<ELEMENT> for 2D planar problems.

Definition at line 49 of file two_layer_spine_mesh.template.cc.

References oomph::TwoLayerSpineMesh< ELEMENT >::build_two_layer_mesh(), oomph::TwoLayerSpineMesh< ELEMENT >::H1, oomph::TwoLayerSpineMesh< ELEMENT >::H2, oomph::TwoLayerSpineMesh< ELEMENT >::Ny1, and oomph::TwoLayerSpineMesh< ELEMENT >::Ny2.

◆ TwoLayerSpineMesh() [2/3]

template<class ELEMENT >
oomph::TwoLayerSpineMesh< ELEMENT >::TwoLayerSpineMesh ( const unsigned &  nx,
const unsigned &  ny1,
const unsigned &  ny2,
const double &  lx,
const double &  h1,
const double &  h2,
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 in bottom and top layer, respectively, axial length and height of top and bottom layers, a boolean flag to make the mesh periodic in the x-direction, and pointer to timestepper (defaults to Steady timestepper)

Constuctor for spine 2D mesh: Pass number of elements in x-direction, number of elements in y-direction in bottom and top layer, respectively, axial length and height of top and bottom layers, a boolean flag to make the mesh periodic in the x-direction, and pointer to timestepper (defaults to Static timestepper).

The mesh contains two layers of elements (of type ELEMENT; e.g SpineElement<QCrouzeixRaviartElement<2>) and an interfacial layer of corresponding Spine interface elements of type INTERFACE_ELEMENT, e.g. SpineLineFluidInterfaceElement<ELEMENT> for 2D planar problems.

Definition at line 98 of file two_layer_spine_mesh.template.cc.

References oomph::TwoLayerSpineMesh< ELEMENT >::build_two_layer_mesh(), oomph::TwoLayerSpineMesh< ELEMENT >::H1, oomph::TwoLayerSpineMesh< ELEMENT >::H2, oomph::TwoLayerSpineMesh< ELEMENT >::Ny1, and oomph::TwoLayerSpineMesh< ELEMENT >::Ny2.

◆ TwoLayerSpineMesh() [3/3]

template<class ELEMENT >
oomph::TwoLayerSpineMesh< ELEMENT >::TwoLayerSpineMesh ( const unsigned &  nx,
const unsigned &  ny1,
const unsigned &  ny2,
const double &  lx,
const double &  h1,
const double &  h2,
const bool &  periodic_in_x,
const bool &  build_mesh,
TimeStepper *  time_stepper_pt = &Mesh::Default_TimeStepper 
)

Constructor: Pass number of elements in x-direction, number of elements in y-direction in bottom and top layer, respectively, axial length and height of top and bottom layers, a boolean flag to make the mesh periodic in the x-direction, a boolean flag to specify whether or not to call the "build_two_layer_mesh" function, and pointer to timestepper (defaults to Steady timestepper)

Constuctor for spine 2D mesh: Pass number of elements in x-direction, number of elements in y-direction in bottom and top layer, respectively, axial length and height of top and bottom layers, a boolean flag to make the mesh periodic in the x-direction, a boolean flag to specify whether or not to call the "build_two_layer_mesh" function, and pointer to timestepper (defaults to Static timestepper).

The mesh contains two layers of elements (of type ELEMENT; e.g SpineElement<QCrouzeixRaviartElement<2>) and an interfacial layer of corresponding Spine interface elements of type INTERFACE_ELEMENT, e.g. SpineLineFluidInterfaceElement<ELEMENT> for 2D planar problems.

Definition at line 156 of file two_layer_spine_mesh.template.cc.

References oomph::RectangularQuadMesh< ELEMENT >::build_mesh(), oomph::TwoLayerSpineMesh< ELEMENT >::build_two_layer_mesh(), oomph::TwoLayerSpineMesh< ELEMENT >::H1, oomph::TwoLayerSpineMesh< ELEMENT >::H2, oomph::TwoLayerSpineMesh< ELEMENT >::Ny1, and oomph::TwoLayerSpineMesh< ELEMENT >::Ny2.

Member Function Documentation

◆ build_two_layer_mesh()

template<class ELEMENT >
void oomph::TwoLayerSpineMesh< ELEMENT >::build_two_layer_mesh ( TimeStepper *  time_stepper_pt)
protectedvirtual

Helper function to actually build the two-layer spine mesh (called from various constructors)

Helper function that actually builds the two-layer spine mesh based on the parameters set in the various constructors.

Definition at line 267 of file two_layer_spine_mesh.template.cc.

References oomph::RectangularQuadMesh< ELEMENT >::build_mesh().

Referenced by oomph::TwoLayerSpineMesh< ELEMENT >::TwoLayerSpineMesh().

◆ interface_lower_boundary_element_pt()

template<class ELEMENT >
FiniteElement*& oomph::TwoLayerSpineMesh< ELEMENT >::interface_lower_boundary_element_pt ( const unsigned long &  i)
inline

Access functions for pointers to elements in bottom layer.

Definition at line 141 of file two_layer_spine_mesh.template.h.

References oomph::TwoLayerSpineMesh< ELEMENT >::Interface_lower_boundary_element_pt.

◆ interface_lower_face_index_at_boundary()

template<class ELEMENT >
int oomph::TwoLayerSpineMesh< ELEMENT >::interface_lower_face_index_at_boundary ( const unsigned &  e)
inline

Index of the face of the elements next to the interface in the lower region (always 2)

Definition at line 167 of file two_layer_spine_mesh.template.h.

◆ interface_upper_boundary_element_pt()

template<class ELEMENT >
FiniteElement*& oomph::TwoLayerSpineMesh< ELEMENT >::interface_upper_boundary_element_pt ( const unsigned long &  i)
inline

Access functions for pointers to elements in upper layer.

Definition at line 135 of file two_layer_spine_mesh.template.h.

References oomph::TwoLayerSpineMesh< ELEMENT >::Interface_upper_boundary_element_pt.

◆ interface_upper_face_index_at_boundary()

template<class ELEMENT >
int oomph::TwoLayerSpineMesh< ELEMENT >::interface_upper_face_index_at_boundary ( const unsigned &  e)
inline

Index of the face of the elements next to the interface in the upper region (always -2)

Definition at line 160 of file two_layer_spine_mesh.template.h.

◆ lower_layer_element_pt()

template<class ELEMENT >
FiniteElement*& oomph::TwoLayerSpineMesh< ELEMENT >::lower_layer_element_pt ( const unsigned long &  i)
inline

Access functions for pointers to elements in bottom layer.

Definition at line 117 of file two_layer_spine_mesh.template.h.

References oomph::TwoLayerSpineMesh< ELEMENT >::Lower_layer_element_pt.

◆ ninterface_lower()

template<class ELEMENT >
unsigned long oomph::TwoLayerSpineMesh< ELEMENT >::ninterface_lower ( ) const
inline

Number of elements in top layer.

Definition at line 153 of file two_layer_spine_mesh.template.h.

References oomph::TwoLayerSpineMesh< ELEMENT >::Interface_lower_boundary_element_pt.

◆ ninterface_upper()

template<class ELEMENT >
unsigned long oomph::TwoLayerSpineMesh< ELEMENT >::ninterface_upper ( ) const
inline

Number of elements in upper layer.

Definition at line 147 of file two_layer_spine_mesh.template.h.

References oomph::TwoLayerSpineMesh< ELEMENT >::Interface_upper_boundary_element_pt.

◆ nlower()

template<class ELEMENT >
unsigned long oomph::TwoLayerSpineMesh< ELEMENT >::nlower ( ) const
inline

Number of elements in top layer.

Definition at line 129 of file two_layer_spine_mesh.template.h.

References oomph::TwoLayerSpineMesh< ELEMENT >::Lower_layer_element_pt.

◆ nupper()

template<class ELEMENT >
unsigned long oomph::TwoLayerSpineMesh< ELEMENT >::nupper ( ) const
inline

Number of elements in upper layer.

Definition at line 123 of file two_layer_spine_mesh.template.h.

References oomph::TwoLayerSpineMesh< ELEMENT >::Upper_layer_element_pt.

◆ spine_node_update()

template<class ELEMENT >
void oomph::TwoLayerSpineMesh< ELEMENT >::spine_node_update ( SpineNode *  spine_node_pt)
inline

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.

Definition at line 175 of file two_layer_spine_mesh.template.h.

References oomph::TwoLayerSpineMesh< ELEMENT >::spine_node_update_lower(), and oomph::TwoLayerSpineMesh< ELEMENT >::spine_node_update_upper().

◆ spine_node_update_lower()

template<class ELEMENT >
void oomph::TwoLayerSpineMesh< ELEMENT >::spine_node_update_lower ( SpineNode *  spine_node_pt)
inlineprotected

Update function for the lower part of the domain.

Definition at line 242 of file two_layer_spine_mesh.template.h.

References oomph::RectangularQuadMesh< ELEMENT >::Ymin.

Referenced by oomph::TwoLayerSpineMesh< ELEMENT >::spine_node_update().

◆ spine_node_update_upper()

template<class ELEMENT >
void oomph::TwoLayerSpineMesh< ELEMENT >::spine_node_update_upper ( SpineNode *  spine_node_pt)
inlineprotected

◆ upper_layer_element_pt()

template<class ELEMENT >
FiniteElement*& oomph::TwoLayerSpineMesh< ELEMENT >::upper_layer_element_pt ( const unsigned long &  i)
inline

Access functions for pointers to elements in upper layer.

Definition at line 111 of file two_layer_spine_mesh.template.h.

References oomph::TwoLayerSpineMesh< ELEMENT >::Upper_layer_element_pt.

◆ x_spacing_function()

template<class ELEMENT >
double oomph::TwoLayerSpineMesh< ELEMENT >::x_spacing_function ( unsigned  xelement,
unsigned  xnode,
unsigned  yelement,
unsigned  ynode 
)
protectedvirtual

The spacing function for the x co-ordinates with two regions.

The spacing function for the x co-ordinate, which is the same as the default function.

Reimplemented from oomph::RectangularQuadMesh< ELEMENT >.

Definition at line 211 of file two_layer_spine_mesh.template.cc.

◆ y_spacing_function()

template<class ELEMENT >
double oomph::TwoLayerSpineMesh< ELEMENT >::y_spacing_function ( unsigned  xelement,
unsigned  xnode,
unsigned  yelement,
unsigned  ynode 
)
protectedvirtual

The spacing function for the y co-ordinates with three regions in each fluid.

The spacing function for the y co-ordinates, which splits the region into two regions (1 and 2), according to the heights H1 and H2, with Ny1 and Ny2 elements respectively.

Reimplemented from oomph::RectangularQuadMesh< ELEMENT >.

Definition at line 227 of file two_layer_spine_mesh.template.cc.

Member Data Documentation

◆ H1

template<class ELEMENT >
double oomph::TwoLayerSpineMesh< ELEMENT >::H1
protected

Height of the lower layer.

Definition at line 207 of file two_layer_spine_mesh.template.h.

Referenced by oomph::TwoLayerSpineMesh< ELEMENT >::TwoLayerSpineMesh().

◆ H2

template<class ELEMENT >
double oomph::TwoLayerSpineMesh< ELEMENT >::H2
protected

Height of the upper layer.

Definition at line 210 of file two_layer_spine_mesh.template.h.

Referenced by oomph::TwoLayerSpineMesh< ELEMENT >::TwoLayerSpineMesh().

◆ Interface_lower_boundary_element_pt

template<class ELEMENT >
Vector<FiniteElement*> oomph::TwoLayerSpineMesh< ELEMENT >::Interface_lower_boundary_element_pt
protected

Vector of pointers to the elements adjacent to the interface on the lower layer.

Definition at line 220 of file two_layer_spine_mesh.template.h.

Referenced by oomph::TwoLayerSpineMesh< ELEMENT >::interface_lower_boundary_element_pt(), and oomph::TwoLayerSpineMesh< ELEMENT >::ninterface_lower().

◆ Interface_upper_boundary_element_pt

template<class ELEMENT >
Vector<FiniteElement*> oomph::TwoLayerSpineMesh< ELEMENT >::Interface_upper_boundary_element_pt
protected

Vector of pointers to the element adjacent to the interface on the upper layer.

Definition at line 224 of file two_layer_spine_mesh.template.h.

Referenced by oomph::TwoLayerSpineMesh< ELEMENT >::interface_upper_boundary_element_pt(), and oomph::TwoLayerSpineMesh< ELEMENT >::ninterface_upper().

◆ Lower_layer_element_pt

template<class ELEMENT >
Vector<FiniteElement*> oomph::TwoLayerSpineMesh< ELEMENT >::Lower_layer_element_pt
protected

Vector of pointers to element in the upper layer.

Definition at line 213 of file two_layer_spine_mesh.template.h.

Referenced by oomph::TwoLayerSpineMesh< ELEMENT >::lower_layer_element_pt(), and oomph::TwoLayerSpineMesh< ELEMENT >::nlower().

◆ Ny1

template<class ELEMENT >
unsigned oomph::TwoLayerSpineMesh< ELEMENT >::Ny1
protected

Number of elements in lower layer.

Definition at line 201 of file two_layer_spine_mesh.template.h.

Referenced by oomph::TwoLayerSpineMesh< ELEMENT >::TwoLayerSpineMesh().

◆ Ny2

template<class ELEMENT >
unsigned oomph::TwoLayerSpineMesh< ELEMENT >::Ny2
protected

Number of elements in upper layer.

Definition at line 204 of file two_layer_spine_mesh.template.h.

Referenced by oomph::TwoLayerSpineMesh< ELEMENT >::TwoLayerSpineMesh().

◆ Upper_layer_element_pt

template<class ELEMENT >
Vector<FiniteElement*> oomph::TwoLayerSpineMesh< ELEMENT >::Upper_layer_element_pt
protected

Vector of pointers to element in the lower layer.

Definition at line 216 of file two_layer_spine_mesh.template.h.

Referenced by oomph::TwoLayerSpineMesh< ELEMENT >::nupper(), and oomph::TwoLayerSpineMesh< ELEMENT >::upper_layer_element_pt().


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