fluid_traction_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-2023 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 // Header file for elements that are used to integrate fluid tractions
27 // This includes the guts (i.e. equations) because we want to inline them
28 // for faster operation, although it slows down the compilation!
29 
30 #ifndef OOMPH_FLUID_TRACTION_ELEMENTS_HEADER
31 #define OOMPH_FLUID_TRACTION_ELEMENTS_HEADER
32 
33 // Config header generated by autoconfig
34 #ifdef HAVE_CONFIG_H
35 #include <oomph-lib-config.h>
36 #endif
37 
38 
39 // OOMPH-LIB headers
40 #include "../generic/Qelements.h"
41 #include "../generic/Telements.h"
42 
43 namespace oomph
44 {
45  //======================================================================
46  /// A class for elements that allow the imposition of an applied traction
47  /// to the Navier--Stokes equations
48  /// The geometrical information can be read from the FaceGeometry<ELEMENT>
49  /// class and and thus, we can be generic enough without the need to have
50  /// a separate equations class
51  //======================================================================
52  template<class ELEMENT>
53  class NavierStokesTractionElement : public virtual FaceGeometry<ELEMENT>,
54  public virtual FaceElement
55  {
56  private:
57  /// Pointer to an imposed traction function
58  void (*Traction_fct_pt)(const double& time,
59  const Vector<double>& x,
60  const Vector<double>& n,
61  Vector<double>& result);
62 
63  protected:
64  /// The "global" intrinsic coordinate of the element when
65  /// viewed as part of a geometric object should be given by
66  /// the FaceElement representation, by default
67  double zeta_nodal(const unsigned& n,
68  const unsigned& k,
69  const unsigned& i) const
70  {
71  return FaceElement::zeta_nodal(n, k, i);
72  }
73 
74  /// Access function that returns the local equation numbers
75  /// for velocity components.
76  /// u_local_eqn(n,i) = local equation number or < 0 if pinned.
77  /// The default is to asssume that n is the local node number
78  /// and the i-th velocity component is the i-th unknown stored at the node.
79  virtual inline int u_local_eqn(const unsigned& n, const unsigned& i)
80  {
81  return nodal_local_eqn(n, i);
82  }
83 
84  /// Function to compute the shape and test functions and to return
85  /// the Jacobian of mapping
86  inline double shape_and_test_at_knot(const unsigned& ipt,
87  Shape& psi,
88  Shape& test) const
89  {
90  // Find number of nodes
91  unsigned n_node = nnode();
92  // Calculate the shape functions
93  shape_at_knot(ipt, psi);
94  // Set the test functions to be the same as the shape functions
95  for (unsigned i = 0; i < n_node; i++)
96  {
97  test[i] = psi[i];
98  }
99  // Return the value of the jacobian
100  return J_eulerian_at_knot(ipt);
101  }
102 
103 
104  /// Function to calculate the traction applied to the fluid
105  void get_traction(const double& time,
106  const Vector<double>& x,
107  const Vector<double>& n,
108  Vector<double>& result)
109  {
110  // If the function pointer is zero return zero
111  if (Traction_fct_pt == 0)
112  {
113  // Loop over dimensions and set body forces to zero
114  for (unsigned i = 0; i < Dim; i++)
115  {
116  result[i] = 0.0;
117  }
118  }
119  // Otherwise call the function
120  else
121  {
122  (*Traction_fct_pt)(time, x, n, result);
123  }
124  }
125 
126 
127  /// This function returns the residuals for the
128  /// traction function.
129  /// flag=1(or 0): do (or don't) compute the Jacobian as well.
131  Vector<double>& residuals, DenseMatrix<double>& jacobian, unsigned flag);
132 
133  /// The highest dimension of the problem
134  unsigned Dim;
135 
136  public:
137  /// Constructor, which takes a "bulk" element and the value of the index
138  /// and its limit
140  FiniteElement* const& element_pt,
141  const int& face_index,
142  const bool& called_from_refineable_constructor = false)
143  : FaceGeometry<ELEMENT>(), FaceElement()
144  {
145  // Attach the geometrical information to the element. N.B. This function
146  // also assigns nbulk_value from the required_nvalue of the bulk element
147  element_pt->build_face_element(face_index, this);
148 
149 #ifdef PARANOID
150  {
151  // Check that the element is not a refineable 3d element
152  if (!called_from_refineable_constructor)
153  {
154  // If it's three-d
155  if (element_pt->dim() == 3)
156  {
157  // Is it refineable
158  RefineableElement* ref_el_pt =
159  dynamic_cast<RefineableElement*>(element_pt);
160  if (ref_el_pt != 0)
161  {
162  if (this->has_hanging_nodes())
163  {
164  throw OomphLibError("This flux element will not work correctly "
165  "if nodes are hanging\n",
166  OOMPH_CURRENT_FUNCTION,
167  OOMPH_EXCEPTION_LOCATION);
168  }
169  }
170  }
171  }
172  }
173 #endif
174 
175  // Set the body force function pointer to zero
176  Traction_fct_pt = 0;
177 
178  // Set the dimension from the dimension of the first node
179  Dim = this->node_pt(0)->ndim();
180  }
181 
182  /// Destructor should not delete anything
184 
185  // Access function for the imposed traction pointer
186  void (*&traction_fct_pt())(const double& t,
187  const Vector<double>& x,
188  const Vector<double>& n,
189  Vector<double>& result)
190  {
191  return Traction_fct_pt;
192  }
193 
194  /// This function returns just the residuals
196  {
197  // Call the generic residuals function with flag set to 0
198  // using a dummy matrix argument
200  residuals, GeneralisedElement::Dummy_matrix, 0);
201  }
202 
203  /// This function returns the residuals and the jacobian
205  DenseMatrix<double>& jacobian)
206  {
207  // Call the generic routine with the flag set to 1
209  residuals, jacobian, 1);
210  }
211 
212  /// Overload the output function
213  void output(std::ostream& outfile)
214  {
215  FiniteElement::output(outfile);
216  }
217 
218  /// Output function: x,y,[z],u,v,[w],p in tecplot format
219  void output(std::ostream& outfile, const unsigned& nplot)
220  {
221  FiniteElement::output(outfile, nplot);
222  }
223  };
224 
225 
226  /// ////////////////////////////////////////////////////////////////////
227  /// ////////////////////////////////////////////////////////////////////
228  /// ////////////////////////////////////////////////////////////////////
229 
230 
231  //============================================================================
232  /// Function that returns the residuals for the imposed traction Navier_Stokes
233  /// equations
234  //============================================================================
235  template<class ELEMENT>
238  Vector<double>& residuals, DenseMatrix<double>& jacobian, unsigned flag)
239  {
240  // Find out how many nodes there are
241  unsigned n_node = nnode();
242 
243  // Get continuous time from timestepper of first node
244  double time = node_pt(0)->time_stepper_pt()->time_pt()->time();
245 
246  // Set up memory for the shape and test functions
247  Shape psif(n_node), testf(n_node);
248 
249  // Set the value of n_intpt
250  unsigned n_intpt = integral_pt()->nweight();
251 
252  // Integers to store local equation numbers
253  int local_eqn = 0;
254 
255  // Loop over the integration points
256  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
257  {
258  // Get the integral weight
259  double w = integral_pt()->weight(ipt);
260 
261  // Find the shape and test functions and return the Jacobian
262  // of the mapping
263  double J = shape_and_test_at_knot(ipt, psif, testf);
264 
265  // Premultiply the weights and the Jacobian
266  double W = w * J;
267 
268  // Need to find position to feed into Traction function
269  Vector<double> interpolated_x(Dim);
270 
271  // Initialise to zero
272  for (unsigned i = 0; i < Dim; i++)
273  {
274  interpolated_x[i] = 0.0;
275  }
276 
277  // Calculate velocities and derivatives
278  for (unsigned l = 0; l < n_node; l++)
279  {
280  // Loop over velocity components
281  for (unsigned i = 0; i < Dim; i++)
282  {
283  interpolated_x[i] += nodal_position(l, i) * psif[l];
284  }
285  }
286 
287  // Get the outer unit normal
288  Vector<double> interpolated_n(Dim);
289  outer_unit_normal(ipt, interpolated_n);
290 
291  // Get the user-defined traction terms
292  Vector<double> traction(Dim);
293  get_traction(time, interpolated_x, interpolated_n, traction);
294 
295  // Now add to the appropriate equations
296 
297  // Loop over the test functions
298  for (unsigned l = 0; l < n_node; l++)
299  {
300  // Loop over the velocity components
301  for (unsigned i = 0; i < Dim; i++)
302  {
303  local_eqn = u_local_eqn(l, i);
304  /*IF it's not a boundary condition*/
305  if (local_eqn >= 0)
306  {
307  // Add the user-defined traction terms
308  residuals[local_eqn] += traction[i] * testf[l] * W;
309 
310  // Assuming the the traction DOES NOT depend upon velocities
311  // or pressures, the jacobian is always zero, so no jacobian
312  // terms are required
313  }
314  } // End of loop over dimension
315  } // End of loop over shape functions
316  }
317  }
318 
319 
320  /// //////////////////////////////////////////////////////////////////////
321  /// //////////////////////////////////////////////////////////////////////
322  /// //////////////////////////////////////////////////////////////////////
323 
324 
325  //======================================================================
326  /// A class for elements that allow the imposition of an applied traction
327  /// to the Navier--Stokes equations
328  /// The geometrical information can be read from the FaceGeometry<ELEMENT>
329  /// class and and thus, we can be generic enough without the need to have
330  /// a separate equations class.
331  ///
332  /// THIS IS THE REFINEABLE VERSION.
333  //======================================================================
334  template<class ELEMENT>
336  : public virtual NavierStokesTractionElement<ELEMENT>,
338  {
339  public:
340  /// Constructor, which takes a "bulk" element and the face index
342  const int& face_index)
343  : // we're calling this from the constructor of the refineable version.
344  NavierStokesTractionElement<ELEMENT>(element_pt, face_index, true)
345  {
346  }
347 
348  /// Destructor should not delete anything
350 
351 
352  /// Number of continuously interpolated values are the
353  /// same as those in the bulk element.
354  unsigned ncont_interpolated_values() const
355  {
356  return dynamic_cast<ELEMENT*>(this->bulk_element_pt())
358  }
359 
360  /// This function returns just the residuals
362  {
363  // Call the generic residuals function using a dummy matrix argument
365  residuals, GeneralisedElement::Dummy_matrix, 0);
366  }
367 
368  /// This function returns the residuals and the Jacobian
370  DenseMatrix<double>& jacobian)
371  {
372  // Call the generic routine
374  residuals, jacobian, 1);
375  }
376 
377 
378  protected:
379  /// This function returns the residuals for the
380  /// traction function.
381  /// flag=1(or 0): do (or don't) compute the Jacobian as well.
383  Vector<double>& residuals, DenseMatrix<double>& jacobian, unsigned flag);
384  };
385 
386 
387  /// ////////////////////////////////////////////////////////////////////
388  /// ////////////////////////////////////////////////////////////////////
389  /// ////////////////////////////////////////////////////////////////////
390 
391 
392  //============================================================================
393  /// Function that returns the residuals for the imposed traction Navier_Stokes
394  /// equations
395  //============================================================================
396  template<class ELEMENT>
399  Vector<double>& residuals, DenseMatrix<double>& jacobian, unsigned flag)
400  {
401  // Get the indices at which the velocity components are stored
402  unsigned u_nodal_index[this->Dim];
403  for (unsigned i = 0; i < this->Dim; i++)
404  {
405  u_nodal_index[i] =
406  dynamic_cast<ELEMENT*>(this->bulk_element_pt())->u_index_nst(i);
407  }
408 
409  // Find out how many nodes there are
410  unsigned n_node = nnode();
411 
412  // Get continuous time from timestepper of first node
413  double time = node_pt(0)->time_stepper_pt()->time_pt()->time();
414 
415  // Set up memory for the shape and test functions
416  Shape psif(n_node), testf(n_node);
417 
418  // Set the value of n_intpt
419  unsigned n_intpt = integral_pt()->nweight();
420 
421  // Integers to store local equation numbers
422  int local_eqn = 0;
423 
424  // Loop over the integration points
425  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
426  {
427  // Get the integral weight
428  double w = integral_pt()->weight(ipt);
429 
430  // Find the shape and test functions and return the Jacobian
431  // of the mapping
432  double J = this->shape_and_test_at_knot(ipt, psif, testf);
433 
434  // Premultiply the weights and the Jacobian
435  double W = w * J;
436 
437  // Need to find position to feed into Traction function
438  Vector<double> interpolated_x(this->Dim);
439 
440  // Initialise to zero
441  for (unsigned i = 0; i < this->Dim; i++)
442  {
443  interpolated_x[i] = 0.0;
444  }
445 
446  // Calculate velocities and derivatives
447  for (unsigned l = 0; l < n_node; l++)
448  {
449  // Loop over velocity components
450  for (unsigned i = 0; i < this->Dim; i++)
451  {
452  interpolated_x[i] += nodal_position(l, i) * psif[l];
453  }
454  }
455 
456  // Get the outer unit normal
457  Vector<double> interpolated_n(this->Dim);
458  this->outer_unit_normal(ipt, interpolated_n);
459 
460  // Get the user-defined traction terms
461  Vector<double> traction(this->Dim);
462  this->get_traction(time, interpolated_x, interpolated_n, traction);
463 
464  // Now add to the appropriate equations
465 
466  // Number of master nodes and storage for the weight of the shape function
467  unsigned n_master = 1;
468  double hang_weight = 1.0;
469 
470  // Pointer to hang info object
471  HangInfo* hang_info_pt = 0;
472 
473  // Loop over the nodes for the test functions/equations
474  //----------------------------------------------------
475  for (unsigned l = 0; l < n_node; l++)
476  {
477  // Local boolean to indicate whether the node is hanging
478  bool is_node_hanging = this->node_pt(l)->is_hanging();
479 
480  // If the node is hanging
481  if (is_node_hanging)
482  {
483  hang_info_pt = this->node_pt(l)->hanging_pt();
484 
485  // Read out number of master nodes from hanging data
486  n_master = hang_info_pt->nmaster();
487  }
488  // Otherwise the node is its own master
489  else
490  {
491  n_master = 1;
492  }
493 
494  // Loop over the master nodes
495  for (unsigned m = 0; m < n_master; m++)
496  {
497  // Loop over velocity components for equations
498  for (unsigned i = 0; i < this->Dim; i++)
499  {
500  // Get the equation number
501  // If the node is hanging
502  if (is_node_hanging)
503  {
504  // Get the equation number from the master node
505  local_eqn = this->local_hang_eqn(hang_info_pt->master_node_pt(m),
506  u_nodal_index[i]);
507  // Get the hang weight from the master node
508  hang_weight = hang_info_pt->master_weight(m);
509  }
510  // If the node is not hanging
511  else
512  {
513  // Local equation number
514  local_eqn = this->nodal_local_eqn(l, u_nodal_index[i]);
515 
516  // Node contributes with full weight
517  hang_weight = 1.0;
518  }
519 
520  // If it's not a boundary condition...
521  if (local_eqn >= 0)
522  {
523  /* //Loop over the test functions */
524  /* for(unsigned l=0;l<n_node;l++) */
525  /* { */
526  /* //Loop over the velocity components */
527  /* for(unsigned i=0;i<Dim;i++) */
528  /* { */
529  /* local_eqn = u_local_eqn(l,i); */
530  /* /\*IF it's not a boundary condition*\/ */
531  /* if(local_eqn >= 0) */
532  /* { */
533 
534 
535  // Add the user-defined traction terms
536  residuals[local_eqn] += traction[i] * testf[l] * W * hang_weight;
537 
538  // Assuming the the traction DOES NOT depend upon velocities
539  // or pressures, the jacobian is always zero, so no jacobian
540  // terms are required
541  }
542  }
543  } // End of loop over dimension
544  } // End of loop over shape functions
545  }
546  }
547 
548 
549 } // namespace oomph
550 
551 #endif
cstr elem_len * i
Definition: cfortran.h:603
char t
Definition: cfortran.h:568
FaceElements are elements that coincide with the faces of higher-dimensional "bulk" elements....
Definition: elements.h:4338
int & face_index()
Index of the face (a number that uniquely identifies the face in the element)
Definition: elements.h:4626
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
In a FaceElement, the "global" intrinsic coordinate of the element along the boundary,...
Definition: elements.h:4497
FiniteElement *& bulk_element_pt()
Pointer to higher-dimensional "bulk" element.
Definition: elements.h:4735
double J_eulerian_at_knot(const unsigned &ipt) const
Return the Jacobian of the mapping from local to global coordinates at the ipt-th integration point O...
Definition: elements.cc:5328
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
Definition: elements.h:4998
A general Finite Element class.
Definition: elements.h:1313
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
virtual void output(std::ostream &outfile)
Output the element data — typically the values at the nodes in a format suitable for post-processing.
Definition: elements.h:3050
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:1432
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
Definition: elements.h:2611
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
virtual void build_face_element(const int &face_index, FaceElement *face_element_pt)
Function for building a lower dimensional FaceElement on the specified face of the FiniteElement....
Definition: elements.cc:5132
virtual void shape_at_knot(const unsigned &ipt, Shape &psi) const
Return the geometric shape function at the ipt-th integration point.
Definition: elements.cc:3220
bool has_hanging_nodes() const
Return boolean to indicate if any of the element's nodes are geometrically hanging.
Definition: elements.h:2470
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
Definition: elements.h:227
Class that contains data for hanging nodes.
Definition: nodes.h:742
double const & master_weight(const unsigned &i) const
Return weight for dofs on i-th master node.
Definition: nodes.h:808
Node *const & master_node_pt(const unsigned &i) const
Return a pointer to the i-th master node.
Definition: nodes.h:791
unsigned nmaster() const
Return the number of master nodes.
Definition: nodes.h:785
A class for elements that allow the imposition of an applied traction to the Navier–Stokes equations ...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
This function returns just the residuals.
void output(std::ostream &outfile)
Overload the output function.
void(* Traction_fct_pt)(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &result)
Pointer to an imposed traction function.
double shape_and_test_at_knot(const unsigned &ipt, Shape &psi, Shape &test) const
Function to compute the shape and test functions and to return the Jacobian of mapping.
virtual int u_local_eqn(const unsigned &n, const unsigned &i)
Access function that returns the local equation numbers for velocity components. u_local_eqn(n,...
void(*&)(const double &t, const Vector< double > &x, const Vector< double > &n, Vector< double > &result) traction_fct_pt()
void fill_in_generic_residual_contribution_fluid_traction(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
This function returns the residuals for the traction function. flag=1(or 0): do (or don't) compute th...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
This function returns the residuals and the jacobian.
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
The "global" intrinsic coordinate of the element when viewed as part of a geometric object should be ...
unsigned Dim
The highest dimension of the problem.
void output(std::ostream &outfile, const unsigned &nplot)
Output function: x,y,[z],u,v,[w],p in tecplot format.
~NavierStokesTractionElement()
Destructor should not delete anything.
void get_traction(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &result)
Function to calculate the traction applied to the fluid.
NavierStokesTractionElement(FiniteElement *const &element_pt, const int &face_index, const bool &called_from_refineable_constructor=false)
Constructor, which takes a "bulk" element and the value of the index and its limit.
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
Definition: nodes.h:1054
/////////////////////////////////////////////////////////////////// /////////////////////////////////...
An OomphLibError object which should be thrown when an run-time error is encountered....
RefineableElements are FiniteElements that may be subdivided into children to provide a better local ...
////////////////////////////////////////////////////////////////////// //////////////////////////////...
void refineable_fill_in_generic_residual_contribution_fluid_traction(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
This function returns the residuals for the traction function. flag=1(or 0): do (or don't) compute th...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
This function returns the residuals and the Jacobian.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
This function returns just the residuals.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values are the same as those in the bulk element.
RefineableNavierStokesTractionElement(FiniteElement *const &element_pt, const int &face_index)
Constructor, which takes a "bulk" element and the face index.
~RefineableNavierStokesTractionElement()
Destructor should not delete anything.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
unsigned Dim
Dimension of zeta tuples (set by get_dim_helper) – needed because we store the scalar coordinates in ...
Definition: multi_domain.cc:60
//////////////////////////////////////////////////////////////////// ////////////////////////////////...