oomph::TwoLayerPerturbedSpineMesh< ELEMENT > Class Template Reference

#include <two_layer_perturbed_spine_mesh.template.h>

+ Inheritance diagram for oomph::TwoLayerPerturbedSpineMesh< ELEMENT >:

Public Member Functions

 TwoLayerPerturbedSpineMesh (const unsigned &nx, const unsigned &ny1, const unsigned &ny2, const double &lx, const double &h1, const double &h2, SpineMesh *&base_mesh_pt, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
 TwoLayerPerturbedSpineMesh (const unsigned &nx, const unsigned &ny1, const unsigned &ny2, const double &lx, const double &h1, const double &h2, const bool &periodic_in_x, SpineMesh *&base_mesh_pt, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
 TwoLayerPerturbedSpineMesh (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, SpineMesh *&base_mesh_pt, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
void element_reorder ()
 
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 lower layer. More...
 
unsigned long nupper () const
 Number of elements in upper layer. More...
 
unsigned long nlower () const
 Number of elements in top layer. More...
 
void perturbed_spine_node_update (PerturbedSpineNode *spine_node_pt)
 
const unsignedny1 () const
 Access function for number of elements in lower layer. More...
 
const unsignedny2 () const
 Access function for number of elements in upper layer. More...
 
const void set_perturbation_to_nodal_positions_indices (const unsigned &cosine_index, const unsigned &sine_index)
 
- 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)
 
 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)
 
 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)
 
 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)
 
const unsignednx () const
 Return number of elements in x direction. More...
 
const unsignedny () 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...
 
- Public Member Functions inherited from oomph::QuadMeshBase
 QuadMeshBase ()
 Constructor (empty) More...
 
 QuadMeshBase (const QuadMeshBase &node)=delete
 Broken copy constructor. More...
 
void operator= (const QuadMeshBase &)=delete
 Broken assignment operator. More...
 
virtual ~QuadMeshBase ()
 Destructor (empty) More...
 
void setup_boundary_element_info ()
 
void setup_boundary_element_info (std::ostream &outfile)
 
- Public Member Functions inherited from oomph::Mesh
 Mesh ()
 Default constructor. More...
 
 Mesh (const Vector< Mesh * > &sub_mesh_pt)
 
void merge_meshes (const Vector< Mesh * > &sub_mesh_pt)
 
virtual void reset_boundary_element_info (Vector< unsigned > &ntmp_boundary_elements, Vector< Vector< unsigned >> &ntmp_boundary_elements_in_region, Vector< FiniteElement * > &deleted_elements)
 Virtual function to perform the reset boundary elements info rutines. More...
 
template<class BULK_ELEMENT >
void doc_boundary_coordinates (const unsigned &b, std::ofstream &the_file)
 
virtual void scale_mesh (const double &factor)
 
 Mesh (const Mesh &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const Mesh &)=delete
 Broken assignment operator. More...
 
virtual ~Mesh ()
 Virtual Destructor to clean up all memory. More...
 
void flush_element_and_node_storage ()
 
void flush_element_storage ()
 
void flush_node_storage ()
 
Node *& node_pt (const unsigned long &n)
 Return pointer to global node n. More...
 
Nodenode_pt (const unsigned long &n) const
 Return pointer to global node n (const version) More...
 
GeneralisedElement *& element_pt (const unsigned long &e)
 Return pointer to element e. More...
 
GeneralisedElementelement_pt (const unsigned long &e) const
 Return pointer to element e (const version) More...
 
const Vector< GeneralisedElement * > & element_pt () const
 Return reference to the Vector of elements. More...
 
Vector< GeneralisedElement * > & element_pt ()
 Return reference to the Vector of elements. More...
 
FiniteElementfinite_element_pt (const unsigned &e) const
 
Node *& boundary_node_pt (const unsigned &b, const unsigned &n)
 Return pointer to node n on boundary b. More...
 
Nodeboundary_node_pt (const unsigned &b, const unsigned &n) const
 Return pointer to node n on boundary b. More...
 
void set_nboundary (const unsigned &nbound)
 Set the number of boundaries in the mesh. More...
 
void remove_boundary_nodes ()
 Clear all pointers to boundary nodes. More...
 
void remove_boundary_nodes (const unsigned &b)
 
void remove_boundary_node (const unsigned &b, Node *const &node_pt)
 Remove a node from the boundary b. More...
 
void add_boundary_node (const unsigned &b, Node *const &node_pt)
 Add a (pointer to) a node to the b-th boundary. More...
 
void copy_boundary_node_data_from_nodes ()
 
bool boundary_coordinate_exists (const unsigned &i) const
 Indicate whether the i-th boundary has an intrinsic coordinate. More...
 
unsigned long nelement () const
 Return number of elements in the mesh. More...
 
unsigned long nnode () const
 Return number of nodes in the mesh. More...
 
unsigned ndof_types () const
 Return number of dof types in mesh. More...
 
unsigned elemental_dimension () const
 Return number of elemental dimension in mesh. More...
 
unsigned nodal_dimension () const
 Return number of nodal dimension in mesh. More...
 
void add_node_pt (Node *const &node_pt)
 Add a (pointer to a) node to the mesh. More...
 
void add_element_pt (GeneralisedElement *const &element_pt)
 Add a (pointer to) an element to the mesh. More...
 
virtual void reorder_nodes (const bool &use_old_ordering=true)
 
virtual void get_node_reordering (Vector< Node * > &reordering, const bool &use_old_ordering=true) const
 
template<class BULK_ELEMENT , template< class > class FACE_ELEMENT>
void build_face_mesh (const unsigned &b, Mesh *const &face_mesh_pt)
 
unsigned self_test ()
 Self-test: Check elements and nodes. Return 0 for OK. More...
 
void max_and_min_element_size (double &max_size, double &min_size)
 
double total_size ()
 
void check_inverted_elements (bool &mesh_has_inverted_elements, std::ofstream &inverted_element_file)
 
void check_inverted_elements (bool &mesh_has_inverted_elements)
 
