oomph::WomersleyMesh< WOMERSLEY_ELEMENT > Class Template Reference

#include <womersley_elements.h>

+ Inheritance diagram for oomph::WomersleyMesh< WOMERSLEY_ELEMENT >:

Public Member Functions

 WomersleyMesh (Mesh *n_st_outflow_mesh_pt, TimeStepper *time_stepper_pt, const unsigned &fixed_coordinate, const unsigned &w_index)
 
- 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)
 
virtual void scale_mesh (const double &factor)
 
 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...
 

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::TemplateFreeWomersleyMeshBase
static bool Suppress_warning_about_unpinned_nst_dofs
 Static bool to suppress warning. 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)
 
- 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 WOMERSLEY_ELEMENT>
class oomph::WomersleyMesh< WOMERSLEY_ELEMENT >

///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// Mesh of Womersley elements whose topology, nodal position etc. matches that of a given mesh of face elements in the outflow cross-section of a full Navier-Stokes mesh.

Constructor & Destructor Documentation

◆ WomersleyMesh()

template<class WOMERSLEY_ELEMENT >
oomph::WomersleyMesh< WOMERSLEY_ELEMENT >::WomersleyMesh ( Mesh n_st_outflow_mesh_pt,
TimeStepper time_stepper_pt,
const unsigned fixed_coordinate,
const unsigned w_index 
)
inline

Constructor: Pass pointer to mesh of face elements in the outflow cross-section of a full Navier-Stokes mesh, the timestepper to be used for the Womersley elements, the coordinate (in the higher-dimensional Navier-Stokes mesh) that is constant in the outflow cross-section and the velocity component (in terms of the nodal index) that represents the outflow component – the latter is used to automatically apply the boundary conditions for the Womersley problem.

How many elements and nodes are there in the original mesh?

Create new elements

