two_layer_spine_mesh.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-2025 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// The mesh
30#include "generic/spines.h"
32
33namespace oomph
34{
35 //======================================================================
36 /// Two-layer spine mesh class derived from standard 2D mesh.
37 /// The mesh contains two layers of spinified fluid elements (of type ELEMENT;
38 /// e.g SpineElement<QCrouzeixRaviartElement<2>).
39 ///
40 /// This mesh paritions the elements into those above and below a
41 /// notional interface and relabels boundaries so that the mesh is as follows
42 /// 3
43 /// ---------------------------------------
44 /// | |
45 /// 4 | | 2
46 /// | 6 |
47 /// ---------------------------------------
48 /// | |
49 /// 5 | | 1
50 /// | |
51 /// ---------------------------------------
52 /// 0
53 /// Update information for the nodes in response to changes in spine length
54 /// is given, but additional equations must be specified in order to
55 /// completely specify the problem.
56 //======================================================================
57 template<class ELEMENT>
58 class TwoLayerSpineMesh : public RectangularQuadMesh<ELEMENT>,
59 public SpineMesh
60 {
61 public:
62 /// Constructor: Pass number of elements in x-direction, number of
63 /// elements in y-direction in bottom and top layer, respectively,
64 /// axial length and height of top and bottom layers and pointer
65 /// to timestepper (defaults to Steady timestepper)
67 const unsigned& nx,
68 const unsigned& ny1,
69 const unsigned& ny2,
70 const double& lx,
71 const double& h1,
72 const double& h2,
73 TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper);
74
75
76 /// Constructor: Pass number of elements in x-direction, number of
77 /// elements in y-direction in bottom and top layer, respectively,
78 /// axial length and height of top and bottom layers, a boolean
79 /// flag to make the mesh periodic in the x-direction, and pointer
80 /// to timestepper (defaults to Steady timestepper)
82 const unsigned& nx,
83 const unsigned& ny1,
84 const unsigned& ny2,
85 const double& lx,
86 const double& h1,
87 const double& h2,
88 const bool& periodic_in_x,
89 TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper);
90
91
92 /// Constructor: Pass number of elements in x-direction, number of
93 /// elements in y-direction in bottom and top layer, respectively,
94 /// axial length and height of top and bottom layers, a boolean
95 /// flag to make the mesh periodic in the x-direction, a boolean flag to
96 /// specify whether or not to call the "build_two_layer_mesh" function,
97 /// and pointer to timestepper (defaults to Steady timestepper)
99 const unsigned& nx,
100 const unsigned& ny1,
101 const unsigned& ny2,
102 const double& lx,
103 const double& h1,
104 const double& h2,
105 const bool& periodic_in_x,
106 const bool& build_mesh,
107 TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper);
108
109 /// Access functions for pointers to elements in upper layer
110 FiniteElement*& upper_layer_element_pt(const unsigned long& i)
111 {
113 }
114
115 /// Access functions for pointers to elements in bottom layer
116 FiniteElement*& lower_layer_element_pt(const unsigned long& i)
117 {
119 }
120
121 /// Number of elements in upper layer
122 unsigned long nupper() const
123 {
124 return Upper_layer_element_pt.size();
125 }
126
127 /// Number of elements in top layer
128 unsigned long nlower() const
129 {
130 return Lower_layer_element_pt.size();
131 }
132
133 /// Access functions for pointers to elements in upper layer
138
139 /// Access functions for pointers to elements in bottom layer
144
145 /// Number of elements in upper layer
146 unsigned long ninterface_upper() const
147 {
149 }
150
151 /// Number of elements in top layer
152 unsigned long ninterface_lower() const
153 {
155 }
156
157 /// Index of the face of the elements next to the interface
158 /// in the upper region (always -2)
160 {
161 return -2;
162 }
163
164 /// Index of the face of the elements next to the interface in
165 /// the lower region (always 2)
167 {
168 return 2;
169 }
170
171 /// General node update function implements pure virtual function
172 /// defined in SpineMesh base class and performs specific update
173 /// actions, depending on the node update fct id stored for each node.
175 {
176 unsigned id = spine_node_pt->node_update_fct_id();
177 switch (id)
178 {
179 case 0:
181 break;
182
183 case 1:
185 break;
186
187 default:
188 std::ostringstream error_message;
189 error_message << "Unknown id passed to spine_node_update " << id
190 << std::endl;
191 throw OomphLibError(error_message.str(),
194 }
195 }
196
197 protected:
198 /// Number of elements in lower layer
199 unsigned Ny1;
200
201 /// Number of elements in upper layer
202 unsigned Ny2;
203
204 /// Height of the lower layer
205 double H1;
206
207 /// Height of the upper layer
208 double H2;
209
210 /// Vector of pointers to element in the upper layer
212
213 /// Vector of pointers to element in the lower layer
215
216 /// Vector of pointers to the elements adjacent to the interface
217 /// on the lower layer
219
220 /// Vector of pointers to the element adjacent to the interface
221 /// on the upper layer
223
224
225 /// The spacing function for the x co-ordinates with two
226 /// regions.
227 double x_spacing_function(unsigned xelement,
228 unsigned xnode,
229 unsigned yelement,
230 unsigned ynode);
231
232 /// The spacing function for the y co-ordinates with three
233 /// regions in each fluid.
234 double y_spacing_function(unsigned xelement,
235 unsigned xnode,
236 unsigned yelement,
237 unsigned ynode);
238
239 /// Update function for the lower part of the domain
241 {
242 // Get fraction along the spine
243 double W = spine_node_pt->fraction();
244 // Get spine height
245 double H = spine_node_pt->h();
246 // Set the value of y
247 spine_node_pt->x(1) = this->Ymin + W * H;
248 }
249
250
251 /// Update function for the upper part of the domain
253 {
254 // Get fraction alon the spine
255 double W = spine_node_pt->fraction();
256
257 // Get spine height
258 double H = spine_node_pt->h();
259
260 // Set the value of y
261 spine_node_pt->x(1) =
262 (this->Ymin + H) + W * (this->Ymax - (this->Ymin + H));
263 }
264
265 /// Helper function to actually build the two-layer spine mesh
266 /// (called from various constructors)
268 };
269
270} // namespace oomph
271
273#endif
Collapsible channel mesh with MacroElement-based node update. The collapsible segment is represented ...
RectangularQuadMesh is a two-dimensional mesh of Quad elements with Nx elements in the "x" (horizonal...
const unsigned & nx() const
Return number of elements in x direction.
double Ymax
Maximum value of y coordinate.
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.
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.
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)
FiniteElement *& interface_lower_boundary_element_pt(const unsigned long &i)
Access functions for pointers to elements in bottom layer.
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.
FiniteElement *& lower_layer_element_pt(const unsigned long &i)
Access functions for pointers to elements in bottom 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.
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 *& upper_layer_element_pt(const unsigned long &i)
Access functions for pointers to elements in upper layer.