unsigned check_for_repeated_nodes (const double &epsilon=1.0e-12)
 
Vector< Node * > prune_dead_nodes ()
 
unsigned nboundary () const
 Return number of boundaries. More...
 
unsigned long nboundary_node (const unsigned &ibound) const
 Return number of nodes on a particular boundary. More...
 
FiniteElementboundary_element_pt (const unsigned &b, const unsigned &e) const
 Return pointer to e-th finite element on boundary b. More...
 
Nodeget_some_non_boundary_node () const
 
unsigned nboundary_element (const unsigned &b) const
 Return number of finite elements that are adjacent to boundary b. More...
 
int face_index_at_boundary (const unsigned &b, const unsigned &e) const
 
void dump (const std::string &dump_file_name, const bool &use_old_ordering=true) const
 Dump the data in the mesh into a file for restart. More...
 
void output_paraview (std::ofstream &file_out, const unsigned &nplot) const
 
void output_fct_paraview (std::ofstream &file_out, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const
 
void output_fct_paraview (std::ofstream &file_out, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) const
 
void output (std::ostream &outfile)
 Output for all elements. More...
 
void output (std::ostream &outfile, const unsigned &n_plot)
 Output at f(n_plot) points in each element. More...
 
void output (FILE *file_pt)
 Output for all elements (C-style output) More...
 
void output (FILE *file_pt, const unsigned &nplot)
 Output at f(n_plot) points in each element (C-style output) More...
 
void output (const std::string &output_filename)
 Output for all elements. More...
 
void output (const std::string &output_filename, const unsigned &n_plot)
 Output at f(n_plot) points in each element. More...
 
void output_fct (std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt)
 Output a given Vector function at f(n_plot) points in each element. More...
 
void output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt)
 
void output_boundaries (std::ostream &outfile)
 Output the nodes on the boundaries (into separate tecplot zones) More...
 
void output_boundaries (const std::string &output_filename)
 
void assign_initial_values_impulsive ()
 Assign initial values for an impulsive start. More...
 
void shift_time_values ()
 
void calculate_predictions ()
 
void set_nodal_and_elemental_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
virtual void set_mesh_level_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
void set_consistent_pinned_values_for_continuation (ContinuationStorageScheme *const &continuation_stepper_pt)
 Set consistent values for pinned data in continuation. More...
 
bool does_pointer_correspond_to_mesh_data (double *const &parameter_pt)
 Does the double pointer correspond to any mesh data. More...
 
void set_nodal_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 Set the timestepper associated with the nodal data in the mesh. More...
 
void set_elemental_internal_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
virtual void compute_norm (double &norm)
 
virtual void compute_norm (Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 Returns the norm of the error and that of the exact solution. More...
 
virtual void compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
bool is_mesh_distributed () const
 Boolean to indicate if Mesh has been distributed. More...
 
OomphCommunicatorcommunicator_pt () const
 
void delete_all_external_storage ()
 Wipe the storage for all externally-based elements. More...
 
- Public Member Functions inherited from oomph::PerturbedSpineMesh
virtual ~PerturbedSpineMesh ()
 Destructor to clean up the memory allocated to the perturbed spines. More...
 
PerturbedSpine *const & perturbed_spine_pt (const unsigned long &i) const
 Return the i-th perturbed spine in the mesh. More...
 
unsigned long nspine () const
 Return the number of perturbed spines in the mesh. More...
 
void add_perturbed_spine_pt (PerturbedSpine *const &spine_pt)
 Add a perturbed spine to the mesh. More...
 
PerturbedSpineNodenode_pt (const unsigned long &n)
 Return a pointer to the n-th global PerturbedSpineNode. More...
 
PerturbedSpineNodeelement_node_pt (const unsigned long &e, const unsigned &n)
 
unsigned long assign_global_eqn_numbers (Vector< double * > &Dof_pt)
 Assign equation numbers for perturbed spines. More...
 
void node_update (const bool &update_all_solid_nodes=false)
 
void dump (std::ofstream &dump_file, const bool &use_old_ordering=true) const
 Overload the dump function so that the spine data is dumped. More...
 
void read (std::ifstream &restart_file)
 Overload the read function so that the spine data is also read. More...
 

Protected Member Functions

double x_spacing_function (unsigned xelement, unsigned xnode, unsigned yelement, unsigned ynode)
 
double y_spacing_function (unsigned xelement, unsigned xnode, unsigned yelement, unsigned ynode)
 
void spine_node_update_lower (PerturbedSpineNode *spine_node_pt)
 Update function for the lower part of the domain. More...
 
void spine_node_update_upper (PerturbedSpineNode *spine_node_pt)
 Update function for the upper part of the domain. More...
 
virtual void build_two_layer_mesh (TimeStepper *time_stepper_pt)
 
- 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)
 
- Protected Member Functions inherited from oomph::Mesh
unsigned long assign_global_eqn_numbers (Vector< double * > &Dof_pt)
 Assign (global) equation numbers to the nodes. More...
 
void describe_dofs (std::ostream &out, const std::string &current_string) const
 
void describe_local_dofs (std::ostream &out, const std::string &current_string) const
 
void assign_local_eqn_numbers (const bool &store_local_dof_pt)
 Assign local equation numbers in all elements. More...
 
void convert_to_boundary_node (Node *&node_pt, const Vector< FiniteElement * > &finite_element_pt)
 
void convert_to_boundary_node (Node *&node_pt)
 

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...
 
SpineMeshBase_mesh_pt
 Pointer to corresponding mesh of base state problem. More...
 
int YC_index
 
int YS_index
 
Vector< FiniteElement * > Lower_layer_element_pt
 Vector of pointers to element in the lower layer. More...
 
Vector< FiniteElement * > Upper_layer_element_pt
 Vector of pointers to element in 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
 
- Protected Attributes inherited from oomph::Mesh
Vector< Vector< Node * > > Boundary_node_pt
 
bool Lookup_for_elements_next_boundary_is_setup
 
