Tdarcy_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//====================================================================
30 #include "Tdarcy_elements.h"
31 
32 namespace oomph
33 {
34 
35 
36 //////////////////////////////////////////////////////////////
37 //////////////////////////////////////////////////////////////
38 /// Lowest-order Raviart-Thomas based Darcy equation element
39 //////////////////////////////////////////////////////////////
40 //////////////////////////////////////////////////////////////
41 
42 
43 //===========================================================================
44  /// \short Constructor. Order 0 elements have 1 pressure dof and no internal
45  /// velocity dofs
46 //===========================================================================
47  template<>
49  TElement<2,3>(), DarcyEquations<2>(), Sign_edge(3,1)
50  {
52  }
53 
54 //===========================================================================
55  /// Destructor
56 //===========================================================================
57  template<>
59  {
60  }
61 
62 
63 //===========================================================================
64  /// Return the number of edge basis functions for q
65 //===========================================================================
66  template<>
68  {
69  return 3;
70  }
71 
72 
73 
74 //===========================================================================
75  /// Return the number of internal basis functions for q
76 //===========================================================================
77  template<>
79  {
80  return 0;
81  }
82 
83 //===========================================================================
84  /// Compute the local form of the q basis at local coordinate s
85 //===========================================================================
86  template<>
88  Shape &q_basis) const
89  {
90  // RT_0 basis functions
91  q_basis(0,0) = Sign_edge[0]*std::sqrt(2)*s[0];
92  q_basis(0,1) = Sign_edge[0]*std::sqrt(2)*s[1];
93 
94  q_basis(1,0) = Sign_edge[1]*(s[0]-1);
95  q_basis(1,1) = Sign_edge[1]*s[1];
96 
97  q_basis(2,0) = Sign_edge[2]*s[0];
98  q_basis(2,1) = Sign_edge[2]*(s[1]-1);
99  }
100 
101 //===========================================================================
102  /// Compute the local form of the q basis and dbasis/ds at local coordinate s
103 //===========================================================================
104  template<>
106  const Vector<double> &s,
107  Shape &div_q_basis_ds) const
108  {
109  div_q_basis_ds(0) = Sign_edge[0]*2*std::sqrt(2);
110  div_q_basis_ds(1) = Sign_edge[1]*2;
111  div_q_basis_ds(2) = Sign_edge[2]*2;
112 
113  // Scale the basis by the ratio of the length of the edge to the length of
114  // the corresponding edge on the reference element
115  // hierher explain please
116  scale_basis(div_q_basis_ds);
117  }
118 
119 //===========================================================================
120  /// Return the number of flux interpolation points along each edge
121  /// of the element
122 //===========================================================================
123  template<>
125  {
126  return 1;
127  }
128 
129 //===========================================================================
130  /// Return the equation number of the n-th pressure degree of freedom
131 //===========================================================================
132  template<>
133  int TRaviartThomasDarcyElement<0>::p_local_eqn(const unsigned &n) const
134  {
136  }
137 
138 //===========================================================================
139  /// Return the nth pressure value
140 //===========================================================================
141  template<>
142  double TRaviartThomasDarcyElement<0>::p_value(const unsigned &n) const
143  {
144  return this->internal_data_pt(P_internal_data_index)->value(n);
145  }
146 
147 //===========================================================================
148  /// Return the total number of pressure basis functions
149 //===========================================================================
150  template<>
152  {
153  return 1;
154  }
155 
156 //===========================================================================
157  /// Compute the pressure basis
158 //===========================================================================
159  template<>
161  Shape &p_basis) const
162  {
163  p_basis(0) = 1.0;
164  }
165 
166 
167 
168 //===========================================================================
169 /// Recovery order for Z2 error estimator
170 //===========================================================================
171  template<>
173  {
174  return 2; // Need to experiment with this...
175  }
176 
177 //===========================================================================
178  /// The number of values stored at each node: A single flux dof is
179  /// stored at each midside node
180 //===========================================================================
181  template<>
183  {0,0,0,1,1,1};
184 
185 
186 //===========================================================================
187  /// Face index associated with edge flux degree of freedom
188 //===========================================================================
189  template<>
191  {2,0,1};
192 
193 
194 //===========================================================================
195  /// Conversion scheme from an edge degree of freedom to the node it's stored at
196 /// Fluxes are stored at mid-side nodes
197 //===========================================================================
198  template<>
199  const unsigned TRaviartThomasDarcyElement<0>::Q_edge_conv[3] = {3,4,5};
200 
201 
202 //===========================================================================
203  /// The points along each edge where the fluxes are taken to be
204 // hierher explain please
205 //===========================================================================
206  template<>
208  {0.5};
209 
210 
211 
212 
213 //////////////////////////////////////////////////////////////
214 //////////////////////////////////////////////////////////////
215 /// Second-order Raviart-Thomas based Darcy equation element
216 //////////////////////////////////////////////////////////////
217 //////////////////////////////////////////////////////////////
218 
219 
220 //===========================================================================
221  /// \short Constructor. Order 1 elements have 3 pressure dofs and 2 internal
222  /// velocity dofs
223 //===========================================================================
224  template<>
226  TElement<2,3>(), DarcyEquations<2>(), Sign_edge(3,1)
227  {
228  // RT_1 elements have 2 internal degrees of freedom for q, and 3 for p
231  }
232 
233 //===========================================================================
234  /// Destructor
235 //===========================================================================
236  template<>
238  {
239  }
240 
241 
242 //===========================================================================
243  /// Return the number of edge basis functions for q
244 //===========================================================================
245  template<>
247  {
248  return 6;
249  }
250 
251 
252 //===========================================================================
253  /// Return the number of internal basis functions for q
254 //===========================================================================
255  template<>
257  {
258  return 2;
259  }
260 
261 //===========================================================================
262  /// Compute the local form of the q basis at local coordinate s
263 //===========================================================================
264  template<>
266  Shape &q_basis) const
267  {
268  // RT_1 basis functions
269 
272  double g1=g1_vect[0];
273  double g2=g2_vect[0];
274  q_basis(0,0) = Sign_edge[0]*std::sqrt(2.0)*s[0]*(s[1]-g2)/(g1-g2);
275  q_basis(0,1) = Sign_edge[0]*std::sqrt(2.0)*s[1]*(s[1]-g2)/(g1-g2);
276 
277  q_basis(1,0) = Sign_edge[0]*std::sqrt(2.0)*s[0]*(s[1]-g1)/(g2-g1);
278  q_basis(1,1) = Sign_edge[0]*std::sqrt(2.0)*s[1]*(s[1]-g1)/(g2-g1);
279 
280 
281  g1_vect=edge_flux_interpolation_point(1,0);
282  g2_vect=edge_flux_interpolation_point(1,1);
283  g1=g1_vect[0];
284  g2=g2_vect[0];
285  q_basis(2,0) = Sign_edge[1]*(s[0]-1.0)*(s[1]-g1)/(g2-g1);
286  q_basis(2,1) = Sign_edge[1]*s[1]*(s[1]-g1)/(g2-g1);
287 
288  q_basis(3,0) = Sign_edge[1]*(s[0]-1.0)*(s[1]-g2)/(g1-g2);
289  q_basis(3,1) = Sign_edge[1]*s[1]*(s[1]-g2)/(g1-g2);
290 
291 
292  g1_vect=edge_flux_interpolation_point(2,0);
293  g2_vect=edge_flux_interpolation_point(2,1);
294  g1=g1_vect[0];
295  g2=g2_vect[0];
296  q_basis(4,0) = Sign_edge[2]*s[0]*(s[0]-g2)/(g1-g2);
297  q_basis(4,1) = Sign_edge[2]*(s[1]-1.0)*(s[0]-g2)/(g1-g2);
298 
299  q_basis(5,0) = Sign_edge[2]*s[0]*(s[0]-g1)/(g2-g1);
300  q_basis(5,1) = Sign_edge[2]*(s[1]-1.0)*(s[0]-g1)/(g2-g1);
301 
302  q_basis(6,0) = s[1]*s[0];
303  q_basis(6,1) = s[1]*(s[1]-1.0);
304 
305  q_basis(7,0) = s[0]*(s[0]-1.0);
306  q_basis(7,1) = s[0]*s[1];
307  }
308 
309 
310 //===========================================================================
311  /// Compute the local form of the q basis and dbasis/ds at local coordinate s
312 //===========================================================================
313  template<>
315  const Vector<double> &s,
316  Shape &div_q_basis_ds) const
317  {
318 
321  double g1=g1_vect[0];
322  double g2=g2_vect[0];
323  div_q_basis_ds(0) = Sign_edge[0]*std::sqrt(2.0)*(3.0*s[1]-2.0*g2)/(g1-g2);
324  div_q_basis_ds(1) = Sign_edge[0]*std::sqrt(2.0)*(2.0*g1-3.0*s[1])/(g1-g2);
325 
326 
327  g1_vect=edge_flux_interpolation_point(1,0);
328  g2_vect=edge_flux_interpolation_point(1,1);
329  g1=g1_vect[0];
330  g2=g2_vect[0];
331  div_q_basis_ds(2) = Sign_edge[1]*(2.0*g1-3.0*s[1])/(g1-g2);
332  div_q_basis_ds(3) = Sign_edge[1]*(3.0*s[1]-2.0*g2)/(g1-g2);
333 
334 
335 
336  g1_vect=edge_flux_interpolation_point(2,0);
337  g2_vect=edge_flux_interpolation_point(2,1);
338  g1=g1_vect[0];
339  g2=g2_vect[0];
340  div_q_basis_ds(4) = Sign_edge[2]*(3.0*s[0]-2.0*g2)/(g1-g2);
341  div_q_basis_ds(5) = Sign_edge[2]*(2.0*g1-3.0*s[0])/(g1-g2);
342 
343  div_q_basis_ds(6) = 3.0*s[1]-1.0;
344  div_q_basis_ds(7) = 3.0*s[0]-1.0;
345 
346  // Scale the basis by the ratio of the length of the edge to the length of
347  // the corresponding edge on the reference element
348  scale_basis(div_q_basis_ds);
349  }
350 
351 //===========================================================================
352 /// Return the number of flux_interpolation points along each edge of
353 /// the element
354 //===========================================================================
355  template<>
357  {
358  return 2;
359  }
360 
361 //===========================================================================
362  /// Return the equation number of the n-th pressure degree of freedom
363 //===========================================================================
364  template<>
365  int TRaviartThomasDarcyElement<1>::p_local_eqn(const unsigned &n) const
366  {
368  }
369 
370 //===========================================================================
371  /// Return the nth pressure value
372 //===========================================================================
373  template<>
374  double TRaviartThomasDarcyElement<1>::p_value(const unsigned &n) const
375  {
376  return this->internal_data_pt(P_internal_data_index)->value(n);
377  }
378 
379 //===========================================================================
380  /// Return the total number of pressure basis functions
381 //===========================================================================
382  template<>
384  {
385  return 3;
386  }
387 
388 //===========================================================================
389  /// Compute the pressure basis
390 //===========================================================================
391  template<>
393  Shape &p_basis) const
394  {
395  p_basis(0) = 1.0;
396  p_basis(1) = s[0];
397  p_basis(2) = s[1];
398  }
399 
400 
401 //===========================================================================
402 /// Recovery order for Z2 error estimator
403 //===========================================================================
404  template<>
406  {
407  return 2; // Need to experiment with this...
408  }
409 
410 //===========================================================================
411  /// The number of values stored at each node. Two flux values are stored
412  /// at each midside node.
413 //===========================================================================
414  template<>
416  {0,0,0,2,2,2};
417 
418 
419 //===========================================================================
420  /// Face index associated with edge flux degree of freedom
421 //===========================================================================
422  template<>
424  {2,2,0,0,1,1};
425 
426 //===========================================================================
427  /// Conversion scheme from an edge degree of freedom to the node it's stored at
428 /// Fluxes are stored at midside nodes
429 //===========================================================================
430  template<>
431  const unsigned TRaviartThomasDarcyElement<1>::Q_edge_conv[3] = {3,4,5};
432 
433 
434 //===========================================================================
435  /// The points along each edge where the fluxes are taken to be
436 //===========================================================================
437  template<>
439  {0.5-std::sqrt(3.0)/6.0,
440  0.5+std::sqrt(3.0)/6.0};
441 
442 
443 
444 //===========================================================================
445  // Force building of templates
446 //===========================================================================
447  template class TRaviartThomasDarcyElement<0>;
448  template class TRaviartThomasDarcyElement<1>;
449 
450 } // End of oomph namespace
451 
double value(const unsigned &i) const
Return i-th stored value. This function is not virtual so that it can be inlined. This means that if ...
Definition: nodes.h:291
Class implementing the generic maths of the Darcy equations using Raviart-Thomas elements with both e...
void scale_basis(Shape &basis) const
Scale the edge basis to allow arbitrary edge mappings.
unsigned nrecovery_order()
Recovery order for Z2 error estimator.
unsigned np_basis() const
Return the total number of pressure basis functions.
TRaviartThomasDarcyElement()
Constructor.
unsigned Q_internal_data_index
The internal data index where the internal q degrees of freedom are stored.
unsigned nedge_flux_interpolation_point() const
Return the number of flux interpolation points along each edge of the element.
unsigned nq_basis_internal() const
Return the number of internal basis functions for q.
static char t char * s
Definition: cfortran.h:572
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
Definition: elements.h:623
Vector< double > edge_flux_interpolation_point(const unsigned &edge, const unsigned &n) const
A class that represents a collection of data; each Data object may contain many different individual ...
Definition: nodes.h:89
unsigned P_internal_data_index
The internal data index where the p degrees of freedom are stored.
std::vector< short > Sign_edge
Unit normal signs associated with each edge to ensure inter-element continuity of the flux...
int p_local_eqn(const unsigned &n) const
Return the equation number of the n-th pressure degree of freedom.
double p_value(const unsigned &n) const
Return the nth pressure value.
unsigned nq_basis_edge() const
Return the number of edge basis functions for q.
void get_p_basis(const Vector< double > &s, Shape &p_basis) const
Compute the pressure basis.
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
void get_div_q_basis_local(const Vector< double > &s, Shape &div_q_basis_ds) const
Return 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
Return the local form of the q basis at local coordinate s.
int internal_local_eqn(const unsigned &i, const unsigned &j) const
Return the local equation number corresponding to the j-th value stored at the i-th internal data...
Definition: elements.h:268