horizontal_single_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_HORIZONTAL_SINGLE_LAYER_SPINE_MESH_HEADER
27 #define OOMPH_HORIZONTAL_SINGLE_LAYER_SPINE_MESH_HEADER
28 
29 // oomph-lib includes
30 #include "../generic/spines.h"
32 
33 // Created by Francisco
34 
35 namespace oomph
36 {
37  //======================================================================
43  //======================================================================
44  template<class ELEMENT>
46  public SpineMesh
47  {
48  public:
53  const unsigned& nx,
54  const unsigned& ny,
55  const double& lx,
56  const double& h,
57  TimeStepper* time_stepper_pt = &Mesh::Default_TimeStepper);
58 
59 
63  virtual void spine_node_update(SpineNode* spine_node_pt)
64  {
65  // Get fraction along the spine
66  double W = spine_node_pt->fraction();
67  // Get spine height
68  double H = spine_node_pt->h();
69  // Set the value of y
70  // spine_node_pt->x(0) = this->Xmin + W*H;
71  spine_node_pt->x(0) = W * H;
72  }
73 
74 
75  protected:
79  TimeStepper* time_stepper_pt);
80  };
81 
82 } // namespace oomph
83 
84 #endif
MatrixXf H
Definition: HessenbergDecomposition_matrixH.cpp:4
Definition: horizontal_single_layer_spine_mesh.template.h:47
virtual void spine_node_update(SpineNode *spine_node_pt)
Definition: horizontal_single_layer_spine_mesh.template.h:63
HorizontalSingleLayerSpineMesh(const unsigned &nx, const unsigned &ny, const double &lx, const double &h, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Definition: horizontal_single_layer_spine_mesh.template.cc:47
virtual void build_horizontal_single_layer_mesh(TimeStepper *time_stepper_pt)
Definition: horizontal_single_layer_spine_mesh.template.cc:75
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: rectangular_quadmesh.template.h:59
const unsigned & ny() const
Return number of elements in y direction.
Definition: rectangular_quadmesh.template.h:231
const unsigned & nx() const
Return number of elements in x direction.
Definition: rectangular_quadmesh.template.h:224
Definition: spines.h:613
Definition: spines.h:328
double & h()
Access function to spine height.
Definition: spines.h:397
double & fraction()
Set reference to fraction along spine.
Definition: spines.h:378
Definition: timesteppers.h:231
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