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-2022 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  //======================================================================
57  //======================================================================
58  template<class ELEMENT>
59  class TwoLayerSpineMesh : public RectangularQuadMesh<ELEMENT>,
60  public SpineMesh
61  {
62  public:
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 
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 
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 
111  FiniteElement*& upper_layer_element_pt(const unsigned long& i)
112  {
113  return Upper_layer_element_pt[i];
114  }
115 
117  FiniteElement*& lower_layer_element_pt(const unsigned long& i)
118  {
119  return Lower_layer_element_pt[i];
120  }
121 
123  unsigned long nupper() const
124  {
125  return Upper_layer_element_pt.size();
126  }
127 
129  unsigned long nlower() const
130  {
131  return Lower_layer_element_pt.size();
132  }
133 
136  {
138  }
139 
142  {
144  }
145 
147  unsigned long ninterface_upper() const
148  {
150  }
151 
153  unsigned long ninterface_lower() const
154  {
156  }
157 
161  {
162  return -2;
163  }
164 
168  {
169  return 2;
170  }
171 
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(),
195  }
196  }
197 
198 
199  protected:
201  unsigned Ny1;
202 
204  unsigned Ny2;
205 
207  double H1;
208 
210  double H2;
211 
214 
217 
221 
225 
226 
229  double x_spacing_function(unsigned xelement,
230  unsigned xnode,
231  unsigned yelement,
232  unsigned ynode);
233 
236  double y_spacing_function(unsigned xelement,
237  unsigned xnode,
238  unsigned yelement,
239  unsigned ynode);
240 
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 
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 
269  virtual void build_two_layer_mesh(TimeStepper* time_stepper_pt);
270  };
271 
272 } // namespace oomph
273 
274 #endif
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Array< double, 1, 3 > e(1./3., 0.5, 2.)
MatrixXf H
Definition: HessenbergDecomposition_matrixH.cpp:4
Definition: elements.h:1313
static Steady< 0 > Default_TimeStepper
The Steady Timestepper.
Definition: mesh.h:75
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition: nodes.h:1060
Definition: oomph_definitions.h:222
Definition: rectangular_quadmesh.template.h:59
double Ymax
Maximum value of y coordinate.
Definition: rectangular_quadmesh.template.h:77
const unsigned & nx() const
Return number of elements in x direction.
Definition: rectangular_quadmesh.template.h:224
double Ymin
Minimum value of y coordinate.
Definition: rectangular_quadmesh.template.h:75
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
Definition: spines.h:613
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
Definition: two_layer_spine_mesh.template.h:61
Vector< FiniteElement * > Upper_layer_element_pt
Vector of pointers to element in the lower layer.
Definition: two_layer_spine_mesh.template.h:216
double H2
Height of the upper layer.
Definition: two_layer_spine_mesh.template.h:210
FiniteElement *& upper_layer_element_pt(const unsigned long &i)
Access functions for pointers to elements in upper layer.
Definition: two_layer_spine_mesh.template.h:111
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)
Definition: two_layer_spine_mesh.template.cc:49
Vector< FiniteElement * > Lower_layer_element_pt
Vector of pointers to element in the upper layer.
Definition: two_layer_spine_mesh.template.h:213
Vector< FiniteElement * > Interface_upper_boundary_element_pt
Definition: two_layer_spine_mesh.template.h:224
double y_spacing_function(unsigned xelement, unsigned xnode, unsigned yelement, unsigned ynode)
Definition: two_layer_spine_mesh.template.cc:227
int interface_upper_face_index_at_boundary(const unsigned &e)
Definition: two_layer_spine_mesh.template.h:160
int interface_lower_face_index_at_boundary(const unsigned &e)
Definition: two_layer_spine_mesh.template.h:167
unsigned Ny2
Number of elements in upper layer.
Definition: two_layer_spine_mesh.template.h:204
unsigned long ninterface_lower() const
Number of elements in top layer.
Definition: two_layer_spine_mesh.template.h:153
double x_spacing_function(unsigned xelement, unsigned xnode, unsigned yelement, unsigned ynode)
Definition: two_layer_spine_mesh.template.cc:211
unsigned Ny1
Number of elements in lower layer.
Definition: two_layer_spine_mesh.template.h:201
FiniteElement *& lower_layer_element_pt(const unsigned long &i)
Access functions for pointers to elements in bottom layer.
Definition: two_layer_spine_mesh.template.h:117
void spine_node_update_lower(SpineNode *spine_node_pt)
Update function for the lower part of the domain.
Definition: two_layer_spine_mesh.template.h:242
unsigned long nupper() const
Number of elements in upper layer.
Definition: two_layer_spine_mesh.template.h:123
FiniteElement *& interface_upper_boundary_element_pt(const unsigned long &i)
Access functions for pointers to elements in upper layer.
Definition: two_layer_spine_mesh.template.h:135
Vector< FiniteElement * > Interface_lower_boundary_element_pt
Definition: two_layer_spine_mesh.template.h:220
unsigned long nlower() const
Number of elements in top layer.
Definition: two_layer_spine_mesh.template.h:129
double H1
Height of the lower layer.
Definition: two_layer_spine_mesh.template.h:207
void spine_node_update(SpineNode *spine_node_pt)
Definition: two_layer_spine_mesh.template.h:175
virtual void build_two_layer_mesh(TimeStepper *time_stepper_pt)
Definition: two_layer_spine_mesh.template.cc:267
unsigned long ninterface_upper() const
Number of elements in upper layer.
Definition: two_layer_spine_mesh.template.h:147
void spine_node_update_upper(SpineNode *spine_node_pt)
Update function for the upper part of the domain.
Definition: two_layer_spine_mesh.template.h:254
FiniteElement *& interface_lower_boundary_element_pt(const unsigned long &i)
Access functions for pointers to elements in bottom layer.
Definition: two_layer_spine_mesh.template.h:141
const double lx
Definition: ConstraintElementsUnitTest.cpp:33
@ W
Definition: quadtree.h:63
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86