Vector< Vector< FiniteElement * > > Boundary_element_pt
 
Vector< Vector< int > > Face_index_at_boundary
 
Vector< Node * > Node_pt
 Vector of pointers to nodes. More...
 
Vector< GeneralisedElement * > Element_pt
 Vector of pointers to generalised elements. More...
 
std::vector< boolBoundary_coordinate_exists
 
- Protected Attributes inherited from oomph::PerturbedSpineMesh
Vector< PerturbedSpine * > PerturbedSpine_pt
 A PerturbedSpine mesh contains a Vector of pointers to perturbed spines. More...
 

Static Private Attributes

static int Perturbation_to_nodal_position_indices_not_set_up
 

Additional Inherited Members

- Public Types inherited from oomph::Mesh
typedef void(FiniteElement::* SteadyExactSolutionFctPt) (const Vector< double > &x, Vector< double > &soln)
 
typedef void(FiniteElement::* UnsteadyExactSolutionFctPt) (const double &time, const Vector< double > &x, Vector< double > &soln)
 
- Static Public Attributes inherited from oomph::Mesh
static Steady< 0 > Default_TimeStepper
 The Steady Timestepper. More...
 
static bool Suppress_warning_about_empty_mesh_level_time_stepper_function
 Static boolean flag to control warning about mesh level timesteppers. More...
 

Detailed Description

template<class ELEMENT>
class oomph::TwoLayerPerturbedSpineMesh< 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 PerturbedSpineElement<QCrouzeixRaviartElement<2>) and an intermediate interface layer of corresponding PerturbedSpine interface elements, of type INTERFACE_ELEMENT, e.g. PerturbedSpineLineFluidInterfaceElement<ELEMENT> for 2D planar problems.

Constructor & Destructor Documentation

◆ TwoLayerPerturbedSpineMesh() [1/3]

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

Constructor: Pass in:

  • number of elements in x-direction
  • number of elements in y-direction in bottom layer
  • number of elements in y-direction in top layer
  • axial length
  • height of bottom layer
  • height of top layer
  • pointer to the base mesh (where we compute the base flow)
  • pointer to timestepper (defaults to Steady timestepper)

Constructor for 2D perturbed spine mesh: Pass in:

  • number of elements in x-direction
  • number of elements in y-direction in bottom layer
  • number of elements in y-direction in top layer
  • axial length
  • height of bottom layer
  • height of top layer
  • pointer to the base mesh (where we compute the base flow)
  • pointer to timestepper (defaults to Steady timestepper)

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

We set YC_index and YS_index, which store the indices at which the cosine and sine parts (respectively) of the perturbation to the nodal y-position are stored in the bulk element, to a "magic" number which indicates that the indices have not been set up yet.

72  :
73  RectangularQuadMesh<ELEMENT>(nx,ny1+ny2,0.0,lx,0.0,h1+h2,false,false,
74  time_stepper_pt),
75  Base_mesh_pt(base_mesh_pt),
78  {
79  // We've called the "generic" constructor for the RectangularQuadMesh
80  // which doesn't do much...
81  // Now set up the parameters that characterise the mesh geometry
82  // then call the build function
83 
84  // Number of elements in bottom and top layers
85  Ny1 = ny1;
86  Ny2 = ny2;
87 
88  // Set height of upper and lower layers
89  H1 = h1;
90  H2 = h2;
91 
92  // Now build the mesh:
93  build_two_layer_mesh(time_stepper_pt);
94 
95  }
Definition: rectangular_quadmesh.template.h:59
const unsigned & nx() const
Return number of elements in x direction.
Definition: rectangular_quadmesh.template.h:224
virtual void build_two_layer_mesh(TimeStepper *time_stepper_pt)
Definition: two_layer_perturbed_spine_mesh.template.cc:296
unsigned Ny1
Number of elements in lower layer.
Definition: two_layer_perturbed_spine_mesh.template.h:187
const unsigned & ny2() const
Access function for number of elements in upper layer.
Definition: two_layer_perturbed_spine_mesh.template.h:171
int YC_index
Definition: two_layer_perturbed_spine_mesh.template.h:209
const unsigned & ny1() const
Access function for number of elements in lower layer.
Definition: two_layer_perturbed_spine_mesh.template.h:168
int YS_index
Definition: two_layer_perturbed_spine_mesh.template.h:215
SpineMesh * Base_mesh_pt
Pointer to corresponding mesh of base state problem.
Definition: two_layer_perturbed_spine_mesh.template.h:203
unsigned Ny2
Number of elements in upper layer.
Definition: two_layer_perturbed_spine_mesh.template.h:190
static int Perturbation_to_nodal_position_indices_not_set_up
Definition: two_layer_perturbed_spine_mesh.template.h:348
double H1
Height of the lower layer.
Definition: two_layer_perturbed_spine_mesh.template.h:193
double H2
Height of the upper layer.
Definition: two_layer_perturbed_spine_mesh.template.h:196
const double lx
Definition: ConstraintElementsUnitTest.cpp:33

◆ TwoLayerPerturbedSpineMesh() [2/3]

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

Constructor: Pass in:

  • number of elements in x-direction
  • number of elements in y-direction in bottom layer
  • number of elements in y-direction in top layer
  • axial length
  • height of bottom layer
  • height of top layer
  • boolean flag to make the mesh periodic in the x-direction
  • pointer to the base mesh (where we compute the base flow)
  • pointer to timestepper (defaults to Steady timestepper)

Constructor for 2D perturbed spine mesh: Pass in:

  • number of elements in x-direction
  • number of elements in y-direction in bottom layer
  • number of elements in y-direction in top layer
  • axial length
  • height of bottom layer
  • height of top layer
  • boolean flag to make the mesh periodic in the x-direction
  • pointer to the base mesh (where we compute the base flow)
  • pointer to timestepper (defaults to Steady timestepper)

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

We set YC_index and YS_index, which store the indices at which the cosine and sine parts (respectively) of the perturbation to the nodal y-position are stored in the bulk element, to a "magic" number which indicates that the indices have not been set up yet.

