Tdarcy_elements.h
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-2025 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//====================================================================
26#ifndef OOMPH_TRAVIART_THOMAS_DARCY_HEADER
27#define OOMPH_TRAVIART_THOMAS_DARCY_HEADER
28
29// Config header
30#ifdef HAVE_CONFIG_H
31#include <oomph-lib-config.h>
32#endif
33
34#include "darcy_elements.h"
35#include "generic/Telements.h"
36
37namespace oomph
38{
39 //================================================================
40 /// Element which solves the Darcy equations using TElements.
41 /// Geometrically the element is always a six noded triangle.
42 /// We use the mid-side nodes to store edge-based flux degrees of
43 /// freedom and internal data for the discontinuous pressure
44 /// and internal flux dofs.
45 //================================================================
46 template<unsigned ORDER>
48 public DarcyEquations<2>
49 {
50 private:
51 /// Face index associated with edge flux degree of freedom
52 static const unsigned Face_index_of_edge_flux[];
53
54 /// Conversion scheme from an edge degree of freedom to the node it's
55 /// stored at
56 static const unsigned Q_edge_conv[];
57
58 /// The points along each edge where the fluxes are interpolated
59 static const double Flux_interpolation_point[];
60
61 /// The internal data index where the internal q degrees of freedom are
62 /// stored
64
65 /// The internal data index where the p degrees of freedom are stored
67
68 /// Unit normal signs associated with each edge to ensure
69 /// inter-element continuity of the flux
70 std::vector<short> Sign_edge;
71
72 public:
73 /// Constructor
75
76 /// Destructor
78
79 /// Number of values required at node n
80 unsigned required_nvalue(const unsigned& n) const
81 {
82 return Initial_Nvalue[n];
83 }
84
85 /// Return the face index associated with specified edge
86 unsigned face_index_of_edge(const unsigned& j) const
87 {
88 return (j + 2) % 3;
89 }
90
91 /// Compute the face element coordinates of the nth flux
92 /// interpolation point along specified edge
94 const unsigned& edge, const unsigned& n, Vector<double>& s) const
95 {
96 // Get the location of the n-th flux interpolation point along
97 // the edge in terms of the distance along the edge itself
100
101 // Convert the edge number to the number of the mid-edge node along that
102 // edge
103 unsigned node_number = Q_edge_conv[edge];
104
105 // The edge basis functions are defined in a clockwise manner, so we have
106 // to effectively "flip" some coordinates
107 switch (node_number)
108 {
109 case 3:
111 break;
112 case 4:
113 s[0] = 1.0 - flux_interpolation_point[0];
114 break;
115 case 5:
117 break;
118 }
119 }
120
121 /// Return the face index associated with edge flux degree of freedom
122 unsigned face_index_of_q_edge_basis_fct(const unsigned& j) const
123 {
125 }
126
127 /// Return the equation nunmber of the n-th edge (flux) degree of freedom
128 int q_edge_local_eqn(const unsigned& n) const
129 {
131 }
132
133 /// Return the equation number of the n-th internal degree of freedom
134 int q_internal_local_eqn(const unsigned& n) const
135 {
137 }
138
139 /// Return vector of pointers to the Data objects that store the
140 /// edge flux values
142 {
143 // It's the mid-side nodes:
145 data_pt[0] = node_pt(3);
146 data_pt[1] = node_pt(4);
147 data_pt[2] = node_pt(5);
148 return data_pt;
149 }
150
151 /// Return pointer to the Data object that stores the internal flux values
153 {
154 return this->internal_data_pt(Q_internal_data_index);
155 }
156
157 /// Return the index of the internal data where the q_internal
158 /// degrees of freedom are stored
159 unsigned q_internal_index() const
160 {
162 }
163
164 /// Return the nodal index at which the nth edge unknown is stored
165 unsigned q_edge_index(const unsigned& n) const
166 {
167 return n % (ORDER + 1);
168 }
169
170 /// Return the local node number of the node where the nth edge
171 /// unknown is stored
172 unsigned q_edge_node_number(const unsigned& n) const
173 {
174 return Q_edge_conv[n / (ORDER + 1)];
175 }
176
177 /// Get pointer to node that stores the edge flux dofs for specified edge
178 Node* edge_flux_node_pt(const unsigned& edge)
179 {
180 return node_pt(Q_edge_conv[edge]);
181 }
182
183 /// Return the values of the edge (flux) degree of freedom
184 double q_edge(const unsigned& n) const
185 {
187 }
188
189 /// Return the values of the internal degree of freedom
190 double q_internal(const unsigned& n) const
191 {
192 return this->internal_data_pt(q_internal_index())->value(n);
193 }
194
195 /// Set the values of the edge (flux) degree of freedom
196 void set_q_edge(const unsigned& n, const double& value)
197 {
199 }
200
201 /// Set the values of the internal degree of freedom
202 void set_q_internal(const unsigned& n, const double& value)
203 {
204 this->internal_data_pt(q_internal_index())->set_value(n, value);
205 }
206
207 /// Return the number of edge basis functions for q
208 unsigned nq_basis_edge() const;
209
210 /// Return the number of internal basis functions for q
211 unsigned nq_basis_internal() const;
212
213 /// Return the local form of the q basis at local coordinate s
215
216 /// Return the local form of the q basis and dbasis/ds at local coordinate s
218 Shape& div_q_basis_ds) const;
219
220 /// Return the number of flux interpolation points along each
221 /// edge of the element
223
224 /// Return the local coordinate of the nth flux interpolation
225 /// point along an edge
227 const unsigned& n) const
228 {
230 coord[0] = (1.0 - sign_edge(edge)) / 2.0 +
232 return coord;
233 }
234
235
236 /// Compute the global coordinates of the flux interpolation
237 /// point associated with the j-th edge-based q basis fct
239 Vector<double>& x) const
240 {
241 unsigned n = j % nedge_flux_interpolation_point();
242 unsigned edge = (j - n) / nedge_flux_interpolation_point();
244 }
245
246
247 /// Compute the global coordinates of the nth flux interpolation
248 /// point along an edge
250 const unsigned& n,
251 Vector<double>& x) const
252 {
253 // Get the location of the n-th flux interpolation point along
254 // the edge in terms of the distance along the edge itself
257
258 // Convert the edge number to the number of the mid-edge node along that
259 // edge
260 unsigned node_number = Q_edge_conv[edge];
261
262 // Storage for the local coords of the flux_interpolation point
264
265 // The edge basis functions are defined in a clockwise manner, so we have
266 // to effectively "flip" the coordinates along edges 0 and 1 to match this
267 switch (node_number)
268 {
269 case 3:
272 break;
273 case 4:
274 s_flux_interpolation[0] = 0.0;
276 break;
277 case 5:
279 s_flux_interpolation[1] = 0.0;
280 break;
281 }
282
283 // Calculate the global coordinates from the local ones
285 }
286
287 /// Pin the nth internal q value
288 void pin_q_internal_value(const unsigned& n)
289 {
291 }
292
293 /// Return the equation number of the n-th pressure degree of freedom
294 int p_local_eqn(const unsigned& n) const;
295
296 /// Return the nth pressure value
297 double p_value(const unsigned& n) const;
298
299 /// Return the total number of pressure basis functions
300 unsigned np_basis() const;
301
302 /// Compute the pressure basis
304
305 /// Pin the nth pressure value
306 void pin_p_value(const unsigned& n)
307 {
308 this->internal_data_pt(P_internal_data_index)->pin(n);
309 }
310
311 /// Return pointer to the Data object that stores the pressure values
313 {
314 return this->internal_data_pt(P_internal_data_index);
315 }
316
317 /// Set the nth pressure value
318 void set_p_value(const unsigned& n, const double& value)
319 {
320 this->internal_data_pt(P_internal_data_index)->set_value(n, value);
321 }
322
323 /// Scale the edge basis to allow arbitrary edge mappings
324 // hierher explain please
326 {
327 // Storage for the lengths of the edges of the element
328 Vector<double> length(3, 0.0);
329
330 // Temporary storage for the vertex positions
331 double x0, y0, x1, y1;
332
333 // loop over the edges of the element and calculate their lengths (in x-y
334 // space)
335 for (unsigned i = 0; i < 3; i++)
336 {
337 x0 = this->node_pt(i)->x(0);
338 y0 = this->node_pt(i)->x(1);
339 x1 = this->node_pt((i + 1) % 3)->x(0);
340 y1 = this->node_pt((i + 1) % 3)->x(1);
341
342 length[i] = std::sqrt(std::pow(y1 - y0, 2) + std::pow(x1 - x0, 2));
343 }
344
345 // lengths of the sides of the reference element (in the same order as the
346 // basis functions)
347 const double ref_length[3] = {std::sqrt(2.0), 1, 1};
348
349 // get the number of basis functions associated with the edges
350 unsigned n_q_basis_edge = nq_basis_edge();
351
352 // rescale the edge basis functions to allow arbitrary edge mappings from
353 // element to ref. element
354 const unsigned n_index2 = basis.nindex2();
355 for (unsigned i = 0; i < n_index2; i++)
356 {
357 for (unsigned l = 0; l < n_q_basis_edge; l++)
358 {
359 basis(l, i) *=
360 (length[l / (ORDER + 1)] / ref_length[l / (ORDER + 1)]);
361 }
362 }
363 }
364
365 /// Accessor for the unit normal sign of edge n (const version)
366 const short& sign_edge(const unsigned& n) const
367 {
368 return Sign_edge[n];
369 }
370
371 /// Accessor for the unit normal sign of edge n
372 short& sign_edge(const unsigned& n)
373 {
374 return Sign_edge[n];
375 }
376
377 /// Output with default number of plot points
378 void output(std::ostream& outfile)
379 {
381 }
382
383 /// Output FE representation of soln: x,y,u1,u2,div_q,p at
384 /// Nplot^DIM plot points
385 void output(std::ostream& outfile, const unsigned& Nplot)
386 {
388 }
389
390 /// Number of vertex nodes in the element
391 unsigned nvertex_node() const
392 {
394 }
395
396 /// Pointer to the j-th vertex node in the element
397 Node* vertex_node_pt(const unsigned& j) const
398 {
400 }
401
402
403 /// Recovery order for Z2 error estimator
404 unsigned nrecovery_order();
405
406 protected:
407 /// Return the geometric basis, and the q, p and divergence basis
408 /// functions and test functions at local coordinate s
410 Shape& psi,
411 Shape& q_basis,
412 Shape& q_test,
413 Shape& p_basis,
414 Shape& p_test,
416 Shape& div_q_test_ds) const
417 {
418 const unsigned n_q_basis = this->nq_basis();
419
422 this->get_p_basis(s, p_basis);
424
425 double J = this->transform_basis(s, q_basis_local, psi, q_basis);
426
427 q_test = q_basis;
428 p_test = p_basis;
430
431 return J;
432 }
433
434 /// Return the geometric basis, and the q, p and divergence basis
435 /// functions and test functions at integration point ipt
436 double shape_basis_test_local_at_knot(const unsigned& ipt,
437 Shape& psi,
438 Shape& q_basis,
439 Shape& q_test,
440 Shape& p_basis,
441 Shape& p_test,
443 Shape& div_q_test_ds) const
444 {
445 Vector<double> s(2);
446 for (unsigned i = 0; i < 2; i++)
447 {
448 s[i] = this->integral_pt()->knot(ipt, i);
449 }
450
452 psi,
453 q_basis,
454 q_test,
455 p_basis,
456 p_test,
459 }
460
461 /// The number of values stored at each node
462 static const unsigned Initial_Nvalue[];
463 };
464
465
466 //============================================================
467 /// Face geometry for TRaviartThomasDarcyElement<0>
468 //============================================================
469 template<>
471 : public virtual TElement<1, 3>
472 {
473 public:
474 /// Constructor: Call constructor of base
475 FaceGeometry() : TElement<1, 3>() {}
476 };
477
478
479 //============================================================
480 /// Face geometry for TRaviartThomasDarcyElement<1>
481 //============================================================
482 template<>
484 : public virtual TElement<1, 3>
485 {
486 public:
487 /// Constructor: Call constructor of base class
488 FaceGeometry() : TElement<1, 3>() {}
489 };
490
491
492} // namespace oomph
493
494#endif
static char t char * s
Definition cfortran.h:568
cstr elem_len * i
Definition cfortran.h:603
Class implementing the generic maths of the Darcy equations using Raviart-Thomas elements with both e...
double transform_basis(const Vector< double > &s, const Shape &q_basis_local, Shape &psi, Shape &q_basis) const
Performs a div-conserving transformation of the vector basis functions from the reference element to ...
unsigned nq_basis() const
Return the total number of computational basis functions for q.
void output(std::ostream &outfile)
Output with default number of plot points.
A class that represents a collection of data; each Data object may contain many different individual ...
Definition nodes.h:86
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition nodes.h:385
void set_value(const unsigned &i, const double &value_)
Set the i-th stored data value to specified value. The only reason that we require an explicit set fu...
Definition nodes.h:271
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:293
FaceGeometry()
Constructor: Call constructor of base.
FaceGeometry()
Constructor: Call constructor of base class.
FaceGeometry class definition: This policy class is used to allow construction of face elements that ...
Definition elements.h:5002
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition elements.h:1967
double nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n. Produces suitably interpolated values for hanging nodes...
Definition elements.h:2597
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
Definition elements.cc:3992
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
Definition elements.h:1436
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition elements.h:2179
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
Definition elements.h:605
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:267
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
Definition nodes.h:906
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition nodes.h:1060
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition shape.h:76
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
General TElement class.
Definition Telements.h:1208
Element which solves the Darcy equations using TElements. Geometrically the element is always a six n...
void edge_flux_interpolation_point_global(const unsigned &j, Vector< double > &x) const
Compute the global coordinates of the flux interpolation point associated with the j-th edge-based q ...
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.
unsigned face_index_of_edge(const unsigned &j) const
Return the face index associated with specified edge.
const short & sign_edge(const unsigned &n) const
Accessor for the unit normal sign of edge n (const version)
void face_local_coordinate_of_flux_interpolation_point(const unsigned &edge, const unsigned &n, Vector< double > &s) const
Compute the face element coordinates of the nth flux interpolation point along specified edge.
void pin_q_internal_value(const unsigned &n)
Pin the nth internal q value.
int p_local_eqn(const unsigned &n) const
Return the equation number of the n-th pressure degree of freedom.
unsigned required_nvalue(const unsigned &n) const
Number of values required at node n.
unsigned nvertex_node() const
Number of vertex nodes in the element.
void edge_flux_interpolation_point_global(const unsigned &edge, const unsigned &n, Vector< double > &x) const
Compute the global coordinates of the nth flux interpolation point along an edge.
Vector< Data * > q_edge_data_pt() const
Return vector of pointers to the Data objects that store the edge flux values.
double p_value(const unsigned &n) const
Return the nth pressure value.
unsigned q_edge_index(const unsigned &n) const
Return the nodal index at which the nth edge unknown is stored.
Data * q_internal_data_pt() const
Return pointer to the Data object that stores the internal flux values.
void output(std::ostream &outfile)
Output with default number of plot points.
double q_internal(const unsigned &n) const
Return the values of the internal degree of freedom.
void set_q_internal(const unsigned &n, const double &value)
Set the values of the internal degree of freedom.
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.
static const unsigned Initial_Nvalue[]
The number of values stored at each node.
unsigned np_basis() const
Return the total number of pressure basis functions.
static const unsigned Q_edge_conv[]
Conversion scheme from an edge degree of freedom to the node it's stored at.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
unsigned face_index_of_q_edge_basis_fct(const unsigned &j) const
Return the face index associated with edge flux degree of freedom.
unsigned q_edge_node_number(const unsigned &n) const
Return the local node number of the node where the nth edge unknown is stored.
double q_edge(const unsigned &n) const
Return the values of the edge (flux) degree of freedom.
static const unsigned Face_index_of_edge_flux[]
Face index associated with edge flux degree of freedom.
short & sign_edge(const unsigned &n)
Accessor for the unit normal sign of edge n.
std::vector< short > Sign_edge
Unit normal signs associated with each edge to ensure inter-element continuity of the flux.
void set_p_value(const unsigned &n, const double &value)
Set the nth pressure value.
void output(std::ostream &outfile, const unsigned &Nplot)
Output FE representation of soln: x,y,u1,u2,div_q,p at Nplot^DIM plot points.
Data * p_data_pt() const
Return pointer to the Data object that stores the pressure values.
int q_internal_local_eqn(const unsigned &n) const
Return the equation number of the n-th internal degree of freedom.
static const double Flux_interpolation_point[]
The points along each edge where the fluxes are interpolated.
unsigned q_internal_index() const
Return the index of the internal data where the q_internal degrees of freedom are stored.
void scale_basis(Shape &basis) const
Scale the edge basis to allow arbitrary edge mappings.
void get_p_basis(const Vector< double > &s, Shape &p_basis) const
Compute the pressure basis.
unsigned P_internal_data_index
The internal data index where the p degrees of freedom are stored.
double shape_basis_test_local(const Vector< double > &s, Shape &psi, Shape &q_basis, Shape &q_test, Shape &p_basis, Shape &p_test, Shape &div_q_basis_ds, Shape &div_q_test_ds) const
Return the geometric basis, and the q, p and divergence basis functions and test functions at local c...
int q_edge_local_eqn(const unsigned &n) const
Return the equation nunmber of the n-th edge (flux) degree of freedom.
void set_q_edge(const unsigned &n, const double &value)
Set the values of the edge (flux) degree of freedom.
unsigned nrecovery_order()
Recovery order for Z2 error estimator.
Node * edge_flux_node_pt(const unsigned &edge)
Get pointer to node that stores the edge flux dofs for specified edge.
void pin_p_value(const unsigned &n)
Pin the nth pressure value.
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.
unsigned nq_basis_edge() const
Return the number of edge basis functions for q.
double shape_basis_test_local_at_knot(const unsigned &ipt, Shape &psi, Shape &q_basis, Shape &q_test, Shape &p_basis, Shape &p_test, Shape &div_q_basis_ds, Shape &div_q_test_ds) const
Return the geometric basis, and the q, p and divergence basis functions and test functions at integra...
Vector< double > edge_flux_interpolation_point(const unsigned &edge, const unsigned &n) const
Return the local coordinate of the nth flux interpolation point along an edge.
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).