oomph::OneDMesh< ELEMENT > Class Template Reference

#include <one_d_mesh.template.h>

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

Public Member Functions

 OneDMesh (const unsigned &n_element, const double &length, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
 OneDMesh (const unsigned &n_element, const double &xmin, const double &xmax, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
- Public Member Functions inherited from oomph::LineMeshBase
 LineMeshBase ()
 Constructor (empty) More...
 
 LineMeshBase (const LineMeshBase &node)=delete
 Broken copy constructor. More...
 
void operator= (const LineMeshBase &)=delete
 Broken assignment operator. More...
 
virtual ~LineMeshBase ()
 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 Member Functions

void check_1d () const
 
void build_mesh (TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 Generic mesh constuction routine, called by all constructors. 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

double Xmin
 Minimum coordinate. More...
 
double Xmax
 Maximum coordinate. More...
 
double Length
 Length of the domain. More...
 
unsigned N
 Number of elements. 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...
 

Detailed Description

template<class ELEMENT>
class oomph::OneDMesh< ELEMENT >

1D mesh consisting of N one-dimensional elements from the QElement family.

\[ x \in [Xmin,Xmax] \]

The mesh has two boundaries:

  • Boundary 0 is at \(x=Xmin\).
  • Boundary 1 is at \(x=Xmax\).

There is one node on each of these boundaries.

Constructor & Destructor Documentation

◆ OneDMesh() [1/2]

template<class ELEMENT >
oomph::OneDMesh< ELEMENT >::OneDMesh ( const unsigned n_element,
const double length,
TimeStepper time_stepper_pt = &Mesh::Default_TimeStepper 
)
inline

Constructor: Pass number of elements, n_element, length of domain, length, and pointer to timestepper (defaults to a Steady timestepper so we don't need to specify one in problems without time-dependence).

61  : Xmin(0.0), Xmax(length), N(n_element)
62  {
63  check_1d();
64 
65  build_mesh(time_stepper_pt);
66  }
void build_mesh(TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Generic mesh constuction routine, called by all constructors.
Definition: one_d_mesh.template.cc:39
unsigned N
Number of elements.
Definition: one_d_mesh.template.h:110
void check_1d() const
Definition: one_d_mesh.template.h:85
double Xmin
Minimum coordinate.
Definition: one_d_mesh.template.h:101
double Xmax
Maximum coordinate.
Definition: one_d_mesh.template.h:104

References oomph::OneDMesh< ELEMENT >::build_mesh(), and oomph::OneDMesh< ELEMENT >::check_1d().

◆ OneDMesh() [2/2]

template<class ELEMENT >
oomph::OneDMesh< ELEMENT >::OneDMesh ( const unsigned n_element,
const double xmin,
const double xmax,
TimeStepper time_stepper_pt = &Mesh::Default_TimeStepper 
)
inline

Constructor: Pass number of elements, n_element, minimum coordinate, xmin, maximum coordinate, xmax, and a pointer to a timestepper.

75  : Xmin(xmin), Xmax(xmax), N(n_element)
76  {
77  check_1d();
78 
79  build_mesh(time_stepper_pt);
80  }

References oomph::OneDMesh< ELEMENT >::build_mesh(), and oomph::OneDMesh< ELEMENT >::check_1d().

Member Function Documentation

◆ build_mesh()

template<class ELEMENT >
void oomph::OneDMesh< ELEMENT >::build_mesh ( TimeStepper time_stepper_pt = &Mesh::Default_TimeStepper)
protected

Generic mesh constuction routine, called by all constructors.

The generic mesh construction routine — this contains all the hard work and is called by all constructors

40  {
41  // Set the length of the domain
42  Length = Xmax - Xmin;
43 
44  // Set the number of boundaries -- there are 2 boundaries in a 1D mesh
45  set_nboundary(2);
46 
47  // Allocate storage for the pointers to the elements
48  Element_pt.resize(N);
49 
50  // Allocate memory for the first element
51  Element_pt[0] = new ELEMENT;
52 
53  // Read out the number of nodes in the element (the member function
54  // nnode_1d() is implemented in QElement)
55  const unsigned n_node =
56  dynamic_cast<ELEMENT*>(finite_element_pt(0))->nnode_1d();
57 
58  // We can now allocate storage for the pointers to the nodes in the mesh
59  Node_pt.resize(1 + (n_node - 1) * N);
60 
61  // Initialise the node counter
62  unsigned node_count = 0;
63 
64  // Initialise the minimum x coordinate in the mesh
65  const double xinit = Xmin;
66 
67  // Calculate the length of the element
68  const double el_length = Length / double(N);
69 
70  // Allocate storage for the local coordinate in the element
71  Vector<double> s_fraction;
72 
73  // If the number of elements is 1, the first element is also the
74  // last element
75  if (N == 1)
76  {
77  // Set the first node
78  // ------------------
79 
80  // Allocate memory for the node, using the element's own construct_node
81  // function -- only the element knows what type of nodes it needs!
82  Node_pt[node_count] =
83  finite_element_pt(0)->construct_boundary_node(0, time_stepper_pt);
84 
85  // Set the position of the node
86  node_pt(node_count)->x(0) = xinit;
87 
88  // Add the node to the boundary 0
89  add_boundary_node(0, Node_pt[node_count]);
90 
91  // Increment the counter for the nodes
92  node_count++;
93 
94  // Now build central nodes (ignore last one which needs special
95  // ------------------------------------------------------------
96  // treatment because it's on the boundary)
97  // ---------------------------------------
98  for (unsigned jnod = 1; jnod < (n_node - 1); jnod++)
99  {
100  // Allocate memory for nodes, as before
101  Node_pt[node_count] =
102  finite_element_pt(0)->construct_node(jnod, time_stepper_pt);
103 
104  // Get the local coordinate of the node
105  finite_element_pt(0)->local_fraction_of_node(jnod, s_fraction);
106 
107  // Set the position of the node (linear mapping)
108  node_pt(node_count)->x(0) = xinit + el_length * s_fraction[0];
109 
110  // Increment the node counter
111  node_count++;
112  }
113 
114  // New final node
115  // --------------
116 
117  // Allocate memory for the node, as before
119  n_node - 1, time_stepper_pt);
120 
121  // Set the position of the node
122  node_pt(node_count)->x(0) = xinit + Length;
123 
124  // Add the node to the boundary 1
125  add_boundary_node(1, Node_pt[node_count]);
126 
127  // Increment the node counter
128  node_count++;
129  }
130 
131  // Otherwise, i.e. if there is more than one element, build all elements
132  else
133  {
134  // -------------
135  // FIRST ELEMENT
136  // -------------
137 
138  // Set the first node
139  // ------------------
140 
141  // Allocate memory for the node, using the element's own construct_node
142  // function -- only the element knows what type of nodes it needs!
143  Node_pt[node_count] =
144  finite_element_pt(0)->construct_boundary_node(0, time_stepper_pt);
145 
146  // Set the position of the node
147  node_pt(node_count)->x(0) = xinit;
148 
149  // Add the node to the boundary 0
150  add_boundary_node(0, Node_pt[node_count]);
151 
152  // Increment the counter for the nodes
153  node_count++;
154 
155  // Now build the other nodes in the first element
156  // ----------------------------------------------
157 
158  // Loop over the other nodes in the first element
159  for (unsigned jnod = 1; jnod < n_node; jnod++)
160  {
161  // Allocate memory for the nodes
162  Node_pt[node_count] =
163  finite_element_pt(0)->construct_node(jnod, time_stepper_pt);
164 
165  // Get the local coordinate of the node
166  finite_element_pt(0)->local_fraction_of_node(jnod, s_fraction);
167 
168  // Set the position of the node (linear mapping)
169  node_pt(node_count)->x(0) = xinit + el_length * s_fraction[0];
170 
171  // Increment the node counter
172  node_count++;
173  }
174 
175  // ----------------
176  // CENTRAL ELEMENTS
177  // ----------------
178 
179  // Loop over central elements in mesh
180  for (unsigned e = 1; e < (N - 1); e++)
181  {
182  // Allocate memory for the new element
183  Element_pt[e] = new ELEMENT;
184 
185  // The first node is the same as the last node in the neighbouring
186  // element on the left
188  finite_element_pt(e - 1)->node_pt((n_node - 1));
189 
190  // Loop over the remaining nodes in the element
191  for (unsigned jnod = 1; jnod < n_node; jnod++)
192  {
193  // Allocate memory for the nodes, as before
194  Node_pt[node_count] =
195  finite_element_pt(e)->construct_node(jnod, time_stepper_pt);
196 
197  // Get the local coordinate of the nodes
198  finite_element_pt(e)->local_fraction_of_node(jnod, s_fraction);
199 
200  // Set the position of the node (linear mapping)
201  node_pt(node_count)->x(0) = xinit + el_length * (e + s_fraction[0]);
202 
203  // Increment the node counter
204  node_count++;
205  }
206  } // End of loop over central elements
207 
208 
209  // FINAL ELEMENT
210  //--------------
211 
212  // Allocate memory for element
213  Element_pt[N - 1] = new ELEMENT;
214 
215  // The first node is the same as the last node in the neighbouring
216  // element on the left
217  finite_element_pt(N - 1)->node_pt(0) =
218  finite_element_pt(N - 2)->node_pt(n_node - 1);
219 
220  // New central nodes (ignore last one which needs special treatment
221  // because it's on the boundary)
222  for (unsigned jnod = 1; jnod < (n_node - 1); jnod++)
223  {
224  // Allocate memory for nodes, as before
225  Node_pt[node_count] =
226  finite_element_pt(N - 1)->construct_node(jnod, time_stepper_pt);
227 
228  // Get the local coordinate of the node
229  finite_element_pt(N - 1)->local_fraction_of_node(jnod, s_fraction);
230 
231  // Set the position of the node
232  node_pt(node_count)->x(0) = xinit + el_length * (N - 1 + s_fraction[0]);
233 
234  // Increment the node counter
235  node_count++;
236  }
237 
238  // New final node
239  // --------------
240 
241  // Allocate memory for the node, as before
243  n_node - 1, time_stepper_pt);
244 
245  // Set the position of the node
246  node_pt(node_count)->x(0) = xinit + Length;
247 
248  // Add the node to the boundary 1
249  add_boundary_node(1, Node_pt[node_count]);
250 
251  // Increment the node counter
252  node_count++;
253  }
254  }
Array< double, 1, 3 > e(1./3., 0.5, 2.)
virtual Node * construct_node(const unsigned &n)
Definition: elements.h:2509
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
virtual void local_fraction_of_node(const unsigned &j, Vector< double > &s_fraction)
Definition: elements.cc:3191
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
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
Node *& node_pt(const unsigned long &n)
Return pointer to global node n.
Definition: mesh.h:436
Vector< GeneralisedElement * > Element_pt
Vector of pointers to generalised elements.
Definition: mesh.h:186
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition: nodes.h:1060
double Length
Length of the domain.
Definition: one_d_mesh.template.h:107

References e(), Global_Physical_Variables::Length, and N.

Referenced by oomph::OneDMesh< ELEMENT >::OneDMesh().

◆ check_1d()

template<class ELEMENT >
void oomph::OneDMesh< ELEMENT >::check_1d ( ) const
inlineprotected

Mesh can only be built with 1D elements (but can be either T or Q so can't use the normal assert_geometric_element function.

86  {
87 #ifdef PARANOID
88  FiniteElement* el_pt = new ELEMENT;
89  if (el_pt->dim() != 1)
90  {
91  std::string err = "OneDMesh is only for 1D elements";
92  throw OomphLibError(
94  }
95  delete el_pt;
96  el_pt = 0;
97 #endif
98  }
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References oomph::FiniteElement::dim(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::Global_string_for_annotation::string().

Referenced by oomph::OneDMesh< ELEMENT >::OneDMesh().

Member Data Documentation

◆ Length

template<class ELEMENT >
double oomph::OneDMesh< ELEMENT >::Length
protected

Length of the domain.

◆ N

template<class ELEMENT >
unsigned oomph::OneDMesh< ELEMENT >::N
protected

Number of elements.

◆ Xmax

template<class ELEMENT >
double oomph::OneDMesh< ELEMENT >::Xmax
protected

Maximum coordinate.

◆ Xmin

template<class ELEMENT >
double oomph::OneDMesh< ELEMENT >::Xmin
protected

Minimum coordinate.


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