oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT > Class Template Reference

#include <specific_node_update_interface_elements.h>

+ Inheritance diagram for oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >:

Public Member Functions

 ElasticUpdateFluidInterfaceElement (FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0)
 
double zeta_nodal (const unsigned &n, const unsigned &k, const unsigned &i) const
 
doublelagrange (const unsigned &n)
 Return the lagrange multiplier at local node n. More...
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 Fill in contribution to residuals and Jacobian. More...
 
void output (std::ostream &outfile)
 Overload the output function. More...
 
void output (std::ostream &outfile, const unsigned &n_plot)
 Output the element. More...
 
void output (FILE *file_pt)
 Overload the C-style output function. More...
 
void output (FILE *file_pt, const unsigned &n_plot)
 C-style Output function. More...
 
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)
 
virtual FluidInterfaceBoundingElementmake_bounding_element (const int &face_index)
 
- Public Member Functions inherited from oomph::Hijacked< FaceGeometry< ELEMENT > >
 Hijacked ()
 Constructor, call the constructors of the base elements. More...
 
 Hijacked (FiniteElement *const &element_pt, const int &face_index)
 Constructor used for hijacking face elements. More...
 
 Hijacked (FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0)
 
Datahijack_internal_value (const unsigned &n, const unsigned &i, const bool &return_data=true)
 
Datahijack_external_value (const unsigned &n, const unsigned &i, const bool &return_data=true)
 
Datahijack_nodal_value (const unsigned &n, const unsigned &i, const bool &return_data=true)
 
Datahijack_nodal_position_value (const unsigned &n, const unsigned &i, const bool &return_data=true)
 
Datahijack_nodal_spine_value (const unsigned &n, const unsigned &i, const bool &return_data=true)
 
void assign_local_eqn_numbers (const bool &store_local_dof_pt)
 
void get_residuals (Vector< double > &residuals)
 
void get_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
- Public Member Functions inherited from oomph::HijackedElementBase
 HijackedElementBase ()
 
virtual ~HijackedElementBase ()
 Destructor, destroy the storage for the equation numbers. More...
 
void unhijack_all_data ()
 
const doubleresidual_multiplier () const
 Return the value of the residual multiplier. More...
 
double *& residual_multiplier_pt ()
 Return the pointer to the residual multiplier. More...
 

Protected Member Functions

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)
 
- Protected Member Functions inherited from oomph::HijackedElementBase
void hijack_global_eqn (long *const &global_eqn_pt)
 
void unhijack_global_eqn (long *const &global_eqn_pt)
 

Private Member Functions

unsigned lagrange_index (const unsigned &n)
 
int kinematic_local_eqn (const unsigned &n)
 
void hijack_kinematic_conditions (const Vector< unsigned > &bulk_node_number)
 

Private Attributes

Vector< unsignedLagrange_index
 

Additional Inherited Members

- Protected Attributes inherited from oomph::HijackedElementBase
std::set< long * > * Hijacked_global_eqn_number_pt
 
Vector< int > * Hijacked_local_eqn_number_pt
 
doubleResidual_multiplier_pt
 
- Static Protected Attributes inherited from oomph::HijackedElementBase
static double Default_residual_multiplier = 0.0
 

Detailed Description

template<class EQUATION_CLASS, class DERIVATIVE_CLASS, class ELEMENT>
class oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >

//////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////// Generic Elastic node update interface template class that can be combined with a given surface equations class and surface derivative class to provide a concrete implementation of any surface element that uses elastic node updates.

Constructor & Destructor Documentation

◆ ElasticUpdateFluidInterfaceElement()

template<class EQUATION_CLASS , class DERIVATIVE_CLASS , class ELEMENT >
oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::ElasticUpdateFluidInterfaceElement ( FiniteElement *const &  element_pt,
const int face_index,
const unsigned id = 0 
)
inline

Constructor, pass a pointer to the bulk element and the face index of the bulk element to which the element is to be attached to. The optional identifier can be used to distinguish the additional nodal value (Lagr mult) created by this element from those created by other FaceElements.

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
780  FluidInterfaceAdditionalValues<EQUATION_CLASS>*
781  interface_additional_values_pt =
782  new FluidInterfaceAdditionalValues<EQUATION_CLASS>();
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  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Vector< unsigned > Lagrange_index
Definition: specific_node_update_interface_elements.h:681
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References oomph::FiniteElement::build_face_element(), i, oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::Lagrange_index, n, oomph::FluidInterfaceAdditionalValues< ELEMENT >::nadditional_values(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::FluidInterfaceAdditionalValues< ELEMENT >::setup_equation_indices().

Member Function Documentation

◆ add_additional_residual_contributions_interface()

template<class EQUATION_CLASS , class DERIVATIVE_CLASS , class ELEMENT >
void oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::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 
)
inline

Helper function to calculate the additional contributions to be added at each integration point. This deals with Lagrange multiplier contribution, as well as any additional contributions by the other equations.

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  }
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
double & lagrange(const unsigned &n)
Return the lagrange multiplier at local node n.
Definition: specific_node_update_interface_elements.h:830
int kinematic_local_eqn(const unsigned &n)
Definition: specific_node_update_interface_elements.h:692
RealScalar s
Definition: level1_cplx_impl.h:130
@ W
Definition: quadtree.h:63

