specific_node_update_interface_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 specific (two-dimensional) fluid free surface elements
27 
28 // Include guards, to prevent multiple includes
29 #ifndef OOMPH_SPECIFIC_NODE_UPDATE_INTERFACE_ELEMENTS_HEADER
30 #define OOMPH_SPECIFIC_NODE_UPDATE_INTERFACE_ELEMENTS_HEADER
31 
32 // Config header generated by autoconfig
33 #ifdef HAVE_CONFIG_H
34 #include <oomph-lib-config.h>
35 #endif
36 
37 // OOMPH-LIB headers
38 #include "../generic/Qelements.h"
39 #include "../generic/spines.h"
40 #include "../generic/hijacked_elements.h"
41 #include "interface_elements.h"
42 
43 namespace oomph
44 {
45  //=======================================================================
46  /// This policy class is used to associate specific bounding
47  /// elements with specific FluidInterface elements. It must be filled
48  /// in for every class that uses the SpineUpdateFluidInterface<...>
49  /// or ElasticUpdateFluidInterface<....> generic template classes.
50  /// Examples for our default Line, Axisymmetric and Surface types
51  /// are included below
52  //=======================================================================
53  template<class ELEMENT>
55  {
56  };
57 
58  //======================================================================
59  /// This policy class is used to allow additional values to be
60  /// added to the nodes from new surface equations, for examples of
61  /// usage see the SurfactantTransportFluidInterfaceElements.
62  /// The use of this class avoids issues with calling virtual
63  /// functions in constructors and avoids having a global look-up
64  /// able, although it functions in much the same way.
65  /// Typically, this will only be filled in by "expert users" and
66  /// is only required if you want to write generic surface-element
67  /// classes. Specific classes can always be overloaded on a
68  /// case-by-case basis.
69  //=====================================================================
70  template<class ELEMENT>
72  {
73  private:
74  // Issue an error message if this base class is called
76  {
77  std::ostringstream error_message;
79  << "The generic FluidInterfaceAdditionalValues<ELEMENT> class has "
80  "been\n"
81  << "called. This should never happen. If you are creating your own\n"
82  << "surface equations you must overload the policy class to specify\n"
83  << "how many additional values are required at the surface nodes.\n"
84  << "For an example, see "
85  "src/fluid_interface/surfactant_transport_elements.h\n"
86  << std::endl;
87  throw OomphLibError(
88  error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
89  }
90 
91  public:
92  /// Empty constructor
94 
95  /// Specific interface that states how many additional values are
96  /// required for the n-th node. Default is zero, but issue error_message.
97  inline unsigned nadditional_values(const unsigned& n)
98  {
99  error_message();
100  return 0;
101  }
102 
103  /// Specify any additional index setup information that is required;
104  /// i.e. the look-up schemes for the additional values.
105  /// Default is empty with error message
106  inline void setup_equation_indices(ELEMENT* const& element_pt,
107  const unsigned& id)
108  {
109  error_message();
110  }
111  };
112 
113 
114  //======================================================================
115  /// Specific policy class for the FluidInterfaceElemetnts,
116  /// which do not require any additional values at the nodes.
117  //=====================================================================
118  template<>
120  {
121  public:
122  /// Empty constructor
124 
125  /// Specific interface that states how many additional values are
126  /// required for the n-th node. No additional values
127  inline unsigned nadditional_values(const unsigned& n)
128  {
129  return 0;
130  }
131 
132  /// Specify any additional index setup information that is required;
133  /// i.e. the look-up schemes for the additional values.
134  /// Empty
135  inline void setup_equation_indices(FluidInterfaceElement* const& element_pt,
136  const unsigned& id)
137  {
138  }
139  };
140 
141 
142  //-------------SPINE NODE UPDATE CLASSES-------------------------------
143  //---------------------------------------------------------------------
144 
145  //======================================================================
146  /// Generic Spine node update interface template class that can
147  /// be combined with a given surface equations class and surface
148  /// derivative class to provide
149  /// a concrete implementation of any surface element that uses spines.
150  //======================================================================
151  template<class EQUATION_CLASS, class DERIVATIVE_CLASS, class ELEMENT>
153  : public virtual Hijacked<SpineElement<FaceGeometry<ELEMENT>>>,
154  public virtual EQUATION_CLASS,
155  public DERIVATIVE_CLASS
156  {
157  private:
158  /// In spine elements, the kinematic condition is the equation
159  /// used to determine the unknown spine heights. Overload the
160  /// function accordingly
161  int kinematic_local_eqn(const unsigned& n)
162  {
163  return this->spine_local_eqn(n);
164  }
165 
166  /// Hijacking the kinematic condition corresponds to hijacking the
167  /// variables associated with the spine heights.
168  void hijack_kinematic_conditions(const Vector<unsigned>& bulk_node_number)
169  {
170  // Loop over all the node numbers that are passed in
171  for (Vector<unsigned>::const_iterator it = bulk_node_number.begin();
172  it != bulk_node_number.end();
173  ++it)
174  {
175  // Hijack the spine heights. (and delete the returned data object)
176  delete this->hijack_nodal_spine_value(*it, 0);
177  }
178  }
179 
180  protected:
181  /// Fill in the specific surface derivative calculations
182  /// by calling the appropriate class function
184  const Shape& psi,
185  const DShape& dpsids,
186  const DenseMatrix<double>& interpolated_t,
187  const Vector<double>& interpolated_x,
188  DShape& surface_gradient,
189  DShape& surface_divergence)
190  {
191  return DERIVATIVE_CLASS::compute_surface_derivatives(psi,
192  dpsids,
193  interpolated_t,
194  interpolated_x,
195  surface_gradient,
196  surface_divergence);
197  }
198 
199 
200  public:
201  /// Constructor, the arguments are a pointer to the "bulk" element
202  /// and the index of the face to be created
203  SpineUpdateFluidInterfaceElement(FiniteElement* const& element_pt,
204  const int& face_index,
205  const unsigned& id = 0)
206  : Hijacked<SpineElement<FaceGeometry<ELEMENT>>>(),
207  EQUATION_CLASS(),
208  DERIVATIVE_CLASS()
209  {
210  // Attach the geometrical information to the element, by
211  // making the face element from the bulk element
212  element_pt->build_face_element(face_index, this);
213 
214 #ifdef PARANOID
215  // Is it refineable
216  RefineableElement* ref_el_pt =
217  dynamic_cast<RefineableElement*>(element_pt);
218  if (ref_el_pt != 0)
219  {
220  if (this->has_hanging_nodes())
221  {
222  throw OomphLibError("This interface element will not work correctly "
223  "if nodes are hanging\n",
224  OOMPH_CURRENT_FUNCTION,
225  OOMPH_EXCEPTION_LOCATION);
226  }
227  }
228 #endif
229 
230  // Find the index at which the velocity unknowns are stored
231  // from the bulk element and allocate the local storage
232  ELEMENT* cast_element_pt = dynamic_cast<ELEMENT*>(element_pt);
233  // Find number of momentum equation required
234  const unsigned n_u_index = cast_element_pt->n_u_nst();
235  this->U_index_interface.resize(n_u_index);
236  for (unsigned i = 0; i < n_u_index; i++)
237  {
238  this->U_index_interface[i] = cast_element_pt->u_index_nst(i);
239  }
240 
241  // Add any additional values required by the equations class
242  unsigned n_node_face = this->nnode();
243  // Create an instance of the policy class that determines
244  // how many additional values are required
246  interface_additional_values_pt =
248  // Do we need to add any values
249  // By default the spines don't need to so we can save
250  // some effort here
251  bool add_values = false;
252  // Storage for the number of additional values required
253  // for each node
254  Vector<unsigned> additional_data_values(n_node_face);
255  for (unsigned i = 0; i < n_node_face; ++i)
256  {
257  // Read out additional values for each node
258  additional_data_values[i] =
259  interface_additional_values_pt->nadditional_values(i);
260  // If any of these are non-zero it will set the boolean to true
261  add_values |= additional_data_values[i];
262  }
263 
264  // If storage is needed allocate it and call
265  // the associated local index setup function
266  if (add_values)
267  {
268  this->add_additional_values(additional_data_values, id);
269  // Call any local setup from the interface policy class
270  interface_additional_values_pt->setup_equation_indices(this, id);
271  }
272 
273  // Delete the policy class, it's work is done.
274  delete interface_additional_values_pt;
275  interface_additional_values_pt = 0;
276  }
277 
278  /// Calculate the contribution to the residuals and the jacobian
279  void fill_in_contribution_to_jacobian(Vector<double>& residuals,
280  DenseMatrix<double>& jacobian)
281  {
282  // Call the generic routine with the flag set to 1
283  EQUATION_CLASS::fill_in_generic_residual_contribution_interface(
284  residuals, jacobian, 1);
285 
286  // Call the generic routine to handle the shape derivatives
287  this->fill_in_jacobian_from_geometric_data(jacobian);
288  }
289 
290  /// \short
291  /// Helper function to calculate the additional contributions
292  /// These are those filled in by the particular equations
294  Vector<double>& residuals,
295  DenseMatrix<double>& jacobian,
296  const unsigned& flag,
297  const Shape& psif,
298  const DShape& dpsifds,
299  const DShape& dpsifdS,
300  const DShape& dpsifdS_div,
301  const Vector<double>& s,
302  const Vector<double>& interpolated_x,
303  const Vector<double>& interpolated_n,
304  const double& W,
305  const double& J)
306  {
307  EQUATION_CLASS::add_additional_residual_contributions_interface(
308  residuals,
309  jacobian,
310  flag,
311  psif,
312  dpsifds,
313  dpsifdS,
314  dpsifdS_div,
315  s,
316  interpolated_x,
317  interpolated_n,
318  W,
319  J);
320  }
321 
322  /// Overload the output function
323  void output(std::ostream& outfile)
324  {
325  EQUATION_CLASS::output(outfile);
326  }
327 
328  /// Output the element
329  void output(std::ostream& outfile, const unsigned& n_plot)
330  {
331  EQUATION_CLASS::output(outfile, n_plot);
332  }
333 
334  /// Overload the C-style output function
335  void output(FILE* file_pt)
336  {
337  EQUATION_CLASS::output(file_pt);
338  }
339 
340  /// C-style Output function
341  void output(FILE* file_pt, const unsigned& n_plot)
342  {
343  EQUATION_CLASS::output(file_pt, n_plot);
344  }
345 
346 
347  /// Create an "bounding" element of the type
348  /// specified by the BoundingElementType policy class
349  /// Here, this allows
350  /// the application of a contact angle boundary condition on the
351  /// the specified face.
353  const int& face_index)
354  {
355  // Create a temporary pointer to the appropriate FaceElement read our from
356  // our policy class
358  DERIVATIVE_CLASS,
359  ELEMENT>>*
360  face_el_pt = new BoundingElementType<
361  SpineUpdateFluidInterfaceElement<EQUATION_CLASS,
362  DERIVATIVE_CLASS,
363  ELEMENT>>;
364 
365  // Attach the geometrical information to the new element
366  this->build_face_element(face_index, face_el_pt);
367 
368  // Set the index at which the velocity nodes are stored
369  face_el_pt->u_index_interface_boundary() = this->U_index_interface;
370 
371  // Set the value of the nbulk_value, the node is not resized
372  // in this bounding element,
373  // so it will just be the actual nvalue here
374  // There is some ambiguity about what this means (however)
375  // We are interpreting it to mean the number of
376  // values in this FaceElement before creating the new
377  // bounding element.
378  const unsigned n_node_bounding = face_el_pt->nnode();
379  for (unsigned n = 0; n < n_node_bounding; n++)
380  {
381  face_el_pt->nbulk_value(n) = face_el_pt->node_pt(n)->nvalue();
382  }
383 
384  // Set of unique geometric data that is used to update the bulk,
385  // but is not used to update the face
386  std::set<Data*> unique_additional_geom_data;
387 
388  // Get all the geometric data for this (bulk) element
389  this->assemble_set_of_all_geometric_data(unique_additional_geom_data);
390 
391  // Now assemble the set of geometric data for the face element
392  std::set<Data*> unique_face_geom_data_pt;
393  face_el_pt->assemble_set_of_all_geometric_data(unique_face_geom_data_pt);
394 
395  // Erase the face geometric data from the additional data
396  for (std::set<Data*>::iterator it = unique_face_geom_data_pt.begin();
397  it != unique_face_geom_data_pt.end();
398  ++it)
399  {
400  unique_additional_geom_data.erase(*it);
401  }
402 
403  // Finally add all unique additional data as external data
404  for (std::set<Data*>::iterator it = unique_additional_geom_data.begin();
405  it != unique_additional_geom_data.end();
406  ++it)
407  {
408  face_el_pt->add_external_data(*it);
409  }
410 
411  // Return the value of the pointer
412  return face_el_pt;
413  }
414  };
415 
416 
417  //=====================================================================
418  /// Spine version of the PointFluidInterfaceBoundingElement
419  //=====================================================================
420  template<class ELEMENT>
422  : public Hijacked<SpineElement<FaceGeometry<FaceGeometry<ELEMENT>>>>,
424 
425  {
426  public:
427  /// Constructor
429  : Hijacked<SpineElement<FaceGeometry<FaceGeometry<ELEMENT>>>>(),
431  {
432  }
433 
434  /// Overload the output function
435  void output(std::ostream& outfile)
436  {
437  FiniteElement::output(outfile);
438  }
439 
440  /// Output the element
441  void output(std::ostream& outfile, const unsigned& n_plot)
442  {
443  FluidInterfaceBoundingElement::output(outfile, n_plot);
444  }
445 
446  /// Overload the C-style output function
447  void output(FILE* file_pt)
448  {
449  FiniteElement::output(file_pt);
450  }
451 
452  /// C-style Output function
453  void output(FILE* file_pt, const unsigned& n_plot)
454  {
455  FluidInterfaceBoundingElement::output(file_pt, n_plot);
456  }
457 
458  /// Calculate the elemental residual vector and the Jacobian
459  void fill_in_contribution_to_jacobian(Vector<double>& residuals,
460  DenseMatrix<double>& jacobian)
461  {
462  // Call the generic routine with the flag set to 1
464  residuals, jacobian, 1);
465  // Call generic FD routine for the external data
466  this->fill_in_jacobian_from_external_by_fd(jacobian);
467 
468  // Call the generic routine to handle the spine variables
469  this->fill_in_jacobian_from_geometric_data(jacobian);
470  }
471 
472  /// Return local equation number associated with the kinematic
473  /// constraint for local node n
474  int kinematic_local_eqn(const unsigned& n)
475  {
476  return this->spine_local_eqn(n);
477  }
478  };
479 
480 
481  //=========================================================================
482  /// Spine version of the LineFluidInterfaceBoundingElement
483  //========================================================================
484  template<class ELEMENT>
486  : public Hijacked<SpineElement<FaceGeometry<FaceGeometry<ELEMENT>>>>,
488 
489  {
490  public:
491  /// Constructor
493  : Hijacked<SpineElement<FaceGeometry<FaceGeometry<ELEMENT>>>>(),
495  {
496  }
497 
498  /// Overload the output function
499  void output(std::ostream& outfile)
500  {
501  FiniteElement::output(outfile);
502  }
503 
504  /// Output the element
505  void output(std::ostream& outfile, const unsigned& n_plot)
506  {
507  FluidInterfaceBoundingElement::output(outfile, n_plot);
508  }
509 
510  /// Overload the C-style output function
511  void output(FILE* file_pt)
512  {
513  FiniteElement::output(file_pt);
514  }
515 
516  /// C-style Output function
517  void output(FILE* file_pt, const unsigned& n_plot)
518  {
519  FluidInterfaceBoundingElement::output(file_pt, n_plot);
520  }
521 
522 
523  /// Calculate the jacobian
524  void fill_in_contribution_to_jacobian(Vector<double>& residuals,
525  DenseMatrix<double>& jacobian)
526  {
527  // Call the generic routine with the flag set to 1
529  residuals, jacobian, 1);
530  // Call generic FD routine for the external data
531  this->fill_in_jacobian_from_external_by_fd(jacobian);
532 
533  // Call the generic routine to handle the spine variables
534  this->fill_in_jacobian_from_geometric_data(jacobian);
535  }
536 
537 
538  /// Local eqn number of the kinematic bc associated with local node n
539  int kinematic_local_eqn(const unsigned& n)
540  {
541  // Kinematic bc is always associated with the n-th spine height
542  return this->spine_local_eqn(n);
543  }
544  };
545 
546 
547  //============GEOMETRIC SPECIALISATIONS============================
548 
549 
550  // Specialise the spine update template class to concrete 1D case
551  template<class ELEMENT>
553  : public SpineUpdateFluidInterfaceElement<FluidInterfaceElement,
554  LineDerivatives,
555  ELEMENT>
556  {
557  public:
558  SpineLineFluidInterfaceElement(FiniteElement* const& element_pt,
559  const int& face_index)
562  ELEMENT>(element_pt, face_index)
563  {
564  }
565  };
566 
567 
568  // Define the BoundingElement type associated with the 1D surface element
569  template<class ELEMENT>
573  ELEMENT>>
575  {
576  public:
579  ELEMENT>>()
581  {
582  }
583  };
584 
585 
586  // Specialise Spine update case to concrete axisymmetric case
587  template<class ELEMENT>
589  : public SpineUpdateFluidInterfaceElement<FluidInterfaceElement,
590  AxisymmetricDerivatives,
591  ELEMENT>
592  {
593  public:
594  SpineAxisymmetricFluidInterfaceElement(FiniteElement* const& element_pt,
595  const int& face_index)
598  ELEMENT>(element_pt, face_index)
599  {
600  }
601  };
602 
603  // Define the bounding element associated with the axisymmetric fluid
604  // interface element
605  template<class ELEMENT>
609  ELEMENT>>
611  {
612  public:
616  ELEMENT>>()
618  {
619  }
620  };
621 
622  // Specialise Spine update case to concrete 2D case
623  template<class ELEMENT>
625  : public SpineUpdateFluidInterfaceElement<FluidInterfaceElement,
626  SurfaceDerivatives,
627  ELEMENT>
628  {
629  public:
630  SpineSurfaceFluidInterfaceElement(FiniteElement* const& element_pt,
631  const int& face_index)
634  ELEMENT>(element_pt, face_index)
635  {
636  }
637  };
638 
639  // Define the bounding element type for the 2D surface
640  template<class ELEMENT>
644  ELEMENT>>
646  {
647  public:
650  ELEMENT>>()
652  {
653  }
654  };
655 
656 
657  /// ////////////////////////////////////////////////////////////////////
658  /// ////////////////////////////////////////////////////////////////////
659  /// ////////////////////////////////////////////////////////////////////
660 
661  //-------------ELASTIC NODE UPDATE CLASSES-------------------------------
662  //---------------------------------------------------------------------
663 
664 
665  //=======================================================================
666  /// Generic Elastic node update interface template class that can
667  /// be combined with a given surface equations class and surface derivative
668  /// class to provide a concrete implementation of any surface element
669  /// that uses elastic node updates.
670  //======================================================================
671  template<class EQUATION_CLASS, class DERIVATIVE_CLASS, class ELEMENT>
673  : public virtual Hijacked<FaceGeometry<ELEMENT>>,
674  public EQUATION_CLASS,
675  public DERIVATIVE_CLASS
676  {
677  private:
678  /// Storage for the location of the Lagrange multiplier
679  /// (If other additional values have been added we need
680  /// to add the Lagrange multiplier at the end)
681  Vector<unsigned> Lagrange_index;
682 
683  /// Return the index at which the lagrange multiplier is
684  /// stored at the n-th node
685  inline unsigned lagrange_index(const unsigned& n)
686  {
687  return this->Lagrange_index[n];
688  }
689 
690  /// Equation number of the kinematic BC associated with node j.
691  /// (This is the equation for the Lagrange multiplier)
692  inline int kinematic_local_eqn(const unsigned& n)
693  {
694  // Get the index of the nodal value associated with Lagrange multiplier
695  return this->nodal_local_eqn(n, this->lagrange_index(n));
696  }
697 
698  /// Hijacking the kinematic condition corresponds to hijacking the
699  /// variables associated with the Lagrange multipliers that are assigned
700  /// on construction of this element.
701  void hijack_kinematic_conditions(const Vector<unsigned>& bulk_node_number)
702  {
703  // Loop over all the nodes that are passed in
704  for (Vector<unsigned>::const_iterator it = bulk_node_number.begin();
705  it != bulk_node_number.end();
706  ++it)
707  {
708  // Hijack the appropriate value and delete the returned Node
709  delete this->hijack_nodal_value(*it, this->lagrange_index(*it));
710  }
711  }
712 
713  protected:
714  /// Fill in the specific surface derivative calculations
715  /// by calling the appropriate function from the derivative class
717  const Shape& psi,
718  const DShape& dpsids,
719  const DenseMatrix<double>& interpolated_t,
720  const Vector<double>& interpolated_x,
721  DShape& surface_gradient,
722  DShape& surface_divergence)
723  {
724  return DERIVATIVE_CLASS::compute_surface_derivatives(psi,
725  dpsids,
726  interpolated_t,
727  interpolated_x,
728  surface_gradient,
729  surface_divergence);
730  }
731 
732 
733  public:
734  /// Constructor, pass a pointer to the bulk element and the face
735  /// index of the bulk element to which the element is to be attached to.
736  /// The optional identifier can be used
737  /// to distinguish the additional nodal value (Lagr mult) created by
738  /// this element from those created by other FaceElements.
739  ElasticUpdateFluidInterfaceElement(FiniteElement* const& element_pt,
740  const int& face_index,
741  const unsigned& id = 0)
742  : FaceGeometry<ELEMENT>(), EQUATION_CLASS(), DERIVATIVE_CLASS()
743  {
744  // Attach the geometrical information to the element
745  // This function also assigned nbulk_value from required_nvalue of the
746  // bulk element
747  element_pt->build_face_element(face_index, this);
748 
749 #ifdef PARANOID
750  // Is it refineable
751  RefineableElement* ref_el_pt =
752  dynamic_cast<RefineableElement*>(element_pt);
753  if (ref_el_pt != 0)
754  {
755  if (this->has_hanging_nodes())
756  {
757  throw OomphLibError(
758  "This flux element will not work correctly if nodes are hanging\n",
759  OOMPH_CURRENT_FUNCTION,
760  OOMPH_EXCEPTION_LOCATION);
761  }
762  }
763 #endif
764 
765  // Find the index at which the velocity unknowns are stored
766  // from the bulk element and resize the local storage scheme
767  ELEMENT* cast_element_pt = dynamic_cast<ELEMENT*>(element_pt);
768  const unsigned n_u_index = cast_element_pt->n_u_nst();
769  this->U_index_interface.resize(n_u_index);
770  for (unsigned i = 0; i < n_u_index; i++)
771  {
772  this->U_index_interface[i] = cast_element_pt->u_index_nst(i);
773  }
774 
775  // Read out the number of nodes on the face
776  unsigned n_node_face = this->nnode();
777 
778  // Create an instance of the policy class that determines
779  // how many additional values are required
781  interface_additional_values_pt =
783 
784  // Set the additional data values in the face
785  // There is always also one additional values at each node --- the
786  // Lagrange multiplier
787  Vector<unsigned> additional_data_values(n_node_face);
788  for (unsigned n = 0; n < n_node_face; n++)
789  {
790  // Now add one to the addtional values at every single node
791  additional_data_values[n] =
792  interface_additional_values_pt->nadditional_values(n) + 1;
793  }
794 
795  // Now add storage for Lagrange multipliers and set the map containing
796  // the position of the first entry of this face element's
797  // additional values.
798  this->add_additional_values(additional_data_values, id);
799 
800  // Now I can just store the lagrange index offset to give the storage
801  // location of the nodes
802  Lagrange_index.resize(n_node_face);
803  for (unsigned n = 0; n < n_node_face; ++n)
804  {
805  Lagrange_index[n] =
806  additional_data_values[n] - 1 +
807  dynamic_cast<BoundaryNodeBase*>(this->node_pt(n))
808  ->index_of_first_value_assigned_by_face_element(id);
809  }
810 
811  // Call any local setup from the interface policy class
812  interface_additional_values_pt->setup_equation_indices(this, id);
813 
814  // Can now delete the policy class
815  delete interface_additional_values_pt;
816  interface_additional_values_pt = 0;
817  }
818 
819  /// The "global" intrinsic coordinate of the element when
820  /// viewed as part of a geometric object should be given by
821  /// the FaceElement representation, by default
822  double zeta_nodal(const unsigned& n,
823  const unsigned& k,
824  const unsigned& i) const
825  {
826  return FaceElement::zeta_nodal(n, k, i);
827  }
828 
829  /// Return the lagrange multiplier at local node n
830  double& lagrange(const unsigned& n)
831  {
832  return *this->node_pt(n)->value_pt(this->lagrange_index(n));
833  }
834 
835 
836  /// Fill in contribution to residuals and Jacobian
837  void fill_in_contribution_to_jacobian(Vector<double>& residuals,
838  DenseMatrix<double>& jacobian)
839  {
840  // Call the generic routine with the flag set to 1
841  EQUATION_CLASS::fill_in_generic_residual_contribution_interface(
842  residuals, jacobian, 1);
843 
844  // Call the generic finite difference routine for the solid variables
845  this->fill_in_jacobian_from_solid_position_by_fd(jacobian);
846  }
847 
848  /// Overload the output function
849  void output(std::ostream& outfile)
850  {
851  EQUATION_CLASS::output(outfile);
852  }
853 
854  /// Output the element
855  void output(std::ostream& outfile, const unsigned& n_plot)
856  {
857  EQUATION_CLASS::output(outfile, n_plot);
858  }
859 
860  /// Overload the C-style output function
861  void output(FILE* file_pt)
862  {
863  EQUATION_CLASS::output(file_pt);
864  }
865 
866  /// C-style Output function
867  void output(FILE* file_pt, const unsigned& n_plot)
868  {
869  EQUATION_CLASS::output(file_pt, n_plot);
870  }
871 
872 
873  /// Helper function to calculate the additional contributions
874  /// to be added at each integration point. This deals with
875  /// Lagrange multiplier contribution, as well as any
876  /// additional contributions by the other equations.
878  Vector<double>& residuals,
879  DenseMatrix<double>& jacobian,
880  const unsigned& flag,
881  const Shape& psif,
882  const DShape& dpsifds,
883  const DShape& dpsifdS,
884  const DShape& dpsifdS_div,
885  const Vector<double>& s,
886  const Vector<double>& interpolated_x,
887  const Vector<double>& interpolated_n,
888  const double& W,
889  const double& J)
890  {
891  EQUATION_CLASS::add_additional_residual_contributions_interface(
892  residuals,
893  jacobian,
894  flag,
895  psif,
896  dpsifds,
897  dpsifdS,
898  dpsifdS_div,
899  s,
900  interpolated_x,
901  interpolated_n,
902  W,
903  J);
904 
905  // Assemble Lagrange multiplier by loop over the shape functions
906  const unsigned n_node = this->nnode();
907  // Read out the dimension of the node
908  const unsigned nodal_dimension = this->nodal_dimension();
909 
910  double interpolated_lagrange = 0.0;
911  for (unsigned l = 0; l < n_node; l++)
912  {
913  // Note same shape functions used for lagrange multiplier field
914  interpolated_lagrange += lagrange(l) * psif(l);
915  }
916 
917  int local_eqn = 0, local_unknown = 0;
918 
919  // Loop over the shape functions to assemble contributions
920  for (unsigned l = 0; l < n_node; l++)
921  {
922  // Loop over the directions
923  for (unsigned i = 0; i < nodal_dimension; i++)
924  {
925  // Now using the same shape functions for the elastic equations,
926  // so we can stay in the loop
927  local_eqn = this->position_local_eqn(l, 0, i);
928  if (local_eqn >= 0)
929  {
930  // Add in the Lagrange multiplier contribution
931  residuals[local_eqn] -=
932  interpolated_lagrange * interpolated_n[i] * psif(l) * J * W;
933 
934  // Do the Jacobian calculation
935  if (flag)
936  {
937  // Loop over the nodes
938  for (unsigned l2 = 0; l2 < n_node; l2++)
939  {
940  // Dependence on solid positions will be handled by FDs
941  // That leaves the Lagrange multiplier contribution
942  local_unknown = this->kinematic_local_eqn(l2);
943  if (local_unknown >= 0)
944  {
945  jacobian(local_eqn, local_unknown) -=
946  psif(l2) * interpolated_n[i] * psif(l) * J * W;
947  }
948  }
949  } // End of Jacobian calculation
950  }
951  }
952  } // End of loop over shape functions
953  }
954 
955 
956  /// Create an "bounding" element (here actually a 2D line element
957  /// of type ElasticLineFluidInterfaceBoundingElement<ELEMENT> that allows
958  /// the application of a contact angle boundary condition on the
959  /// the specified face.
961  const int& face_index)
962  {
963  // Create a temporary pointer to the appropriate FaceElement
965  DERIVATIVE_CLASS,
966  ELEMENT>>*
967  face_el_pt = new BoundingElementType<
968  ElasticUpdateFluidInterfaceElement<EQUATION_CLASS,
969  DERIVATIVE_CLASS,
970  ELEMENT>>;
971 
972  // Attach the geometrical information to the new element
973  this->build_face_element(face_index, face_el_pt);
974 
975  // Set the index at which the velocity nodes are stored
976  face_el_pt->u_index_interface_boundary() = this->U_index_interface;
977 
978  // Set the value of the nbulk_value, the node is not resized
979  // in this bounding element,
980  // so it will just be the actual nvalue here
981  // There is some ambiguity about what this means (however)
982  const unsigned n_node_bounding = face_el_pt->nnode();
983  // Storage for lagrange multiplier index at the face nodes
984  Vector<unsigned> local_lagrange_index(n_node_bounding);
985  for (unsigned n = 0; n < n_node_bounding; n++)
986  {
987  face_el_pt->nbulk_value(n) = face_el_pt->node_pt(n)->nvalue();
988  // At the moment the assumption is that it is stored at all nodes, but
989  // that is consistent with the assumption in this element
990  local_lagrange_index[n] =
991  this->Lagrange_index[face_el_pt->bulk_node_number(n)];
992  }
993 
994  // Pass the ID and offset down
995  face_el_pt->set_lagrange_index(local_lagrange_index);
996 
997  // Find the nodes
998  std::set<SolidNode*> set_of_solid_nodes;
999  const unsigned n_node = this->nnode();
1000  for (unsigned n = 0; n < n_node; n++)
1001  {
1002  set_of_solid_nodes.insert(static_cast<SolidNode*>(this->node_pt(n)));
1003  }
1004 
1005  // Delete the nodes from the face
1006  // n_node = face_el_pt->nnode();
1007  for (unsigned n = 0; n < n_node_bounding; n++)
1008  {
1009  // Set the value of the nbulk_value, from the present element
1010  face_el_pt->nbulk_value(n) =
1011  this->nbulk_value(face_el_pt->bulk_node_number(n));
1012 
1013  // Now delete the nodes from the set
1014  set_of_solid_nodes.erase(
1015  static_cast<SolidNode*>(face_el_pt->node_pt(n)));
1016  }
1017 
1018  // Now add these as external data
1019  for (std::set<SolidNode*>::iterator it = set_of_solid_nodes.begin();
1020  it != set_of_solid_nodes.end();
1021  ++it)
1022  {
1023  face_el_pt->add_external_data((*it)->variable_position_pt());
1024  }
1025 
1026 
1027  // Return the value of the pointer
1028  return face_el_pt;
1029  }
1030  };
1031 
1032 
1033  //=========================================================================
1034  /// Pseudo-elasticity version of the PointFluidInterfaceBoundingElement
1035  //========================================================================
1036  template<class ELEMENT>
1038  : public FaceGeometry<FaceGeometry<ELEMENT>>,
1040  public virtual SolidFiniteElement
1041 
1042  {
1043  private:
1044  /// Short Storage for the index of the Lagrange multiplier at the chosen
1045  /// nodes
1046  Vector<unsigned> Lagrange_index;
1047 
1048  public:
1049  /// Set the Id and offset
1050  void set_lagrange_index(const Vector<unsigned>& lagrange_index)
1051  {
1052  Lagrange_index = lagrange_index;
1053  }
1054 
1055 
1056  /// Specify the value of nodal zeta from the face geometry
1057  /// The "global" intrinsic coordinate of the element when
1058  /// viewed as part of a geometric object should be given by
1059  /// the FaceElement representation, by default
1060  double zeta_nodal(const unsigned& n,
1061  const unsigned& k,
1062  const unsigned& i) const
1063  {
1064  return FaceElement::zeta_nodal(n, k, i);
1065  }
1066 
1067 
1068  /// Constructor
1070  : FaceGeometry<FaceGeometry<ELEMENT>>(),
1072  {
1073  }
1074 
1075  /// Overload the output function
1076  void output(std::ostream& outfile)
1077  {
1078  FiniteElement::output(outfile);
1079  }
1080 
1081  /// Output the element
1082  void output(std::ostream& outfile, const unsigned& n_plot)
1083  {
1084  FluidInterfaceBoundingElement::output(outfile, n_plot);
1085  }
1086 
1087  /// Overload the C-style output function
1088  void output(FILE* file_pt)
1089  {
1090  FiniteElement::output(file_pt);
1091  }
1092 
1093  /// C-style Output function
1094  void output(FILE* file_pt, const unsigned& n_plot)
1095  {
1096  FluidInterfaceBoundingElement::output(file_pt, n_plot);
1097  }
1098 
1099  /// Calculate the element's residual vector and Jacobian
1100  void fill_in_contribution_to_jacobian(Vector<double>& residuals,
1101  DenseMatrix<double>& jacobian)
1102  {
1103  // Call the generic routine with the flag set to 1
1105  residuals, jacobian, 1);
1106  // Call the generic FD routine to get external data
1107  this->fill_in_jacobian_from_external_by_fd(jacobian);
1108 
1109  // Call the generic finite difference routine to handle the solid
1110  // variables
1111  this->fill_in_jacobian_from_solid_position_by_fd(jacobian);
1112  }
1113 
1114  /// Set the kinematic local equation
1115  inline int kinematic_local_eqn(const unsigned& n)
1116  {
1117  return this->nodal_local_eqn(n, this->Lagrange_index[n]);
1118  }
1119  };
1120 
1121 
1122  //=========================================================================
1123  /// Pseudo-elasticity version of the LineFluidInterfaceBoundingElement
1124  //========================================================================
1125  template<class ELEMENT>
1127  : public FaceGeometry<FaceGeometry<ELEMENT>>,
1129  public virtual SolidFiniteElement
1130 
1131  {
1132  /// Short Storage for the index of Lagrange multiplier
1133  Vector<unsigned> Lagrange_index;
1134 
1135  public:
1136  /// Set the Id
1137  void set_lagrange_index(const Vector<unsigned>& lagrange_index)
1138  {
1139  Lagrange_index = lagrange_index;
1140  }
1141 
1142  /// Constructor
1144  : FaceGeometry<FaceGeometry<ELEMENT>>(),
1146  {
1147  }
1148 
1149  /// Specify the value of nodal zeta from the face geometry
1150  /// The "global" intrinsic coordinate of the element when
1151  /// viewed as part of a geometric object should be given by
1152  /// the FaceElement representation, by default
1153  double zeta_nodal(const unsigned& n,
1154  const unsigned& k,
1155  const unsigned& i) const
1156  {
1157  return FaceElement::zeta_nodal(n, k, i);
1158  }
1159 
1160 
1161  /// Overload the output function
1162  void output(std::ostream& outfile)
1163  {
1164  FiniteElement::output(outfile);
1165  }
1166 
1167  /// Output the element
1168  void output(std::ostream& outfile, const unsigned& n_plot)
1169  {
1170  FluidInterfaceBoundingElement::output(outfile, n_plot);
1171  }
1172 
1173  /// Overload the C-style output function
1174  void output(FILE* file_pt)
1175  {
1176  FiniteElement::output(file_pt);
1177  }
1178 
1179  /// C-style Output function
1180  void output(FILE* file_pt, const unsigned& n_plot)
1181  {
1182  FluidInterfaceBoundingElement::output(file_pt, n_plot);
1183  }
1184 
1185  /// Calculate the elemental residual vector and Jacobian
1186  void fill_in_contribution_to_jacobian(Vector<double>& residuals,
1187  DenseMatrix<double>& jacobian)
1188  {
1189  // Call the generic routine with the flag set to 1
1191  residuals, jacobian, 1);
1192 
1193  // Call the generic FD routine to get externals
1194  this->fill_in_jacobian_from_external_by_fd(jacobian);
1195 
1196  // Call the generic finite difference routine to handle the solid
1197  // variables
1198  this->fill_in_jacobian_from_solid_position_by_fd(jacobian);
1199  }
1200 
1201  /// Local eqn number of kinematic bc associated with local node n
1202  int kinematic_local_eqn(const unsigned& n)
1203  {
1204  // Read out the kinematic constraint from the Id which is passed down
1205  // from the constructing element
1206  return this->nodal_local_eqn(n, this->Lagrange_index[n]);
1207  }
1208  };
1209 
1210 
1211  //==================GEOMETRIC SPECIALISATIONS==========================
1212 
1213 
1214  /// Specialise the elastic update template class to concrete 1D case
1215  template<class ELEMENT>
1217  : public ElasticUpdateFluidInterfaceElement<FluidInterfaceElement,
1218  LineDerivatives,
1219  ELEMENT>
1220  {
1221  public:
1222  ElasticLineFluidInterfaceElement(FiniteElement* const& element_pt,
1223  const int& face_index,
1224  const unsigned& id = 0)
1227  ELEMENT>(element_pt, face_index, id)
1228  {
1229  }
1230  };
1231 
1232  /// Define the BoundingElement type associated with the 1D surface element
1233  template<class ELEMENT>
1237  ELEMENT>>
1239  {
1240  public:
1244  ELEMENT>>()
1246  {
1247  }
1248  };
1249 
1250 
1251  /// Specialise the Elastic update case to axisymmetric equations
1252  template<class ELEMENT>
1254  : public ElasticUpdateFluidInterfaceElement<FluidInterfaceElement,
1255  AxisymmetricDerivatives,
1256  ELEMENT>
1257  {
1258  public:
1259  ElasticAxisymmetricFluidInterfaceElement(FiniteElement* const& element_pt,
1260  const int& face_index,
1261  const unsigned& id = 0)
1264  ELEMENT>(element_pt, face_index, id)
1265  {
1266  }
1267  };
1268 
1269  // Define the bounding element associated with the axsymmetric elastic fluid
1270  // interface element
1271  template<class ELEMENT>
1275  ELEMENT>>
1277  {
1278  public:
1282  ELEMENT>>()
1284  {
1285  }
1286  };
1287 
1288 
1289  /// Specialise Elastic update case to the concrete 2D case
1290  template<class ELEMENT>
1292  : public ElasticUpdateFluidInterfaceElement<FluidInterfaceElement,
1293  SurfaceDerivatives,
1294  ELEMENT>
1295  {
1296  public:
1297  ElasticSurfaceFluidInterfaceElement(FiniteElement* const& element_pt,
1298  const int& face_index,
1299  const unsigned& id = 0)
1302  ELEMENT>(element_pt, face_index, id)
1303  {
1304  }
1305  };
1306 
1307 
1308  // Define the bounding element associated with the 2D surface elements
1309  template<class ELEMENT>
1313  ELEMENT>>
1315  {
1316  public:
1320  ELEMENT>>()
1322  {
1323  }
1324  };
1325 
1326 
1327 } // namespace oomph
1328 
1329 #endif
Class that establishes the surface derivative functions for AxisymmetricInterfaceElements....
This policy class is used to associate specific bounding elements with specific FluidInterface elemen...
Specialise the Elastic update case to axisymmetric equations.
ElasticAxisymmetricFluidInterfaceElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0)
Pseudo-elasticity version of the LineFluidInterfaceBoundingElement.
void output(FILE *file_pt)
Overload the C-style output function.
int kinematic_local_eqn(const unsigned &n)
Local eqn number of kinematic bc associated with local node n.
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Specify the value of nodal zeta from the face geometry The "global" intrinsic coordinate of the eleme...
void output(std::ostream &outfile)
Overload the output function.
void set_lagrange_index(const Vector< unsigned > &lagrange_index)
Set the Id.
void output(FILE *file_pt, const unsigned &n_plot)
C-style Output function.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Calculate the elemental residual vector and Jacobian.
void output(std::ostream &outfile, const unsigned &n_plot)
Output the element.
Vector< unsigned > Lagrange_index
Short Storage for the index of Lagrange multiplier.
Specialise the elastic update template class to concrete 1D case.
ElasticLineFluidInterfaceElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0)
Pseudo-elasticity version of the PointFluidInterfaceBoundingElement.
void output(std::ostream &outfile, const unsigned &n_plot)
Output the element.
void output(std::ostream &outfile)
Overload the output function.
void output(FILE *file_pt, const unsigned &n_plot)
C-style Output function.
void output(FILE *file_pt)
Overload the C-style output function.
Vector< unsigned > Lagrange_index
Short Storage for the index of the Lagrange multiplier at the chosen nodes.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Calculate the element's residual vector and Jacobian.
void set_lagrange_index(const Vector< unsigned > &lagrange_index)
Set the Id and offset.
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Specify the value of nodal zeta from the face geometry The "global" intrinsic coordinate of the eleme...
int kinematic_local_eqn(const unsigned &n)
Set the kinematic local equation.
Specialise Elastic update case to the concrete 2D case.
ElasticSurfaceFluidInterfaceElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0)
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
void output(FILE *file_pt, const unsigned &n_plot)
C-style Output function.
void add_additional_residual_contributions_interface(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag, const Shape &psif, const DShape &dpsifds, const DShape &dpsifdS, const DShape &dpsifdS_div, const Vector< double > &s, const Vector< double > &interpolated_x, const Vector< double > &interpolated_n, const double &W, const double &J)
Helper function to calculate the additional contributions to be added at each integration point....
unsigned lagrange_index(const unsigned &n)
Return the index at which the lagrange multiplier is stored at the n-th node.
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 ...
void output(std::ostream &outfile)
Overload the output function.
ElasticUpdateFluidInterfaceElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0)
Constructor, pass a pointer to the bulk element and the face index of the bulk element to which the e...
Vector< unsigned > Lagrange_index
Storage for the location of the Lagrange multiplier (If other additional values have been added we ne...
double & lagrange(const unsigned &n)
Return the lagrange multiplier at local node n.
void output(FILE *file_pt)
Overload the C-style output function.
virtual FluidInterfaceBoundingElement * make_bounding_element(const int &face_index)
Create an "bounding" element (here actually a 2D line element of type ElasticLineFluidInterfaceBoundi...
int kinematic_local_eqn(const unsigned &n)
Equation number of the kinematic BC associated with node j. (This is the equation for the Lagrange mu...
void output(std::ostream &outfile, const unsigned &n_plot)
Output the element.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution to residuals and Jacobian.
void hijack_kinematic_conditions(const Vector< unsigned > &bulk_node_number)
Hijacking the kinematic condition corresponds to hijacking the variables associated with the Lagrange...
double compute_surface_derivatives(const Shape &psi, const DShape &dpsids, const DenseMatrix< double > &interpolated_t, const Vector< double > &interpolated_x, DShape &surface_gradient, DShape &surface_divergence)
Fill in the specific surface derivative calculations by calling the appropriate function from the der...
Specific policy class for the FluidInterfaceElemetnts, which do not require any additional values at ...
void setup_equation_indices(FluidInterfaceElement *const &element_pt, const unsigned &id)
Specify any additional index setup information that is required; i.e. the look-up schemes for the add...
unsigned nadditional_values(const unsigned &n)
Specific interface that states how many additional values are required for the n-th node....
This policy class is used to allow additional values to be added to the nodes from new surface equati...
unsigned nadditional_values(const unsigned &n)
Specific interface that states how many additional values are required for the n-th node....
void setup_equation_indices(ELEMENT *const &element_pt, const unsigned &id)
Specify any additional index setup information that is required; i.e. the look-up schemes for the add...
Base class for elements at the boundary of free surfaces or interfaces, used typically to impose cont...
void output(std::ostream &outfile)
Overload the output function.
Base class establishing common interfaces and functions for all Navier-Stokes-like fluid interface el...
Class that establishes the surface derivative functions for LineElements. These are defined in a sepa...
Specialisation of the interface boundary constraint to a line.
void fill_in_generic_residual_contribution_interface_boundary(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Overload the helper function to calculate the residuals and (if flag==true) the Jacobian – this funct...
Specialisation of the interface boundary constraint to a point.
void fill_in_generic_residual_contribution_interface_boundary(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Overload the helper function to calculate the residuals and (if flag==1) the Jacobian – this function...
SpineAxisymmetricFluidInterfaceElement(FiniteElement *const &element_pt, const int &face_index)
Spine version of the LineFluidInterfaceBoundingElement.
int kinematic_local_eqn(const unsigned &n)
Local eqn number of the kinematic bc associated with local node n.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Calculate the jacobian.
void output(std::ostream &outfile)
Overload the output function.
void output(std::ostream &outfile, const unsigned &n_plot)
Output the element.
void output(FILE *file_pt, const unsigned &n_plot)
C-style Output function.
void output(FILE *file_pt)
Overload the C-style output function.
SpineLineFluidInterfaceElement(FiniteElement *const &element_pt, const int &face_index)
Spine version of the PointFluidInterfaceBoundingElement.
void output(std::ostream &outfile, const unsigned &n_plot)
Output the element.
void output(FILE *file_pt)
Overload the C-style output function.
void output(FILE *file_pt, const unsigned &n_plot)
C-style Output function.
int kinematic_local_eqn(const unsigned &n)
Return local equation number associated with the kinematic constraint for local node n.
void output(std::ostream &outfile)
Overload the output function.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Calculate the elemental residual vector and the Jacobian.
SpineSurfaceFluidInterfaceElement(FiniteElement *const &element_pt, const int &face_index)
Generic Spine node update interface template class that can be combined with a given surface equation...
SpineUpdateFluidInterfaceElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0)
Constructor, the arguments are a pointer to the "bulk" element and the index of the face to be create...
void add_additional_residual_contributions_interface(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag, const Shape &psif, const DShape &dpsifds, const DShape &dpsifdS, const DShape &dpsifdS_div, const Vector< double > &s, const Vector< double > &interpolated_x, const Vector< double > &interpolated_n, const double &W, const double &J)
Helper function to calculate the additional contributions These are those filled in by the particular...
void output(std::ostream &outfile, const unsigned &n_plot)
Output the element.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Calculate the contribution to the residuals and the jacobian.
double compute_surface_derivatives(const Shape &psi, const DShape &dpsids, const DenseMatrix< double > &interpolated_t, const Vector< double > &interpolated_x, DShape &surface_gradient, DShape &surface_divergence)
Fill in the specific surface derivative calculations by calling the appropriate class function.
virtual FluidInterfaceBoundingElement * make_bounding_element(const int &face_index)
Create an "bounding" element of the type specified by the BoundingElementType policy class Here,...
int kinematic_local_eqn(const unsigned &n)
In spine elements, the kinematic condition is the equation used to determine the unknown spine height...
void output(FILE *file_pt)
Overload the C-style output function.
void hijack_kinematic_conditions(const Vector< unsigned > &bulk_node_number)
Hijacking the kinematic condition corresponds to hijacking the variables associated with the spine he...
void output(std::ostream &outfile)
Overload the output function.
void output(FILE *file_pt, const unsigned &n_plot)
C-style Output function.
Class that establishes the surface derivative functions for SurfaceInterfaceElements (2D surfaces in ...