navier_stokes_flux_control_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 // Include guard to prevent multiple inclusions of the header
27 #ifndef OOMPH_NAVIER_STOKES_FLUX_CONTROL_ELEMENTS
28 #define OOMPH_NAVIER_STOKES_FLUX_CONTROL_ELEMENTS
29 
30 // Config header generated by autoconfig
31 #ifdef HAVE_CONFIG_H
32 #include <oomph-lib-config.h>
33 #endif
34 
35 // OOMPH-LIB headers
36 #include "../generic/nodes.h"
37 #include "../navier_stokes/navier_stokes_surface_power_elements.h"
38 
39 namespace oomph
40 {
41  /// ///////////////////////////////////////////////////////////////////////
42  /// ///////////////////////////////////////////////////////////////////////
43  /// ///////////////////////////////////////////////////////////////////////
44 
45 
46  //======================================================================
47  /// A template free base class for an element to imposes an applied
48  /// boundary pressure to the Navier-Stokes equations in order to
49  /// control a volume flux when used in conjunction with a
50  /// NetFluxControlElement or
51  /// NetFluxControlElementForWomersleyPressureControl).
52  //======================================================================
54  : public virtual GeneralisedElement
55  {
56  public:
57  /// Empty constructor
59 
60  /// Empty virtual destructor
62 
63  /// Pure virtual function to calculate integral of the volume flux
64  virtual double get_volume_flux() = 0;
65 
66  /// Function adds to the external data the Data object whose
67  /// single value is the pressure applied by the element
68  void add_pressure_data(Data* pressure_data_pt)
69  {
70  Pressure_data_id = add_external_data(pressure_data_pt);
71  }
72 
73  protected:
74  /// Access function gives id of external Data object whose
75  /// single value is the pressure applied by the element
76  unsigned& pressure_data_id()
77  {
78  return Pressure_data_id;
79  }
80 
81  private:
82  /// Id of external Data object whose single value is the
83  /// pressure applied by the elements
84  unsigned Pressure_data_id;
85  };
86 
87 
88  /// ///////////////////////////////////////////////////////////////////////
89  /// ///////////////////////////////////////////////////////////////////////
90  /// ///////////////////////////////////////////////////////////////////////
91 
92 
93  //======================================================================
94  /// A class for an element that controls the net fluid flux across a
95  /// boundary by the imposition of an unknown applied pressure to the
96  /// Navier-Stokes equations. This element is used with a mesh of
97  /// NavierStokesFluxControlElement elements which are attached
98  /// to the boundary.
99  /// Note: fill_in_contribution_to_jacobian() does not calculate
100  /// Jacobian contributions for this element as they are calculated by
101  /// NavierStokesFluxControlElement::fill_in_contribution_to_jacobian(...)
102  //======================================================================
104  {
105  public:
106  /// Constructor takes a mesh of
107  /// TemplateFreeNavierStokesFluxControlElementBase elements
108  /// that impose the pressure to control the flux, plus a pointer to
109  /// the double which contains the desired flux value
110  NetFluxControlElement(Mesh* flux_control_mesh_pt,
111  double* prescribed_flux_value_pt)
112  : Flux_control_mesh_pt(flux_control_mesh_pt),
113  Prescribed_flux_value_pt(prescribed_flux_value_pt)
114  {
115  // Construct Pressure_data_pt
116  Pressure_data_pt = new Data(1);
117 
118  // Add the new Data to internal Data for this element
120 
121  // There's no need to add the external data for this element since
122  // this elements Jacobian contributions are calculated by the
123  // NavierStokesFluxControlElements
124 
125  // Loop over elements in the Flux_control_mesh to add this element's
126  // Data to the external Data in the elements in the flux control mesh
127  unsigned n_el = Flux_control_mesh_pt->nelement();
128  for (unsigned e = 0; e < n_el; e++)
129  {
130  // Get pointer to the element
132 
133  // Perform cast to TemplateFreeNavierStokesFluxControlElementBase
134  // pointer
137 
138  flux_el_pt->add_pressure_data(Pressure_data_pt);
139  }
140 
141  // Default value for Dof_number_for_unknown, indiating that it's
142  // uninitialised
143  Dof_number_for_unknown = UINT_MAX;
144  }
145 
146 
147  /// Empty Destructor - Data gets deleted automatically
149 
150  /// Broken copy constructor
152 
153  /// Broken assignment operator
154  // Commented out broken assignment operator because this can lead to a
155  // conflict warning when used in the virtual inheritence hierarchy.
156  // Essentially the compiler doesn't realise that two separate
157  // implementations of the broken function are the same and so, quite
158  // rightly, it shouts.
159  /*void operator=(const NetFluxControlElement&) = delete;*/
160 
161  /// Spatial dimension of the problem
162  unsigned dim() const
163  {
164  return Dim;
165  }
166 
167  /// Function to return a pointer to the Data object whose
168  /// single value is the pressure applied by the
169  /// NavierStokesFluxControlElement elements
171  {
172  return Pressure_data_pt;
173  }
174 
175 
176  /// Add the element's contribution to its residual vector:
177  /// i.e. the flux constraint.
179  {
180  // Call the generic routine
182  }
183 
184  /// This function returns the residuals, but adds nothing to the
185  /// Jacobian as this element's Jacobian contributions are calculated by
186  /// the NavierStokesFluxControlElements which impose the traction
187  /// used to control the flux.
189  DenseMatrix<double>& jacobian)
190  {
191  // Call the generic routine
193  }
194 
195 
196  /// The number of "DOF types" that degrees of freedom in this element
197  /// are sub-divided into - it's set to Dof_number_for_unknown+1
198  /// because it's expected this element is added to a fluid mesh
199  /// containing navier stokes elements
200  unsigned ndof_types() const
201  {
202 #ifdef PARANOID
203  if (Dof_number_for_unknown == UINT_MAX)
204  {
205  std::ostringstream error_message;
206  error_message << "Dof_number_for_unknown hasn't been set yet!\n"
207  << "Please do so using the dof_number_for_unknown()\n"
208  << "access function\n";
209  throw OomphLibError(error_message.str(),
210  OOMPH_CURRENT_FUNCTION,
211  OOMPH_EXCEPTION_LOCATION);
212  }
213 #endif
214  return Dof_number_for_unknown + 1;
215  }
216 
217  /// Function to set / get the nodal value of the "DOF type" to which
218  /// the degree of freedom in this element (the pressure that enforces
219  /// the required volume flux!) is added to.
220  /// This should be set to the Navier-Stokes pressure DOF type
221  /// (usually the dimension of the problem, for example, in 3D, the DOF types
222  /// for single-physics Navier-Stokes elements are usually
223  /// labelled 0, 1, 2, 3 for u, v and w
224  /// velocities and pressure respectively. It is important to note that this
225  /// is dimension dependent, so should not be hard coded in!! In
226  /// particularly, this should not simply be set to the dimension of the
227  /// problem if there is further splitting of the velocity DOF types) if this
228  /// element is added to a fluid mesh containing Navier-Stokes elements.
230  {
231  return Dof_number_for_unknown;
232  }
233 
234  /// Create a list of pairs for all unknowns in this element,
235  /// so that the first entry in each pair contains the global equation
236  /// number of the unknown, while the second one contains the number
237  /// of the "DOF type" that this unknown is associated with.
238  /// (Function can obviously only be called if the equation numbering
239  /// scheme has been set up.) The single degree of freedom is given the
240  /// DOF type number of Dof_number_for_unknown since it's expected this
241  /// unknown is added to the Navier-Stokes pressure DOF block (it is also
242  /// assumed that the user has set the Dof_number_for_unknown variable to
243  /// the velocity DOF type using the function dof_number_for_unknown()).
245  std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list) const
246  {
247 #ifdef PARANOID
248  if (Dof_number_for_unknown == UINT_MAX)
249  {
250  std::ostringstream error_message;
251  error_message << "Dof_number_for_unknown hasn't been set yet!\n"
252  << "Please do so using the dof_number_for_unknown()\n"
253  << "access function\n";
254  throw OomphLibError(error_message.str(),
255  OOMPH_CURRENT_FUNCTION,
256  OOMPH_EXCEPTION_LOCATION);
257  }
258 #endif
259 
260  // pair to store dof lookup prior to being added to list
261  std::pair<unsigned, unsigned> dof_lookup;
262 
263  dof_lookup.first = this->eqn_number(0);
264  dof_lookup.second = Dof_number_for_unknown;
265 
266  // add to list
267  dof_lookup_list.push_front(dof_lookup);
268  }
269 
270  protected:
271  /// This function returns the residuals for the
272  /// flux control master element.
274  Vector<double>& residuals)
275  {
276  // Initialise volume flux
277  double volume_flux = 0.0;
278 
279  // Loop over elements in Flux_control_mesh_pt and calculate flux
280  unsigned n_el = Flux_control_mesh_pt->nelement();
281  for (unsigned e = 0; e < n_el; e++)
282  {
283  // Get a pointer to the element
285 
286  // Cast to NavierStokesFluxControlElement
287  TemplateFreeNavierStokesFluxControlElementBase* flux_control_el_pt = 0;
288  flux_control_el_pt =
290 
291 #ifdef PARANOID
292  if (flux_control_el_pt == 0)
293  {
294  throw OomphLibError("Element must be used with a mesh of "
295  "NavierStokesFluxControlElements",
296  OOMPH_CURRENT_FUNCTION,
297  OOMPH_EXCEPTION_LOCATION);
298  }
299 #endif
300 
301  // Add the elemental volume flux
302  volume_flux += flux_control_el_pt->get_volume_flux();
303  }
304 
305  residuals[0] += *Prescribed_flux_value_pt - volume_flux;
306  }
307 
308 
309  private:
310  /// Data object whose single value is the pressure
311  /// applied by the elements in the Flux_control_mesh_pt
313 
314  /// Mesh of elements which impose the pressure which controls
315  /// the net flux
317 
318  /// Pointer to the value that stores the prescribed flux
320 
321  /// The id number of the "DOF type" to which the degree
322  /// of freedom in this element is added to. This should be set to the
323  /// number id of the Navier-Stokes pressure DOF block (which is dimension
324  /// dependent!) if this element is added to a fluid mesh
325  /// containing navier stokes elements
327 
328  /// spatial dim of NS system
329  unsigned Dim;
330  };
331 
332 
333  /// ///////////////////////////////////////////////////////////////////////
334  /// ///////////////////////////////////////////////////////////////////////
335  /// ///////////////////////////////////////////////////////////////////////
336 
337 
338  //======================================================================
339  /// A class of element to impose an applied boundary pressure to
340  /// Navier-Stokes elements to control to control a volume flux. A mesh of
341  /// these elements are used in conjunction with a NetFluxControlElement.
342  /// The template arguement ELEMENT is a Navier-Stokes "bulk" element.
343  ///
344  /// Note: This element calculates Jacobian contributions for both itself
345  /// and also for the NetFluxControlElement with respect to its unknowns.
346  //======================================================================
347  template<class ELEMENT>
350  public virtual NavierStokesSurfacePowerElement<ELEMENT>
351  {
352  public:
353  /// Constructor, which takes a "bulk" element and face index
355  FiniteElement* const& element_pt,
356  const int& face_index,
357  const bool& called_from_refineable_constructor = false)
358  : NavierStokesSurfacePowerElement<ELEMENT>(element_pt, face_index)
359  {
360 #ifdef PARANOID
361  {
362  // Check that the element is not a refineable 3d element
363  if (!called_from_refineable_constructor)
364  {
365  ELEMENT* elem_pt = new ELEMENT;
366  // If it's three-d
367  if (elem_pt->dim() == 3)
368  {
369  // Is it refineable
370  if (dynamic_cast<RefineableElement*>(elem_pt))
371  {
372  // Throw Error
373  std::ostringstream error_message;
374  error_message
375  << "This element does not work properly with refineable bulk \n"
376  << "elements in 3D. Please use the refineable version\n"
377  << "instead.\n";
378  throw OomphLibError(error_message.str(),
379  OOMPH_CURRENT_FUNCTION,
380  OOMPH_EXCEPTION_LOCATION);
381  }
382  }
383  }
384  }
385 #endif
386 
387  // Set the dimension from the dimension of the first node (since Dim is
388  // private in the parent class)
389  Dim = this->node_pt(0)->ndim();
390  }
391 
392  /// Destructor should not delete anything
394 
395  /// This function returns just the residuals
397  {
398  // Call the generic residuals function using a dummy matrix argument
400  residuals, GeneralisedElement::Dummy_matrix, 0);
401  }
402 
403  /// This function returns the residuals and the Jacobian
404  /// including the Jacobian contribution from the flux control
405  /// master element with respect to dof in this
406  /// element
408  DenseMatrix<double>& jacobian)
409  {
410  // Call the generic routine
412  residuals, jacobian, 1);
413  }
414 
415  /// Function to get the integral of the volume flux
417  {
419  }
420 
421  protected:
422  /// Access function that returns the local equation numbers
423  /// for velocity components.
424  /// u_local_eqn(n,i) = local equation number or < 0 if pinned.
425  /// The default is to asssume that n is the local node number
426  /// and the i-th velocity component is the i-th unknown stored at the node.
427  virtual inline int u_local_eqn(const unsigned& n, const unsigned& i)
428  {
429  return this->nodal_local_eqn(n, i);
430  }
431 
432  /// Function to compute the shape and test functions and to return
433  /// the Jacobian of mapping
434  inline double shape_and_test_at_knot(const unsigned& ipt,
435  Shape& psi,
436  Shape& test) const
437  {
438  // Find number of nodes
439  unsigned n_node = this->nnode();
440  // Calculate the shape functions
441  this->shape_at_knot(ipt, psi);
442  // Set the test functions to be the same as the shape functions
443  for (unsigned i = 0; i < n_node; i++)
444  {
445  test[i] = psi[i];
446  }
447  // Return the value of the jacobian
448  return this->J_eulerian_at_knot(ipt);
449  }
450 
451 
452  /// This function returns the residuals for the traction function
453  /// flag=1(or 0): do (or don't) compute the Jacobian as well.
454  /// This function also calculates the Jacobian contribution for the
455  /// NetFluxControlElement
457  Vector<double>& residuals, DenseMatrix<double>& jacobian, unsigned flag)
458  {
459  // Find out how many nodes there are
460  unsigned n_node = this->nnode();
461 
462  // Set up memory for the shape and test functions
463  Shape psif(n_node), testf(n_node);
464 
465  // Set the value of n_intpt
466  unsigned n_intpt = this->integral_pt()->nweight();
467 
468  // Integers to store local equation numbers
469  int local_eqn = 0;
470 
471  // Get the pressure at the outflow
472  double pressure = this->external_data_pt(pressure_data_id())->value(0);
473 
474  // Loop over the integration points
475  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
476  {
477  // Get the integral weight
478  double w = this->integral_pt()->weight(ipt);
479 
480  // Find the shape and test functions and return the Jacobian
481  // of the mapping
482  double J = shape_and_test_at_knot(ipt, psif, testf);
483 
484  // Premultiply the weights and the Jacobian
485  double W = w * J;
486 
487  // Get the outer unit normal
488  Vector<double> unit_normal(Dim);
489  this->outer_unit_normal(ipt, unit_normal);
490 
491  // Calculate the traction
492  Vector<double> traction(Dim);
493  for (unsigned i = 0; i < Dim; i++)
494  {
495  traction[i] = -pressure * unit_normal[i];
496  }
497 
498  // Loop over the test functions
499  for (unsigned l = 0; l < n_node; l++)
500  {
501  // Loop over the velocity components
502  for (unsigned i = 0; i < Dim; i++)
503  {
504  local_eqn = u_local_eqn(l, i);
505 
506  /*IF it's not a boundary condition*/
507  if (local_eqn >= 0)
508  {
509  // Add the user-defined traction terms
510  residuals[local_eqn] += traction[i] * testf[l] * W;
511 
512  // Calculate the Jacobian if required. It is assumed
513  // that traction DOES NOT depend upon velocities
514  // or pressures in the Navier Stokes elements, but
515  // depend in the Data value which holds the
516  // pressure.
517  if (flag)
518  {
519  // Get equation number of the pressure data unknown
520  int local_unknown =
522 
523  // IF it's not a boundary condition
524  if (local_unknown >= 0)
525  {
526  // Add to Jacobian for this element
527  double jac_contribution = -unit_normal[i] * testf[l] * W;
528  jacobian(local_eqn, local_unknown) += jac_contribution;
529 
530  // Add to Jacobian for master element
531  jacobian(local_unknown, local_eqn) += jac_contribution;
532  }
533  }
534  }
535  } // End of loop over dimension
536  } // End of loop over shape functions
537  }
538  }
539 
540  protected:
541  /// The highest dimension of the problem
542  unsigned Dim;
543  };
544 
545 
546  /// //////////////////////////////////////////////////////////////////////
547  /// //////////////////////////////////////////////////////////////////////
548  /// //////////////////////////////////////////////////////////////////////
549 
550 
551  //======================================================================
552  /// A class of element to impose an applied boundary pressure to
553  /// Navier-Stokes elements to control to control a volume flux. A mesh of
554  /// these elements are used in conjunction with a NetFluxControlElement.
555  /// The template arguement ELEMENT is a Navier-Stokes "bulk" element.
556  ///
557  /// Note: This element calculates Jacobian contributions for both itself
558  /// and also for the NetFluxControlElement with respect to its unknowns.
559  ///
560  /// THIS IS THE REFINEABLE VERSION.
561  //======================================================================
562  template<class ELEMENT>
564  : public virtual NavierStokesFluxControlElement<ELEMENT>,
566  {
567  public:
568  /// Constructor, which takes a "bulk" element and the face index
570  const int& face_index)
571  : NavierStokesSurfacePowerElement<ELEMENT>(element_pt, face_index),
572  // we're calling this from the constructor of the refineable version.
573  NavierStokesFluxControlElement<ELEMENT>(element_pt, face_index, true)
574  {
575  }
576 
577  /// Destructor should not delete anything
579 
580 
581  /// Number of continuously interpolated values are the
582  /// same as those in the bulk element.
583  unsigned ncont_interpolated_values() const
584  {
585  return dynamic_cast<ELEMENT*>(this->bulk_element_pt())
587  }
588 
589  /// This function returns just the residuals
591  {
592  // Call the generic residuals function using a dummy matrix argument
594  residuals, GeneralisedElement::Dummy_matrix, 0);
595  }
596 
597  /// This function returns the residuals and the Jacobian
598  /// including the Jacobian contribution from the flux control
599  /// master element with respect to dof in this
600  /// element
602  DenseMatrix<double>& jacobian)
603  {
604  // Call the generic routine
606  residuals, jacobian, 1);
607  }
608 
609  protected:
610  /// This function returns the residuals for the traction function
611  /// flag=1(or 0): do (or don't) compute the Jacobian as well.
612  /// This function also calculates the Jacobian contribution for the
613  /// NetFluxControlElement
615  Vector<double>& residuals, DenseMatrix<double>& jacobian, unsigned flag)
616  {
617  // Get the indices at which the velocity components are stored
618  unsigned u_nodal_index[this->Dim];
619  for (unsigned i = 0; i < this->Dim; i++)
620  {
621  u_nodal_index[i] =
622  dynamic_cast<ELEMENT*>(this->bulk_element_pt())->u_index_nst(i);
623  }
624 
625  // Pointer to hang info object
626  HangInfo* hang_info_pt = 0;
627 
628  // Find out how many nodes there are
629  unsigned n_node = this->nnode();
630 
631  // Set up memory for the shape and test functions
632  Shape psif(n_node), testf(n_node);
633 
634  // Set the value of n_intpt
635  unsigned n_intpt = this->integral_pt()->nweight();
636 
637  // Integers to store local equation numbers
638  int local_eqn = 0;
639 
640  // Get the pressure at the outflow
641  double pressure =
642  this->external_data_pt(this->pressure_data_id())->value(0);
643 
644  // Loop over the integration points
645  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
646  {
647  // Get the integral weight
648  double w = this->integral_pt()->weight(ipt);
649 
650  // Find the shape and test functions and return the Jacobian
651  // of the mapping
652  double J = this->shape_and_test_at_knot(ipt, psif, testf);
653 
654  // Premultiply the weights and the Jacobian
655  double W = w * J;
656 
657  // Get the outer unit normal
658  Vector<double> unit_normal(this->Dim);
659  this->outer_unit_normal(ipt, unit_normal);
660 
661  // Calculate the traction
662  Vector<double> traction(this->Dim);
663  for (unsigned i = 0; i < this->Dim; i++)
664  {
665  traction[i] = -pressure * unit_normal[i];
666  }
667 
668 
669  // Number of master nodes and storage for the weight of the shape
670  // function
671  unsigned n_master = 1;
672  double hang_weight = 1.0;
673 
674  // Loop over the nodes for the test functions/equations
675  //----------------------------------------------------
676  for (unsigned l = 0; l < n_node; l++)
677  {
678  // Local boolean to indicate whether the node is hanging
679  bool is_node_hanging = this->node_pt(l)->is_hanging();
680 
681  // If the node is hanging
682  if (is_node_hanging)
683  {
684  hang_info_pt = this->node_pt(l)->hanging_pt();
685 
686  // Read out number of master nodes from hanging data
687  n_master = hang_info_pt->nmaster();
688  }
689  // Otherwise the node is its own master
690  else
691  {
692  n_master = 1;
693  }
694 
695  // Loop over the master nodes
696  for (unsigned m = 0; m < n_master; m++)
697  {
698  // Loop over velocity components for equations
699  for (unsigned i = 0; i < this->Dim; i++)
700  {
701  // Get the equation number
702  // If the node is hanging
703  if (is_node_hanging)
704  {
705  // Get the equation number from the master node
706  local_eqn = this->local_hang_eqn(
707  hang_info_pt->master_node_pt(m), u_nodal_index[i]);
708  // Get the hang weight from the master node
709  hang_weight = hang_info_pt->master_weight(m);
710  }
711  // If the node is not hanging
712  else
713  {
714  // Local equation number
715  local_eqn = this->nodal_local_eqn(l, u_nodal_index[i]);
716 
717  // Node contributes with full weight
718  hang_weight = 1.0;
719  }
720 
721  // If it's not a boundary condition...
722  if (local_eqn >= 0)
723  {
724  // Add the user-defined traction terms
725  residuals[local_eqn] +=
726  traction[i] * testf[l] * W * hang_weight;
727 
728  // Calculate the Jacobian if required. It is assumed
729  // that traction DOES NOT depend upon velocities
730  // or pressures in the Navier Stokes elements, but
731  // depend in the Data value which holds the
732  // pressure.
733  if (flag)
734  {
735  // Get equation number of the pressure data unknown
736  int local_unknown =
737  this->external_local_eqn(this->pressure_data_id(), 0);
738 
739  // IF it's not a boundary condition
740  if (local_unknown >= 0)
741  {
742  // Add to Jacobian for this element
743  double jac_contribution =
744  -unit_normal[i] * testf[l] * W * hang_weight;
745  jacobian(local_eqn, local_unknown) += jac_contribution;
746 
747  // Add to Jacobian for master element
748  jacobian(local_unknown, local_eqn) += jac_contribution;
749  }
750  }
751  }
752  } // End of loop over dimension
753  } // End of loop over master nodes
754  } // End of loop over nodes
755  }
756  }
757  };
758 
759 
760 } // namespace oomph
761 
762 #endif
e
Definition: cfortran.h:571
cstr elem_len * i
Definition: cfortran.h:603
A class that represents a collection of data; each Data object may contain many different individual ...
Definition: nodes.h:86
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
int & face_index()
Index of the face (a number that uniquely identifies the face in the element)
Definition: elements.h:4630
void outer_unit_normal(const Vector< double > &s, Vector< double > &unit_normal) const
Compute outer unit normal at the specified local coordinate.
Definition: elements.cc:6036
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
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
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 nnode() const
Return the number of nodes.
Definition: elements.h:2214
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1967
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
A Generalised Element class.
Definition: elements.h:73
unsigned add_internal_data(Data *const &data_pt, const bool &fd=true)
Add a (pointer to an) internal data object to the element and return the index required to obtain it ...
Definition: elements.cc:67
unsigned long eqn_number(const unsigned &ieqn_local) const
Return the global equation number corresponding to the ieqn_local-th local equation number.
Definition: elements.h:708
Data *& external_data_pt(const unsigned &i)
Return a pointer to i-th external data object.
Definition: elements.h:659
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
int external_local_eqn(const unsigned &i, const unsigned &j)
Return the local equation number corresponding to the j-th value stored at the i-th external data.
Definition: elements.h:311
unsigned add_external_data(Data *const &data_pt, const bool &fd=true)
Add a (pointer to an) external data object to the element and return its index (i....
Definition: elements.cc:312
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
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
A general mesh class.
Definition: mesh.h:67
GeneralisedElement *& element_pt(const unsigned long &e)
Return pointer to element e.
Definition: mesh.h:448
unsigned long nelement() const
Return number of elements in the mesh.
Definition: mesh.h:590
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
double get_volume_flux()
Function to get the integral of the volume flux.
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.
~NavierStokesFluxControlElement()
Destructor should not delete anything.
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 the...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
This function returns the residuals and the Jacobian including the Jacobian contribution from the flu...
unsigned Dim
The highest dimension of the problem.
NavierStokesFluxControlElement(FiniteElement *const &element_pt, const int &face_index, const bool &called_from_refineable_constructor=false)
Constructor, which takes a "bulk" element and face index.
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.
A class of elements that allow the determination of the power input and various other fluxes over the...
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector: i.e. the flux constraint.
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into - it's set to ...
void fill_in_generic_residual_contribution_flux_control(Vector< double > &residuals)
This function returns the residuals for the flux control master element.
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contain...
Mesh * Flux_control_mesh_pt
Mesh of elements which impose the pressure which controls the net flux.
unsigned & dof_number_for_unknown()
Function to set / get the nodal value of the "DOF type" to which the degree of freedom in this elemen...
~NetFluxControlElement()
Empty Destructor - Data gets deleted automatically.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
This function returns the residuals, but adds nothing to the Jacobian as this element's Jacobian cont...
unsigned Dof_number_for_unknown
The id number of the "DOF type" to which the degree of freedom in this element is added to....
NetFluxControlElement(const NetFluxControlElement &dummy)=delete
Broken copy constructor.
NetFluxControlElement(Mesh *flux_control_mesh_pt, double *prescribed_flux_value_pt)
Constructor takes a mesh of TemplateFreeNavierStokesFluxControlElementBase elements that impose the p...
Data * pressure_data_pt() const
Function to return a pointer to the Data object whose single value is the pressure applied by the Nav...
unsigned dim() const
Broken assignment operator.
Data * Pressure_data_pt
Data object whose single value is the pressure applied by the elements in the Flux_control_mesh_pt.
double * Prescribed_flux_value_pt
Pointer to the value that stores the prescribed flux.
HangInfo *const & hanging_pt() const
Return pointer to hanging node data (this refers to the geometric hanging node status) (const version...
Definition: nodes.h:1228
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
Definition: nodes.h:1054
bool is_hanging() const
Test whether the node is geometrically hanging.
Definition: nodes.h:1285
/////////////////////////////////////////////////////////////////// /////////////////////////////////...
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 ...
int local_hang_eqn(Node *const &node_pt, const unsigned &i)
Access function that returns the local equation number for the hanging node variables (values stored ...
////////////////////////////////////////////////////////////////////// //////////////////////////////...
unsigned ncont_interpolated_values() const
Number of continuously interpolated values are the same as those in the bulk element.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
This function returns the residuals and the Jacobian including the Jacobian contribution from the flu...
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 the...
~RefineableNavierStokesFluxControlElement()
Destructor should not delete anything.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
This function returns just the residuals.
RefineableNavierStokesFluxControlElement(FiniteElement *const &element_pt, const int &face_index)
Constructor, which takes a "bulk" element and the face index.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
virtual double get_volume_flux()=0
Pure virtual function to calculate integral of the volume flux.
void add_pressure_data(Data *pressure_data_pt)
Function adds to the external data the Data object whose single value is the pressure applied by the ...
unsigned & pressure_data_id()
Access function gives id of external Data object whose single value is the pressure applied by the el...
unsigned Pressure_data_id
Id of external Data object whose single value is the pressure applied by the elements.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...