spacetime_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-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 // 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_SPACETIME_FLUID_TRACTION_ELEMENTS_HEADER
31 #define OOMPH_SPACETIME_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 // OOMPH-LIB headers
39 #include "generic/Qelements.h"
40 #include "generic/Telements.h"
41 
42 namespace oomph
43 {
44  //======================================================================
45  /// A class for elements that allow the imposition of an applied traction
46  /// to the Navier--Stokes equations
47  /// The geometrical information can be read from the FaceGeometry<ELEMENT>
48  /// class and and thus, we can be generic enough without the need to have
49  /// a separate equations class
50  //======================================================================
51  template<class ELEMENT>
52  class NavierStokesSpaceTimeTractionElement
53  : 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  public:
64  /// Constructor, which takes a "bulk" element and the value of the index
65  /// and its limit
67  FiniteElement* const& element_pt,
68  const int& face_index,
69  const bool& called_from_refineable_constructor = false)
70  : FaceGeometry<ELEMENT>(), FaceElement()
71  {
72  // Attach the geometrical information to the element. N.B. This function
73  // also assigns nbulk_value from the required_nvalue of the bulk element
74  element_pt->build_face_element(face_index, this);
75 
76 #ifdef PARANOID
77  {
78  // Check that the element is not a refineable 3D element
79  if (!called_from_refineable_constructor)
80  {
81  // If it's 3D
82  if (element_pt->dim() == 3)
83  {
84  // Is it refineable
85  RefineableElement* ref_el_pt =
86  dynamic_cast<RefineableElement*>(element_pt);
87 
88  // If the upcast was successful
89  if (ref_el_pt != 0)
90  {
91  // If there are any hanging nodes
92  if (this->has_hanging_nodes())
93  {
94  // Create an output stream
95  std::ostringstream error_message_stream;
96 
97  // Create an error message
98  error_message_stream
99  << "This flux element will not work correctly "
100  << "if nodes are hanging" << std::endl;
101 
102  // Throw the error message
103  throw OomphLibError(error_message_stream.str(),
104  OOMPH_CURRENT_FUNCTION,
105  OOMPH_EXCEPTION_LOCATION);
106  }
107  } // if (ref_el_pt!=0)
108  } // if (element_pt->dim()==3)
109  } // if (!called_from_refineable_constructor)
110  }
111 #endif
112 
113  // Set the body force function pointer to zero
114  Traction_fct_pt = 0;
115 
116  // Set the dimension from the dimension of the first node
117  Dim = this->node_pt(0)->ndim();
118  } // End of NavierStokesSpaceTimeTractionElement
119 
120 
121  /// Destructor should not delete anything
123 
124 
125  /// Access function for the imposed traction pointer
126  void (*&traction_fct_pt())(const double& t,
127  const Vector<double>& x,
128  const Vector<double>& n,
129  Vector<double>& result)
130  {
131  // Return the function pointer
132  return Traction_fct_pt;
133  } // End of traction_fct_pt
134 
135 
136  /// This function returns just the residuals
138  {
139  // Call the generic residuals function with flag set to 0
140  // using a dummy matrix argument
142  residuals, GeneralisedElement::Dummy_matrix, 0);
143  }
144 
145 
146  /// This function returns the residuals and the jacobian
148  DenseMatrix<double>& jacobian)
149  {
150  // Call the generic routine with the flag set to 1
152  residuals, jacobian, 1);
153  }
154 
155 
156  /// Overload the output function
157  void output(std::ostream& outfile)
158  {
159  FiniteElement::output(outfile);
160  }
161 
162 
163  /// Output function: x,y,[z],u,v,[w],p in tecplot format
164  void output(std::ostream& outfile, const unsigned& nplot)
165  {
166  FiniteElement::output(outfile, nplot);
167  }
168 
169  protected:
170  /// This function returns the residuals for the traction function.
171  /// flag=1(or 0): do (or don't) compute the Jacobian as well.
173  Vector<double>& residuals,
174  DenseMatrix<double>& jacobian,
175  const unsigned& flag);
176 
177 
178  /// The "global" intrinsic coordinate of the element when
179  /// viewed as part of a geometric object should be given by
180  /// the FaceSpaceTimeElement representation, by default
181  double zeta_nodal(const unsigned& n,
182  const unsigned& k,
183  const unsigned& i) const
184  {
185  return FaceElement::zeta_nodal(n, k, i);
186  } // End of zeta_nodal
187 
188 
189  /// Access function that returns the local equation numbers
190  /// for velocity components.
191  /// u_local_eqn(n,i)=local equation number or < 0 if pinned.
192  /// The default is to asssume that n is the local node number
193  /// and the i-th velocity component is the i-th unknown stored at the node.
194  virtual inline int u_local_eqn(const unsigned& n, const unsigned& i)
195  {
196  return nodal_local_eqn(n, i);
197  }
198 
199  /// Function to compute the shape and test functions and to return
200  /// the Jacobian of mapping
201  inline double shape_and_test_at_knot(const unsigned& ipt,
202  Shape& psi,
203  Shape& test) const
204  {
205  // Calculate the shape functions
206  shape_at_knot(ipt, psi);
207 
208  // The Galerkin method uses the same test functions as shape functions
209  test = psi;
210 
211  // Return the value of the jacobian
212  return J_eulerian_at_knot(ipt);
213  } // End of shape_and_test_at_knot
214 
215 
216  /// Function to calculate the traction applied to the fluid
217  void get_traction(const double& time,
218  const Vector<double>& x,
219  const Vector<double>& n,
220  Vector<double>& result)
221  {
222  // If the function pointer is zero return zero
223  if (Traction_fct_pt == 0)
224  {
225  // Loop over the (spatial) dimensions
226  for (unsigned i = 0; i < Dim - 1; i++)
227  {
228  // Set the i-th traction component to zero
229  result[i] = 0.0;
230  }
231  }
232  // If the traction function pointer has been set
233  else
234  {
235  // Call the function pointer
236  (*Traction_fct_pt)(time, x, n, result);
237  } // if (Traction_fct_pt==0)
238  } // End of get_traction
239 
240 
241  /// The highest dimension of the problem
242  unsigned Dim;
243  };
244 
245  /// ////////////////////////////////////////////////////////////////////
246  /// ////////////////////////////////////////////////////////////////////
247  /// ////////////////////////////////////////////////////////////////////
248 
249  //=========================================================================
250  /// Function that returns the residuals for the imposed traction
251  /// Navier-Stokes equations
252  //=========================================================================
253  template<class ELEMENT>
256  Vector<double>& residuals,
257  DenseMatrix<double>& jacobian,
258  const unsigned& flag)
259  {
260  // Find out how many nodes there are
261  unsigned n_node = nnode();
262 
263  // Set up memory for the shape functions
264  Shape psif(n_node);
265 
266  // Set up memory for the test functions
267  Shape testf(n_node);
268 
269  // Number of integration points
270  unsigned n_intpt = integral_pt()->nweight();
271 
272  // Integers to store local equation numbers
273  int local_eqn = 0;
274 
275  // Loop over the integration points
276  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
277  {
278  // Get the integral weight
279  double w = integral_pt()->weight(ipt);
280 
281  // Find the shape and test functions and return the Jacobian
282  // of the mapping
283  double J = shape_and_test_at_knot(ipt, psif, testf);
284 
285  // Premultiply the weights and the Jacobian
286  double W = w * J;
287 
288  // Get continuous time
289  double interpolated_t = 0.0;
290 
291  // Need to find position to feed into traction function
292  Vector<double> interpolated_x(Dim - 1, 0.0);
293 
294  // Allocate space for the outer unit normal
295  Vector<double> interpolated_spacetime_n(Dim, 0.0);
296 
297  // Allocate space for the outer unit normal (in the given time slice)
298  Vector<double> interpolated_n(Dim - 1, 0.0);
299 
300  // Allocate space for the traction vector
301  Vector<double> traction(Dim - 1, 0.0);
302 
303  // Calculate velocities and derivatives
304  for (unsigned l = 0; l < n_node; l++)
305  {
306  // Update the interpolated time value
307  interpolated_t += raw_nodal_position(l, Dim - 1) * psif(l);
308 
309  // Loop over velocity components
310  for (unsigned i = 0; i < Dim - 1; i++)
311  {
312  // Update the i-th interpolated coordinate value
313  interpolated_x[i] += nodal_position(l, i) * psif(l);
314  }
315  } // for (unsigned l=0;l<n_node;l++)
316 
317  // Calculate the outer unit normal
318  outer_unit_normal(ipt, interpolated_spacetime_n);
319 
320 #ifdef PARANOID
321  // Sanity check: Make sure the normal component in the time-direction is
322  // zero. If it isn't then something has gone wrong...
323  if (std::fabs(interpolated_spacetime_n[Dim - 1]) > 1.0e-15)
324  {
325  // Create an output stream
326  std::ostringstream error_message_stream;
327 
328  // Create the error message
329  error_message_stream
330  << "The component of the outer unit normal in the "
331  << "time-direction\nshould be zero but the actual "
332  << "value is: " << interpolated_spacetime_n[Dim - 1] << std::endl;
333 
334  // Throw the error message
335  throw OomphLibError(error_message_stream.str(),
336  OOMPH_CURRENT_FUNCTION,
337  OOMPH_EXCEPTION_LOCATION);
338  }
339 #endif
340 
341  // Loop over spatial coordinates
342  for (unsigned i = 0; i < Dim - 1; i++)
343  {
344  // Copy the i-th normal component over
345  interpolated_n[i] = interpolated_spacetime_n[i];
346  }
347 
348  // Get the user-defined traction
349  get_traction(interpolated_t, interpolated_x, interpolated_n, traction);
350 
351  // Now add to the appropriate equations:
352 
353  // Loop over the test functions
354  for (unsigned l = 0; l < n_node; l++)
355  {
356  // Loop over the velocity components
357  for (unsigned i = 0; i < Dim - 1; i++)
358  {
359  // Get the local equation number associated with this unknown and node
360  local_eqn = u_local_eqn(l, i);
361 
362  // If it's not a boundary condition
363  if (local_eqn >= 0)
364  {
365  // Add the user-defined traction terms
366  residuals[local_eqn] += traction[i] * testf[l] * W;
367 
368  // Assuming the traction DOES NOT depend upon velocities
369  // or pressures, the jacobian is always zero, so no jacobian
370  // terms are required
371  }
372  } // for (unsigned i=0;i<Dim;i++)
373  } // for (unsigned l=0;l<n_node;l++)
374  } // for (unsigned ipt=0;ipt<n_intpt;ipt++)
375  } // End of fill_in_generic_residual_contribution_fluid_traction
376 
377  /// //////////////////////////////////////////////////////////////////////
378  /// //////////////////////////////////////////////////////////////////////
379  /// //////////////////////////////////////////////////////////////////////
380 
381  //======================================================================
382  /// A class for elements that allow the imposition of an applied traction
383  /// to the Navier--Stokes equations
384  /// The geometrical information can be read from the FaceGeometry<ELEMENT>
385  /// class and and thus, we can be generic enough without the need to have
386  /// a separate equations class.
387  ///
388  /// THIS IS THE REFINEABLE VERSION.
389  //======================================================================
390  template<class ELEMENT>
391  class RefineableNavierStokesSpaceTimeTractionElement
392  : public virtual NavierStokesSpaceTimeTractionElement<ELEMENT>,
393  public virtual NonRefineableElementWithHangingNodes
394  {
395  public:
396  /// Constructor, which takes a "bulk" element and the face index
398  FiniteElement* const& element_pt, const int& face_index)
399  : // we're calling this from the constructor of the refineable version.
401  element_pt, face_index, true)
402  {
403  }
404 
405 
406  /// Destructor should not delete anything
408 
409 
410  /// Number of continuously interpolated values are the
411  /// same as those in the bulk element.
412  unsigned ncont_interpolated_values() const
413  {
414  return dynamic_cast<ELEMENT*>(this->bulk_element_pt())
416  }
417 
418  /// This function returns just the residuals
420  {
421  // Call the generic residuals function using a dummy matrix argument
423  residuals, GeneralisedElement::Dummy_matrix, 0);
424  }
425 
426  /// This function returns the residuals and the Jacobian
428  DenseMatrix<double>& jacobian)
429  {
430  // Call the generic routine
432  residuals, jacobian, 1);
433  }
434 
435 
436  protected:
437  /// This function returns the residuals for the traction function.
438  /// flag=1(or 0): do (or don't) compute the Jacobian as well.
440  Vector<double>& residuals,
441  DenseMatrix<double>& jacobian,
442  const unsigned& flag);
443  };
444 
445  /// ////////////////////////////////////////////////////////////////////
446  /// ////////////////////////////////////////////////////////////////////
447  /// ////////////////////////////////////////////////////////////////////
448 
449  //============================================================================
450  /// Function that returns the residuals for the imposed traction Navier_Stokes
451  /// equations
452  //============================================================================
453  template<class ELEMENT>
456  Vector<double>& residuals,
457  DenseMatrix<double>& jacobian,
458  const unsigned& flag)
459  {
460  // Get the indices at which the velocity components are stored
461  unsigned u_nodal_index[this->Dim - 1];
462 
463  // Loop over the velocity components
464  for (unsigned i = 0; i < this->Dim - 1; i++)
465  {
466  // Get the i-th velocity component index
467  u_nodal_index[i] =
468  dynamic_cast<ELEMENT*>(this->bulk_element_pt())->u_index_nst(i);
469  }
470 
471  // Find out how many nodes there are
472  unsigned n_node = nnode();
473 
474  // Set up memory for the shape functions
475  Shape psif(n_node);
476 
477  // Set up memory for the test functions
478  Shape testf(n_node);
479 
480  // Number of integration points
481  unsigned n_intpt = integral_pt()->nweight();
482 
483  // Integers to store local equation numbers
484  int local_eqn = 0;
485 
486  // Loop over the integration points
487  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
488  {
489  // Get the integral weight
490  double w = integral_pt()->weight(ipt);
491 
492  // Find the shape and test functions and return the Jacobian
493  // of the mapping
494  double J = this->shape_and_test_at_knot(ipt, psif, testf);
495 
496  // Premultiply the weights and the Jacobian
497  double W = w * J;
498 
499  // Get continuous time
500  double interpolated_t = 0.0;
501 
502  // Need to find position to feed into traction function
503  Vector<double> interpolated_x(this->Dim - 1, 0.0);
504 
505  // Allocate space for the outer unit normal
506  Vector<double> interpolated_spacetime_n(this->Dim, 0.0);
507 
508  // Allocate space for the outer unit normal (in the given time slice)
509  Vector<double> interpolated_n(this->Dim - 1, 0.0);
510 
511  // Allocate space for the traction vector
512  Vector<double> traction(this->Dim - 1, 0.0);
513 
514  // Calculate velocities and derivatives
515  for (unsigned l = 0; l < n_node; l++)
516  {
517  // Update the interpolated time value
518  interpolated_t += raw_nodal_position(l, this->Dim - 1) * psif(l);
519 
520  // Loop over velocity components
521  for (unsigned i = 0; i < this->Dim - 1; i++)
522  {
523  // Update the i-th interpolated coordinate value
524  interpolated_x[i] += this->nodal_position(l, i) * psif(l);
525  }
526  } // for (unsigned l=0;l<n_node;l++)
527 
528  // Calculate the outer unit normal
529  this->outer_unit_normal(ipt, interpolated_spacetime_n);
530 
531 #ifdef PARANOID
532  // Sanity check: Make sure the normal component in the time-direction is
533  // zero. If it isn't then something has gone wrong...
534  if (std::fabs(interpolated_n[this->Dim - 1]) > 1.0e-15)
535  {
536  // Create an output stream
537  std::ostringstream error_message_stream;
538 
539  // Create the error message
540  error_message_stream << "The component of the outer unit normal in the "
541  << "time-direction\nshould be zero but the actual "
542  << "value is: " << interpolated_n[this->Dim - 1]
543  << std::endl;
544 
545  // Throw the error message
546  throw OomphLibError(error_message_stream.str(),
547  OOMPH_CURRENT_FUNCTION,
548  OOMPH_EXCEPTION_LOCATION);
549  }
550 #endif
551 
552  // Loop over spatial coordinates
553  for (unsigned i = 0; i < this->Dim - 1; i++)
554  {
555  // Copy the i-th normal component over
556  interpolated_n[i] = interpolated_spacetime_n[i];
557  }
558 
559  // Get the user-defined traction
560  this->get_traction(
561  interpolated_t, interpolated_x, interpolated_n, traction);
562 
563  // Now add to the appropriate equations:
564 
565  // Number of master nodes and storage for the weight of the shape function
566  unsigned n_master = 1;
567 
568  // The weight of the shape function
569  double hang_weight = 1.0;
570 
571  // Pointer to hang info object
572  HangInfo* hang_info_pt = 0;
573 
574  // Loop over the nodes for the test functions/equations
575  //-----------------------------------------------------
576  for (unsigned l = 0; l < n_node; l++)
577  {
578  // Local boolean to indicate whether the node is hanging
579  bool is_node_hanging = this->node_pt(l)->is_hanging();
580 
581  // If the node is hanging
582  if (is_node_hanging)
583  {
584  // Get the HangInfo pointer associated with this hanging node
585  hang_info_pt = this->node_pt(l)->hanging_pt();
586 
587  // Read out number of master nodes from hanging data
588  n_master = hang_info_pt->nmaster();
589  }
590  // Otherwise the node is its own master
591  else
592  {
593  // If it's not hanging then it is its (one and only) master node
594  n_master = 1;
595  }
596 
597  // Loop over the master nodes
598  for (unsigned m = 0; m < n_master; m++)
599  {
600  // Loop over velocity components for equations
601  for (unsigned i = 0; i < this->Dim - 1; i++)
602  {
603  // Get the equation number:
604 
605  // If the node is hanging
606  if (is_node_hanging)
607  {
608  // Get the equation number from the master node
609  local_eqn = this->local_hang_eqn(hang_info_pt->master_node_pt(m),
610  u_nodal_index[i]);
611 
612  // Get the hang weight from the master node
613  hang_weight = hang_info_pt->master_weight(m);
614  }
615  // If the node is not hanging
616  else
617  {
618  // Local equation number
619  local_eqn = this->nodal_local_eqn(l, u_nodal_index[i]);
620 
621  // Node contributes with full weight
622  hang_weight = 1.0;
623  }
624 
625  // If it's not a boundary condition...
626  if (local_eqn >= 0)
627  {
628  // Add the user-defined traction terms
629  residuals[local_eqn] += traction[i] * testf[l] * W * hang_weight;
630 
631  // Assuming the the traction DOES NOT depend upon velocities
632  // or pressures, the jacobian is always zero, so no jacobian
633  // terms are required.
634  }
635  } // for (unsigned i=0;i<this->Dim-1;i++)
636  } // for (unsigned m=0;m<n_master;m++)
637  } // for (unsigned l=0;l<n_node;l++)
638  } // for (unsigned ipt=0;ipt<n_intpt;ipt++)
639  } // End of refineable_fill_in_generic_residual_contribution_fluid_traction
640 } // End of namespace oomph
641 
642 #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:4342
int & face_index()
Index of the face (a number that uniquely identifies the face in the element)
Definition: elements.h:4630
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:4501
FiniteElement *& bulk_element_pt()
Pointer to higher-dimensional "bulk" element.
Definition: elements.h:4739
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:5358
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
Definition: elements.h:5002
A general Finite Element class.
Definition: elements.h:1317
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2179
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:3054
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
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
Definition: elements.h:2615
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:5162
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:3250
bool has_hanging_nodes() const
Return boolean to indicate if any of the element's nodes are geometrically hanging.
Definition: elements.h:2474
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
A class for elements that allow the imposition of an applied traction to the Navier–Stokes equations ...
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 fill_in_contribution_to_residuals(Vector< double > &residuals)
This function returns just the residuals.
NavierStokesSpaceTimeTractionElement(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.
void output(std::ostream &outfile, const unsigned &nplot)
Output function: x,y,[z],u,v,[w],p in tecplot format.
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.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
This function returns the residuals and the jacobian.
void output(std::ostream &outfile)
Overload the output function.
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 ...
~NavierStokesSpaceTimeTractionElement()
Destructor should not delete anything.
void fill_in_generic_residual_contribution_fluid_traction(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
This function returns the residuals for the traction function. flag=1(or 0): do (or don't) compute th...
void(* Traction_fct_pt)(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &result)
Pointer to an imposed traction function.
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.
void(*&)(const double &t, const Vector< double > &x, const Vector< double > &n, Vector< double > &result) traction_fct_pt()
Access function for the imposed traction pointer.
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 fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
This function returns the residuals and the Jacobian.
RefineableNavierStokesSpaceTimeTractionElement(FiniteElement *const &element_pt, const int &face_index)
Constructor, which takes a "bulk" element and the face index.
void refineable_fill_in_generic_residual_contribution_fluid_traction(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
This function returns the residuals for the traction function. flag=1(or 0): do (or don't) compute th...
~RefineableNavierStokesSpaceTimeTractionElement()
Destructor should not delete anything.
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.
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
//////////////////////////////////////////////////////////////////// ////////////////////////////////...