Taxisym_poroelasticity_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-2024 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  /// ///////////////////////////////////////////////////////////////////////////
31  /// ///////////////////////////////////////////////////////////////////////////
32  /// Lowest-order Raviart-Thomas based Darcy/axisym lin elast equation element
33  /// ///////////////////////////////////////////////////////////////////////////
34  /// ///////////////////////////////////////////////////////////////////////////
35 
36  //==========================================================================
37  /// Constructor. Order 0 elements have 1 pressure dof and no internal
38  /// flux dofs
39  //==========================================================================
40  template<>
42  : TElement<2, 3>(), AxisymmetricPoroelasticityEquations(), Sign_edge(3, 1)
43  {
45  }
46 
47  //==========================================================================
48  /// Destructor
49  //==========================================================================
50  template<>
52  {
53  }
54 
55  //==========================================================================
56  /// Return the number of edge basis functions for flux q
57  //==========================================================================
58  template<>
60  {
61  return 3;
62  }
63 
64  //==========================================================================
65  /// Return the number of internal basis functions for flux q
66  //==========================================================================
67  template<>
69  {
70  return 0;
71  }
72 
73  //==========================================================================
74  /// Compute the local form of the q basis at local coordinate s
75  //==========================================================================
76  template<>
78  const Vector<double>& s, Shape& q_basis) const
79  {
80  // RT_0 basis functions
81  q_basis(0, 0) = Sign_edge[0] * std::sqrt(2) * s[0];
82  q_basis(0, 1) = Sign_edge[0] * std::sqrt(2) * s[1];
83 
84  q_basis(1, 0) = Sign_edge[1] * (s[0] - 1);
85  q_basis(1, 1) = Sign_edge[1] * s[1];
86 
87  q_basis(2, 0) = Sign_edge[2] * s[0];
88  q_basis(2, 1) = Sign_edge[2] * (s[1] - 1);
89  }
90 
91  //===========================================================================
92  /// Compute the local form of the q basis and dbasis/ds at local coordinate s
93  //==========================================================================
94  template<>
96  const Vector<double>& s, Shape& div_q_basis_ds) const
97  {
98  div_q_basis_ds(0) = Sign_edge[0] * 2 * std::sqrt(2);
99  div_q_basis_ds(1) = Sign_edge[1] * 2;
100  div_q_basis_ds(2) = Sign_edge[2] * 2;
101 
102  // Scale the basis by the ratio of the length of the edge to the length of
103  // the corresponding edge on the reference element
104  scale_basis(div_q_basis_ds);
105  }
106 
107  //==========================================================================
108  /// Return the number of flux_interpolation points along each
109  /// edge of the element
110  //==========================================================================
111  template<>
113  0>::nedge_flux_interpolation_point() const
114  {
115  return 1;
116  }
117 
118  //==========================================================================
119  /// Return the total number of pressure basis functions
120  //==========================================================================
121  template<>
123  {
124  return 1;
125  }
126 
127  //==========================================================================
128  /// Return the pressure basis
129  //==========================================================================
130  template<>
132  const Vector<double>& s, Shape& p_basis) const
133  {
134  p_basis(0) = 1.0;
135  }
136 
137  //==========================================================================
138  /// The number of values stored at each node
139  //==========================================================================
140  template<>
142  2, 2, 2, 3, 3, 3};
143 
144 
145  //===========================================================================
146  /// Face index associated with edge flux degree of freedom
147  //===========================================================================
148  template<>
149  const unsigned
151  2, 0, 1};
152 
153  //==========================================================================
154  /// Conversion scheme from an edge degree of freedom to the node it's stored
155  /// at
156  //==========================================================================
157  template<>
159  3, 4, 5};
160 
161 
162  //==========================================================================
163  /// The points along each edge where the fluxes are taken to be
164  //==========================================================================
165  template<>
166  const double
168 
169 
170  /// ////////////////////////////////////////////////////////////////////////
171  /// ////////////////////////////////////////////////////////////////////////
172  /// Second-orderRaviart-Thomas based Darcy/axisym lin elast equation element
173  /// ////////////////////////////////////////////////////////////////////////
174  /// ////////////////////////////////////////////////////////////////////////
175 
176  //==========================================================================
177  /// Constructor. Order 1 elements have 3 pressure dofs and 2 internal
178  /// velocity dofs
179  //==========================================================================
180  template<>
182  : TElement<2, 3>(), AxisymmetricPoroelasticityEquations(), Sign_edge(3, 1)
183  {
184  // RT_1 elements have 2 internal degrees of freedom for u, and 3 for p
187  }
188 
189  //==========================================================================
190  /// Destructor
191  //==========================================================================
192  template<>
194  {
195  }
196 
197  //==========================================================================
198  /// Return the number of edge basis functions for flux q
199  //==========================================================================
200  template<>
202  {
203  return 6;
204  }
205 
206  //==========================================================================
207  /// Return the number of internal basis functions for flux q
208  //==========================================================================
209  template<>
211  {
212  return 2;
213  }
214 
215  //==========================================================================
216  /// Returns the local form of the q basis at local coordinate s
217  //==========================================================================
218  template<>
220  const Vector<double>& s, Shape& q_basis) const
221  {
222  // RT_1 basis functions
223  Vector<double> g1_vect, g2_vect;
224  g1_vect = edge_flux_interpolation_point(0, 0);
225  g2_vect = edge_flux_interpolation_point(0, 1);
226  double g1 = g1_vect[0];
227  double g2 = g2_vect[0];
228 
229  q_basis(0, 0) =
230  Sign_edge[0] * std::sqrt(2.0) * s[0] * (s[1] - g2) / (g1 - g2);
231  q_basis(0, 1) =
232  Sign_edge[0] * std::sqrt(2.0) * s[1] * (s[1] - g2) / (g1 - g2);
233 
234  q_basis(1, 0) =
235  Sign_edge[0] * std::sqrt(2.0) * s[0] * (s[1] - g1) / (g2 - g1);
236  q_basis(1, 1) =
237  Sign_edge[0] * std::sqrt(2.0) * s[1] * (s[1] - g1) / (g2 - g1);
238 
239  g1_vect = edge_flux_interpolation_point(1, 0);
240  g2_vect = edge_flux_interpolation_point(1, 1);
241  g1 = g1_vect[0];
242  g2 = g2_vect[0];
243  q_basis(2, 0) = Sign_edge[1] * (s[0] - 1.0) * (s[1] - g1) / (g2 - g1);
244  q_basis(2, 1) = Sign_edge[1] * s[1] * (s[1] - g1) / (g2 - g1);
245 
246  q_basis(3, 0) = Sign_edge[1] * (s[0] - 1.0) * (s[1] - g2) / (g1 - g2);
247  q_basis(3, 1) = Sign_edge[1] * s[1] * (s[1] - g2) / (g1 - g2);
248 
249  g1_vect = edge_flux_interpolation_point(2, 0);
250  g2_vect = edge_flux_interpolation_point(2, 1);
251 
252  g1 = g1_vect[0];
253  g2 = g2_vect[0];
254  q_basis(4, 0) = Sign_edge[2] * s[0] * (s[0] - g2) / (g1 - g2);
255  q_basis(4, 1) = Sign_edge[2] * (s[1] - 1.0) * (s[0] - g2) / (g1 - g2);
256 
257  q_basis(5, 0) = Sign_edge[2] * s[0] * (s[0] - g1) / (g2 - g1);
258  q_basis(5, 1) = Sign_edge[2] * (s[1] - 1.0) * (s[0] - g1) / (g2 - g1);
259 
260  q_basis(6, 0) = s[1] * s[0];
261  q_basis(6, 1) = s[1] * (s[1] - 1.0);
262 
263  q_basis(7, 0) = s[0] * (s[0] - 1.0);
264  q_basis(7, 1) = s[0] * s[1];
265  }
266 
267  //==========================================================================
268  /// Returns the local form of the q basis and dbasis/ds at local coordinate s
269  //==========================================================================
270  template<>
272  const Vector<double>& s, Shape& div_q_basis_ds) const
273  {
274  Vector<double> g1_vect, g2_vect;
275  g1_vect = edge_flux_interpolation_point(0, 0);
276  g2_vect = edge_flux_interpolation_point(0, 1);
277  double g1 = g1_vect[0];
278  double g2 = g2_vect[0];
279  div_q_basis_ds(0) =
280  Sign_edge[0] * std::sqrt(2.0) * (3.0 * s[1] - 2.0 * g2) / (g1 - g2);
281  div_q_basis_ds(1) =
282  Sign_edge[0] * std::sqrt(2.0) * (2.0 * g1 - 3.0 * s[1]) / (g1 - g2);
283 
284  g1_vect = edge_flux_interpolation_point(1, 0);
285  g2_vect = edge_flux_interpolation_point(1, 1);
286  g1 = g1_vect[0];
287  g2 = g2_vect[0];
288  div_q_basis_ds(2) = Sign_edge[1] * (2.0 * g1 - 3.0 * s[1]) / (g1 - g2);
289  div_q_basis_ds(3) = Sign_edge[1] * (3.0 * s[1] - 2.0 * g2) / (g1 - g2);
290 
291  g1_vect = edge_flux_interpolation_point(2, 0);
292  g2_vect = edge_flux_interpolation_point(2, 1);
293  g1 = g1_vect[0];
294  g2 = g2_vect[0];
295  div_q_basis_ds(4) = Sign_edge[2] * (3.0 * s[0] - 2.0 * g2) / (g1 - g2);
296  div_q_basis_ds(5) = Sign_edge[2] * (2.0 * g1 - 3.0 * s[0]) / (g1 - g2);
297 
298  div_q_basis_ds(6) = 3.0 * s[1] - 1.0;
299  div_q_basis_ds(7) = 3.0 * s[0] - 1.0;
300 
301  // Scale the basis by the ratio of the length of the edge to the length of
302  // the corresponding edge on the reference element
303  scale_basis(div_q_basis_ds);
304  }
305 
306  //==========================================================================
307  /// Returns the number of flux_interpolation points along each edge of
308  /// the element
309  //==========================================================================
310  template<>
312  1>::nedge_flux_interpolation_point() const
313  {
314  return 2;
315  }
316 
317  //==========================================================================
318  /// Return the total number of pressure basis functions
319  //==========================================================================
320  template<>
322  {
323  return 3;
324  }
325 
326  //==========================================================================
327  /// Return the pressure basis
328  //==========================================================================
329  template<>
331  const Vector<double>& s, Shape& p_basis) const
332  {
333  p_basis(0) = 1.0;
334  p_basis(1) = s[0];
335  p_basis(2) = s[1];
336  }
337 
338  //==========================================================================
339  /// The number of values stored at each node
340  //==========================================================================
341  template<>
343  2, 2, 2, 4, 4, 4};
344 
345 
346  //===========================================================================
347  /// Face index associated with edge flux degree of freedom
348  //===========================================================================
349  template<>
350  const unsigned
352  2, 0, 1};
353 
354 
355  //==========================================================================
356  /// Conversion scheme from an edge degree of freedom to the node it's stored
357  /// at
358  //==========================================================================
359  template<>
361  3, 4, 5};
362 
363  //==========================================================================
364  /// The points along each edge where the fluxes are taken to be
365  //==========================================================================
366  template<>
367  const double
369  0.5 - std::sqrt(3.0) / 6.0, 0.5 + std::sqrt(3.0) / 6.0};
370 
371 
372  //==========================================================================
373  // Force building of templates
374  //==========================================================================
377 
378 } // namespace oomph
static char t char * s
Definition: cfortran.h:568
Class implementing the generic maths of the axisym poroelasticity equations: axisym linear elasticity...
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:67
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
================================================================= Element which solves the Darcy/line...
void get_p_basis(const Vector< double > &s, Shape &p_basis) const
Return the pressure basis.
unsigned Q_internal_data_index
The internal data index where the internal q degrees of freedom are stored.
unsigned P_internal_data_index
The internal data index where the p degrees of freedom are stored.
unsigned nq_basis_edge() const
Return the number of edge basis functions for flux q.
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.
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 nq_basis_internal() const
Return the number of internal basis functions for flux q.
unsigned np_basis() const
Return the total number of pressure basis functions.
General TElement class.
Definition: Telements.h:1208
//////////////////////////////////////////////////////////////////// ////////////////////////////////...