133  :
134  RectangularQuadMesh<ELEMENT >(nx,ny1+ny2,0.0,lx,0.0,h1+h2,
135  periodic_in_x,false,time_stepper_pt),
136  Base_mesh_pt(base_mesh_pt),
139  {
140  // We've called the "generic" constructor for the RectangularQuadMesh
141  // which doesn't do much...
142  // Now set up the parameters that characterise the mesh geometry
143  // then call the constructor
144 
145  // Number of elements in bottom and top layers
146  Ny1 = ny1;
147  Ny2 = ny2;
148 
149  // Set height of upper and lower layers
150  H1 = h1;
151  H2 = h2;
152 
153  // Now build the mesh:
154  build_two_layer_mesh(time_stepper_pt);
155 
156  }

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

◆ TwoLayerPerturbedSpineMesh() [3/3]

template<class ELEMENT >
oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::TwoLayerPerturbedSpineMesh ( 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,
SpineMesh *&  base_mesh_pt,
TimeStepper time_stepper_pt = &Mesh::Default_TimeStepper 
)

Constructor: Pass in:

  • number of elements in x-direction
  • number of elements in y-direction in bottom layer
  • number of elements in y-direction in top layer
  • axial length
  • height of bottom layer
  • height of top layer
  • boolean flag to make the mesh periodic in the x-direction
  • boolean flag to specify whether or not to call the "build_two_layer_mesh" function
  • pointer to the base mesh (where we compute the base flow)
  • pointer to timestepper (defaults to Steady timestepper)

Constructor for 2D perturbed spine mesh: Pass in:

  • number of elements in x-direction
  • number of elements in y-direction in bottom layer
  • number of elements in y-direction in top layer
  • axial length
  • height of bottom layer
  • height of top layer
  • boolean flag to make the mesh periodic in the x-direction
  • boolean flag to specify whether or not to call the "build_two_layer_mesh" function
  • pointer to the base mesh (where we compute the base flow)
  • pointer to timestepper (defaults to Steady timestepper)

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

We set YC_index and YS_index, which store the indices at which the cosine and sine parts (respectively) of the perturbation to the nodal y-position are stored in the bulk element, to a "magic" number which indicates that the indices have not been set up yet.

