Taxisym_poroelasticity_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-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//====================================================================
26 #ifndef OOMPH_TAXISYM_POROELASTICITY_ELEMENTS_HEADER
27 #define OOMPH_TAXISYM_POROELASTICITY_ELEMENTS_HEADER
28 
29 // Config header generated by autoconfig
30 #ifdef HAVE_CONFIG_H
31 #include <oomph-lib-config.h>
32 #endif
33 
35 #include "../generic/Telements.h"
36 
37 namespace oomph
38 {
39  /// =================================================================
40  /// Element which solves the Darcy/linear elasticity equations
41  /// using TElements
42  /// Geometrically the element is always a six noded triangle.
43  /// We use the mid-side nodes to store edge-based flux degrees of
44  /// freedom and internal data for the discontinuous pressure
45  /// and internal flux dofs.
46  /// =================================================================
47  template<unsigned ORDER>
49  : public TElement<2, 3>,
51  {
52  private:
53  /// The number of values stored at each node
54  static const unsigned Initial_Nvalue[];
55 
56  /// Face index associated with edge flux degree of freedom
57  static const unsigned Face_index_of_edge_flux[];
58 
59  /// Conversion scheme from an edge degree of freedom to the node
60  /// it's stored at
61  static const unsigned Q_edge_conv[];
62 
63  /// The points along each edge where the fluxes are taken to be
64  static const double Flux_interpolation_point[];
65 
66  /// The internal data index where the internal q degrees of freedom are
67  /// stored
69 
70  /// The internal data index where the p degrees of freedom are stored
72 
73  /// Unit normal signs associated with each edge to ensure
74  /// inter-element continuity of the flux
75  std::vector<short> Sign_edge;
76 
77  public:
78  /// Constructor
80 
81  /// Destructor
83 
84  /// Number of values required at node n
85  unsigned required_nvalue(const unsigned& n) const
86  {
87  return Initial_Nvalue[n];
88  }
89 
90  /// Return the face index associated with specified edge
91  unsigned face_index_of_edge(const unsigned& j) const
92  {
93  return (j + 2) % 3;
94  }
95 
96  /// Compute the face element coordinates of the nth flux
97  /// interpolation point along specified edge
99  const unsigned& edge, const unsigned& n, Vector<double>& s) const
100  {
101  // Get the location of the n-th flux interpolation point along
102  // the edge in terms of the distance along the edge itself
103  Vector<double> flux_interpolation_point =
105 
106  // Convert the edge number to the number of the mid-edge node along that
107  // edge
108  unsigned node_number = Q_edge_conv[edge];
109 
110  // The edge basis functions are defined in a clockwise manner, so we have
111  // to effectively "flip" some coordinates
112  switch (node_number)
113  {
114  case 3:
115  s[0] = flux_interpolation_point[0];
116  break;
117  case 4:
118  s[0] = 1.0 - flux_interpolation_point[0];
119  break;
120  case 5:
121  s[0] = flux_interpolation_point[0];
122  break;
123  }
124  }
125 
126  /// Return the face index associated with j-th edge flux degree of freedom
127  unsigned face_index_of_q_edge_basis_fct(const unsigned& j) const
128  {
129  return Face_index_of_edge_flux[j];
130  }
131 
132  /// Return the nodal index of the j-th solid displacement unknown
133  /// [0: r; 1: z]
134  unsigned u_index_axisym_poroelasticity(const unsigned& j) const
135  {
136 #ifdef RANGE_CHECKING
137  if (j >= 2)
138  {
139  std::ostringstream error_message;
140  error_message << "Range Error: j " << j << " is not in the range (0,1)";
141  throw OomphLibError(error_message.str(),
142  OOMPH_CURRENT_FUNCTION,
143  OOMPH_EXCEPTION_LOCATION);
144  }
145 #endif
146  return j;
147  }
148 
149 
150  /// Return the equation number of the j-th edge (flux) degree of freedom
151  int q_edge_local_eqn(const unsigned& j) const
152  {
153 #ifdef RANGE_CHECKING
154  if (j >= nq_basis_edge())
155  {
156  std::ostringstream error_message;
157  error_message << "Range Error: j " << j << " is not in the range (0,"
158  << nq_basis_edge() - 1 << ")";
159  throw OomphLibError(error_message.str(),
160  OOMPH_CURRENT_FUNCTION,
161  OOMPH_EXCEPTION_LOCATION);
162  }
163 #endif
164  return this->nodal_local_eqn(q_edge_node_number(j), q_edge_index(j));
165  }
166 
167  /// Return the equation number of the j-th internal degree of freedom
168  int q_internal_local_eqn(const unsigned& j) const
169  {
170 #ifdef RANGE_CHECKING
171  if (j >= nq_basis_internal())
172  {
173  std::ostringstream error_message;
174  error_message << "Range Error: j " << j << " is not in the range (0,"
175  << nq_basis_internal() - 1 << ")";
176  throw OomphLibError(error_message.str(),
177  OOMPH_CURRENT_FUNCTION,
178  OOMPH_EXCEPTION_LOCATION);
179  }
180 #endif
182  }
183 
184  /// Return vector of pointers to the Data objects that store the
185  /// edge flux values
187  {
188  // It's the mid-side nodes:
189  Vector<Data*> data_pt(3);
190  data_pt[0] = node_pt(3);
191  data_pt[1] = node_pt(4);
192  data_pt[2] = node_pt(5);
193  return data_pt;
194  }
195 
196  /// Return pointer to the Data object that stores the internal flux values
198  {
199  return this->internal_data_pt(Q_internal_data_index);
200  }
201 
202  /// Return the index of the internal data where the q_internal
203  /// degrees of freedom are stored
204  unsigned q_internal_index() const
205  {
206  return Q_internal_data_index;
207  }
208 
209  /// Return the nodal index at which the jth edge unknown is stored
210  unsigned q_edge_index(const unsigned& j) const
211  {
212 #ifdef RANGE_CHECKING
213  if (j >= (nq_basis_edge()))
214  {
215  std::ostringstream error_message;
216  error_message << "Range Error: j " << j << " is not in the range (0,"
217  << nq_basis_edge() - 1 << ")";
218  throw OomphLibError(error_message.str(),
219  OOMPH_CURRENT_FUNCTION,
220  OOMPH_EXCEPTION_LOCATION);
221  }
222 #endif
223  return j % (ORDER + 1) + 2;
224  }
225 
226  /// Return the number of the node where the jth edge unknown is stored
227  unsigned q_edge_node_number(const unsigned& j) const
228  {
229 #ifdef RANGE_CHECKING
230  if (j >= (nq_basis_edge()))
231  {
232  std::ostringstream error_message;
233  error_message << "Range Error: j " << j << " is not in the range (0,"
234  << nq_basis_edge() - 1 << ")";
235  throw OomphLibError(error_message.str(),
236  OOMPH_CURRENT_FUNCTION,
237  OOMPH_EXCEPTION_LOCATION);
238  }
239 #endif
240  return Q_edge_conv[j / (ORDER + 1)];
241  }
242 
243 
244  /// Get pointer to node that stores the edge flux dofs for specified edge
245  Node* edge_flux_node_pt(const unsigned& edge)
246  {
247  return node_pt(Q_edge_conv[edge]);
248  }
249 
250  /// Return the values of the j-th edge (flux) degree of freedom
251  double q_edge(const unsigned& j) const
252  {
253 #ifdef RANGE_CHECKING
254  if (j >= (nq_basis_edge()))
255  {
256  std::ostringstream error_message;
257  error_message << "Range Error: j " << j << " is not in the range (0,"
258  << nq_basis_edge() - 1 << ")";
259  throw OomphLibError(error_message.str(),
260  OOMPH_CURRENT_FUNCTION,
261  OOMPH_EXCEPTION_LOCATION);
262  }
263 #endif
265  }
266 
267  /// Return the values of the j-th edge (flux) degree of
268  /// freedom at time history level t
269  double q_edge(const unsigned& t, const unsigned& j) const
270  {
271 #ifdef RANGE_CHECKING
272  if (j >= (nq_basis_edge()))
273  {
274  std::ostringstream error_message;
275  error_message << "Range Error: j " << j << " is not in the range (0,"
276  << nq_basis_edge() - 1 << ")";
277  throw OomphLibError(error_message.str(),
278  OOMPH_CURRENT_FUNCTION,
279  OOMPH_EXCEPTION_LOCATION);
280  }
281 #endif
283  }
284 
285  /// Return the values of the internal degree of freedom
286  double q_internal(const unsigned& j) const
287  {
288 #ifdef RANGE_CHECKING
289  if (j >= nq_basis_internal())
290  {
291  std::ostringstream error_message;
292  error_message << "Range Error: j " << j << " is not in the range (0,"
293  << nq_basis_internal() - 1 << ")";
294  throw OomphLibError(error_message.str(),
295  OOMPH_CURRENT_FUNCTION,
296  OOMPH_EXCEPTION_LOCATION);
297  }
298 #endif
299  return this->internal_data_pt(q_internal_index())->value(j);
300  }
301 
302  /// Return the value of the j-th internal degree of freedom at
303  /// time history level t
304  double q_internal(const unsigned& t, const unsigned& j) const
305  {
306 #ifdef RANGE_CHECKING
307  if (j >= nq_basis_internal())
308  {
309  std::ostringstream error_message;
310  error_message << "Range Error: j " << j << " is not in the range (0,"
311  << nq_basis_internal() - 1 << ")";
312  throw OomphLibError(error_message.str(),
313  OOMPH_CURRENT_FUNCTION,
314  OOMPH_EXCEPTION_LOCATION);
315  }
316 #endif
317  return this->internal_data_pt(q_internal_index())->value(t, j);
318  }
319 
320  /// Pin the j-th edge (flux) degree of freedom and set it to specified value
321  void pin_q_edge_value(const unsigned& j, const double& value)
322  {
325  }
326 
327  /// Set the values of the j-th edge (flux) degree of freedom
328  void set_q_edge(const unsigned& j, const double& value)
329  {
331  }
332 
333  /// Set the values of the j-th edge (flux) degree of freedom at
334  /// time history level t
335  void set_q_edge(const unsigned& t, const unsigned& j, const double& value)
336  {
338  }
339 
340 
341  /// Set the values of the j-th internal degree of freedom
342  void set_q_internal(const unsigned& j, const double& value)
343  {
344  this->internal_data_pt(q_internal_index())->set_value(j, value);
345  }
346 
347  /// Set the values of the j-th internal degree of freedom at
348  /// time history level t
349  void set_q_internal(const unsigned& t,
350  const unsigned& j,
351  const double& value)
352  {
353  this->internal_data_pt(q_internal_index())->set_value(t, j, value);
354  }
355 
356  /// Return the number of edge basis functions for flux q
357  unsigned nq_basis_edge() const;
358 
359  /// Return the number of internal basis functions for flux q
360  unsigned nq_basis_internal() const;
361 
362  /// Returns the local form of the q basis at local coordinate s
363  void get_q_basis_local(const Vector<double>& s, Shape& q_basis) const;
364 
365  /// Returns the local form of the q basis and dbasis/ds at local coordinate
366  /// s
368  Shape& div_q_basis_ds) const;
369 
370  /// Returns the number of flux_interpolation points along each
371  /// edge of the element
373 
374  /// Returns the local coordinate of the jth flux_interpolation point
375  /// along specified edge
377  const unsigned& j) const
378  {
379 #ifdef RANGE_CHECKING
380  if (edge >= 3)
381  {
382  std::ostringstream error_message;
383  error_message << "Range Error: edge " << edge
384  << " is not in the range (0,2)";
385  throw OomphLibError(error_message.str(),
386  OOMPH_CURRENT_FUNCTION,
387  OOMPH_EXCEPTION_LOCATION);
388  }
390  {
391  std::ostringstream error_message;
392  error_message << "Range Error: j " << j << " is not in the range (0,"
393  << nedge_flux_interpolation_point() - 1 << ")";
394  throw OomphLibError(error_message.str(),
395  OOMPH_CURRENT_FUNCTION,
396  OOMPH_EXCEPTION_LOCATION);
397  }
398 #endif
399  Vector<double> coord(1);
400  coord[0] = (1.0 - sign_edge(edge)) / 2.0 +
402  return coord;
403  }
404 
405  /// Compute the global coordinates of the jth flux_interpolation
406  /// point along specified edge
407  void edge_flux_interpolation_point_global(const unsigned& edge,
408  const unsigned& j,
409  Vector<double>& x) const
410  {
411 #ifdef RANGE_CHECKING
412  if (edge >= 3)
413  {
414  std::ostringstream error_message;
415  error_message << "Range Error: edge " << edge
416  << " is not in the range (0,2)";
417  throw OomphLibError(error_message.str(),
418  OOMPH_CURRENT_FUNCTION,
419  OOMPH_EXCEPTION_LOCATION);
420  }
422  {
423  std::ostringstream error_message;
424  error_message << "Range Error: j " << j << " is not in the range (0,"
425  << nedge_flux_interpolation_point() - 1 << ")";
426  throw OomphLibError(error_message.str(),
427  OOMPH_CURRENT_FUNCTION,
428  OOMPH_EXCEPTION_LOCATION);
429  }
430 #endif
431 
432  // Get the location of the n-th flux_interpolation point along the
433  // edge in terms of the distance along the edge itself
434  Vector<double> flux_interpolation_point =
436 
437  // Convert the edge number to the number of the mid-edge node along that
438  // edge
439  unsigned node_number = Q_edge_conv[edge];
440 
441  // Storage for the local coords of the flux_interpolation point
442  Vector<double> s_flux_interpolation(2, 0);
443 
444  // The edge basis functions are defined in a clockwise manner, so we have
445  // to effectively "flip" the coordinates along edges 0 and 1 to match this
446  switch (node_number)
447  {
448  case 3:
449  s_flux_interpolation[0] = 1.0 - flux_interpolation_point[0];
450  s_flux_interpolation[1] = flux_interpolation_point[0];
451  break;
452  case 4:
453  s_flux_interpolation[0] = 0.0;
454  s_flux_interpolation[1] = 1.0 - flux_interpolation_point[0];
455  break;
456  case 5:
457  s_flux_interpolation[0] = flux_interpolation_point[0];
458  s_flux_interpolation[1] = 0.0;
459  break;
460  }
461 
462  // Calculate the global coordinates from the local ones
463  interpolated_x(s_flux_interpolation, x);
464  }
465 
466  /// Pin the jth internal q value and set it to specified value
467  void pin_q_internal_value(const unsigned& j, const double& q)
468  {
469 #ifdef RANGE_CHECKING
470  if (j >= nq_basis_internal())
471  {
472  std::ostringstream error_message;
473  error_message << "Range Error: j " << j << " is not in the range (0,"
474  << nq_basis_internal() - 1 << ")";
475  throw OomphLibError(error_message.str(),
476  OOMPH_CURRENT_FUNCTION,
477  OOMPH_EXCEPTION_LOCATION);
478  }
479 #endif
482  }
483 
484  /// Return the equation number of the j-th pressure degree of freedom
485  int p_local_eqn(const unsigned& j) const
486  {
487 #ifdef RANGE_CHECKING
488  if (j >= np_basis())
489  {
490  std::ostringstream error_message;
491  error_message << "Range Error: j " << j << " is not in the range (0,"
492  << np_basis() - 1 << ")";
493  throw OomphLibError(error_message.str(),
494  OOMPH_CURRENT_FUNCTION,
495  OOMPH_EXCEPTION_LOCATION);
496  }
497 #endif
498  return this->internal_local_eqn(P_internal_data_index, j);
499  }
500 
501  /// Return the jth pressure value
502  double p_value(const unsigned& j) const
503  {
504 #ifdef RANGE_CHECKING
505  if (j >= np_basis())
506  {
507  std::ostringstream error_message;
508  error_message << "Range Error: j " << j << " is not in the range (0,"
509  << np_basis() - 1 << ")";
510  throw OomphLibError(error_message.str(),
511  OOMPH_CURRENT_FUNCTION,
512  OOMPH_EXCEPTION_LOCATION);
513  }
514 #endif
515  return this->internal_data_pt(P_internal_data_index)->value(j);
516  }
517 
518  /// Return the total number of pressure basis functions
519  unsigned np_basis() const;
520 
521  /// Return the pressure basis
522  void get_p_basis(const Vector<double>& s, Shape& p_basis) const;
523 
524  /// Pin the jth pressure value and set to specified value
525  void pin_p_value(const unsigned& j, const double& p)
526  {
527 #ifdef RANGE_CHECKING
528  if (j >= np_basis())
529  {
530  std::ostringstream error_message;
531  error_message << "Range Error: j " << j << " is not in the range (0,"
532  << np_basis() - 1 << ")";
533  throw OomphLibError(error_message.str(),
534  OOMPH_CURRENT_FUNCTION,
535  OOMPH_EXCEPTION_LOCATION);
536  }
537 #endif
538  this->internal_data_pt(P_internal_data_index)->pin(j);
539  this->internal_data_pt(P_internal_data_index)->set_value(j, p);
540  }
541 
542  /// Return pointer to the Data object that stores the pressure values
543  Data* p_data_pt() const
544  {
545  return this->internal_data_pt(P_internal_data_index);
546  }
547 
548  /// Set the jth pressure value
549  void set_p_value(const unsigned& j, const double& value)
550  {
551  this->internal_data_pt(P_internal_data_index)->set_value(j, value);
552  }
553 
554  /// Scale the edge basis to allow arbitrary edge mappings
555  void scale_basis(Shape& basis) const
556  {
557  // Storage for the lengths of the edges of the element
558  Vector<double> length(3, 0.0);
559 
560  // Temporary storage for the vertex positions
561  double x0, y0, x1, y1;
562 
563  // loop over the edges of the element and calculate their lengths (in x-y
564  // space)
565  for (unsigned i = 0; i < 3; i++)
566  {
567  x0 = this->node_pt(i)->x(0);
568  y0 = this->node_pt(i)->x(1);
569  x1 = this->node_pt((i + 1) % 3)->x(0);
570  y1 = this->node_pt((i + 1) % 3)->x(1);
571 
572  length[i] = std::sqrt(std::pow(y1 - y0, 2) + std::pow(x1 - x0, 2));
573  }
574 
575  // lengths of the sides of the reference element (in the same order as the
576  // basis functions)
577  const double ref_length[3] = {std::sqrt(2.0), 1, 1};
578 
579  // get the number of basis functions associated with the edges
580  unsigned n_q_basis_edge = nq_basis_edge();
581 
582  // rescale the edge basis functions to allow arbitrary edge mappings from
583  // element to ref. element
584  const unsigned n_index2 = basis.nindex2();
585  for (unsigned i = 0; i < n_index2; i++)
586  {
587  for (unsigned l = 0; l < n_q_basis_edge; l++)
588  {
589  basis(l, i) *=
590  (length[l / (ORDER + 1)] / ref_length[l / (ORDER + 1)]);
591  }
592  }
593  }
594 
595  /// Accessor for the unit normal sign of edge n (const version)
596  const short& sign_edge(const unsigned& n) const
597  {
598  return Sign_edge[n];
599  }
600 
601  /// Accessor for the unit normal sign of edge n
602  short& sign_edge(const unsigned& n)
603  {
604  return Sign_edge[n];
605  }
606 
607  /// Output with default number of plot points
608  void output(std::ostream& outfile)
609  {
611  }
612 
613  /// Output FE representation of soln: x,y,u1,u2,div_q,p at
614  /// Nplot^DIM plot points
615  void output(std::ostream& outfile, const unsigned& Nplot)
616  {
618  }
619 
620 
621  /// Number of vertex nodes in the element
622  unsigned nvertex_node() const
623  {
625  }
626 
627  /// Pointer to the j-th vertex node in the element
628  Node* vertex_node_pt(const unsigned& j) const
629  {
631  }
632 
633  /// Recovery order for Z2 error estimator
634  unsigned nrecovery_order()
635  {
636  return 2; // need to experiment with this...
637  }
638 
639  protected:
640  /// Returns the geometric basis, and the u, p and divergence basis
641  /// functions and test functions at local coordinate s
643  Shape& psi,
644  DShape& dpsi,
645  Shape& u_basis,
646  Shape& u_test,
647  DShape& du_basis_dx,
648  DShape& du_test_dx,
649  Shape& q_basis,
650  Shape& q_test,
651  Shape& p_basis,
652  Shape& p_test,
653  Shape& div_q_basis_ds,
654  Shape& div_q_test_ds) const
655  {
656  const unsigned n_q_basis = this->nq_basis();
657 
658  Shape q_basis_local(n_q_basis, 2);
659  this->get_q_basis_local(s, q_basis_local);
660  this->get_p_basis(s, p_basis);
661  this->get_div_q_basis_local(s, div_q_basis_ds);
662 
663  double J = this->transform_basis(s, q_basis_local, psi, dpsi, q_basis);
664 
665  // u_basis consists of the normal Lagrangian shape functions
666  u_basis = psi;
667  du_basis_dx = dpsi;
668 
669  u_test = psi;
670  du_test_dx = dpsi;
671 
672  q_test = q_basis;
673  p_test = p_basis;
674  div_q_test_ds = div_q_basis_ds;
675 
676  return J;
677  }
678 
679  /// Returns the geometric basis, and the u, p and divergence basis
680  /// functions and test functions at integration point ipt
681  double shape_basis_test_local_at_knot(const unsigned& ipt,
682  Shape& psi,
683  DShape& dpsi,
684  Shape& u_basis,
685  Shape& u_test,
686  DShape& du_basis_dx,
687  DShape& du_test_dx,
688  Shape& q_basis,
689  Shape& q_test,
690  Shape& p_basis,
691  Shape& p_test,
692  Shape& div_q_basis_ds,
693  Shape& div_q_test_ds) const
694  {
695  Vector<double> s(2);
696  for (unsigned i = 0; i < 2; i++)
697  {
698  s[i] = this->integral_pt()->knot(ipt, i);
699  }
700 
701  return shape_basis_test_local(s,
702  psi,
703  dpsi,
704  u_basis,
705  u_test,
706  du_basis_dx,
707  du_test_dx,
708  q_basis,
709  q_test,
710  p_basis,
711  p_test,
712  div_q_basis_ds,
713  div_q_test_ds);
714  }
715  };
716 
717 
718  //================================================================
719  /// Face geometry for TAxisymmetricPoroelasticityElement<0>
720  //================================================================
721  template<>
723  : public virtual TElement<1, 3>
724  {
725  public:
726  /// Constructor: Call constructor of base
727  FaceGeometry() : TElement<1, 3>() {}
728  };
729 
730 
731  //================================================================
732  /// Face geometry for TAxisymmetricPoroelasticityElement<1>
733  //================================================================
734  template<>
736  : public virtual TElement<1, 3>
737  {
738  public:
739  /// Constructor: Call constructor of base class
740  FaceGeometry() : TElement<1, 3>() {}
741  };
742 
743 
744 } // namespace oomph
745 
746 #endif
static char t char * s
Definition: cfortran.h:568
cstr elem_len * i
Definition: cfortran.h:603
char t
Definition: cfortran.h:568
Class implementing the generic maths of the axisym poroelasticity equations: axisym linear elasticity...
double transform_basis(const Vector< double > &s, const Shape &q_basis_local, Shape &psi, DShape &dpsi, Shape &q_basis) const
Performs a div-conserving transformation of the vector basis functions from the reference element to ...
virtual 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 for the derivatives of shape functions The class design is essentially the same as Shape,...
Definition: shape.h:278
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
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
Definition: elements.h:5002
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2179
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
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1967
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
Definition: elements.h:622
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
An OomphLibError object which should be thrown when an run-time error is encountered....
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
unsigned nindex2() const
Return the range of index 2 of the shape function object.
Definition: shape.h:265
================================================================= Element which solves the Darcy/line...
void set_q_internal(const unsigned &t, const unsigned &j, const double &value)
Set the values of the j-th internal degree of freedom at time history level t.
void edge_flux_interpolation_point_global(const unsigned &edge, const unsigned &j, Vector< double > &x) const
Compute the global coordinates of the jth flux_interpolation point along specified edge.
unsigned nvertex_node() const
Number of vertex nodes in the element.
double q_internal(const unsigned &t, const unsigned &j) const
Return the value of the j-th internal degree of freedom at time history level t.
void get_p_basis(const Vector< double > &s, Shape &p_basis) const
Return the pressure basis.
double q_edge(const unsigned &j) const
Return the values of the j-th edge (flux) degree of freedom.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
std::vector< short > Sign_edge
Unit normal signs associated with each edge to ensure inter-element continuity of the flux.
Node * edge_flux_node_pt(const unsigned &edge)
Get pointer to node that stores the edge flux dofs for specified edge.
unsigned Q_internal_data_index
The internal data index where the internal q degrees of freedom are stored.
Data * q_internal_data_pt() const
Return pointer to the Data object that stores the internal flux values.
unsigned q_edge_index(const unsigned &j) const
Return the nodal index at which the jth edge unknown is stored.
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.
short & sign_edge(const unsigned &n)
Accessor for the unit normal sign of edge n.
unsigned P_internal_data_index
The internal data index where the p degrees of freedom are stored.
unsigned required_nvalue(const unsigned &n) const
Number of values required at node n.
double p_value(const unsigned &j) const
Return the jth pressure value.
static const unsigned Initial_Nvalue[]
The number of values stored at each node.
unsigned q_internal_index() const
Return the index of the internal data where the q_internal degrees of freedom are stored.
unsigned nq_basis_edge() const
Return the number of edge basis functions for flux q.
unsigned face_index_of_edge(const unsigned &j) const
Return the face index associated with specified edge.
unsigned q_edge_node_number(const unsigned &j) const
Return the number of the node where the jth edge unknown is stored.
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 set_q_internal(const unsigned &j, const double &value)
Set the values of the j-th internal degree of freedom.
void set_p_value(const unsigned &j, const double &value)
Set the jth pressure value.
Data * p_data_pt() const
Return pointer to the Data object that stores the pressure values.
static const double Flux_interpolation_point[]
The points along each edge where the fluxes are taken to be.
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.
void pin_p_value(const unsigned &j, const double &p)
Pin the jth pressure value and set to specified value.
void pin_q_internal_value(const unsigned &j, const double &q)
Pin the jth internal q value and set it to specified value.
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 set_q_edge(const unsigned &t, const unsigned &j, const double &value)
Set the values of the j-th edge (flux) degree of freedom at time history level t.
const short & sign_edge(const unsigned &n) const
Accessor for the unit normal sign of edge n (const version)
int p_local_eqn(const unsigned &j) const
Return the equation number of the j-th pressure degree of freedom.
double shape_basis_test_local(const Vector< double > &s, Shape &psi, DShape &dpsi, Shape &u_basis, Shape &u_test, DShape &du_basis_dx, DShape &du_test_dx, Shape &q_basis, Shape &q_test, Shape &p_basis, Shape &p_test, Shape &div_q_basis_ds, Shape &div_q_test_ds) const
Returns the geometric basis, and the u, p and divergence basis functions and test functions at local ...
void set_q_edge(const unsigned &j, const double &value)
Set the values of the j-th edge (flux) degree of freedom.
void pin_q_edge_value(const unsigned &j, const double &value)
Pin the j-th edge (flux) degree of freedom and set it to specified value.
void output(std::ostream &outfile)
Output with default number of plot points.
unsigned nedge_flux_interpolation_point() const
Returns the number of flux_interpolation points along each edge of the element.
int q_edge_local_eqn(const unsigned &j) const
Return the equation number of the j-th edge (flux) degree of freedom.
unsigned face_index_of_q_edge_basis_fct(const unsigned &j) const
Return the face index associated with j-th edge flux degree of freedom.
double shape_basis_test_local_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsi, Shape &u_basis, Shape &u_test, DShape &du_basis_dx, DShape &du_test_dx, Shape &q_basis, Shape &q_test, Shape &p_basis, Shape &p_test, Shape &div_q_basis_ds, Shape &div_q_test_ds) const
Returns the geometric basis, and the u, p and divergence basis functions and test functions at integr...
void scale_basis(Shape &basis) const
Scale the edge basis to allow arbitrary edge mappings.
unsigned nq_basis_internal() const
Return the number of internal basis functions for flux q.
unsigned u_index_axisym_poroelasticity(const unsigned &j) const
Return the nodal index of the j-th solid displacement unknown [0: r; 1: z].
int q_internal_local_eqn(const unsigned &j) const
Return the equation number of the j-th internal degree of freedom.
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.
static const unsigned Face_index_of_edge_flux[]
Face index associated with edge flux degree of freedom.
double q_edge(const unsigned &t, const unsigned &j) const
Return the values of the j-th edge (flux) degree of freedom at time history level t.
Vector< Data * > q_edge_data_pt() const
Return vector of pointers to the Data objects that store the edge flux values.
double q_internal(const unsigned &j) const
Return the values of the internal degree of freedom.
unsigned nrecovery_order()
Recovery order for Z2 error estimator.
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.
General TElement class.
Definition: Telements.h:1208
//////////////////////////////////////////////////////////////////// ////////////////////////////////...