oomph::ChannelWithLeafletMesh< ELEMENT > Class Template Reference

Channel with leaflet mesh. More...

#include <channel_with_leaflet_mesh.template.h>

+ Inheritance diagram for oomph::ChannelWithLeafletMesh< ELEMENT >:

Public Member Functions

 ChannelWithLeafletMesh (GeomObject *leaflet_pt, const double &lleft, const double &lright, const double &hleaflet, const double &htot, const unsigned &nleft, const unsigned &nright, const unsigned &ny1, const unsigned &ny2, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
virtual ~ChannelWithLeafletMesh ()
 Destructor : empty. More...
 
ChannelWithLeafletDomaindomain_pt ()
 Access function to domain. More...
 
- Public Member Functions inherited from oomph::SimpleRectangularQuadMesh< ELEMENT >
 SimpleRectangularQuadMesh (const unsigned &Nx, const unsigned &Ny, const double &Lx, const double &Ly, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
const unsignednx () const
 Access function for number of elements in x directions. More...
 
const unsignedny () const
 Access function for number of elements in y directions. More...
 
- 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)
 
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...
 

Protected Attributes

ChannelWithLeafletDomainDomain_pt
 Pointer to domain. More...
 
GeomObjectLeaflet_pt
 Pointer to GeomObject that represents the leaflet. More...
 
- 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
 

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...
 
- 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)
 

Detailed Description

template<class ELEMENT>
class oomph::ChannelWithLeafletMesh< ELEMENT >

Channel with leaflet mesh.

Constructor & Destructor Documentation

◆ ChannelWithLeafletMesh()

template<class ELEMENT >
oomph::ChannelWithLeafletMesh< ELEMENT >::ChannelWithLeafletMesh ( GeomObject leaflet_pt,
const double lleft,
const double lright,
const double hleaflet,
const double htot,
const unsigned nleft,
const unsigned nright,
const unsigned ny1,
const unsigned ny2,
TimeStepper time_stepper_pt = &Mesh::Default_TimeStepper 
)

Constructor: Pass pointer to GeomObject that represents the leaflet, the length of the domain to left and right of the leaflet, the height of the leaflet and the overall height of the channel, the number of element columns to the left and right of the leaflet, the number of rows of elements from the bottom of the channel to the end of the leaflet, the number of rows of elements above the end of the leaflet. Timestepper defaults to Steady default Timestepper defined in the Mesh base class