1620  {
1622  unsigned nelem = n_st_outflow_mesh_pt->nelement();
1623 
1624  // Navier-Stokes outflow mesh may not have any nodes stored (it usually
1625  // just acts as a container for traction elements) -->
1626  // Count number of distinct Navier-Stokes nodes by adding
1627  // the elements' nodes to a set
1628  std::set<Node*> n_st_nodes;
1629  for (unsigned e = 0; e < nelem; e++)
1630  {
1631  FiniteElement* el_pt = n_st_outflow_mesh_pt->finite_element_pt(e);
1632  unsigned nnod_el = el_pt->nnode();
1633  for (unsigned j = 0; j < nnod_el; j++)
1634  {
1635  n_st_nodes.insert(el_pt->node_pt(j));
1636 
1637  // Careful: It there are hanging nodes this won't work!
1638  if (el_pt->node_pt(j)->is_hanging())
1639  {
1640  throw OomphLibError(
1641  "Cannot build WomersleyMesh from mesh with hanging nodes!",
1644  }
1645  }
1646  }
1647 
1648  // Extract size then wipe
1649  unsigned nnode_n_st = n_st_nodes.size();
1650  n_st_nodes.clear();
1651 
1652  // Create enough storage
1653  Node_pt.resize(nnode_n_st);
1654 
1656  for (unsigned e = 0; e < nelem; e++)
1657  {
1658  add_element_pt(new WOMERSLEY_ELEMENT);
1659 #ifdef PARANOID
1660  if (finite_element_pt(e)->nnode() !=
1661  n_st_outflow_mesh_pt->finite_element_pt(e)->nnode())
1662  {
1663  throw OomphLibError(
1664  "Number of nodes in existing and new elements don't match",
1667  }
1668 #endif
1669  }
1670 
1671  // Map to record which Navier-Stokes nodes have been visited (default
1672  // return is false)
1673  std::map<Node*, bool> n_st_node_done;
1674 
1675  // Map to store the Womersley node that corresponds to a
1676  // Navier Stokes node
1677  std::map<Node*, Node*> equivalent_womersley_node_pt;
1678 
1679  // Initialise count of newly created nodes
1680  unsigned node_count = 0;
1681 
1682 
1683  // This is awkward do diagnose: We're assuming that
1684  // the boundary conditions have been applied for the
1685  // underlying Navier-Stokes problem before calling
1686  // this function (otherwise it's really tricky to
1687  // apply the right boundary conditions here), but it's
1688  // hard to police. Issue definite (but suppressable)
1689  // warning if nothing has been pinned at all.
1690  unsigned n_pinned_nodes = 0;
1691 
1692  // Loop over nst and womersley elements in tandem to sort out
1693  // which new nodes are required
1694  for (unsigned e = 0; e < nelem; e++)
1695  {
1696  FiniteElement* n_st_el_pt = n_st_outflow_mesh_pt->finite_element_pt(e);
1697  unsigned nnod_el = n_st_el_pt->nnode();
1698  for (unsigned j = 0; j < nnod_el; j++)
1699  {
1700  // Has the Navier Stokes node been done yet?
1701  Node* n_st_node_pt = n_st_el_pt->node_pt(j);
1702 
1703  // Hasn't been done: Create new node in Womersley element
1704  if (!n_st_node_done[n_st_node_pt])
1705  {
1706  // Create a new node in the Womersley element
1707  Node* new_node_pt =
1708  finite_element_pt(e)->construct_node(j, time_stepper_pt);
1709 
1710  // Add newly created node
1711  Node_pt[node_count] = new_node_pt;
1712  node_count++;
1713 
1714 
1715  // Set coordinates
1716  unsigned dim = n_st_node_pt->ndim();
1717  unsigned icount = 0;
1718  for (unsigned i = 0; i < dim; i++)
1719  {
1720  if (i != fixed_coordinate)
1721  {
1722  new_node_pt->x(icount) = n_st_node_pt->x(i);
1723  icount++;
1724  }
1725  }
1726 
1727  // Set pin status
1728  if (n_st_node_pt->is_pinned(w_index))
1729  {
1730  new_node_pt->pin(0);
1731  n_pinned_nodes++;
1732  }
1733  else
1734  {
1735  // shouldn't need this, but...
1736  new_node_pt->unpin(0);
1737  }
1738 
1739  // Record which Womersley node the
1740  // Navier Stokes node is associated with
1741  equivalent_womersley_node_pt[n_st_node_pt] = new_node_pt;
1742 
1743  // Now the Navier-Stokes node has been done
1744  n_st_node_done[n_st_node_pt] = true;
1745  }
1746  // The node has already been done -- set pointer to existing
1747  // Womersley node
1748  else
1749  {
1751  equivalent_womersley_node_pt[n_st_node_pt];
1752  }
1753  }
1754 
1755  bool passed = true;
1757  if (!passed)
1758  {
1759  // Reverse the nodes
1760  unsigned nnod = finite_element_pt(e)->nnode();
1761  Vector<Node*> orig_nod_pt(nnod);
1762  for (unsigned j = 0; j < nnod; j++)
1763  {
1764  orig_nod_pt[j] = finite_element_pt(e)->node_pt(j);
1765  }
1766  for (unsigned j = 0; j < nnod; j++)
1767  {
1768  finite_element_pt(e)->node_pt(j) = orig_nod_pt[nnod - j - 1];
1769  }
1770  bool passed = true;
1772  if (!passed)
1773  {
1774  throw OomphLibError("Element remains inverted even after reversing "
1775  "the local node numbers",
1778  }
1779  }
1780  }
1781 
1782 
1783 #ifdef PARANOID
1785  {
1786  if (n_pinned_nodes == 0)
1787  {
1788  std::stringstream bla;
1789  bla << "Boundary conditions must be applied in Navier-Stokes\n"
1790  << "problem before attaching impedance elements.\n"
1791  << "Note: This warning can be suppressed by setting the\n"
1792  << "global static boolean\n\n"
1793  << " "
1794  "TemplateFreeWomersleyMeshBase::Suppress_warning_about_"
1795  "unpinned_nst_dofs\n\n"
1796  << "to true\n";
1797  OomphLibWarning(
1799  }
1800  }
1801 #endif
1802 
1803 #ifdef PARANOID
1804  if (nnode() != nnode_n_st)
1805  {
1806  throw OomphLibError(
1807  "Number of nodes in the new mesh don't match that in the old one",
1810  }
1811 #endif
1812  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Array< double, 1, 3 > e(1./3., 0.5, 2.)
virtual Node * construct_node(const unsigned &n)
Definition: elements.h:2509
void check_J_eulerian_at_knots(bool &passed) const
Definition: elements.cc:4237
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
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
unsigned long nnode() const
Return number of nodes in the mesh.
Definition: mesh.h:596
void add_element_pt(GeneralisedElement *const &element_pt)
Add a (pointer to) an element to the mesh.
Definition: mesh.h:617
static bool Suppress_warning_about_unpinned_nst_dofs
Static bool to suppress warning.
Definition: womersley_elements.h:1589
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References oomph::Mesh::add_element_pt(), oomph::FiniteElement::check_J_eulerian_at_knots(), oomph::FiniteElement::construct_node(), e(), oomph::Mesh::finite_element_pt(), i, oomph::Node::is_hanging(), oomph::Data::is_pinned(), j, oomph::Node::ndim(), oomph::Mesh::nelement(), oomph::FiniteElement::nnode(), oomph::Mesh::nnode(), oomph::FiniteElement::node_pt(), oomph::Mesh::Node_pt, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::Data::pin(), oomph::TemplateFreeWomersleyMeshBase::Suppress_warning_about_unpinned_nst_dofs, oomph::Data::unpin(), and oomph::Node::x().


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