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

Public Member Functions

CylinderAndInterfaceDomaindomain_pt ()
 
 CylinderAndInterfaceMesh (const double &length, const double &height, TimeStepper *time_stepper_pt)
 
- 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)
 
- 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 setup_boundary_element_info ()
 
virtual void setup_boundary_element_info (std::ostream &outfile)
 
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...
 

Protected Attributes

CylinderAndInterfaceDomainDomain_pt
 
- 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
 

Private Attributes

double Height
 

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::SolidMesh
static SolidICProblem Solid_IC_problem
 
- 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...
 
- 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)
 

Constructor & Destructor Documentation

◆ CylinderAndInterfaceMesh()

template<class ELEMENT >
CylinderAndInterfaceMesh< ELEMENT >::CylinderAndInterfaceMesh ( const double length,
const double height,
TimeStepper time_stepper_pt 
)
inline
453  : Height(height)
454  {
455  //Create the domain
457 
458  //Initialise the node counter
459  unsigned node_count=0;
460  //Vectors Used to get data from domains
461  Vector<double> s(2), r(2);
462 
463  //Setup temporary storage for the Node
464  Vector<Node *> Tmp_node_pt;
465 
466  //Now blindly loop over the macro elements and associate and finite
467  //element with each
468  unsigned Nmacro_element = Domain_pt->nmacro_element();
469  for(unsigned e=0;e<Nmacro_element;e++)
470  {
471  //Create the FiniteElement and add to the Element_pt Vector
472  Element_pt.push_back(new ELEMENT);
473 
474  //Read out the number of linear points in the element
475  unsigned Np =
476  dynamic_cast<ELEMENT*>(finite_element_pt(e))->nnode_1d();
477 
478  //Loop over nodes in the column
479  for(unsigned l1=0;l1<Np;l1++)
480  {
481  //Loop over the nodes in the row
482  for(unsigned l2=0;l2<Np;l2++)
483  {
484  //Allocate the memory for the node
485  Tmp_node_pt.push_back(finite_element_pt(e)->
486  construct_node(l1*Np+l2,time_stepper_pt));
487 
488  //Read out the position of the node from the macro element
489  s[0] = -1.0 + 2.0*(double)l2/(double)(Np-1);
490  s[1] = -1.0 + 2.0*(double)l1/(double)(Np-1);
492 
493  //Set the position of the node
494  Tmp_node_pt[node_count]->x(0) = r[0];
495  Tmp_node_pt[node_count]->x(1) = r[1];
496 
497  //Increment the node number
498  node_count++;
499  }
500  }
501  } //End of loop over macro elements
502 
503  //Now the elements have been created, but there will be nodes in
504  //common, need to loop over the common edges and sort it, by reassigning
505  //pointers and the deleting excess nodes
506 
507  //Read out the number of linear points in the element
508  unsigned Np =
509  dynamic_cast<ELEMENT*>(finite_element_pt(0))->nnode_1d();
510 
511  //DelaunayEdge between Elements 0 and 1
512  for(unsigned n=0;n<Np;n++)
513  {
514  //Set the nodes in element 1 to be the same as in element 0
515  finite_element_pt(1)->node_pt(Np*n)
516  = finite_element_pt(0)->node_pt(n*Np+Np-1);
517 
518  //Remove the nodes in element 1 from the temporaray node list
519  delete Tmp_node_pt[Np*Np + Np*n];
520  Tmp_node_pt[Np*Np + Np*n] = 0;
521  }
522 
523  //DelaunayEdge between Elements 1 and 2
524  for(unsigned n=0;n<Np;n++)
525  {
526  //Set the nodes in element 2 to be the same as in element 1
527  finite_element_pt(2)->node_pt(n*Np)
528  = finite_element_pt(1)->node_pt((Np-1)*Np+Np-1-n);
529 
530  //Remove the nodes in element 2 from the temporaray node list
531  delete Tmp_node_pt[2*Np*Np + n*Np];
532  Tmp_node_pt[2*Np*Np + n*Np] = 0;
533  }
534 
535  //DelaunayEdge between Elements 1 and 4
536  for(unsigned n=0;n<Np;n++)
537  {
538  //Set the nodes in element 4 to be the same as in element 1
539  finite_element_pt(4)->node_pt(n*Np)
540  = finite_element_pt(1)->node_pt(n);
541 
542  //Remove the nodes in element 2 from the temporaray node list
543  delete Tmp_node_pt[4*Np*Np + n*Np];
544  Tmp_node_pt[4*Np*Np + n*Np] = 0;
545  }
546 
547  //DelaunayEdge between Element 2 and 3
548  for(unsigned n=0;n<Np;n++)
549  {
550  //Set the nodes in element 3 to be the same as in element 2
551  finite_element_pt(3)->node_pt(Np*(Np-1)+n)
552  = finite_element_pt(2)->node_pt(Np*n+Np-1);
553 
554  //Remove the nodes in element 3 from the temporaray node list
555  delete Tmp_node_pt[3*Np*Np + Np*(Np-1)+n];
556  Tmp_node_pt[3*Np*Np + Np*(Np-1)+n] = 0;
557  }
558 
559 
560  //DelaunayEdge between Element 4 and 3
561  for(unsigned n=0;n<Np;n++)
562  {
563  //Set the nodes in element 3 to be the same as in element 4
565  = finite_element_pt(4)->node_pt(Np*(Np-n-1)+Np-1);
566 
567  //Remove the nodes in element 3 from the temporaray node list
568  delete Tmp_node_pt[3*Np*Np + n];
569  Tmp_node_pt[3*Np*Np + n] = 0;
570  }
571 
572 
573  //DelaunayEdge between Element 3 and 5
574  for(unsigned n=0;n<Np;n++)
575  {
576  //Set the nodes in element 5 to be the same as in element 3
577  finite_element_pt(5)->node_pt(n*Np)
578  = finite_element_pt(3)->node_pt(Np*n+Np-1);
579 
580  //Remove the nodes in element 5 from the temporaray node list
581  delete Tmp_node_pt[5*Np*Np + n*Np];
582  Tmp_node_pt[5*Np*Np + n*Np] = 0;
583  }
584 
585  //Now set the actual true nodes
586  for(unsigned n=0;n<node_count;n++)
587  {
588  if(Tmp_node_pt[n]!=0) {Node_pt.push_back(Tmp_node_pt[n]);}
589  }
590 
591 
592 
593  //Finally set the nodes on the boundaries
594  set_nboundary(5);
595 
596  for(unsigned n=0;n<Np;n++)
597  {
598  //Left hand side
599  Node* temp_node_pt = finite_element_pt(0)->node_pt(n*Np);
600  this->convert_to_boundary_node(temp_node_pt);
601  add_boundary_node(3,temp_node_pt);
602 
603 
604  //Right hand side
605  temp_node_pt = finite_element_pt(5)->node_pt(n*Np+Np-1);
606  this->convert_to_boundary_node(temp_node_pt);
607  add_boundary_node(1,temp_node_pt);
608 
609  //LH part of lower boundary
610  temp_node_pt = finite_element_pt(0)->node_pt(n);
611  this->convert_to_boundary_node(temp_node_pt);
612  add_boundary_node(0,temp_node_pt);
613 
614  //First part of upper boundary
615  temp_node_pt = finite_element_pt(0)->node_pt(Np*(Np-1)+n);
616  this->convert_to_boundary_node(temp_node_pt);
617  add_boundary_node(2,temp_node_pt);
618 
619  //First part of hole boundary
620  temp_node_pt = finite_element_pt(4)->node_pt(Np*(Np-1)+n);
621  this->convert_to_boundary_node(temp_node_pt);
622  add_boundary_node(4,temp_node_pt);
623  }
624 
625  for(unsigned n=1;n<Np;n++)
626  {
627  //Middle of lower boundary
628  Node* temp_node_pt = finite_element_pt(4)->node_pt(n);
629  this->convert_to_boundary_node(temp_node_pt);
630  add_boundary_node(0,temp_node_pt);
631 
632  //Middle of upper boundary
633  temp_node_pt = finite_element_pt(2)->node_pt(Np*(Np-1)+n);
634  this->convert_to_boundary_node(temp_node_pt);
635  add_boundary_node(2,temp_node_pt);
636 
637  //Next part of hole
638  temp_node_pt = finite_element_pt(3)->node_pt(n*Np);
639  this->convert_to_boundary_node(temp_node_pt);
640  add_boundary_node(4,temp_node_pt);
641  }
642 
643  for(unsigned n=1;n<Np;n++)
644  {
645  //Final part of lower boundary
646  Node* temp_node_pt = finite_element_pt(5)->node_pt(n);
647  this->convert_to_boundary_node(temp_node_pt);
648  add_boundary_node(0,temp_node_pt);
649 
650  //Middle of upper boundary
651  temp_node_pt = finite_element_pt(5)->node_pt(Np*(Np-1)+n);
652  this->convert_to_boundary_node(temp_node_pt);
653  add_boundary_node(2,temp_node_pt);
654 
655  //Next part of hole
656  temp_node_pt = finite_element_pt(2)->node_pt(Np-n-1);
657  this->convert_to_boundary_node(temp_node_pt);
658  add_boundary_node(4,temp_node_pt);
659  }
660 
661  for(unsigned n=1;n<Np-1;n++)
662  {
663  //Final part of hole
664  Node* temp_node_pt = finite_element_pt(1)->node_pt(Np*(Np-n-1)+Np-1);
665  this->convert_to_boundary_node(temp_node_pt);
666  add_boundary_node(4,temp_node_pt);
667  }
668 
669  //Now loop over all the nodes and set their Lagrangian coordinates
670  unsigned Nnode = nnode();
671  for(unsigned n=0;n<Nnode;n++)
672  {
673  //Cast node to an elastic node
674  SolidNode* temp_pt = static_cast<SolidNode*>(Node_pt[n]);
675  for(unsigned i=0;i<2;i++)
676  {temp_pt->xi(i) = temp_pt->x(i);}
677  }
678  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Array< double, 1, 3 > e(1./3., 0.5, 2.)
Definition: adaptive_interface.cc:99
double Height
Definition: adaptive_interface.cc:440
CylinderAndInterfaceDomain * Domain_pt
Definition: adaptive_interface.cc:444
unsigned nmacro_element()
Number of macro elements in domain.
Definition: domain.h:123
MacroElement * macro_element_pt(const unsigned &i)
Access to i-th macro element.
Definition: domain.h:116
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
void macro_map(const Vector< double > &s, Vector< double > &r)
The mapping from local to global coordinates at the current time : r(s)
Definition: macro_element.h:126
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
Vector< Node * > Node_pt
Vector of pointers to nodes.
Definition: mesh.h:183
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
unsigned long nnode() const
Return number of nodes in the mesh.
Definition: mesh.h:596
void convert_to_boundary_node(Node *&node_pt, const Vector< FiniteElement * > &finite_element_pt)
Definition: mesh.cc:2590
Vector< GeneralisedElement * > Element_pt
Vector of pointers to generalised elements.
Definition: mesh.h:186
Definition: nodes.h:906
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition: nodes.h:1060
Definition: nodes.h:1686
double & xi(const unsigned &i)
Reference to i-th Lagrangian position.
Definition: nodes.h:1883
RealScalar s
Definition: level1_cplx_impl.h:130
double height(const double &x)
Height of domain.
Definition: simple_spine_channel.cc:429
r
Definition: UniformPSDSelfTest.py:20

References e(), Global_Physical_Variables::height(), i, oomph::Domain::macro_element_pt(), oomph::MacroElement::macro_map(), n, oomph::Domain::nmacro_element(), UniformPSDSelfTest::r, s, oomph::Node::x(), and oomph::SolidNode::xi().

Member Function Documentation

◆ domain_pt()

template<class ELEMENT >
CylinderAndInterfaceDomain* CylinderAndInterfaceMesh< ELEMENT >::domain_pt ( )
inline
449 {return Domain_pt;}

Member Data Documentation

◆ Domain_pt

template<class ELEMENT >
CylinderAndInterfaceDomain* CylinderAndInterfaceMesh< ELEMENT >::Domain_pt
protected

◆ Height

template<class ELEMENT >
double CylinderAndInterfaceMesh< ELEMENT >::Height
private

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