ElasticTwoLayerMesh< ELEMENT > Class Template Reference
+ Inheritance diagram for ElasticTwoLayerMesh< ELEMENT >:

Public Member Functions

 ElasticTwoLayerMesh (const unsigned &nx, const unsigned &ny1, const unsigned &ny2, const double &lx, const double &h1, const double &h2, const bool &periodic_in_x=false, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
FiniteElement *& upper_layer_element_pt (const unsigned long &i)
 Access functions for pointers to elements in upper layer. More...
 
FiniteElement *& lower_layer_element_pt (const unsigned long &i)
 Access functions for pointers to elements in bottom layer. More...
 
unsigned long nupper () const
 Number of elements in upper layer. More...
 
unsigned long nlower () const
 Number of elements in top layer. More...
 
FiniteElement *& interface_upper_boundary_element_pt (const unsigned long &i)
 Access functions for pointers to elements in upper layer. More...
 
FiniteElement *& interface_lower_boundary_element_pt (const unsigned long &i)
 Access functions for pointers to elements in bottom layer. More...
 
unsigned long ninterface_upper () const
 Number of elements in upper layer. More...
 
unsigned long ninterface_lower () const
 Number of elements in top layer. More...
 
int interface_upper_face_index_at_boundary (const unsigned &e)
 
int interface_lower_face_index_at_boundary (const unsigned &e)
 
- Public Member Functions inherited from oomph::ElasticRectangularQuadMesh< ELEMENT >
 ElasticRectangularQuadMesh (const unsigned &nx, const unsigned &ny, const double &lx, const double &ly, const Vector< double > &origin, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
 ElasticRectangularQuadMesh (const unsigned &nx, const unsigned &ny, const double &lx, const double &ly, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
 ElasticRectangularQuadMesh (const unsigned &nx, const unsigned &ny, const double &lx, const double &ly, const bool &periodic_in_x, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
- Public Member Functions inherited from oomph::RectangularQuadMesh< ELEMENT >
 RectangularQuadMesh (const unsigned &nx, const unsigned &ny, const double &lx, const double &ly, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
 RectangularQuadMesh (const unsigned &nx, const unsigned &ny, const double &xmin, const double &xmax, const double &ymin, const double &ymax, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
 RectangularQuadMesh (const unsigned &nx, const unsigned &ny, const double &lx, const double &ly, const bool &periodic_in_x, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
 RectangularQuadMesh (const unsigned &nx, const unsigned &ny, const double &xmin, const double &xmax, const double &ymin, const double &ymax, const bool &periodic_in_x, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
const unsignednx () const
 Return number of elements in x direction. More...
 
const unsignedny () const
 Return number of elements in y direction. More...
 
const double x_min () const
 Return the minimum value of x coordinate. More...
 
const double x_max () const
 Return the maximum value of x coordinate. More...
 
const double y_min () const
 Return the minimum value of y coordinate. More...
 
const double y_max () const
 Return the maximum value of y coordinate. More...
 
virtual void element_reorder ()
 
virtual double x_spacing_function (unsigned xelement, unsigned xnode, unsigned yelement, unsigned ynode)
 
virtual double y_spacing_function (unsigned xelement, unsigned xnode, unsigned yelement, unsigned ynode)
 
- Public Member Functions inherited from oomph::QuadMeshBase
 QuadMeshBase ()
 Constructor (empty) More...
 
 QuadMeshBase (const QuadMeshBase &node)=delete
 Broken copy constructor. More...
 
void operator= (const QuadMeshBase &)=delete
 Broken assignment operator. More...
 
virtual ~QuadMeshBase ()
 Destructor (empty) More...
 
void setup_boundary_element_info ()
 
void setup_boundary_element_info (std::ostream &outfile)
 
- Public Member Functions inherited from oomph::Mesh
 Mesh ()
 Default constructor. More...
 
 Mesh (const Vector< Mesh * > &sub_mesh_pt)
 
void merge_meshes (const Vector< Mesh * > &sub_mesh_pt)
 
virtual void reset_boundary_element_info (Vector< unsigned > &ntmp_boundary_elements, Vector< Vector< unsigned >> &ntmp_boundary_elements_in_region, Vector< FiniteElement * > &deleted_elements)
 Virtual function to perform the reset boundary elements info rutines. More...
 
template<class BULK_ELEMENT >
void doc_boundary_coordinates (const unsigned &b, std::ofstream &the_file)
 
 Mesh (const Mesh &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const Mesh &)=delete
 Broken assignment operator. More...
 
virtual ~Mesh ()
 Virtual Destructor to clean up all memory. More...
 
void flush_element_and_node_storage ()
 
void flush_element_storage ()
 
void flush_node_storage ()
 
Node *& node_pt (const unsigned long &n)
 Return pointer to global node n. More...
 
Nodenode_pt (const unsigned long &n) const
 Return pointer to global node n (const version) More...
 
GeneralisedElement *& element_pt (const unsigned long &e)
 Return pointer to element e. More...
 
GeneralisedElementelement_pt (const unsigned long &e) const
 Return pointer to element e (const version) More...
 
const Vector< GeneralisedElement * > & element_pt () const
 Return reference to the Vector of elements. More...
 
Vector< GeneralisedElement * > & element_pt ()
 Return reference to the Vector of elements. More...
 
FiniteElementfinite_element_pt (const unsigned &e) const
 
Node *& boundary_node_pt (const unsigned &b, const unsigned &n)
 Return pointer to node n on boundary b. More...
 
Nodeboundary_node_pt (const unsigned &b, const unsigned &n) const
 Return pointer to node n on boundary b. More...
 
void set_nboundary (const unsigned &nbound)
 Set the number of boundaries in the mesh. More...
 
void remove_boundary_nodes ()
 Clear all pointers to boundary nodes. More...
 
void remove_boundary_nodes (const unsigned &b)
 
void remove_boundary_node (const unsigned &b, Node *const &node_pt)
 Remove a node from the boundary b. More...
 
void add_boundary_node (const unsigned &b, Node *const &node_pt)
 Add a (pointer to) a node to the b-th boundary. More...
 
void copy_boundary_node_data_from_nodes ()
 
bool boundary_coordinate_exists (const unsigned &i) const
 Indicate whether the i-th boundary has an intrinsic coordinate. More...
 
unsigned long nelement () const
 Return number of elements in the mesh. More...
 
unsigned long nnode () const
 Return number of nodes in the mesh. More...
 
unsigned ndof_types () const
 Return number of dof types in mesh. More...
 
unsigned elemental_dimension () const
 Return number of elemental dimension in mesh. More...
 
unsigned nodal_dimension () const
 Return number of nodal dimension in mesh. More...
 
void add_node_pt (Node *const &node_pt)
 Add a (pointer to a) node to the mesh. More...
 
void add_element_pt (GeneralisedElement *const &element_pt)
 Add a (pointer to) an element to the mesh. More...
 
virtual void node_update (const bool &update_all_solid_nodes=false)
 
virtual void reorder_nodes (const bool &use_old_ordering=true)
 
virtual void get_node_reordering (Vector< Node * > &reordering, const bool &use_old_ordering=true) const
 
template<class BULK_ELEMENT , template< class > class FACE_ELEMENT>
void build_face_mesh (const unsigned &b, Mesh *const &face_mesh_pt)
 
unsigned self_test ()
 Self-test: Check elements and nodes. Return 0 for OK. More...
 
void max_and_min_element_size (double &max_size, double &min_size)
 
double total_size ()
 
void check_inverted_elements (bool &mesh_has_inverted_elements, std::ofstream &inverted_element_file)
 
void check_inverted_elements (bool &mesh_has_inverted_elements)
 
unsigned check_for_repeated_nodes (const double &epsilon=1.0e-12)
 
Vector< Node * > prune_dead_nodes ()
 
unsigned nboundary () const
 Return number of boundaries. More...
 
unsigned long nboundary_node (const unsigned &ibound) const
 Return number of nodes on a particular boundary. More...
 
FiniteElementboundary_element_pt (const unsigned &b, const unsigned &e) const
 Return pointer to e-th finite element on boundary b. More...
 
Nodeget_some_non_boundary_node () const
 
unsigned nboundary_element (const unsigned &b) const
 Return number of finite elements that are adjacent to boundary b. More...
 
int face_index_at_boundary (const unsigned &b, const unsigned &e) const
 
virtual void dump (std::ofstream &dump_file, const bool &use_old_ordering=true) const
 Dump the data in the mesh into a file for restart. More...
 
void dump (const std::string &dump_file_name, const bool &use_old_ordering=true) const
 Dump the data in the mesh into a file for restart. More...
 
virtual void read (std::ifstream &restart_file)
 Read solution from restart file. More...
 
void output_paraview (std::ofstream &file_out, const unsigned &nplot) const
 
void output_fct_paraview (std::ofstream &file_out, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const
 
void output_fct_paraview (std::ofstream &file_out, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) const
 
void output (std::ostream &outfile)
 Output for all elements. More...
 
void output (std::ostream &outfile, const unsigned &n_plot)
 Output at f(n_plot) points in each element. More...
 
void output (FILE *file_pt)
 Output for all elements (C-style output) More...
 
void output (FILE *file_pt, const unsigned &nplot)
 Output at f(n_plot) points in each element (C-style output) More...
 
void output (const std::string &output_filename)
 Output for all elements. More...
 
void output (const std::string &output_filename, const unsigned &n_plot)
 Output at f(n_plot) points in each element. More...
 
void output_fct (std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt)
 Output a given Vector function at f(n_plot) points in each element. More...
 
void output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt)
 
void output_boundaries (std::ostream &outfile)
 Output the nodes on the boundaries (into separate tecplot zones) More...
 
void output_boundaries (const std::string &output_filename)
 
void assign_initial_values_impulsive ()
 Assign initial values for an impulsive start. More...
 
void shift_time_values ()
 
void calculate_predictions ()
 
void set_nodal_and_elemental_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
virtual void set_mesh_level_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
void set_consistent_pinned_values_for_continuation (ContinuationStorageScheme *const &continuation_stepper_pt)
 Set consistent values for pinned data in continuation. More...
 
bool does_pointer_correspond_to_mesh_data (double *const &parameter_pt)
 Does the double pointer correspond to any mesh data. More...
 
void set_nodal_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 Set the timestepper associated with the nodal data in the mesh. More...
 
void set_elemental_internal_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
virtual void compute_norm (double &norm)
 
virtual void compute_norm (Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 Returns the norm of the error and that of the exact solution. More...
 
virtual void compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
bool is_mesh_distributed () const
 Boolean to indicate if Mesh has been distributed. More...
 
OomphCommunicatorcommunicator_pt () const
 
void delete_all_external_storage ()
 Wipe the storage for all externally-based elements. More...
 
- Public Member Functions inherited from oomph::SolidMesh
 SolidMesh ()
 Default constructor. More...
 
 SolidMesh (const SolidMesh &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const SolidMesh &)=delete
 Broken assignment operator. More...
 
 SolidMesh (const Vector< SolidMesh * > &sub_mesh_pt)
 
SolidNodenode_pt (const unsigned long &n)
 Return a pointer to the n-th global SolidNode. More...
 
SolidNodeboundary_node_pt (const unsigned &b, const unsigned &n)
 Return n-th SolidNodes on b-th boundary. More...
 
SolidNodeelement_node_pt (const unsigned long &e, const unsigned &n)
 
void set_lagrangian_nodal_coordinates ()
 
void scale_mesh (const double &factor)
 

Private Attributes

Vector< FiniteElement * > Lower_layer_element_pt
 Vector of pointers to element in the upper layer. More...
 
Vector< FiniteElement * > Upper_layer_element_pt
 Vector of pointers to element in the lower layer. More...
 
Vector< FiniteElement * > Interface_lower_boundary_element_pt
 
Vector< FiniteElement * > Interface_upper_boundary_element_pt
 

Additional Inherited Members

- Public Types inherited from oomph::Mesh
typedef void(FiniteElement::* SteadyExactSolutionFctPt) (const Vector< double > &x, Vector< double > &soln)
 
typedef void(FiniteElement::* UnsteadyExactSolutionFctPt) (const double &time, const Vector< double > &x, Vector< double > &soln)
 
- Static Public Attributes inherited from oomph::Mesh
static Steady< 0 > Default_TimeStepper
 The Steady Timestepper. More...
 
static bool Suppress_warning_about_empty_mesh_level_time_stepper_function
 Static boolean flag to control warning about mesh level timesteppers. More...
 
- Static Public Attributes inherited from oomph::SolidMesh
static SolidICProblem Solid_IC_problem
 
- Protected Member Functions inherited from oomph::RectangularQuadMesh< ELEMENT >
void build_mesh (TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 Generic mesh construction function: contains all the hard work. More...
 
 RectangularQuadMesh (const unsigned &nx, const unsigned &ny, const double &xmin, const double &xmax, const double &ymin, const double &ymax, const bool &periodic_in_x, const bool &build, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
- Protected Member Functions inherited from oomph::Mesh
unsigned long assign_global_eqn_numbers (Vector< double * > &Dof_pt)
 Assign (global) equation numbers to the nodes. More...
 
void describe_dofs (std::ostream &out, const std::string &current_string) const
 
void describe_local_dofs (std::ostream &out, const std::string &current_string) const
 
void assign_local_eqn_numbers (const bool &store_local_dof_pt)
 Assign local equation numbers in all elements. More...
 
void convert_to_boundary_node (Node *&node_pt, const Vector< FiniteElement * > &finite_element_pt)
 
void convert_to_boundary_node (Node *&node_pt)
 
- Protected Attributes inherited from oomph::RectangularQuadMesh< ELEMENT >
unsigned Nx
 Nx: number of elements in x-direction. More...
 
unsigned Ny
 Ny: number of elements in y-direction. More...
 
unsigned Np
 Np: number of (linear) points in the element. More...
 
double Xmin
 Minimum value of x coordinate. More...
 
double Xmax
 Maximum value of x coordinate. More...
 
double Ymin
 Minimum value of y coordinate. More...
 
double Ymax
 Maximum value of y coordinate. More...
 
bool Xperiodic
 
- Protected Attributes inherited from oomph::Mesh
Vector< Vector< Node * > > Boundary_node_pt
 
bool Lookup_for_elements_next_boundary_is_setup
 
Vector< Vector< FiniteElement * > > Boundary_element_pt
 
Vector< Vector< int > > Face_index_at_boundary
 
Vector< Node * > Node_pt
 Vector of pointers to nodes. More...
 
Vector< GeneralisedElement * > Element_pt
 Vector of pointers to generalised elements. More...
 
std::vector< boolBoundary_coordinate_exists
 

Detailed Description

template<class ELEMENT>
class ElasticTwoLayerMesh< ELEMENT >

Two layer mesh which employs a pseudo-solid node-update strategy. This class is essentially a wrapper to an ElasticRectangularQuadMesh, with an additional boundary to represent the interface between the two fluid layers. In addition, the mesh paritions the elements into those above and below the interface and relabels boundaries so that we can impose a volume constraint on the lower or upper fluid.

                            3
          ---------------------------------------
          |                                     |
        4 |                                     | 2
          |                 6                   |
          ---------------------------------------
          |                                     |
        5 |                                     | 1
          |                                     |
          ---------------------------------------
                            0 

Constructor & Destructor Documentation

◆ ElasticTwoLayerMesh()

template<class ELEMENT >
ElasticTwoLayerMesh< ELEMENT >::ElasticTwoLayerMesh ( const unsigned nx,
const unsigned ny1,
const unsigned ny2,
const double lx,
const double h1,
const double h2,
const bool periodic_in_x = false,
TimeStepper time_stepper_pt = &Mesh::Default_TimeStepper 
)
inline

Constructor: Pass number of elements in x-direction, number of elements in y-direction in bottom and top layer, respectively, axial length and height of top and bottom layers, a boolean flag to make the mesh periodic in the x-direction, and pointer to timestepper (defaults to Steady timestepper)

472  :
473  RectangularQuadMesh<ELEMENT>(nx,ny1+ny2,lx,h1+h2,
474  periodic_in_x,time_stepper_pt),
476  periodic_in_x,time_stepper_pt)
477  {
478  //Set up the pointers to elements in the upper and lower fluid
479  //Calculate numbers of nodes in upper and lower regions
480  unsigned long n_lower = nx*ny1;
481  unsigned long n_upper = nx*ny2;
482  //Loop over lower elements and push back onto the array
483  Lower_layer_element_pt.resize(n_lower);
484  for(unsigned e=0;e<n_lower;e++)
485  {
487  }
488  //Loop over upper elements and push back onto the array
489  Upper_layer_element_pt.resize(n_upper);
490  for(unsigned e=0;e<n_upper;e++)
491  {
492  Upper_layer_element_pt[e] = this->finite_element_pt(n_lower + e);
493  }
494  //end of upper and lower fluid element assignment
495 
496  //Set the elements adjacent to the interface on both sides
499  {
500  unsigned count_lower=nx*(ny1-1);
501  unsigned count_upper= count_lower + nx;
502  for(unsigned e=0;e<nx;e++)
503  {
505  this->finite_element_pt(count_lower); ++count_lower;
507  this->finite_element_pt(count_upper); ++count_upper;
508  }
509  } //end of bulk elements next to interface setup
510 
511 
512  // Reset the number of boundaries
513  this->set_nboundary(7);
514 
515  //Relabel the boundary nodes
516  //Storage for the boundary coordinates that will be transferred directly
517  Vector<double> b_coord;
518  {
519  //Store the interface level
520  const double y_interface = h1 - this->Ymin;
521 
522  //Nodes on boundary 3 are now on boundaries 4 and 5
523  unsigned n_boundary_node = this->nboundary_node(3);
524  //Loop over the nodes remove them from boundary 3
525  //and add them to boundary 4 or 5 depending on their y coordinate
526  for(unsigned n=0;n<n_boundary_node;n++)
527  {
528  //Cache pointer to the node
529  Node* const nod_pt = this->boundary_node_pt(3,n);
530  //Get the boundary coordinates if set
531  if(this->boundary_coordinate_exists(3))
532  {
533  b_coord.resize(nod_pt->ncoordinates_on_boundary(3));
534  nod_pt->get_coordinates_on_boundary(3,b_coord);
535  //Indicate that the boundary coordinates are to be set on the
536  //new boundaries
537  this->Boundary_coordinate_exists[4]=true;
538  this->Boundary_coordinate_exists[5]=true;
539  }
540 
541  //Now remove the node from the old boundary
542  nod_pt->remove_from_boundary(3);
543 
544  //Find the height of the node
545  double y = nod_pt->x(1);
546  //If it's above or equal to the interface, it's on boundary 4
547  if(y >= y_interface)
548  {
549  this->add_boundary_node(4,nod_pt);
550  //Add the boundary coordinate if it has been set up
551  if(this->Boundary_coordinate_exists[4])
552  {
553  nod_pt->set_coordinates_on_boundary(4,b_coord);
554  }
555  }
556  //otherwise it's on boundary 5
557  if(y <= y_interface)
558  {
559  this->add_boundary_node(5,nod_pt);
560  //Add the boundary coordinate if it has been set up
561  if(this->Boundary_coordinate_exists[5])
562  {
563  nod_pt->set_coordinates_on_boundary(5,b_coord);
564  }
565  }
566  }
567 
568  //Now clear the boundary node information from boundary 3
569  this->Boundary_node_pt[3].clear();
570 
571  //Relabel the nodes on boundary 2 to be on boundary 3
572  n_boundary_node = this->nboundary_node(2);
573  //Loop over the nodes remove them from boundary 2
574  //and add them to boundary 3
575  for(unsigned n=0;n<n_boundary_node;n++)
576  {
577  //Cache pointer to the node
578  Node* const nod_pt = this->boundary_node_pt(2,n);
579  //Get the boundary coordinates if set
580  if(this->boundary_coordinate_exists(2))
581  {
582  b_coord.resize(nod_pt->ncoordinates_on_boundary(2));
583  nod_pt->get_coordinates_on_boundary(2,b_coord);
584  this->Boundary_coordinate_exists[3]=true;
585  }
586 
587  //Now remove the node from the boundary 2
588  nod_pt->remove_from_boundary(2);
589  //and add to boundary 3
590  this->add_boundary_node(3,nod_pt);
591  if(this->Boundary_coordinate_exists[3])
592  {
593  nod_pt->set_coordinates_on_boundary(3,b_coord);
594  }
595  }
596 
597  //Clear the information from boundary 2
598  this->Boundary_node_pt[2].clear();
599 
600  //Storage for nodes that should be removed from boundary 1
601  std::list<Node*> nodes_to_be_removed;
602 
603  //Nodes on boundary 1 are now on boundaries 1 and 2
604  n_boundary_node = this->nboundary_node(1);
605  //Loop over the nodes remove some of them from boundary 1
606  for(unsigned n=0;n<n_boundary_node;n++)
607  {
608  //Cache pointer to the node
609  Node* const nod_pt = this->boundary_node_pt(1,n);
610 
611  //Find the height of the node
612  double y = nod_pt->x(1);
613  //If it's above or equal to the interface it's on boundary 2
614  if(y >= y_interface)
615  {
616  //Get the boundary coordinates if set
617  if(this->boundary_coordinate_exists(1))
618  {
619  b_coord.resize(nod_pt->ncoordinates_on_boundary(1));
620  nod_pt->get_coordinates_on_boundary(1,b_coord);
621  this->Boundary_coordinate_exists[2]=true;
622  }
623 
624  //Now remove the node from the boundary 1 if above interace
625  if(y > y_interface)
626  {
627  nodes_to_be_removed.push_back(nod_pt);
628  }
629  //Always add to boundary 2
630  this->add_boundary_node(2,nod_pt);
631  //Add the boundary coordinate if it has been set up
632  if(this->Boundary_coordinate_exists[2])
633  {
634  nod_pt->set_coordinates_on_boundary(2,b_coord);
635  }
636  }
637  }
638 
639  //Loop over the nodes that are to be removed from boundary 1 and remove
640  //them
641  for(std::list<Node*>::iterator it=nodes_to_be_removed.begin();
642  it!=nodes_to_be_removed.end();++it)
643  {
644  this->remove_boundary_node(1,*it);
645  }
646  nodes_to_be_removed.clear();
647  } //end of boundary relabelling
648 
649  //Add the nodes to the interface
650 
651  //Read out number of linear points in the element
652  unsigned n_p = dynamic_cast<ELEMENT*>
653  (this->finite_element_pt(0))->nnode_1d();
654 
655  //Add the nodes on the interface to the boundary 6
656  //Storage for boundary coordinates (x-coordinate)
657  b_coord.resize(1);
659  //Starting index of the nodes
660  unsigned n_start=0;
661  for(unsigned e=0;e<nx;e++)
662  {
663  //If we are past the
664  if(e>0) {n_start=1;}
665  //Get the layer of elements just above the interface
666  FiniteElement* el_pt = this->finite_element_pt(nx*ny1+e);
667  //The first n_p nodes lie on the boundary
668  for(unsigned n=n_start;n<n_p;n++)
669  {
670  Node* nod_pt = el_pt->node_pt(n);
671  this->convert_to_boundary_node(nod_pt);
672  this->add_boundary_node(6,nod_pt);
673  b_coord[0] = nod_pt->x(0);
674  nod_pt->set_coordinates_on_boundary(6,b_coord);
675  }
676  }
677 
678  // Set up the boundary element information
680  }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Array< double, 1, 3 > e(1./3., 0.5, 2.)
Vector< FiniteElement * > Lower_layer_element_pt
Vector of pointers to element in the upper layer.
Definition: static_two_layer.cc:725
Vector< FiniteElement * > Interface_lower_boundary_element_pt
Definition: static_two_layer.cc:732
Vector< FiniteElement * > Upper_layer_element_pt
Vector of pointers to element in the lower layer.
Definition: static_two_layer.cc:728
Vector< FiniteElement * > Interface_upper_boundary_element_pt
Definition: static_two_layer.cc:736
Definition: rectangular_quadmesh.template.h:423
Definition: elements.h:1313
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
void add_boundary_node(const unsigned &b, Node *const &node_pt)
Add a (pointer to) a node to the b-th boundary.
Definition: mesh.cc:243
unsigned long nboundary_node(const unsigned &ibound) const
Return number of nodes on a particular boundary.
Definition: mesh.h:833
void remove_boundary_node(const unsigned &b, Node *const &node_pt)
Remove a node from the boundary b.
Definition: mesh.cc:221
std::vector< bool > Boundary_coordinate_exists
Definition: mesh.h:190
FiniteElement * finite_element_pt(const unsigned &e) const
Definition: mesh.h:473
void set_nboundary(const unsigned &nbound)
Set the number of boundaries in the mesh.
Definition: mesh.h:505
virtual void setup_boundary_element_info()
Definition: mesh.h:275
void convert_to_boundary_node(Node *&node_pt, const Vector< FiniteElement * > &finite_element_pt)
Definition: mesh.cc:2590
Vector< Vector< Node * > > Boundary_node_pt
Definition: mesh.h:83
bool boundary_coordinate_exists(const unsigned &i) const
Indicate whether the i-th boundary has an intrinsic coordinate.
Definition: mesh.h:565
Definition: nodes.h:906
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition: nodes.h:1060
virtual unsigned ncoordinates_on_boundary(const unsigned &b)
Definition: nodes.cc:2363
virtual void remove_from_boundary(const unsigned &b)
Definition: nodes.cc:2350
virtual void set_coordinates_on_boundary(const unsigned &b, const unsigned &k, const Vector< double > &boundary_zeta)
Definition: nodes.cc:2394
virtual void get_coordinates_on_boundary(const unsigned &b, const unsigned &k, Vector< double > &boundary_zeta)
Definition: nodes.cc:2379
Definition: rectangular_quadmesh.template.h:59
const unsigned & nx() const
Return number of elements in x direction.
Definition: rectangular_quadmesh.template.h:224
double Ymin
Minimum value of y coordinate.
Definition: rectangular_quadmesh.template.h:75
SolidNode * boundary_node_pt(const unsigned &b, const unsigned &n)
Return n-th SolidNodes on b-th boundary.
Definition: mesh.h:2612
Scalar * y
Definition: level1_cplx_impl.h:128
const double lx
Definition: ConstraintElementsUnitTest.cpp:33

References oomph::Mesh::add_boundary_node(), oomph::Mesh::Boundary_coordinate_exists, oomph::Mesh::boundary_coordinate_exists(), oomph::Mesh::Boundary_node_pt, oomph::SolidMesh::boundary_node_pt(), oomph::Mesh::convert_to_boundary_node(), e(), oomph::Mesh::finite_element_pt(), oomph::Node::get_coordinates_on_boundary(), ElasticTwoLayerMesh< ELEMENT >::Interface_lower_boundary_element_pt, ElasticTwoLayerMesh< ELEMENT >::Interface_upper_boundary_element_pt, ElasticTwoLayerMesh< ELEMENT >::Lower_layer_element_pt, n, oomph::Mesh::nboundary_node(), oomph::Node::ncoordinates_on_boundary(), oomph::FiniteElement::node_pt(), oomph::RectangularQuadMesh< ELEMENT >::nx(), oomph::Mesh::remove_boundary_node(), oomph::Node::remove_from_boundary(), oomph::Node::set_coordinates_on_boundary(), oomph::Mesh::set_nboundary(), oomph::Mesh::setup_boundary_element_info(), ElasticTwoLayerMesh< ELEMENT >::Upper_layer_element_pt, oomph::Node::x(), y, and oomph::RectangularQuadMesh< ELEMENT >::Ymin.

Member Function Documentation

◆ interface_lower_boundary_element_pt()

template<class ELEMENT >
FiniteElement* & ElasticTwoLayerMesh< ELEMENT >::interface_lower_boundary_element_pt ( const unsigned long &  i)
inline

Access functions for pointers to elements in bottom layer.

int i
Definition: BiCGSTAB_step_by_step.cpp:9

References i, and ElasticTwoLayerMesh< ELEMENT >::Interface_lower_boundary_element_pt.

◆ interface_lower_face_index_at_boundary()

template<class ELEMENT >
int ElasticTwoLayerMesh< ELEMENT >::interface_lower_face_index_at_boundary ( const unsigned e)
inline

Index of the face of the elements next to the interface in the lower region (always 2)

720  {return 2;}

◆ interface_upper_boundary_element_pt()

template<class ELEMENT >
FiniteElement* & ElasticTwoLayerMesh< ELEMENT >::interface_upper_boundary_element_pt ( const unsigned long &  i)
inline

Access functions for pointers to elements in upper layer.

References i, and ElasticTwoLayerMesh< ELEMENT >::Interface_upper_boundary_element_pt.

◆ interface_upper_face_index_at_boundary()

template<class ELEMENT >
int ElasticTwoLayerMesh< ELEMENT >::interface_upper_face_index_at_boundary ( const unsigned e)
inline

Index of the face of the elements next to the interface in the upper region (always -2)

715  {return -2;}

◆ lower_layer_element_pt()

template<class ELEMENT >
FiniteElement* & ElasticTwoLayerMesh< ELEMENT >::lower_layer_element_pt ( const unsigned long &  i)
inline

Access functions for pointers to elements in bottom layer.

688  {return Lower_layer_element_pt[i];}

References i, and ElasticTwoLayerMesh< ELEMENT >::Lower_layer_element_pt.

◆ ninterface_lower()

template<class ELEMENT >
unsigned long ElasticTwoLayerMesh< ELEMENT >::ninterface_lower ( ) const
inline

Number of elements in top layer.

710  {return Interface_lower_boundary_element_pt.size();}

References ElasticTwoLayerMesh< ELEMENT >::Interface_lower_boundary_element_pt.

◆ ninterface_upper()

template<class ELEMENT >
unsigned long ElasticTwoLayerMesh< ELEMENT >::ninterface_upper ( ) const
inline

Number of elements in upper layer.

706  {return Interface_upper_boundary_element_pt.size();}

References ElasticTwoLayerMesh< ELEMENT >::Interface_upper_boundary_element_pt.

◆ nlower()

template<class ELEMENT >
unsigned long ElasticTwoLayerMesh< ELEMENT >::nlower ( ) const
inline

Number of elements in top layer.

694 {return Lower_layer_element_pt.size();}

References ElasticTwoLayerMesh< ELEMENT >::Lower_layer_element_pt.

◆ nupper()

template<class ELEMENT >
unsigned long ElasticTwoLayerMesh< ELEMENT >::nupper ( ) const
inline

Number of elements in upper layer.

691 {return Upper_layer_element_pt.size();}

References ElasticTwoLayerMesh< ELEMENT >::Upper_layer_element_pt.

◆ upper_layer_element_pt()

template<class ELEMENT >
FiniteElement* & ElasticTwoLayerMesh< ELEMENT >::upper_layer_element_pt ( const unsigned long &  i)
inline

Access functions for pointers to elements in upper layer.

684  {return Upper_layer_element_pt[i];}

References i, and ElasticTwoLayerMesh< ELEMENT >::Upper_layer_element_pt.

Member Data Documentation

◆ Interface_lower_boundary_element_pt

template<class ELEMENT >
Vector<FiniteElement*> ElasticTwoLayerMesh< ELEMENT >::Interface_lower_boundary_element_pt
private

◆ Interface_upper_boundary_element_pt

template<class ELEMENT >
Vector<FiniteElement *> ElasticTwoLayerMesh< ELEMENT >::Interface_upper_boundary_element_pt
private

◆ Lower_layer_element_pt

template<class ELEMENT >
Vector<FiniteElement *> ElasticTwoLayerMesh< ELEMENT >::Lower_layer_element_pt
private

◆ Upper_layer_element_pt

template<class ELEMENT >
Vector<FiniteElement *> ElasticTwoLayerMesh< ELEMENT >::Upper_layer_element_pt
private

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