two_layer_spine_mesh.template.h
Go to the documentation of this file.
1 // LIC// ====================================================================
2 // LIC// This file forms part of oomph-lib, the object-oriented,
3 // LIC// multi-physics finite-element library, available
4 // LIC// at http://www.oomph-lib.org.
5 // LIC//
6 // LIC// Copyright (C) 2006-2023 Matthias Heil and Andrew Hazel
7 // LIC//
8 // LIC// This library is free software; you can redistribute it and/or
9 // LIC// modify it under the terms of the GNU Lesser General Public
10 // LIC// License as published by the Free Software Foundation; either
11 // LIC// version 2.1 of the License, or (at your option) any later version.
12 // LIC//
13 // LIC// This library is distributed in the hope that it will be useful,
14 // LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // LIC// Lesser General Public License for more details.
17 // LIC//
18 // LIC// You should have received a copy of the GNU Lesser General Public
19 // LIC// License along with this library; if not, write to the Free Software
20 // LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21 // LIC// 02110-1301 USA.
22 // LIC//
23 // LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24 // LIC//
25 // LIC//====================================================================
26 #ifndef OOMPH_TWO_LAYER_SPINE_MESH_HEADER
27 #define OOMPH_TWO_LAYER_SPINE_MESH_HEADER
28 
29 
30 // The mesh
31 #include "../generic/spines.h"
33 
34 namespace oomph
35 {
36  //======================================================================
37  /// Two-layer spine mesh class derived from standard 2D mesh.
38  /// The mesh contains two layers of spinified fluid elements (of type ELEMENT;
39  /// e.g SpineElement<QCrouzeixRaviartElement<2>).
40  ///
41  /// This mesh paritions the elements into those above and below a
42  /// notional interface and relabels boundaries so that the mesh is as follows
43  /// 3
44  /// ---------------------------------------
45  /// | |
46  /// 4 | | 2
47  /// | 6 |
48  /// ---------------------------------------
49  /// | |
50  /// 5 | | 1
51  /// | |
52  /// ---------------------------------------
53  /// 0
54  /// Update information for the nodes in response to changes in spine length
55  /// is given, but additional equations must be specified in order to
56  /// completely specify the problem.
57  //======================================================================
58  template<class ELEMENT>
59  class TwoLayerSpineMesh : public RectangularQuadMesh<ELEMENT>,
60  public SpineMesh
61  {
62  public:
63  /// Constructor: Pass number of elements in x-direction, number of
64  /// elements in y-direction in bottom and top layer, respectively,
65  /// axial length and height of top and bottom layers and pointer
66  /// to timestepper (defaults to Steady timestepper)
68  const unsigned& nx,
69  const unsigned& ny1,
70  const unsigned& ny2,
71  const double& lx,
72  const double& h1,
73  const double& h2,
74  TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper);
75 
76 
77  /// Constructor: Pass number of elements in x-direction, number of
78  /// elements in y-direction in bottom and top layer, respectively,
79  /// axial length and height of top and bottom layers, a boolean
80  /// flag to make the mesh periodic in the x-direction, and pointer
81  /// to timestepper (defaults to Steady timestepper)
83  const unsigned& nx,
84  const unsigned& ny1,
85  const unsigned& ny2,
86  const double& lx,
87  const double& h1,
88  const double& h2,
89  const bool& periodic_in_x,
90  TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper);
91 
92 
93  /// Constructor: Pass number of elements in x-direction, number of
94  /// elements in y-direction in bottom and top layer, respectively,
95  /// axial length and height of top and bottom layers, a boolean
96  /// flag to make the mesh periodic in the x-direction, a boolean flag to
97  /// specify whether or not to call the "build_two_layer_mesh" function,
98  /// and pointer to timestepper (defaults to Steady timestepper)
100  const unsigned& nx,
101  const unsigned& ny1,
102  const unsigned& ny2,
103  const double& lx,
104  const double& h1,
105  const double& h2,
106  const bool& periodic_in_x,
107  const bool& build_mesh,
108  TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper);
109 
110  /// Access functions for pointers to elements in upper layer
111  FiniteElement*& upper_layer_element_pt(const unsigned long& i)
112  {
113  return Upper_layer_element_pt[i];
114  }
115 
116  /// Access functions for pointers to elements in bottom layer
117  FiniteElement*& lower_layer_element_pt(const unsigned long& i)
118  {
119  return Lower_layer_element_pt[i];
120  }
121 
122  /// Number of elements in upper layer
123  unsigned long nupper() const
124  {
125  return Upper_layer_element_pt.size();
126  }
127 
128  /// Number of elements in top layer
129  unsigned long nlower() const
130  {
131  return Lower_layer_element_pt.size();
132  }
133 
134  /// Access functions for pointers to elements in upper layer
136  {
138  }
139 
140  /// Access functions for pointers to elements in bottom layer
142  {
144  }
145 
146  /// Number of elements in upper layer
147  unsigned long ninterface_upper() const
148  {
150  }
151 
152  /// Number of elements in top layer
153  unsigned long ninterface_lower() const
154  {
156  }
157 
158  /// Index of the face of the elements next to the interface
159  /// in the upper region (always -2)
161  {
162  return -2;
163  }
164 
165  /// Index of the face of the elements next to the interface in
166  /// the lower region (always 2)
168  {
169  return 2;
170  }
171 
172  /// General node update function implements pure virtual function
173  /// defined in SpineMesh base class and performs specific update
174  /// actions, depending on the node update fct id stored for each node.
175  void spine_node_update(SpineNode* spine_node_pt)
176  {
177  unsigned id = spine_node_pt->node_update_fct_id();
178  switch (id)
179  {
180  case 0:
181  spine_node_update_lower(spine_node_pt);
182  break;
183 
184  case 1:
185  spine_node_update_upper(spine_node_pt);
186  break;
187 
188  default:
189  std::ostringstream error_message;
190  error_message << "Unknown id passed to spine_node_update " << id
191  << std::endl;
192  throw OomphLibError(error_message.str(),
193  OOMPH_CURRENT_FUNCTION,
194  OOMPH_EXCEPTION_LOCATION);
195  }
196  }
197 
198 
199  protected:
200  /// Number of elements in lower layer
201  unsigned Ny1;
202 
203  /// Number of elements in upper layer
204  unsigned Ny2;
205 
206  /// Height of the lower layer
207  double H1;
208 
209  /// Height of the upper layer
210  double H2;
211 
212  /// Vector of pointers to element in the upper layer
214 
215  /// Vector of pointers to element in the lower layer
217 
218  /// Vector of pointers to the elements adjacent to the interface
219  /// on the lower layer
221 
222  /// Vector of pointers to the element adjacent to the interface
223  /// on the upper layer
225 
226 
227  /// The spacing function for the x co-ordinates with two
228  /// regions.
229  double x_spacing_function(unsigned xelement,
230  unsigned xnode,
231  unsigned yelement,
232  unsigned ynode);
233 
234  /// The spacing function for the y co-ordinates with three
235  /// regions in each fluid.
236  double y_spacing_function(unsigned xelement,
237  unsigned xnode,
238  unsigned yelement,
239  unsigned ynode);
240 
241  /// Update function for the lower part of the domain
242  void spine_node_update_lower(SpineNode* spine_node_pt)
243  {
244  // Get fraction along the spine
245  double W = spine_node_pt->fraction();
246  // Get spine height
247  double H = spine_node_pt->h();
248  // Set the value of y
249  spine_node_pt->x(1) = this->Ymin + W * H;
250  }
251 
252 
253  /// Update function for the upper part of the domain
254  void spine_node_update_upper(SpineNode* spine_node_pt)
255  {
256  // Get fraction alon the spine
257  double W = spine_node_pt->fraction();
258 
259  // Get spine height
260  double H = spine_node_pt->h();
261 
262  // Set the value of y
263  spine_node_pt->x(1) =
264  (this->Ymin + H) + W * (this->Ymax - (this->Ymin + H));
265  }
266 
267  /// Helper function to actually build the two-layer spine mesh
268  /// (called from various constructors)
269  virtual void build_two_layer_mesh(TimeStepper* time_stepper_pt);
270  };
271 
272 } // namespace oomph
273 
274 #endif
e
Definition: cfortran.h:571
cstr elem_len * i
Definition: cfortran.h:603
A general Finite Element class.
Definition: elements.h:1313
static Steady< 0 > Default_TimeStepper
Default Steady Timestepper, to be used in default arguments to Mesh constructors.
Definition: mesh.h:75
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition: nodes.h:1060
An OomphLibError object which should be thrown when an run-time error is encountered....
RectangularQuadMesh is a two-dimensional mesh of Quad elements with Nx elements in the "x" (horizonal...
double Ymax
Maximum value of y coordinate.
const unsigned & nx() const
Return number of elements in x direction.
double Ymin
Minimum value of y coordinate.
void build_mesh(TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Generic mesh construction function: contains all the hard work.
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
Definition: spines.h:613
Class for nodes that live on spines. The assumption is that each Node lies at a fixed fraction on a s...
Definition: spines.h:328
double & h()
Access function to spine height.
Definition: spines.h:397
unsigned & node_update_fct_id()
Access function to ID of node update function (within specific mesh)
Definition: spines.h:384
double & fraction()
Set reference to fraction along spine.
Definition: spines.h:378
////////////////////////////////////////////////////////////////////// //////////////////////////////...
Definition: timesteppers.h:231
Two-layer spine mesh class derived from standard 2D mesh. The mesh contains two layers of spinified f...
Vector< FiniteElement * > Upper_layer_element_pt
Vector of pointers to element in the lower layer.
double H2
Height of the upper layer.
FiniteElement *& upper_layer_element_pt(const unsigned long &i)
Access functions for pointers to elements in upper layer.
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 ...
Vector< FiniteElement * > Lower_layer_element_pt
Vector of pointers to element in the upper layer.
Vector< FiniteElement * > Interface_upper_boundary_element_pt
Vector of pointers to the element adjacent to the interface on the upper layer.
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.
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)
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)
unsigned Ny2
Number of elements in upper layer.
unsigned long ninterface_lower() const
Number of elements in top layer.
double x_spacing_function(unsigned xelement, unsigned xnode, unsigned yelement, unsigned ynode)
The spacing function for the x co-ordinates with two regions.
unsigned Ny1
Number of elements in lower layer.
FiniteElement *& lower_layer_element_pt(const unsigned long &i)
Access functions for pointers to elements in bottom layer.
void spine_node_update_lower(SpineNode *spine_node_pt)
Update function for the lower part of the domain.
unsigned long nupper() const
Number of elements in upper layer.
FiniteElement *& interface_upper_boundary_element_pt(const unsigned long &i)
Access functions for pointers to elements in upper layer.
Vector< FiniteElement * > Interface_lower_boundary_element_pt
Vector of pointers to the elements adjacent to the interface on the lower layer.
unsigned long nlower() const
Number of elements in top layer.
double H1
Height of the lower layer.
void spine_node_update(SpineNode *spine_node_pt)
General node update function implements pure virtual function defined in SpineMesh base class and per...
virtual void build_two_layer_mesh(TimeStepper *time_stepper_pt)
Helper function to actually build the two-layer spine mesh (called from various constructors)
unsigned long ninterface_upper() const
Number of elements in upper layer.
void spine_node_update_upper(SpineNode *spine_node_pt)
Update function for the upper part of the domain.
FiniteElement *& interface_lower_boundary_element_pt(const unsigned long &i)
Access functions for pointers to elements in bottom layer.
A slight extension to the standard template vector class so that we can include "graceful" array rang...
Definition: Vector.h:58
//////////////////////////////////////////////////////////////////// ////////////////////////////////...