59  : SimpleRectangularQuadMesh<ELEMENT>(
60  nright + nleft, ny1 + ny2, lright + lleft, htot, time_stepper_pt)
61  {
62  // Mesh can only be built with 2D Qelements.
63  MeshChecker::assert_geometric_element<QElementGeometricBase, ELEMENT>(2);
64 
65  // Copy pointer to leaflet into private storage
66  Leaflet_pt = leaflet_pt;
67 
68  // Create the ChannelWithLeafletDomain with the wall
69  // represented by the geometric object
70  Domain_pt = new ChannelWithLeafletDomain(
71  leaflet_pt, lleft, lright, hleaflet, htot, nleft, nright, ny1, ny2);
72 
73 
74  // Total number of (macro/finite)elements
75  unsigned nmacro = (ny1 + ny2) * (nleft + nright);
76 
77  // Loop over all elements and set macro element pointer
78  for (unsigned e = 0; e < nmacro; e++)
79  {
80  // Get pointer to finite element
81  FiniteElement* el_pt = this->finite_element_pt(e);
82 
83  // Set pointer to macro element
84  el_pt->set_macro_elem_pt(this->Domain_pt->macro_element_pt(e));
85  }
86 
87  // Update the nodal positions corresponding to their
88  // macro element representations.
89  this->node_update();
90 
91  // Change the numbers of boundaries
92  this->set_nboundary(7);
93 
94  // Remove the nodes of boundary 0
95  this->remove_boundary_nodes(0);
96 
97  // Get the number of nodes along the element edge from first element
98  unsigned nnode_1d = this->finite_element_pt(0)->nnode_1d();
99 
100  // Boundary 0 will be before the wall, boundary 6 after it
101  for (unsigned e = 0; e < nleft; e++)
102  {
103  for (unsigned i = 0; i < nnode_1d; i++)
104  {
105  Node* nod_pt = this->finite_element_pt(e)->node_pt(i);
106  // do not add the last node : it will go to boundary 6
107  if (e != nleft - 1 || i != 2)
108  {
109  this->add_boundary_node(0, nod_pt);
110  }
111  }
112  }
113 
114  for (unsigned e = nleft; e < nleft + nright; e++)
115  {
116  for (unsigned i = 0; i < nnode_1d; i++)
117  {
118  Node* nod_pt = this->finite_element_pt(e)->node_pt(i);
119  this->add_boundary_node(6, nod_pt);
120  }
121  }
122 
123  // Vector of Lagrangian coordinates used as boundary coordinate
124  Vector<double> zeta(1);
125 
126  // Set the wall as a boundary
127  for (unsigned k = 0; k < ny1; k++)
128  {
129  unsigned e = (nleft + nright) * k + nleft - 1;
130  for (unsigned i = 0; i < nnode_1d; i++)
131  {
132  Node* nod_pt =
133  this->finite_element_pt(e)->node_pt((i + 1) * nnode_1d - 1);
134  this->convert_to_boundary_node(nod_pt);
135  this->add_boundary_node(4, nod_pt);
136 
137  // Set coordinates
138  zeta[0] = double(k) * hleaflet / double(ny1) +
139  double(i) * hleaflet / double(ny1) / double(nnode_1d - 1);
140  nod_pt->set_coordinates_on_boundary(4, zeta);
141  }
142  }
143 
144 
145  // Duplicate the nodes of the wall and assign then as a boundary.
146  // This will make one boundary for the east of the elements at the
147  // left of the wall, and one for the west of the elements at the right
148  // of the wall.
149  // This is necessary to use TaylorHoodElements, because it will allow
150  // a discontinuity of the pressure accross the wall.
151  // We separate the lower element from the rest as we do not want to
152  // add the same node twice to the boundary, and the upper element as its
153  // upper node must be the same one than the node of the upper element
154  // at the right of the wall
155 
156  // Lower element
157  unsigned e = nleft - 1;
158 
159  // Allocate storage for newly created node outside
160  // so we can refer to the most recently created one below.
161  Node* nod_pt = 0;
162 
163  // duplicate the east nodes and add them to the 6th boundary
164  // Add the first node to the 0th boundary (horizontal)
165  for (unsigned i = 0; i < nnode_1d; i++)
166  {
167  nod_pt = this->finite_element_pt(e)->construct_boundary_node(
168  (i + 1) * nnode_1d - 1, time_stepper_pt);
169  this->add_boundary_node(5, nod_pt);
170  if (i == 0)
171  {
172  this->add_boundary_node(0, nod_pt);
173  }
174  this->add_node_pt(nod_pt);
175  // Set coordinates
176  zeta[0] = i * hleaflet / double(ny1) / double(nnode_1d - 1);
177  nod_pt->set_coordinates_on_boundary(5, zeta);
178  }
179 
180 
181  // Other elements just at the left of the wall
182  for (unsigned k = 1; k < ny1; k++)
183  {
184  e = (nleft + nright) * k + nleft - 1;
185 
186  // add the upper node of the previous element
187  this->finite_element_pt(e)->node_pt(nnode_1d - 1) = nod_pt;
188  this->add_boundary_node(5, nod_pt);
189  // Set coordinates
190  zeta[0] = k * hleaflet / double(ny1);
191  nod_pt->set_coordinates_on_boundary(5, zeta);
192 
193  // Loop over other nodes on element's eastern edge
194  for (unsigned i = 1; i < nnode_1d; i++)
195  {
196  // Don't duplicate the node at the very top of the "obstacle"
197  if ((k == ny1 - 1) && (i == nnode_1d - 1))
198  {
199  // Get the node but don't do anything else
200  nod_pt = this->finite_element_pt(e)->node_pt(nnode_1d * nnode_1d - 1);
201  }
202  else
203  {
204  // Overwrite the node with a boundary node
205  nod_pt = this->finite_element_pt(e)->construct_boundary_node(
206  (i + 1) * nnode_1d - 1, time_stepper_pt);
207  this->add_node_pt(nod_pt);
208  }
209 
210  // Add node to boundary
211  this->add_boundary_node(5, nod_pt);
212  // Set coordinates
213  zeta[0] = double(k) * hleaflet / double(ny1) +
214  double(i) * hleaflet / double(ny1) / double(nnode_1d - 1);
215  nod_pt->set_coordinates_on_boundary(5, zeta);
216  }
217  }
218 
219  this->node_update();
220 
221  // Re-setup lookup scheme that establishes which elements are located
222  // on the mesh boundaries (doesn't need to be wiped)
224 
225  // We have parametrised boundary 4 and 5
226  this->Boundary_coordinate_exists[4] = true;
227  this->Boundary_coordinate_exists[5] = true;
228 
229  } // end of constructor
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Array< double, 1, 3 > e(1./3., 0.5, 2.)
ChannelWithLeafletDomain * Domain_pt
Pointer to domain.
Definition: channel_with_leaflet_mesh.template.h:89
GeomObject * Leaflet_pt
Pointer to GeomObject that represents the leaflet.
Definition: channel_with_leaflet_mesh.template.h:92
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
virtual unsigned nnode_1d() const
Definition: elements.h:2218
virtual Node * construct_boundary_node(const unsigned &n)
Definition: elements.h:2538
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
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
void remove_boundary_nodes()
Clear all pointers to boundary nodes.
Definition: mesh.cc:204
virtual void node_update(const bool &update_all_solid_nodes=false)
Definition: mesh.cc:287
void convert_to_boundary_node(Node *&node_pt, const Vector< FiniteElement * > &finite_element_pt)
Definition: mesh.cc:2590
void add_node_pt(Node *const &node_pt)
Add a (pointer to a) node to the mesh.
Definition: mesh.h:611
void setup_boundary_element_info()
Definition: quad_mesh.h:73
char char char int int * k
Definition: level2_impl.h:374
EIGEN_STRONG_INLINE const Eigen::CwiseBinaryOp< Eigen::internal::scalar_zeta_op< typename DerivedX::Scalar >, const DerivedX, const DerivedQ > zeta(const Eigen::ArrayBase< DerivedX > &x, const Eigen::ArrayBase< DerivedQ > &q)
Definition: SpecialFunctionsArrayAPI.h:152

