Tporoelasticity_elements.cc
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//====================================================================
27 
28 namespace oomph
29 {
30  /// Lowest-order Raviart-Thomas based Darcy equation element
31 
32  /// Constructor. Order 0 elements have 1 pressure dof and no internal
33  /// velocity dofs
34  template<>
36  : TElement<2, 3>(), PoroelasticityEquations<2>(), Sign_edge(3, 1)
37  {
39  }
40 
41  /// Destructor
42  template<>
44  {
45  }
46 
47  /// Return the total number of computational basis functions for u
48  template<>
50  {
51  return 3;
52  }
53 
54  /// Return the number of edge basis functions for u
55  template<>
57  {
58  return 3;
59  }
60 
61  /// Returns the local form of the q basis at local coordinate s
62  template<>
64  Shape& q_basis) const
65  {
66  // RT_0 basis functions
67  q_basis(0, 0) = Sign_edge[0] * std::sqrt(2) * s[0];
68  q_basis(0, 1) = Sign_edge[0] * std::sqrt(2) * s[1];
69 
70  q_basis(1, 0) = Sign_edge[1] * (s[0] - 1);
71  q_basis(1, 1) = Sign_edge[1] * s[1];
72 
73  q_basis(2, 0) = Sign_edge[2] * s[0];
74  q_basis(2, 1) = Sign_edge[2] * (s[1] - 1);
75  }
76 
77  /// Returns the local form of the q basis and dbasis/ds at local coordinate s
78  template<>
80  const Vector<double>& s, Shape& div_q_basis_ds) const
81  {
82  div_q_basis_ds(0) = Sign_edge[0] * 2 * std::sqrt(2);
83  div_q_basis_ds(1) = Sign_edge[1] * 2;
84  div_q_basis_ds(2) = Sign_edge[2] * 2;
85 
86  // Scale the basis by the ratio of the length of the edge to the length of
87  // the corresponding edge on the reference element
88  scale_basis(div_q_basis_ds);
89  }
90 
91  /// Returns the number of gauss points along each edge of the element
92  template<>
94  {
95  return 1;
96  }
97 
98  /// Return the total number of pressure basis functions
99  template<>
101  {
102  return 1;
103  }
104 
105  /// Return the pressure basis
106  template<>
108  Shape& p_basis) const
109  {
110  p_basis(0) = 1.0;
111  }
112 
113  /// The number of values stored at each node
114  template<>
116  //{0,0,0,1,1,1};
117  {2, 2, 2, 3, 3, 3};
118 
119  /// Conversion scheme from an edge degree of freedom to the node it's stored
120  /// at
121  template<>
122  const unsigned TPoroelasticityElement<0>::Q_edge_conv[3] = {3, 4, 5};
123 
124  /// The points along each edge where the fluxes are taken to be
125  template<>
126  const double TPoroelasticityElement<0>::Gauss_point[1] = {0.5};
127 
128  /// Second-order Raviart-Thomas based Darcy equation element
129 
130  /// Constructor. Order 1 elements have 3 pressure dofs and 2 internal
131  /// velocity dofs
132  template<>
134  : TElement<2, 3>(), PoroelasticityEquations<2>(), Sign_edge(3, 1)
135  {
136  // RT_1 elements have 2 internal degrees of freedom for u, and 3 for p
139  }
140 
141  /// Destructor
142  template<>
144  {
145  }
146 
147  /// Return the total number of computational basis functions for u
148  template<>
150  {
151  return 8;
152  }
153 
154  /// Return the number of edge basis functions for u
155  template<>
157  {
158  return 6;
159  }
160 
161  /// Returns the local form of the q basis at local coordinate s
162  template<>
164  Shape& q_basis) const
165  {
166  // RT_1 basis functions
167 
168  double g1, g2;
169 
170  g1 = edge_gauss_point(0, 0);
171  g2 = edge_gauss_point(0, 1);
172  q_basis(0, 0) =
173  Sign_edge[0] * std::sqrt(2) * s[0] * (s[1] - g2) / (g1 - g2);
174  q_basis(0, 1) =
175  Sign_edge[0] * std::sqrt(2) * s[1] * (s[1] - g2) / (g1 - g2);
176 
177  q_basis(1, 0) =
178  Sign_edge[0] * std::sqrt(2) * s[0] * (s[1] - g1) / (g2 - g1);
179  q_basis(1, 1) =
180  Sign_edge[0] * std::sqrt(2) * s[1] * (s[1] - g1) / (g2 - g1);
181 
182  g1 = edge_gauss_point(1, 0);
183  g2 = edge_gauss_point(1, 1);
184  q_basis(2, 0) = Sign_edge[1] * (s[0] - 1) * (s[1] - g1) / (g2 - g1);
185  q_basis(2, 1) = Sign_edge[1] * s[1] * (s[1] - g1) / (g2 - g1);
186 
187  q_basis(3, 0) = Sign_edge[1] * (s[0] - 1) * (s[1] - g2) / (g1 - g2);
188  q_basis(3, 1) = Sign_edge[1] * s[1] * (s[1] - g2) / (g1 - g2);
189 
190  g1 = edge_gauss_point(2, 0);
191  g2 = edge_gauss_point(2, 1);
192  q_basis(4, 0) = Sign_edge[2] * s[0] * (s[0] - g2) / (g1 - g2);
193  q_basis(4, 1) = Sign_edge[2] * (s[1] - 1) * (s[0] - g2) / (g1 - g2);
194 
195  q_basis(5, 0) = Sign_edge[2] * s[0] * (s[0] - g1) / (g2 - g1);
196  q_basis(5, 1) = Sign_edge[2] * (s[1] - 1) * (s[0] - g1) / (g2 - g1);
197 
198  q_basis(6, 0) = s[1] * s[0];
199  q_basis(6, 1) = s[1] * (s[1] - 1);
200 
201  q_basis(7, 0) = s[0] * (s[0] - 1);
202  q_basis(7, 1) = s[0] * s[1];
203  }
204 
205  /// Returns the local form of the q basis and dbasis/ds at local coordinate s
206  template<>
208  const Vector<double>& s, Shape& div_q_basis_ds) const
209  {
210  double g1, g2;
211 
212  g1 = edge_gauss_point(0, 0);
213  g2 = edge_gauss_point(0, 1);
214  div_q_basis_ds(0) =
215  Sign_edge[0] * std::sqrt(2) * (3 * s[1] - 2 * g2) / (g1 - g2);
216  div_q_basis_ds(1) =
217  Sign_edge[0] * std::sqrt(2) * (2 * g1 - 3 * s[1]) / (g1 - g2);
218 
219  g1 = edge_gauss_point(1, 0);
220  g2 = edge_gauss_point(1, 1);
221  div_q_basis_ds(2) = Sign_edge[1] * (2 * g1 - 3 * s[1]) / (g1 - g2);
222  div_q_basis_ds(3) = Sign_edge[1] * (3 * s[1] - 2 * g2) / (g1 - g2);
223 
224  g1 = edge_gauss_point(2, 0);
225  g2 = edge_gauss_point(2, 1);
226  div_q_basis_ds(4) = Sign_edge[2] * (3 * s[0] - 2 * g2) / (g1 - g2);
227  div_q_basis_ds(5) = Sign_edge[2] * (2 * g1 - 3 * s[0]) / (g1 - g2);
228 
229  div_q_basis_ds(6) = 3 * s[1] - 1;
230  div_q_basis_ds(7) = 3 * s[0] - 1;
231 
232  // Scale the basis by the ratio of the length of the edge to the length of
233  // the corresponding edge on the reference element
234  scale_basis(div_q_basis_ds);
235  }
236 
237  /// Returns the number of gauss points along each edge of the element
238  template<>
240  {
241  return 2;
242  }
243 
244  /// Return the total number of pressure basis functions
245  template<>
247  {
248  return 3;
249  }
250 
251  /// Return the pressure basis
252  template<>
254  Shape& p_basis) const
255  {
256  p_basis(0) = 1.0;
257  p_basis(1) = s[0];
258  p_basis(2) = s[1];
259  }
260 
261  /// The number of values stored at each node
262  template<>
264  2, 2, 2, 4, 4, 4};
265 
266  /// Conversion scheme from an edge degree of freedom to the node it's stored
267  /// at
268  template<>
269  const unsigned TPoroelasticityElement<1>::Q_edge_conv[3] = {3, 4, 5};
270 
271  /// The points along each edge where the fluxes are taken to be
272  template<>
274  0.5 - std::sqrt(3.0) / 6.0, 0.5 + std::sqrt(3.0) / 6.0};
275 
276  // Force building of templates
277  template class TPoroelasticityElement<0>;
278  template class TPoroelasticityElement<1>;
279 
280 } // namespace oomph
static char t char * s
Definition: cfortran.h:568
A class that represents a collection of data; each Data object may contain many different individual ...
Definition: nodes.h:86
unsigned add_internal_data(Data *const &data_pt, const bool &fd=true)
Add a (pointer to an) internal data object to the element and return the index required to obtain it ...
Definition: elements.cc:62
Class implementing the generic maths of the poroelasticity equations: linear elasticity coupled with ...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
General TElement class.
Definition: Telements.h:1208
Element which solves the Darcy equations using TElements.
TPoroelasticityElement()
Constructor.
unsigned Q_internal_data_index
The internal data index where the internal q degrees of freedom are stored.
void get_div_q_basis_local(const Vector< double > &s, Shape &div_q_basis_ds) const
Returns the local form of the q basis and dbasis/ds at local coordinate s.
unsigned nedge_gauss_point() const
Returns the number of gauss points along each edge of the element.
void get_p_basis(const Vector< double > &s, Shape &p_basis) const
Return the pressure basis.
unsigned np_basis() const
Return the total number of pressure basis functions.
unsigned nq_basis_edge() const
Return the number of edge basis functions for u.
unsigned P_internal_data_index
The internal data index where the p degrees of freedom are stored.
void get_q_basis_local(const Vector< double > &s, Shape &q_basis) const
Returns the local form of the q basis at local coordinate s.
~TPoroelasticityElement()
Destructor.
unsigned nq_basis() const
Return the total number of computational basis functions for u.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...