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// Version 1.0; svn revision $LastChangedRevision$
7 //LIC//
8 //LIC// $LastChangedDate$
9 //LIC//
10 //LIC// Copyright (C) 2006-2016 Matthias Heil and Andrew Hazel
11 //LIC//
12 //LIC// This library is free software; you can redistribute it and/or
13 //LIC// modify it under the terms of the GNU Lesser General Public
14 //LIC// License as published by the Free Software Foundation; either
15 //LIC// version 2.1 of the License, or (at your option) any later version.
16 //LIC//
17 //LIC// This library is distributed in the hope that it will be useful,
18 //LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
19 //LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 //LIC// Lesser General Public License for more details.
21 //LIC//
22 //LIC// You should have received a copy of the GNU Lesser General Public
23 //LIC// License along with this library; if not, write to the Free Software
24 //LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
25 //LIC// 02110-1301 USA.
26 //LIC//
27 //LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
28 //LIC//
29 //LIC//====================================================================
31 
32 namespace oomph
33 {
34 
35  //////////////////////////////////////////////////////////////////////////////
36  //////////////////////////////////////////////////////////////////////////////
37  /// Lowest-order Raviart-Thomas based Darcy/axisym lin elast equation element
38  //////////////////////////////////////////////////////////////////////////////
39  //////////////////////////////////////////////////////////////////////////////
40 
41  //==========================================================================
42  /// \short Constructor. Order 0 elements have 1 pressure dof and no internal
43  /// flux dofs
44  //==========================================================================
45  template<>
47  TElement<2,3>(),
49  Sign_edge(3,1)
50  {
52  }
53 
54  //==========================================================================
55  /// Destructor
56  //==========================================================================
57  template<>
59  {
60  }
61 
62  //==========================================================================
63  /// Return the number of edge basis functions for flux q
64  //==========================================================================
65  template<>
67  {
68  return 3;
69  }
70 
71  //==========================================================================
72  /// Return the number of internal basis functions for flux q
73  //==========================================================================
74  template<>
76  {
77  return 0;
78  }
79 
80  //==========================================================================
81  /// Compute the local form of the q basis at local coordinate s
82  //==========================================================================
83  template<>
85  const Vector<double> &s,
86  Shape &q_basis) const
87  {
88  // RT_0 basis functions
89  q_basis(0,0) = Sign_edge[0]*std::sqrt(2)*s[0];
90  q_basis(0,1) = Sign_edge[0]*std::sqrt(2)*s[1];
91 
92  q_basis(1,0) = Sign_edge[1]*(s[0]-1);
93  q_basis(1,1) = Sign_edge[1]*s[1];
94 
95  q_basis(2,0) = Sign_edge[2]*s[0];
96  q_basis(2,1) = Sign_edge[2]*(s[1]-1);
97  }
98 
99  //===========================================================================
100  /// Compute the local form of the q basis and dbasis/ds at local coordinate s
101  //==========================================================================
102  template<>
104  const Vector<double> &s,
105  Shape &div_q_basis_ds) const
106  {
107  div_q_basis_ds(0) = Sign_edge[0]*2*std::sqrt(2);
108  div_q_basis_ds(1) = Sign_edge[1]*2;
109  div_q_basis_ds(2) = Sign_edge[2]*2;
110 
111  // Scale the basis by the ratio of the length of the edge to the length of
112  // the corresponding edge on the reference element
113  scale_basis(div_q_basis_ds);
114  }
115 
116  //==========================================================================
117  /// Return the number of flux_interpolation points along each
118  /// edge of the element
119  //==========================================================================
120  template<>
123  {
124  return 1;
125  }
126 
127  //==========================================================================
128  /// Return the total number of pressure basis functions
129  //==========================================================================
130  template<>
132  {
133  return 1;
134  }
135 
136  //==========================================================================
137  /// Return the pressure basis
138  //==========================================================================
139  template<>
141  const Vector<double> &s,
142  Shape &p_basis) const
143  {
144  p_basis(0) = 1.0;
145  }
146 
147  //==========================================================================
148  /// The number of values stored at each node
149  //==========================================================================
150  template<>
152  {2,2,2,3,3,3};
153 
154 
155 //===========================================================================
156  /// Face index associated with edge flux degree of freedom
157 //===========================================================================
158  template<>
160  Face_index_of_edge_flux[3] = {2,0,1};
161 
162  //==========================================================================
163  /// Conversion scheme from an edge degree of freedom to the node it's stored
164  /// at
165  //==========================================================================
166  template<>
168 
169 
170  //==========================================================================
171  /// The points along each edge where the fluxes are taken to be
172  //==========================================================================
173  template<>
175  Flux_interpolation_point[1] = {0.5};
176 
177 
178 
179 ///////////////////////////////////////////////////////////////////////////
180 ///////////////////////////////////////////////////////////////////////////
181  /// Second-orderRaviart-Thomas based Darcy/axisym lin elast equation element
182 ///////////////////////////////////////////////////////////////////////////
183 ///////////////////////////////////////////////////////////////////////////
184 
185  //==========================================================================
186  /// \short Constructor. Order 1 elements have 3 pressure dofs and 2 internal
187  /// velocity dofs
188  //==========================================================================
189  template<>
191  TElement<2,3>(),
193  Sign_edge(3,1)
194  {
195  // RT_1 elements have 2 internal degrees of freedom for u, and 3 for p
198  }
199 
200  //==========================================================================
201  /// Destructor
202  //==========================================================================
203  template<>
205  {}
206 
207  //==========================================================================
208  /// Return the number of edge basis functions for flux q
209  //==========================================================================
210  template<>
212  {
213  return 6;
214  }
215 
216  //==========================================================================
217  /// Return the number of internal basis functions for flux q
218  //==========================================================================
219  template<>
221  {
222  return 2;
223  }
224 
225  //==========================================================================
226  /// Returns the local form of the q basis at local coordinate s
227  //==========================================================================
228  template<>
230  const Vector<double> &s,
231  Shape &q_basis) const
232  {
233  // RT_1 basis functions
234  Vector<double> g1_vect,g2_vect;
235  g1_vect=edge_flux_interpolation_point(0,0);
236  g2_vect=edge_flux_interpolation_point(0,1);
237  double g1=g1_vect[0];
238  double g2=g2_vect[0];
239 
240  q_basis(0,0) = Sign_edge[0]*std::sqrt(2.0)*s[0]*(s[1]-g2)/(g1-g2);
241  q_basis(0,1) = Sign_edge[0]*std::sqrt(2.0)*s[1]*(s[1]-g2)/(g1-g2);
242 
243  q_basis(1,0) = Sign_edge[0]*std::sqrt(2.0)*s[0]*(s[1]-g1)/(g2-g1);
244  q_basis(1,1) = Sign_edge[0]*std::sqrt(2.0)*s[1]*(s[1]-g1)/(g2-g1);
245 
246  g1_vect=edge_flux_interpolation_point(1,0);
247  g2_vect=edge_flux_interpolation_point(1,1);
248  g1=g1_vect[0];
249  g2=g2_vect[0];
250  q_basis(2,0) = Sign_edge[1]*(s[0]-1.0)*(s[1]-g1)/(g2-g1);
251  q_basis(2,1) = Sign_edge[1]*s[1]*(s[1]-g1)/(g2-g1);
252 
253  q_basis(3,0) = Sign_edge[1]*(s[0]-1.0)*(s[1]-g2)/(g1-g2);
254  q_basis(3,1) = Sign_edge[1]*s[1]*(s[1]-g2)/(g1-g2);
255 
256  g1_vect=edge_flux_interpolation_point(2,0);
257  g2_vect=edge_flux_interpolation_point(2,1);
258 
259  g1=g1_vect[0];
260  g2=g2_vect[0];
261  q_basis(4,0) = Sign_edge[2]*s[0]*(s[0]-g2)/(g1-g2);
262  q_basis(4,1) = Sign_edge[2]*(s[1]-1.0)*(s[0]-g2)/(g1-g2);
263 
264  q_basis(5,0) = Sign_edge[2]*s[0]*(s[0]-g1)/(g2-g1);
265  q_basis(5,1) = Sign_edge[2]*(s[1]-1.0)*(s[0]-g1)/(g2-g1);
266 
267  q_basis(6,0) = s[1]*s[0];
268  q_basis(6,1) = s[1]*(s[1]-1.0);
269 
270  q_basis(7,0) = s[0]*(s[0]-1.0);
271  q_basis(7,1) = s[0]*s[1];
272  }
273 
274  //==========================================================================
275  /// Returns the local form of the q basis and dbasis/ds at local coordinate s
276  //==========================================================================
277  template<>
279  const Vector<double> &s,
280  Shape &div_q_basis_ds) const
281  {
282 
283  Vector<double> g1_vect,g2_vect;
284  g1_vect=edge_flux_interpolation_point(0,0);
285  g2_vect=edge_flux_interpolation_point(0,1);
286  double g1=g1_vect[0];
287  double g2=g2_vect[0];
288  div_q_basis_ds(0) = Sign_edge[0]*std::sqrt(2.0)*(3.0*s[1]-2.0*g2)/(g1-g2);
289  div_q_basis_ds(1) = Sign_edge[0]*std::sqrt(2.0)*(2.0*g1-3.0*s[1])/(g1-g2);
290 
291  g1_vect=edge_flux_interpolation_point(1,0);
292  g2_vect=edge_flux_interpolation_point(1,1);
293  g1=g1_vect[0];
294  g2=g2_vect[0];
295  div_q_basis_ds(2) = Sign_edge[1]*(2.0*g1-3.0*s[1])/(g1-g2);
296  div_q_basis_ds(3) = Sign_edge[1]*(3.0*s[1]-2.0*g2)/(g1-g2);
297 
298  g1_vect=edge_flux_interpolation_point(2,0);
299  g2_vect=edge_flux_interpolation_point(2,1);
300  g1=g1_vect[0];
301  g2=g2_vect[0];
302  div_q_basis_ds(4) = Sign_edge[2]*(3.0*s[0]-2.0*g2)/(g1-g2);
303  div_q_basis_ds(5) = Sign_edge[2]*(2.0*g1-3.0*s[0])/(g1-g2);
304 
305  div_q_basis_ds(6) = 3.0*s[1]-1.0;
306  div_q_basis_ds(7) = 3.0*s[0]-1.0;
307 
308  // Scale the basis by the ratio of the length of the edge to the length of
309  // the corresponding edge on the reference element
310  scale_basis(div_q_basis_ds);
311  }
312 
313  //==========================================================================
314  /// Returns the number of flux_interpolation points along each edge of
315  /// the element
316  //==========================================================================
317  template<>
320  {
321  return 2;
322  }
323 
324  //==========================================================================
325  /// Return the total number of pressure basis functions
326  //==========================================================================
327  template<>
329  {
330  return 3;
331  }
332 
333  //==========================================================================
334  /// Return the pressure basis
335  //==========================================================================
336  template<>
338  const Vector<double> &s,
339  Shape &p_basis) const
340  {
341  p_basis(0) = 1.0;
342  p_basis(1) = s[0];
343  p_basis(2) = s[1];
344  }
345 
346  //==========================================================================
347  /// The number of values stored at each node
348  //==========================================================================
349  template<>
351  {2,2,2,4,4,4};
352 
353 
354 //===========================================================================
355  /// Face index associated with edge flux degree of freedom
356 //===========================================================================
357  template<>
359  Face_index_of_edge_flux[3] = {2,0,1};
360 
361 
362  //==========================================================================
363 /// Conversion scheme from an edge degree of freedom to the node it's stored at
364  //==========================================================================
365  template<>
367 
368  //==========================================================================
369  /// The points along each edge where the fluxes are taken to be
370  //==========================================================================
371  template<>
373  Flux_interpolation_point[2] =
374  {0.5-std::sqrt(3.0)/6.0,
375  0.5+std::sqrt(3.0)/6.0};
376 
377 
378  //==========================================================================
379  // Force building of templates
380  //==========================================================================
383 
384 } // End of oomph namespace
385 
Vector< double > edge_flux_interpolation_point(const unsigned &edge, const unsigned &j) const
Returns the local coordinate of the jth flux_interpolation point along specified edge.
unsigned nq_basis_internal() const
Return the number of internal basis functions for flux q.
unsigned Q_internal_data_index
The internal data index where the internal q degrees of freedom are stored.
unsigned nq_basis_edge() const
Return the number of edge basis functions for flux q.
static char t char * s
Definition: cfortran.h:572
A class that represents a collection of data; each Data object may contain many different individual ...
Definition: nodes.h:89
void get_p_basis(const Vector< double > &s, Shape &p_basis) const
Return the pressure basis.
std::vector< short > Sign_edge
Unit normal signs associated with each edge to ensure inter-element continuity of the flux...
unsigned P_internal_data_index
The internal data index where the p degrees of freedom are stored.
unsigned np_basis() const
Return the total number of pressure basis functions.
Class implementing the generic maths of the axisym poroelasticity equations: axisym linear elasticity...
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:66
unsigned nedge_flux_interpolation_point() const
Returns the number of flux_interpolation points along each edge of the element.
void scale_basis(Shape &basis) const
Scale the edge basis to allow arbitrary edge mappings.
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.
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.