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// 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  /// Lowest-order Raviart-Thomas based Darcy equation element
35 
36  /// \short Constructor. Order 0 elements have 1 pressure dof and no internal
37  /// velocity dofs
38  template<>
40  TElement<2,3>(), PoroelasticityEquations<2>(), Sign_edge(3,1)
41  {
43  }
44 
45  /// Destructor
46  template<>
48  {
49  }
50 
51  /// Return the total number of computational basis functions for u
52  template<>
54  {
55  return 3;
56  }
57 
58  /// Return the number of edge basis functions for u
59  template<>
61  {
62  return 3;
63  }
64 
65  /// Returns the local form of the q basis at local coordinate s
66  template<>
68  Shape &q_basis) const
69  {
70  // RT_0 basis functions
71  q_basis(0,0) = Sign_edge[0]*std::sqrt(2)*s[0];
72  q_basis(0,1) = Sign_edge[0]*std::sqrt(2)*s[1];
73 
74  q_basis(1,0) = Sign_edge[1]*(s[0]-1);
75  q_basis(1,1) = Sign_edge[1]*s[1];
76 
77  q_basis(2,0) = Sign_edge[2]*s[0];
78  q_basis(2,1) = Sign_edge[2]*(s[1]-1);
79  }
80 
81  /// Returns the local form of the q basis and dbasis/ds at local coordinate s
82  template<>
84  const Vector<double> &s,
85  Shape &div_q_basis_ds) const
86  {
87  div_q_basis_ds(0) = Sign_edge[0]*2*std::sqrt(2);
88  div_q_basis_ds(1) = Sign_edge[1]*2;
89  div_q_basis_ds(2) = Sign_edge[2]*2;
90 
91  // Scale the basis by the ratio of the length of the edge to the length of
92  // the corresponding edge on the reference element
93  scale_basis(div_q_basis_ds);
94  }
95 
96  /// Returns the number of gauss points along each edge of the element
97  template<>
99  {
100  return 1;
101  }
102 
103  /// Return the total number of pressure basis functions
104  template<>
106  {
107  return 1;
108  }
109 
110  /// Return the pressure basis
111  template<>
113  Shape &p_basis) const
114  {
115  p_basis(0) = 1.0;
116  }
117 
118  /// The number of values stored at each node
119  template<>
121  //{0,0,0,1,1,1};
122  {2,2,2,3,3,3};
123 
124  /// Conversion scheme from an edge degree of freedom to the node it's stored at
125  template<>
126  const unsigned TPoroelasticityElement<0>::Q_edge_conv[3] = {3,4,5};
127 
128  /// The points along each edge where the fluxes are taken to be
129  template<>
130  const double TPoroelasticityElement<0>::Gauss_point[1] = {0.5};
131 
132  /// Second-order Raviart-Thomas based Darcy equation element
133 
134  /// \short Constructor. Order 1 elements have 3 pressure dofs and 2 internal
135  /// velocity dofs
136  template<>
138  TElement<2,3>(), PoroelasticityEquations<2>(), Sign_edge(3,1)
139  {
140  // RT_1 elements have 2 internal degrees of freedom for u, and 3 for p
143  }
144 
145  /// Destructor
146  template<>
148  {
149  }
150 
151  /// Return the total number of computational basis functions for u
152  template<>
154  {
155  return 8;
156  }
157 
158  /// Return the number of edge basis functions for u
159  template<>
161  {
162  return 6;
163  }
164 
165  /// Returns the local form of the q basis at local coordinate s
166  template<>
168  Shape &q_basis) const
169  {
170  // RT_1 basis functions
171 
172  double g1,g2;
173 
174  g1=edge_gauss_point(0,0);
175  g2=edge_gauss_point(0,1);
176  q_basis(0,0) = Sign_edge[0]*std::sqrt(2)*s[0]*(s[1]-g2)/(g1-g2);
177  q_basis(0,1) = Sign_edge[0]*std::sqrt(2)*s[1]*(s[1]-g2)/(g1-g2);
178 
179  q_basis(1,0) = Sign_edge[0]*std::sqrt(2)*s[0]*(s[1]-g1)/(g2-g1);
180  q_basis(1,1) = 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,
209  Shape &div_q_basis_ds) const
210  {
211  double g1,g2;
212 
213  g1=edge_gauss_point(0,0);
214  g2=edge_gauss_point(0,1);
215  div_q_basis_ds(0) = Sign_edge[0]*std::sqrt(2)*(3*s[1]-2*g2)/(g1-g2);
216  div_q_basis_ds(1) = Sign_edge[0]*std::sqrt(2)*(2*g1-3*s[1])/(g1-g2);
217 
218  g1=edge_gauss_point(1,0);
219  g2=edge_gauss_point(1,1);
220  div_q_basis_ds(2) = Sign_edge[1]*(2*g1-3*s[1])/(g1-g2);
221  div_q_basis_ds(3) = Sign_edge[1]*(3*s[1]-2*g2)/(g1-g2);
222 
223  g1=edge_gauss_point(2,0);
224  g2=edge_gauss_point(2,1);
225  div_q_basis_ds(4) = Sign_edge[2]*(3*s[0]-2*g2)/(g1-g2);
226  div_q_basis_ds(5) = Sign_edge[2]*(2*g1-3*s[0])/(g1-g2);
227 
228  div_q_basis_ds(6) = 3*s[1]-1;
229  div_q_basis_ds(7) = 3*s[0]-1;
230 
231  // Scale the basis by the ratio of the length of the edge to the length of
232  // the corresponding edge on the reference element
233  scale_basis(div_q_basis_ds);
234  }
235 
236  /// Returns the number of gauss points along each edge of the element
237  template<>
239  {
240  return 2;
241  }
242 
243  /// Return the total number of pressure basis functions
244  template<>
246  {
247  return 3;
248  }
249 
250  /// Return the pressure basis
251  template<>
253  Shape &p_basis) const
254  {
255  p_basis(0) = 1.0;
256  p_basis(1) = s[0];
257  p_basis(2) = s[1];
258  }
259 
260  /// The number of values stored at each node
261  template<>
263  {2,2,2,4,4,4};
264 
265  /// Conversion scheme from an edge degree of freedom to the node it's stored at
266  template<>
267  const unsigned TPoroelasticityElement<1>::Q_edge_conv[3] = {3,4,5};
268 
269  /// The points along each edge where the fluxes are taken to be
270  template<>
272  {0.5-std::sqrt(3.0)/6.0,
273  0.5+std::sqrt(3.0)/6.0};
274 
275  // Force building of templates
276  template class TPoroelasticityElement<0>;
277  template class TPoroelasticityElement<1>;
278 
279 } // End of oomph namespace
280 
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.
~TPoroelasticityElement()
Destructor.
TPoroelasticityElement()
Constructor.
void scale_basis(Shape &basis) const
Scale the edge basis to allow arbitrary edge mappings.
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.
Element which solves the Darcy equations using TElements.
unsigned nq_basis_edge() const
Return the number of edge basis functions for u.
static char t char * s
Definition: cfortran.h:572
double edge_gauss_point(const unsigned &edge, const unsigned &n) const
unsigned np_basis() const
Return the total number of pressure basis functions.
A class that represents a collection of data; each Data object may contain many different individual ...
Definition: nodes.h:89
Class implementing the generic maths of the poroelasticity equations: linear elasticity coupled with ...
unsigned nq_basis() const
Return the total number of computational basis functions for u.
std::vector< short > Sign_edge
Unit normal signs associated with each edge to ensure inter-element continuity of the flux...
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 P_internal_data_index
The internal data index where the p degrees of freedom are stored.
unsigned Q_internal_data_index
The internal data index where the internal q degrees of freedom are stored.
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.