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-2022 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  //=======================================================================
52  //=======================================================================
53  template<class ELEMENT>
55  {
56  };
57 
58  //======================================================================
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(
89  }
90 
91  public:
94 
97  inline unsigned nadditional_values(const unsigned& n)
98  {
99  error_message();
100  return 0;
101  }
102 
106  inline void setup_equation_indices(ELEMENT* const& element_pt,
107  const unsigned& id)
108  {
109  error_message();
110  }
111  };
112 
113 
114  //======================================================================
117  //=====================================================================
118  template<>
120  {
121  public:
124 
127  inline unsigned nadditional_values(const unsigned& n)
128  {
129  return 0;
130  }
131 
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  //======================================================================
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:
161  int kinematic_local_eqn(const unsigned& n)
162  {
163  return this->spine_local_eqn(n);
164  }
165 
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:
184  const Shape& psi,
185  const DShape& dpsids,
186  const DenseMatrix<double>& interpolated_t,
188  DShape& surface_gradient,
189  DShape& surface_divergence)
190  {
191  return DERIVATIVE_CLASS::compute_surface_derivatives(psi,
192  dpsids,
193  interpolated_t,
195  surface_gradient,
196  surface_divergence);
197  }
198 
199 
200  public:
204  const int& face_index,
205  const unsigned& id = 0)
206  : Hijacked<SpineElement<FaceGeometry<ELEMENT>>>(),
207  EQUATION_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",
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 
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 
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,
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,
317  interpolated_n,
318  W,
319  J);
320  }
321 
323  void output(std::ostream& outfile)
324  {
325  EQUATION_CLASS::output(outfile);
326  }
327 
329  void output(std::ostream& outfile, const unsigned& n_plot)
330  {
331  EQUATION_CLASS::output(outfile, n_plot);
332  }
333 
335  void output(FILE* file_pt)
336  {
337  EQUATION_CLASS::output(file_pt);
338  }
339 
341  void output(FILE* file_pt, const unsigned& n_plot)
342  {
343  EQUATION_CLASS::output(file_pt, n_plot);
344  }
345 
346 
353  const int& face_index)
354  {
355  // Create a temporary pointer to the appropriate FaceElement read our from
356  // our policy class
359  ELEMENT>>*
360  face_el_pt = new BoundingElementType<
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  //=====================================================================
419  //=====================================================================
420  template<class ELEMENT>
422  : public Hijacked<SpineElement<FaceGeometry<FaceGeometry<ELEMENT>>>>,
424 
425  {
426  public:
431  {
432  }
433 
435  void output(std::ostream& outfile)
436  {
437  FiniteElement::output(outfile);
438  }
439 
441  void output(std::ostream& outfile, const unsigned& n_plot)
442  {
443  FluidInterfaceBoundingElement::output(outfile, n_plot);
444  }
445 
447  void output(FILE* file_pt)
448  {
449  FiniteElement::output(file_pt);
450  }
451 
453  void output(FILE* file_pt, const unsigned& n_plot)
454  {
455  FluidInterfaceBoundingElement::output(file_pt, n_plot);
456  }
457 
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 
474  int kinematic_local_eqn(const unsigned& n)
475  {
476  return this->spine_local_eqn(n);
477  }
478  };
479 
480 
481  //=========================================================================
483  //========================================================================
484  template<class ELEMENT>
486  : public Hijacked<SpineElement<FaceGeometry<FaceGeometry<ELEMENT>>>>,
488 
489  {
490  public:
495  {
496  }
497 
499  void output(std::ostream& outfile)
500  {
501  FiniteElement::output(outfile);
502  }
503 
505  void output(std::ostream& outfile, const unsigned& n_plot)
506  {
507  FluidInterfaceBoundingElement::output(outfile, n_plot);
508  }
509 
511  void output(FILE* file_pt)
512  {
513  FiniteElement::output(file_pt);
514  }
515 
517  void output(FILE* file_pt, const unsigned& n_plot)
518  {
519  FluidInterfaceBoundingElement::output(file_pt, n_plot);
520  }
521 
522 
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 
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:
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:
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:
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 
660 
661  //-------------ELASTIC NODE UPDATE CLASSES-------------------------------
662  //---------------------------------------------------------------------
663 
664 
665  //=======================================================================
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:
682 
685  inline unsigned lagrange_index(const unsigned& n)
686  {
687  return this->Lagrange_index[n];
688  }
689 
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 
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:
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:
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",
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 
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 
830  double& lagrange(const unsigned& n)
831  {
832  return *this->node_pt(n)->value_pt(this->lagrange_index(n));
833  }
834 
835 
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 
849  void output(std::ostream& outfile)
850  {
851  EQUATION_CLASS::output(outfile);
852  }
853 
855  void output(std::ostream& outfile, const unsigned& n_plot)
856  {
857  EQUATION_CLASS::output(outfile, n_plot);
858  }
859 
861  void output(FILE* file_pt)
862  {
863  EQUATION_CLASS::output(file_pt);
864  }
865 
867  void output(FILE* file_pt, const unsigned& n_plot)
868  {
869  EQUATION_CLASS::output(file_pt, n_plot);
870  }
871 
872 
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 
961  const int& face_index)
962  {
963  // Create a temporary pointer to the appropriate FaceElement
966  ELEMENT>>*
967  face_el_pt = new BoundingElementType<
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  //=========================================================================
1035  //========================================================================
1036  template<class ELEMENT>
1038  : public FaceGeometry<FaceGeometry<ELEMENT>>,
1040  public virtual SolidFiniteElement
1041 
1042  {
1043  private:
1047 
1048  public:
1050  void set_lagrange_index(const Vector<unsigned>& lagrange_index)
1051  {
1052  Lagrange_index = lagrange_index;
1053  }
1054 
1055 
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 
1070  : FaceGeometry<FaceGeometry<ELEMENT>>(),
1072  {
1073  }
1074 
1076  void output(std::ostream& outfile)
1077  {
1078  FiniteElement::output(outfile);
1079  }
1080 
1082  void output(std::ostream& outfile, const unsigned& n_plot)
1083  {
1084  FluidInterfaceBoundingElement::output(outfile, n_plot);
1085  }
1086 
1088  void output(FILE* file_pt)
1089  {
1090  FiniteElement::output(file_pt);
1091  }
1092 
1094  void output(FILE* file_pt, const unsigned& n_plot)
1095  {
1096  FluidInterfaceBoundingElement::output(file_pt, n_plot);
1097  }
1098 
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
1112  }
1113 
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  //=========================================================================
1124  //========================================================================
1125  template<class ELEMENT>
1127  : public FaceGeometry<FaceGeometry<ELEMENT>>,
1129  public virtual SolidFiniteElement
1130 
1131  {
1134 
1135  public:
1137  void set_lagrange_index(const Vector<unsigned>& lagrange_index)
1138  {
1139  Lagrange_index = lagrange_index;
1140  }
1141 
1144  : FaceGeometry<FaceGeometry<ELEMENT>>(),
1146  {
1147  }
1148 
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 
1162  void output(std::ostream& outfile)
1163  {
1164  FiniteElement::output(outfile);
1165  }
1166 
1168  void output(std::ostream& outfile, const unsigned& n_plot)
1169  {
1170  FluidInterfaceBoundingElement::output(outfile, n_plot);
1171  }
1172 
1174  void output(FILE* file_pt)
1175  {
1176  FiniteElement::output(file_pt);
1177  }
1178 
1180  void output(FILE* file_pt, const unsigned& n_plot)
1181  {
1182  FluidInterfaceBoundingElement::output(file_pt, n_plot);
1183  }
1184 
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
1199  }
1200 
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 
1215  template<class ELEMENT>
1217  : public ElasticUpdateFluidInterfaceElement<FluidInterfaceElement,
1218  LineDerivatives,
1219  ELEMENT>
1220  {
1221  public:
1223  const int& face_index,
1224  const unsigned& id = 0)
1227  ELEMENT>(element_pt, face_index, id)
1228  {
1229  }
1230  };
1231 
1233  template<class ELEMENT>
1237  ELEMENT>>
1239  {
1240  public:
1244  ELEMENT>>()
1246  {
1247  }
1248  };
1249 
1250 
1252  template<class ELEMENT>
1254  : public ElasticUpdateFluidInterfaceElement<FluidInterfaceElement,
1255  AxisymmetricDerivatives,
1256  ELEMENT>
1257  {
1258  public:
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 
1290  template<class ELEMENT>
1292  : public ElasticUpdateFluidInterfaceElement<FluidInterfaceElement,
1293  SurfaceDerivatives,
1294  ELEMENT>
1295  {
1296  public:
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
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
Definition: interface_elements.h:685
Definition: nodes.h:1996
Definition: specific_node_update_interface_elements.h:55
Definition: shape.h:278
Specialise the Elastic update case to axisymmetric equations.
Definition: specific_node_update_interface_elements.h:1257
ElasticAxisymmetricFluidInterfaceElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0)
Definition: specific_node_update_interface_elements.h:1259
Pseudo-elasticity version of the LineFluidInterfaceBoundingElement.
Definition: specific_node_update_interface_elements.h:1131
void output(FILE *file_pt)
Overload the C-style output function.
Definition: specific_node_update_interface_elements.h:1174
int kinematic_local_eqn(const unsigned &n)
Local eqn number of kinematic bc associated with local node n.
Definition: specific_node_update_interface_elements.h:1202
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: specific_node_update_interface_elements.h:1153
void output(std::ostream &outfile)
Overload the output function.
Definition: specific_node_update_interface_elements.h:1162
void set_lagrange_index(const Vector< unsigned > &lagrange_index)
Set the Id.
Definition: specific_node_update_interface_elements.h:1137
void output(FILE *file_pt, const unsigned &n_plot)
C-style Output function.
Definition: specific_node_update_interface_elements.h:1180
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Calculate the elemental residual vector and Jacobian.
Definition: specific_node_update_interface_elements.h:1186
void output(std::ostream &outfile, const unsigned &n_plot)
Output the element.
Definition: specific_node_update_interface_elements.h:1168
Vector< unsigned > Lagrange_index
Short Storage for the index of Lagrange multiplier.
Definition: specific_node_update_interface_elements.h:1133
ElasticLineFluidInterfaceBoundingElement()
Constructor.
Definition: specific_node_update_interface_elements.h:1143
Specialise the elastic update template class to concrete 1D case.
Definition: specific_node_update_interface_elements.h:1220
ElasticLineFluidInterfaceElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0)
Definition: specific_node_update_interface_elements.h:1222
Pseudo-elasticity version of the PointFluidInterfaceBoundingElement.
Definition: specific_node_update_interface_elements.h:1042
ElasticPointFluidInterfaceBoundingElement()
Constructor.
Definition: specific_node_update_interface_elements.h:1069
void output(std::ostream &outfile, const unsigned &n_plot)
Output the element.
Definition: specific_node_update_interface_elements.h:1082
void output(std::ostream &outfile)
Overload the output function.
Definition: specific_node_update_interface_elements.h:1076
void output(FILE *file_pt, const unsigned &n_plot)
C-style Output function.
Definition: specific_node_update_interface_elements.h:1094
void output(FILE *file_pt)
Overload the C-style output function.
Definition: specific_node_update_interface_elements.h:1088
Vector< unsigned > Lagrange_index
Definition: specific_node_update_interface_elements.h:1046
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Calculate the element's residual vector and Jacobian.
Definition: specific_node_update_interface_elements.h:1100
void set_lagrange_index(const Vector< unsigned > &lagrange_index)
Set the Id and offset.
Definition: specific_node_update_interface_elements.h:1050
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: specific_node_update_interface_elements.h:1060
int kinematic_local_eqn(const unsigned &n)
Set the kinematic local equation.
Definition: specific_node_update_interface_elements.h:1115
Specialise Elastic update case to the concrete 2D case.
Definition: specific_node_update_interface_elements.h:1295
ElasticSurfaceFluidInterfaceElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0)
Definition: specific_node_update_interface_elements.h:1297
Definition: specific_node_update_interface_elements.h:676
void output(FILE *file_pt, const unsigned &n_plot)
C-style Output function.
Definition: specific_node_update_interface_elements.h:867
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)
Definition: specific_node_update_interface_elements.h:877
unsigned lagrange_index(const unsigned &n)
Definition: specific_node_update_interface_elements.h:685
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: specific_node_update_interface_elements.h:822
void output(std::ostream &outfile)
Overload the output function.
Definition: specific_node_update_interface_elements.h:849
ElasticUpdateFluidInterfaceElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0)
Definition: specific_node_update_interface_elements.h:739
Vector< unsigned > Lagrange_index
Definition: specific_node_update_interface_elements.h:681
double & lagrange(const unsigned &n)
Return the lagrange multiplier at local node n.
Definition: specific_node_update_interface_elements.h:830
void output(FILE *file_pt)
Overload the C-style output function.
Definition: specific_node_update_interface_elements.h:861
virtual FluidInterfaceBoundingElement * make_bounding_element(const int &face_index)
Definition: specific_node_update_interface_elements.h:960
int kinematic_local_eqn(const unsigned &n)
Definition: specific_node_update_interface_elements.h:692
void output(std::ostream &outfile, const unsigned &n_plot)
Output the element.
Definition: specific_node_update_interface_elements.h:855
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution to residuals and Jacobian.
Definition: specific_node_update_interface_elements.h:837
void hijack_kinematic_conditions(const Vector< unsigned > &bulk_node_number)
Definition: specific_node_update_interface_elements.h:701
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)
Definition: specific_node_update_interface_elements.h:716
void assemble_set_of_all_geometric_data(std::set< Data * > &unique_geom_data_pt)
Return a set of all geometric data associated with the element.
Definition: element_with_moving_nodes.cc:43
void fill_in_jacobian_from_geometric_data(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: element_with_moving_nodes.cc:323
int & face_index()
Definition: elements.h:4626
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: elements.h:4497
Definition: elements.h:4998
Definition: elements.h:1313
virtual void output(std::ostream &outfile)
Definition: elements.h:3050
virtual void build_face_element(const int &face_index, FaceElement *face_element_pt)
Definition: elements.cc:5132
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
Definition: elements.cc:3962
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
bool has_hanging_nodes() const
Definition: elements.h:2470
Definition: specific_node_update_interface_elements.h:120
void setup_equation_indices(FluidInterfaceElement *const &element_pt, const unsigned &id)
Definition: specific_node_update_interface_elements.h:135
unsigned nadditional_values(const unsigned &n)
Definition: specific_node_update_interface_elements.h:127
Definition: specific_node_update_interface_elements.h:72
unsigned nadditional_values(const unsigned &n)
Definition: specific_node_update_interface_elements.h:97
void setup_equation_indices(ELEMENT *const &element_pt, const unsigned &id)
Definition: specific_node_update_interface_elements.h:106
void error_message()
Definition: specific_node_update_interface_elements.h:75
Definition: interface_elements.h:54
void output(std::ostream &outfile)
Overload the output function.
Definition: interface_elements.h:247
Definition: interface_elements.h:321
void fill_in_jacobian_from_external_by_fd(Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
Definition: elements.cc:1199
Definition: hijacked_elements.h:132
Data * hijack_nodal_spine_value(const unsigned &n, const unsigned &i, const bool &return_data=true)
Definition: hijacked_elements.h:279
Data * hijack_nodal_value(const unsigned &n, const unsigned &i, const bool &return_data=true)
Definition: hijacked_elements.h:214
Definition: interface_elements.h:661
Specialisation of the interface boundary constraint to a line.
Definition: interface_elements.h:297
void fill_in_generic_residual_contribution_interface_boundary(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Add contribution to element's residual vector and Jacobian.
Definition: interface_elements.cc:224
Definition: oomph_definitions.h:222
Specialisation of the interface boundary constraint to a point.
Definition: interface_elements.h:277
void fill_in_generic_residual_contribution_interface_boundary(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Add contribution to element's residual vector and Jacobian.
Definition: interface_elements.cc:86
Definition: refineable_elements.h:97
Definition: shape.h:76
Definition: elements.h:3561
virtual void fill_in_jacobian_from_solid_position_by_fd(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: elements.cc:6985
Definition: nodes.h:1686
Definition: specific_node_update_interface_elements.h:592
SpineAxisymmetricFluidInterfaceElement(FiniteElement *const &element_pt, const int &face_index)
Definition: specific_node_update_interface_elements.h:594
Definition: spines.h:477
int spine_local_eqn(const unsigned &n)
Definition: spines.h:521
Spine version of the LineFluidInterfaceBoundingElement.
Definition: specific_node_update_interface_elements.h:489
int kinematic_local_eqn(const unsigned &n)
Local eqn number of the kinematic bc associated with local node n.
Definition: specific_node_update_interface_elements.h:539
SpineLineFluidInterfaceBoundingElement()
Constructor.
Definition: specific_node_update_interface_elements.h:492
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Calculate the jacobian.
Definition: specific_node_update_interface_elements.h:524
void output(std::ostream &outfile)
Overload the output function.
Definition: specific_node_update_interface_elements.h:499
void output(std::ostream &outfile, const unsigned &n_plot)
Output the element.
Definition: specific_node_update_interface_elements.h:505
void output(FILE *file_pt, const unsigned &n_plot)
C-style Output function.
Definition: specific_node_update_interface_elements.h:517
void output(FILE *file_pt)
Overload the C-style output function.
Definition: specific_node_update_interface_elements.h:511
Definition: specific_node_update_interface_elements.h:556
SpineLineFluidInterfaceElement(FiniteElement *const &element_pt, const int &face_index)
Definition: specific_node_update_interface_elements.h:558
Spine version of the PointFluidInterfaceBoundingElement.
Definition: specific_node_update_interface_elements.h:425
void output(std::ostream &outfile, const unsigned &n_plot)
Output the element.
Definition: specific_node_update_interface_elements.h:441
void output(FILE *file_pt)
Overload the C-style output function.
Definition: specific_node_update_interface_elements.h:447
SpinePointFluidInterfaceBoundingElement()
Constructor.
Definition: specific_node_update_interface_elements.h:428
void output(FILE *file_pt, const unsigned &n_plot)
C-style Output function.
Definition: specific_node_update_interface_elements.h:453
int kinematic_local_eqn(const unsigned &n)
Definition: specific_node_update_interface_elements.h:474
void output(std::ostream &outfile)
Overload the output function.
Definition: specific_node_update_interface_elements.h:435
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Calculate the elemental residual vector and the Jacobian.
Definition: specific_node_update_interface_elements.h:459
Definition: specific_node_update_interface_elements.h:628
SpineSurfaceFluidInterfaceElement(FiniteElement *const &element_pt, const int &face_index)
Definition: specific_node_update_interface_elements.h:630
Definition: specific_node_update_interface_elements.h:156
SpineUpdateFluidInterfaceElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0)
Definition: specific_node_update_interface_elements.h:203
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...
Definition: specific_node_update_interface_elements.h:293
void output(std::ostream &outfile, const unsigned &n_plot)
Output the element.
Definition: specific_node_update_interface_elements.h:329
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Calculate the contribution to the residuals and the jacobian.
Definition: specific_node_update_interface_elements.h:279
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)
Definition: specific_node_update_interface_elements.h:183
virtual FluidInterfaceBoundingElement * make_bounding_element(const int &face_index)
Definition: specific_node_update_interface_elements.h:352
int kinematic_local_eqn(const unsigned &n)
Definition: specific_node_update_interface_elements.h:161
void output(FILE *file_pt)
Overload the C-style output function.
Definition: specific_node_update_interface_elements.h:335
void hijack_kinematic_conditions(const Vector< unsigned > &bulk_node_number)
Definition: specific_node_update_interface_elements.h:168
void output(std::ostream &outfile)
Overload the output function.
Definition: specific_node_update_interface_elements.h:323
void output(FILE *file_pt, const unsigned &n_plot)
C-style Output function.
Definition: specific_node_update_interface_elements.h:341
Definition: interface_elements.h:709
RealScalar s
Definition: level1_cplx_impl.h:130
char char char int int * k
Definition: level2_impl.h:374
@ W
Definition: quadtree.h:63
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
void output(std::ostream &outfile, const unsigned &nplot)
Overload output function.
Definition: overloaded_element_body.h:490