197  :
198  RectangularQuadMesh<ELEMENT >(nx,ny1+ny2,0.0,lx,0.0,h1+h2,
199  periodic_in_x,false,time_stepper_pt),
200  Base_mesh_pt(base_mesh_pt),
203  {
204  // We've called the "generic" constructor for the RectangularQuadMesh
205  // which doesn't do much...
206  // Now set up the parameters that characterise the mesh geometry
207  // then call the constructor
208 
209  // Number of elements in bottom and top layers
210  Ny1 = ny1;
211  Ny2 = ny2;
212 
213  // Set height of upper and lower layers
214  H1 = h1;
215  H2 = h2;
216 
217  // Only build the mesh here if build_mesh=true
218  // This is useful when calling this constructor from a derived class
219  // (such as PatrickAxisym2x6TwoLayerPerturbedSpineMesh) where the mesh
220  // building needs to be called from *its* constructor and this constructor
221  // is only used to pass arguments to the RectangularQuadMesh constructor.
222  if(build_mesh) { build_two_layer_mesh(time_stepper_pt); }
223 
224 }
void build_mesh(TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Generic mesh construction function: contains all the hard work.
Definition: rectangular_quadmesh.template.cc:43

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

Member Function Documentation

◆ build_two_layer_mesh()

template<class ELEMENT >
void oomph::TwoLayerPerturbedSpineMesh< 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

297  {
298 
299  // Build the underlying quad mesh:
301 
302  // Set up the pointers to elements in the upper and lower fluid
303  // ------------------------------------------------------------
304 
305  // Calculate numbers of nodes in upper and lower regions
306  const unsigned long n_lower = this->Nx*Ny1;
307  const unsigned long n_upper = this->Nx*Ny2;
308 
309  // Loop over lower elements and push back
310  Lower_layer_element_pt.reserve(n_lower);
311  for(unsigned e=0;e<n_lower;e++)
312  {
313  Lower_layer_element_pt.push_back(this->finite_element_pt(e));
314  }
315  // Loop over upper elements and push back
316  Upper_layer_element_pt.reserve(n_upper);
317  for(unsigned e=n_lower;e<(n_lower+n_upper);e++)
318  {
319  Upper_layer_element_pt.push_back(this->finite_element_pt(e));
320  }
321 
322  // Allocate memory for the spines and fractions along spines
323  // ---------------------------------------------------------
324 
325  // Read out number of linear points in the element
326  const unsigned n_p =
327  dynamic_cast<ELEMENT*>(finite_element_pt(0))->nnode_1d();
328 
329  // Allocate store for the spines:
330  if(this->Xperiodic) { PerturbedSpine_pt.reserve((n_p-1)*this->Nx); }
331  else { PerturbedSpine_pt.reserve((n_p-1)*this->Nx+1); }
332 
333  // FIRST SPINE
334  // -----------
335 
336  // Element 0
337  // Node 0
338 
339  // Create a perturbed spine and pass it a pointer to the corresponding
340  // base spine
341  PerturbedSpine* new_perturbed_spine_pt =
342  new PerturbedSpine(time_stepper_pt,Base_mesh_pt->spine_pt(0));
343  PerturbedSpine_pt.push_back(new_perturbed_spine_pt);
344 
345  // Get pointer to node
346  PerturbedSpineNode* nod_pt = element_node_pt(0,0);
347 
348  // Set the pointer to the spine
349  nod_pt->perturbed_spine_pt() = new_perturbed_spine_pt;
350  // Set the fraction
351  nod_pt->fraction() = 0.0;
352  // Pointer to the mesh that implements the update fct
353  nod_pt->spine_mesh_pt() = this;
354  // ID of update function within this mesh: 0 = lower; 1 = upper
355  nod_pt->node_update_fct_id() = 0;
356 
357  // Pass the spine a pointer to this node
358  new_perturbed_spine_pt->node_at_bottom_of_spine_pt() = nod_pt;
359 
360  // Loop vertically along the spine
361  // Loop over the elements in fluid 1
362  for(unsigned long i=0;i<Ny1;i++)
363  {
364  // Loop over the vertical nodes, apart from the first
365  for(unsigned l1=1;l1<n_p;l1++)
366  {
367  // Get pointer to node
368  PerturbedSpineNode* nod_pt=element_node_pt(i*this->Nx,l1*n_p);
369  // Set the pointer to the spine
370  nod_pt->perturbed_spine_pt() = new_perturbed_spine_pt;
371  // Set the fraction
372  nod_pt->fraction() = (nod_pt->x(1)-this->Ymin)/(H1);
373  // Pointer to the mesh that implements the update fct
374  nod_pt->spine_mesh_pt() = this;
375  // ID of update function within this mesh: 0 = lower; 1 = upper
376  nod_pt->node_update_fct_id() = 0;
377  }
378  }
379 
380  // Loop over the elements in fluid 2
381  for(unsigned long i=0;i<Ny2;i++)
382  {
383  // Loop over vertical nodes, apart from the first
384  for(unsigned l1=1;l1<n_p;l1++)
385  {
386  // Get pointer to node
387  PerturbedSpineNode* nod_pt=element_node_pt((Ny1+i)*this->Nx,l1*n_p);
388 
389  //Set the pointer to the spine
390  nod_pt->perturbed_spine_pt() = new_perturbed_spine_pt;
391  //Set the fraction
392  nod_pt->fraction() =(nod_pt->x(1)-(this->Ymin+H1))/(H2);
393  // Pointer to the mesh that implements the update fct
394  nod_pt->spine_mesh_pt() = this;
395  // ID of update function within this mesh: 0 = lower; 1 = upper
396  nod_pt->node_update_fct_id() = 1;
397 
398  // Pass the spine a pointer to the "top" node
399  if(i==(Ny2-1) && l1==(n_p-1))
400  {
401  new_perturbed_spine_pt->node_at_top_of_spine_pt() = nod_pt;
402  }
403  }
404  }
405 
406 
407  //LOOP OVER OTHER SPINES
408  //----------------------
409 
410  //Now loop over the elements horizontally
411  for(unsigned long j=0;j<this->Nx;j++)
412  {
413  //Loop over the nodes in the elements horizontally, ignoring
414  //the first column
415 
416  // Last spine needs special treatment in x-periodic meshes:
417  unsigned n_pmax=n_p;
418  if ((this->Xperiodic)&&(j==this->Nx-1)) n_pmax=n_p-1;
419 
420  for(unsigned l2=1;l2<n_pmax;l2++)
421  {
422  // Create a perturbed spine and pass it a pointer to the
423  // corresponding base spine
424  new_perturbed_spine_pt =
425  new PerturbedSpine(time_stepper_pt,
426  Base_mesh_pt->spine_pt(((n_pmax-1)*j)+l2));
427  PerturbedSpine_pt.push_back(new_perturbed_spine_pt);
428 
429  // Get pointer to node
431 
432  //Set the pointer to the spine
433  nod_pt->perturbed_spine_pt() = new_perturbed_spine_pt;
434  //Set the fraction
435  nod_pt->fraction() = 0.0;
436  // Pointer to the mesh that implements the update fct
437  nod_pt->spine_mesh_pt() = this;
438  // ID of update function within this mesh: 0 = lower; 1 = upper
439  nod_pt->node_update_fct_id() = 0;
440 
441  // Pass the spine a pointer to this node
442  new_perturbed_spine_pt->node_at_bottom_of_spine_pt() = nod_pt;
443 
444  //Loop vertically along the spine
445  //Loop over the elements in fluid 1
446  for(unsigned long i=0;i<Ny1;i++)
447  {
448  //Loop over the vertical nodes, apart from the first
449  for(unsigned l1=1;l1<n_p;l1++)
450  {
451  // Get pointer to node
452  PerturbedSpineNode* nod_pt=element_node_pt(i*this->Nx+j,l1*n_p+l2);
453  //Set the pointer to the spine
454  nod_pt->perturbed_spine_pt() = new_perturbed_spine_pt;
455  //Set the fraction
456  nod_pt->fraction() = (nod_pt->x(1)-this->Ymin)/H1;
457  // Pointer to the mesh that implements the update fct
458  nod_pt->spine_mesh_pt() = this;
459  // ID of update function within this mesh: 0 = lower; 1 = upper
460  nod_pt->node_update_fct_id() = 0;
461  }
462  }
463 
464  //Loop over the elements in fluid 2
465  for(unsigned long i=0;i<Ny2;i++)
466  {
467  //Loop over vertical nodes, apart from the first
468  for(unsigned l1=1;l1<n_p;l1++)
469  {
470  // Get pointer to node
471  PerturbedSpineNode* nod_pt=element_node_pt((Ny1+i)*this->Nx+j,l1*n_p+l2);
472 
473  //Set the pointer to the spine
474  nod_pt->perturbed_spine_pt() = new_perturbed_spine_pt;
475  //Set the fraction
476  nod_pt->fraction() = (nod_pt->x(1)-(this->Ymin+H1))/H2;
477  // Pointer to the mesh that implements the update fct
478  nod_pt->spine_mesh_pt() = this;
479  // ID of update function within this mesh: 0 = lower; 1 = upper
480  nod_pt->node_update_fct_id() = 1;
481 
482  // Pass the spine a pointer to the "top" node
483  if(i==(Ny2-1) && l1==(n_p-1))
484  {
485  new_perturbed_spine_pt->node_at_top_of_spine_pt() = nod_pt;
486  }
487  }
488  }
489  }
490  }
491 
492 
493  // Last spine needs special treatment for periodic meshes
494  // because it's the same as the first one...
495  if (this->Xperiodic)
496  {
497  // Last spine is the same as first one...
498  PerturbedSpine* final_perturbed_spine_pt=PerturbedSpine_pt[0];
499 
500  // Get pointer to node
501  PerturbedSpineNode* nod_pt=element_node_pt((this->Nx-1),(n_p-1));
502 
503  //Set the pointer to the spine
504  nod_pt->perturbed_spine_pt() = final_perturbed_spine_pt;
505  //Set the fraction to be the same as for the nodes on the first row
506  nod_pt->fraction() = element_node_pt(0,0)->fraction();
507  // Pointer to the mesh that implements the update fct
508  nod_pt->spine_mesh_pt() = element_node_pt(0,0)->spine_mesh_pt();
509  // ID of update function within this mesh: 0 = lower; 1 = upper
511 
512  //Now loop vertically along the spine
513  for(unsigned i=0;i<(Ny1+Ny2);i++)
514  {
515  //Loop over the vertical nodes, apart from the first
516  for(unsigned l1=1;l1<n_p;l1++)
517  {
518  // Get pointer to node
519  PerturbedSpineNode* nod_pt =
520  element_node_pt(i*this->Nx+(this->Nx-1),l1*n_p+(n_p-1));
521 
522  //Set the pointer to the spine
523  nod_pt->perturbed_spine_pt() = final_perturbed_spine_pt;
524  //Set the fraction to be the same as in first row
525  nod_pt->fraction() = element_node_pt(i*this->Nx,l1*n_p)->fraction();
526  // ID of update function within this mesh: 0 = lower; 1 = upper
527  nod_pt->node_update_fct_id() =
528  element_node_pt(i*this->Nx,l1*n_p)->node_update_fct_id();
529  // Pointer to the mesh that implements the update fct
530  nod_pt->spine_mesh_pt() = element_node_pt(i*this->Nx,l1*n_p)
531  ->spine_mesh_pt();
532  }
533  }
534  }
535 
536 }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Array< double, 1, 3 > e(1./3., 0.5, 2.)
FiniteElement * finite_element_pt(const unsigned &e) const
Definition: mesh.h:473
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition: nodes.h:1060
PerturbedSpineNode * element_node_pt(const unsigned long &e, const unsigned &n)
Definition: perturbed_spines.h:613
Vector< PerturbedSpine * > PerturbedSpine_pt
A PerturbedSpine mesh contains a Vector of pointers to perturbed spines.
Definition: perturbed_spines.h:557
Class for nodes that ‘live’ on perturbed spines.
Definition: perturbed_spines.h:276
unsigned & node_update_fct_id()
Access function to ID of node update function (within specific mesh)
Definition: perturbed_spines.h:320
double & fraction()
Set reference to fraction along spine.
Definition: perturbed_spines.h:317
PerturbedSpine *& perturbed_spine_pt()
Access function to perturbed spine.
Definition: perturbed_spines.h:314
PerturbedSpineMesh *& spine_mesh_pt()
Definition: perturbed_spines.h:324
Definition: perturbed_spines.h:49
PerturbedSpineNode *& node_at_bottom_of_spine_pt()
Access function to SpineNode at bottom of spine.
Definition: perturbed_spines.h:121
PerturbedSpineNode *& node_at_top_of_spine_pt()
Access function to SpineNode at top of spine.
Definition: perturbed_spines.h:134
unsigned Nx
Nx: number of elements in x-direction.
Definition: rectangular_quadmesh.template.h:63
bool Xperiodic
Definition: rectangular_quadmesh.template.h:81
double Ymin
Minimum value of y coordinate.
Definition: rectangular_quadmesh.template.h:75
Spine *& spine_pt(const unsigned long &i)
Return the i-th spine in the mesh.
Definition: spines.h:623
Vector< FiniteElement * > Lower_layer_element_pt
Vector of pointers to element in the lower layer.
Definition: two_layer_perturbed_spine_mesh.template.h:218
Vector< FiniteElement * > Upper_layer_element_pt
Vector of pointers to element in the upper layer.
Definition: two_layer_perturbed_spine_mesh.template.h:221
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References oomph::RectangularQuadMesh< ELEMENT >::build_mesh(), e(), and GlobalParameters::Nx.

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

◆ element_reorder()

template<class ELEMENT >
void oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::element_reorder
virtual

Reorder the elements so we loop over them vertically first (advantageous in "wide" domains if a frontal solver is used).

Reorder the elements, so we loop over them vertically first (advantageous in "wide" domains if a frontal solver is used).

Reimplemented from oomph::RectangularQuadMesh< ELEMENT >.

545 {
546 
547  unsigned Nx = this->Nx;
548  //Find out how many elements there are
549  unsigned long Nelement = nelement();
550  //Find out how many fluid elements there are
551  unsigned long Nfluid = Nx*(Ny1+Ny2);
552  //Create a dummy array of elements
554 
555  //Loop over the elements in horizontal order
556  for(unsigned long j=0;j<Nx;j++)
557  {
558  //Loop over the elements in lower layer vertically
559  for(unsigned long i=0;i<Ny1;i++)
560  {
561  //Push back onto the new stack
562  dummy.push_back(finite_element_pt(Nx*i + j));
563  }
564 
565  //Push back the line element onto the stack
566  dummy.push_back(finite_element_pt(Nfluid+j));
567 
568  //Loop over the elements in upper layer vertically
569  for(unsigned long i=Ny1;i<(Ny1+Ny2);i++)
570  {
571  //Push back onto the new stack
572  dummy.push_back(finite_element_pt(Nx*i + j));
573  }
574  }
575 
576  //Now copy the reordered elements into the element_pt
577  for(unsigned long e=0;e<Nelement;e++)
578  {
579  Element_pt[e] = dummy[e];
580  }
581 
582 }
Vector< GeneralisedElement * > Element_pt
Vector of pointers to generalised elements.
Definition: mesh.h:186
unsigned long nelement() const
Return number of elements in the mesh.
Definition: mesh.h:590

References e(), i, j, and GlobalParameters::Nx.

◆ lower_layer_element_pt()

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

◆ nlower()

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

◆ nupper()

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

◆ ny1()

template<class ELEMENT >
const unsigned& oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::ny1 ( ) const
inline

Access function for number of elements in lower layer.

168 { return Ny1; }

References oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::Ny1.

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

◆ ny2()

template<class ELEMENT >
const unsigned& oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::ny2 ( ) const
inline

Access function for number of elements in upper layer.

171 { return Ny2; }

References oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::Ny2.

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

◆ perturbed_spine_node_update()

template<class ELEMENT >
void oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::perturbed_spine_node_update ( PerturbedSpineNode spine_node_pt)
inlinevirtual

General node update function implements pure virtual function defined in PerturbedSpineMesh base class and performs specific update actions, depending on the node update fct id stored for each node. This function sets the nodal positions to be the same as in the base state problem through the perturbed spine's pointer to its corresponding base spine. Additionally this function updates the perturbations to the nodal positions, which are stored as (pinned) values at the nodes.

Implements oomph::PerturbedSpineMesh.

145  {
146  unsigned id=spine_node_pt->node_update_fct_id();
147  switch(id)
148  {
149  case 0:
150  spine_node_update_lower(spine_node_pt);
151  break;
152 
153  case 1:
154  spine_node_update_upper(spine_node_pt);
155  break;
156 
157  default:
158  std::ostringstream error_message;
159  error_message << "Unknown id passed to perturbed_spine_node_update " << id
160  << std::endl;
161  throw OomphLibError(error_message.str(),
162  "TwoLayerPerturbedSpineMesh::perturbed_spine_node_update()",
164  }
165  }
Definition: oomph_definitions.h:222
void spine_node_update_upper(PerturbedSpineNode *spine_node_pt)
Update function for the upper part of the domain.
Definition: two_layer_perturbed_spine_mesh.template.h:285
void spine_node_update_lower(PerturbedSpineNode *spine_node_pt)
Update function for the lower part of the domain.
Definition: two_layer_perturbed_spine_mesh.template.h:234
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61

References oomph::PerturbedSpineNode::node_update_fct_id(), OOMPH_EXCEPTION_LOCATION, oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::spine_node_update_lower(), and oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::spine_node_update_upper().

◆ set_perturbation_to_nodal_positions_indices()

template<class ELEMENT >
const void oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::set_perturbation_to_nodal_positions_indices ( const unsigned cosine_index,
const unsigned sine_index 
)
inline

Set up the internal data YC_index and YS_index, which store the indices at which the perturbations to the nodal y-position are stored in the bulk element

179  {
180  YC_index = cosine_index;
181  YS_index = sine_index;
182  }

References oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::YC_index, and oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::YS_index.

Referenced by PerturbedStateProblem< BASE_ELEMENT, PERTURBED_ELEMENT >::PerturbedStateProblem().

◆ spine_node_update_lower()

template<class ELEMENT >
void oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::spine_node_update_lower ( PerturbedSpineNode spine_node_pt)
inlineprotected

Update function for the lower part of the domain.

235  {
236  // Get fraction along the spine
237  const double w = spine_node_pt->fraction();
238 
239  // Get base spine height
240  const double h =
241  spine_node_pt->perturbed_spine_pt()->base_spine_pt()->height();
242 
243  // Set the nodal y-position
244  spine_node_pt->x(1) = this->Ymin + w*h;
245 
246  // Get perturbed spine height (cosine and sine parts)
247  const double HC = spine_node_pt->perturbed_spine_pt()->height(0);
248  const double HS = spine_node_pt->perturbed_spine_pt()->height(1);
249 
250  // Set the perturbation (cosine and sine parts) to the nodal
251  // y-position. Note that YC_index and YS_index are the indices at
252  // which the perturbations to the nodal y-position are stored in
253  // the bulk elements
254  if(YC_index>=0)
255  {
256  spine_node_pt->set_value(YC_index,w*HC);
257  }
258  else
259  {
260  std::ostringstream error_message;
261  error_message << "The cosine part of the perturbation to the nodal\n"
262  << "y-position has not been set up." << std::endl;
263  throw OomphLibError(
264  error_message.str(),
265  "TwoLayerPerturbedSpineMesh::spine_node_update_lower(...)",
267  }
268  if(YS_index>=0)
269  {
270  spine_node_pt->set_value(YS_index,w*HS);
271  }
272  else
273  {
274  std::ostringstream error_message;
275  error_message << "The sine part of the perturbation to the nodal\n"
276  << "y-position has not been set up." << std::endl;
277  throw OomphLibError(
278  error_message.str(),
279  "TwoLayerPerturbedSpineMesh::spine_node_update_lower(...)",
281  }
282  }
RowVector3d w
Definition: Matrix_resize_int.cpp:3
void set_value(const unsigned &i, const double &value_)
Definition: nodes.h:271
Spine *& base_spine_pt()
Access function to pointer to base spine.
Definition: perturbed_spines.h:98
double & height(const unsigned &i)
Access function to i-th component of spine "height".
Definition: perturbed_spines.h:101
double & height()
Access function to spine height.
Definition: spines.h:150

References oomph::PerturbedSpine::base_spine_pt(), oomph::PerturbedSpineNode::fraction(), oomph::Spine::height(), oomph::PerturbedSpine::height(), OOMPH_EXCEPTION_LOCATION, oomph::PerturbedSpineNode::perturbed_spine_pt(), oomph::Data::set_value(), w, oomph::Node::x(), oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::YC_index, oomph::RectangularQuadMesh< ELEMENT >::Ymin, and oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::YS_index.

Referenced by oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::perturbed_spine_node_update().

◆ spine_node_update_upper()

template<class ELEMENT >
void oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::spine_node_update_upper ( PerturbedSpineNode spine_node_pt)
inlineprotected

Update function for the upper part of the domain.

286  {
287  // Get fraction along the spine
288  const double w = spine_node_pt->fraction();
289 
290  // Get base spine height
291  const double h =
292  spine_node_pt->perturbed_spine_pt()->base_spine_pt()->height();
293 
294  // Set the nodal y-position
295  spine_node_pt->x(1) = (this->Ymin + h) + w*(this->Ymax - (this->Ymin+h));
296 
297  // Get perturbed spine height (cosine and sine parts)
298  const double HC = spine_node_pt->perturbed_spine_pt()->height(0);
299  const double HS = spine_node_pt->perturbed_spine_pt()->height(1);
300 
301  // Set the perturbation (cosine and sine parts) to the nodal
302  // y-position. Note that YC_index and YS_index are the indices at
303  // which the perturbations to the nodal y-position are stored in
304  // the bulk elements
305  if(YC_index>=0)
306  {
307  // To see analysis of why this is okay, see the pdf file
308  // "perturbed_node_update_procedure.pdf" (emailed to myself)
309  spine_node_pt->set_value(YC_index,HC*(1.0-w));
310  }
311  else
312  {
313  std::ostringstream error_message;
314  error_message << "The cosine part of the perturbation to the nodal\n"
315  << "y-position has not been set up." << std::endl;
316  throw OomphLibError(
317  error_message.str(),
318  "TwoLayerPerturbedSpineMesh::spine_node_update_upper(...)",
320  }
321  if(YS_index>=0)
322  {
323  // To see analysis of why this is okay, see the pdf file
324  // "perturbed_node_update_procedure.pdf" (emailed to myself)
325  spine_node_pt->set_value(YS_index,HS*(1.0-w));
326  }
327  else
328  {
329  std::ostringstream error_message;
330  error_message << "The sine part of the perturbation to the nodal\n"
331  << "y-position has not been set up." << std::endl;
332  throw OomphLibError(
333  error_message.str(),
334  "TwoLayerPerturbedSpineMesh::spine_node_update_upper(...)",
336  }
337  }
double Ymax
Maximum value of y coordinate.
Definition: rectangular_quadmesh.template.h:77

References oomph::PerturbedSpine::base_spine_pt(), oomph::PerturbedSpineNode::fraction(), oomph::Spine::height(), oomph::PerturbedSpine::height(), OOMPH_EXCEPTION_LOCATION, oomph::PerturbedSpineNode::perturbed_spine_pt(), oomph::Data::set_value(), w, oomph::Node::x(), oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::YC_index, oomph::RectangularQuadMesh< ELEMENT >::Ymax, oomph::RectangularQuadMesh< ELEMENT >::Ymin, and oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::YS_index.

Referenced by oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::perturbed_spine_node_update().

◆ upper_layer_element_pt()

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

◆ x_spacing_function()

template<class ELEMENT >
double oomph::TwoLayerPerturbedSpineMesh< 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 >.

235  {
236  // Calculate the values of equal increments in nodal values
237  double xstep = (this->Xmax-this->Xmin)/((this->Np-1)*this->Nx);
238 
239  // Return the appropriate value
240  return (this->Xmin + xstep*((this->Np-1)*xelement + xnode));
241  }
unsigned Np
Np: number of (linear) points in the element.
Definition: rectangular_quadmesh.template.h:67
double Xmax
Maximum value of x coordinate.
Definition: rectangular_quadmesh.template.h:72
double Xmin
Minimum value of x coordinate.
Definition: rectangular_quadmesh.template.h:70

References GlobalParameters::Nx.

◆ y_spacing_function()

template<class ELEMENT >
double oomph::TwoLayerPerturbedSpineMesh< 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 >.

252  {
253  // Set up some spacing parameters
254 
255  // The lower region a starts at Ymin
257 
258  // The upper region a ends at Ymax
260 
261  // Number of nodes per element
262  unsigned n_p = RectangularQuadMesh<ELEMENT >::Np;
263 
264  // The lower region starts at Ymin
265  double y1init = Ymin;
266 
267  // The upper region starts at H1 - Ymin
268  double y2init = H1 - Ymin;
269 
270  // Calculate the space between each node in each region, assuming
271  // uniform spacing
272 
273  // Lower region has a length (H1-Ymin)
274  double y1step = (H1-Ymin)/((n_p-1)*Ny1);
275 
276  // Upper region has a length (Ymax-H1)
277  double y2step = (Ymax-H1)/((n_p-1)*Ny2);
278 
279  // Now return the actual node position (different in the two regions)
280  if(yelement < Ny1)
281  {
282  return (y1init + y1step*((n_p-1)*yelement + ynode));
283  }
284  else
285  {
286  return (y2init + y2step*((n_p-1)*(yelement-Ny1) + ynode));
287  }
288  }

Member Data Documentation

◆ Base_mesh_pt

template<class ELEMENT >
SpineMesh* oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::Base_mesh_pt
protected

Pointer to corresponding mesh of base state problem.

◆ H1

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

◆ H2

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

◆ Lower_layer_element_pt

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

◆ Ny1

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

◆ Ny2

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

◆ Perturbation_to_nodal_position_indices_not_set_up

template<class ELEMENT >
int oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::Perturbation_to_nodal_position_indices_not_set_up
staticprivate

Static "magic" number that indicates that the indices at which the perturbations to the nodal y-positions are stored have not been set up

◆ Upper_layer_element_pt

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

◆ YC_index

template<class ELEMENT >
int oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::YC_index
protected

Index at which the cosine part of the perturbation to the nodal y-position is stored in the bulk element. The mesh needs to know this so that the value of the nodal position perturbation can be updated during the "perturbed spine node update" procedure

Referenced by oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::set_perturbation_to_nodal_positions_indices(), oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::spine_node_update_lower(), and oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::spine_node_update_upper().

◆ YS_index

template<class ELEMENT >
int oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::YS_index
protected

Index at which the sine part of the perturbation to the nodal y-position is stored in the bulk element. The mesh needs to know this so that the value of the nodal position perturbation can be updated during the "perturbed spine node update" procedure

Referenced by oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::set_perturbation_to_nodal_positions_indices(), oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::spine_node_update_lower(), and oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::spine_node_update_upper().


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