References oomph::Mesh::add_boundary_node(), oomph::Mesh::add_node_pt(), oomph::Mesh::Boundary_coordinate_exists, oomph::FiniteElement::construct_boundary_node(), oomph::Mesh::convert_to_boundary_node(), oomph::ChannelWithLeafletMesh< ELEMENT >::Domain_pt, e(), oomph::Mesh::finite_element_pt(), i, k, oomph::ChannelWithLeafletMesh< ELEMENT >::Leaflet_pt, oomph::Domain::macro_element_pt(), oomph::FiniteElement::nnode_1d(), oomph::FiniteElement::node_pt(), oomph::Mesh::node_update(), oomph::Mesh::remove_boundary_nodes(), oomph::Node::set_coordinates_on_boundary(), oomph::FiniteElement::set_macro_elem_pt(), oomph::Mesh::set_nboundary(), oomph::QuadMeshBase::setup_boundary_element_info(), and Eigen::zeta().

◆ ~ChannelWithLeafletMesh()

template<class ELEMENT >
virtual oomph::ChannelWithLeafletMesh< ELEMENT >::~ChannelWithLeafletMesh ( )
inlinevirtual

Destructor : empty.

79 {}

Member Function Documentation

◆ domain_pt()

Member Data Documentation

◆ Domain_pt

◆ Leaflet_pt


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