References i, J, oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::kinematic_local_eqn(), oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::lagrange(), s, and oomph::QuadTreeNames::W.

◆ compute_surface_derivatives()

template<class EQUATION_CLASS , class DERIVATIVE_CLASS , class ELEMENT >
double oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::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 
)
inlineprotected

Fill in the specific surface derivative calculations by calling the appropriate function from the derivative class

723  {
724  return DERIVATIVE_CLASS::compute_surface_derivatives(psi,
725  dpsids,
726  interpolated_t,
727  interpolated_x,
728  surface_gradient,
729  surface_divergence);
730  }

◆ fill_in_contribution_to_jacobian()

template<class EQUATION_CLASS , class DERIVATIVE_CLASS , class ELEMENT >
void oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::fill_in_contribution_to_jacobian ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian 
)
inline

Fill in contribution to residuals and 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  }

◆ hijack_kinematic_conditions()

template<class EQUATION_CLASS , class DERIVATIVE_CLASS , class ELEMENT >
void oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::hijack_kinematic_conditions ( const Vector< unsigned > &  bulk_node_number)
inlineprivate

Hijacking the kinematic condition corresponds to hijacking the variables associated with the Lagrange multipliers that are assigned on construction of this element.

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  }
unsigned lagrange_index(const unsigned &n)
Definition: specific_node_update_interface_elements.h:685
Data * hijack_nodal_value(const unsigned &n, const unsigned &i, const bool &return_data=true)
Definition: hijacked_elements.h:214

References oomph::Hijacked< FaceGeometry< ELEMENT > >::hijack_nodal_value(), and oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::lagrange_index().

◆ kinematic_local_eqn()

template<class EQUATION_CLASS , class DERIVATIVE_CLASS , class ELEMENT >
int oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::kinematic_local_eqn ( const unsigned n)
inlineprivate

Equation number of the kinematic BC associated with node j. (This is the equation for the Lagrange multiplier)

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  }

References oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::lagrange_index(), and n.

Referenced by oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::add_additional_residual_contributions_interface().

◆ lagrange()

template<class EQUATION_CLASS , class DERIVATIVE_CLASS , class ELEMENT >
double& oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::lagrange ( const unsigned n)
inline

◆ lagrange_index()

template<class EQUATION_CLASS , class DERIVATIVE_CLASS , class ELEMENT >
unsigned oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::lagrange_index ( const unsigned n)
inlineprivate

◆ make_bounding_element()

template<class EQUATION_CLASS , class DERIVATIVE_CLASS , class ELEMENT >
virtual FluidInterfaceBoundingElement* oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::make_bounding_element ( const int face_index)
inlinevirtual

Create an "bounding" element (here actually a 2D line element of type ElasticLineFluidInterfaceBoundingElement<ELEMENT> that allows the application of a contact angle boundary condition on the the specified face.

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  }
ElasticUpdateFluidInterfaceElement(FiniteElement *const &element_pt, const int &face_index, const unsigned &id=0)
Definition: specific_node_update_interface_elements.h:739

References n.

◆ output() [1/4]

template<class EQUATION_CLASS , class DERIVATIVE_CLASS , class ELEMENT >
void oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::output ( FILE *  file_pt)
inline

Overload the C-style output function.

862  {
863  EQUATION_CLASS::output(file_pt);
864  }
void output(std::ostream &outfile, const unsigned &nplot)
Overload output function.
Definition: overloaded_element_body.h:490

References output().

◆ output() [2/4]

template<class EQUATION_CLASS , class DERIVATIVE_CLASS , class ELEMENT >
void oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::output ( FILE *  file_pt,
const unsigned n_plot 
)
inline

C-style Output function.

868  {
869  EQUATION_CLASS::output(file_pt, n_plot);
870  }

References output().

◆ output() [3/4]

template<class EQUATION_CLASS , class DERIVATIVE_CLASS , class ELEMENT >
void oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::output ( std::ostream &  outfile)
inline

Overload the output function.

850  {
851  EQUATION_CLASS::output(outfile);
852  }

References output().

◆ output() [4/4]

template<class EQUATION_CLASS , class DERIVATIVE_CLASS , class ELEMENT >
void oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::output ( std::ostream &  outfile,
const unsigned n_plot 
)
inline

Output the element.

856  {
857  EQUATION_CLASS::output(outfile, n_plot);
858  }

References output().

◆ zeta_nodal()

template<class EQUATION_CLASS , class DERIVATIVE_CLASS , class ELEMENT >
double oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::zeta_nodal ( const unsigned n,
const unsigned k,
const unsigned i 
) const
inline

The "global" intrinsic coordinate of the element when viewed as part of a geometric object should be given by the FaceElement representation, by default

825  {
826  return FaceElement::zeta_nodal(n, k, i);
827  }
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: elements.h:4497
char char char int int * k
Definition: level2_impl.h:374

References i, k, n, and oomph::FaceElement::zeta_nodal().

Member Data Documentation

◆ Lagrange_index

template<class EQUATION_CLASS , class DERIVATIVE_CLASS , class ELEMENT >
Vector<unsigned> oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::Lagrange_index
private

Storage for the location of the Lagrange multiplier (If other additional values have been added we need to add the Lagrange multiplier at the end)

Referenced by oomph::ElasticUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::ElasticUpdateFluidInterfaceElement().


The documentation for this class was generated from the following file: