oomph::PMLQuadMeshBase< ELEMENT > Class Template Reference

PML mesh class. Policy class for 2D PML meshes. More...

#include <pml_meshes.h>

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

Public Member Functions

 PMLQuadMeshBase (const unsigned &n_pml_x, const unsigned &n_pml_y, const double &x_pml_start, const double &x_pml_end, const double &y_pml_start, const double &y_pml_end, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 Constructor: Create the underlying rectangular quad mesh. More...
 
void pml_locate_zeta (const Vector< double > &x, FiniteElement *&coarse_mesh_el_pt)
 
- 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)
 
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...
 
- 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 oomph::PMLQuadMeshBase< ELEMENT >

PML mesh class. Policy class for 2D PML meshes.

Constructor & Destructor Documentation

◆ PMLQuadMeshBase()

template<class ELEMENT >
oomph::PMLQuadMeshBase< ELEMENT >::PMLQuadMeshBase ( const unsigned n_pml_x,
const unsigned n_pml_y,
const double x_pml_start,
const double x_pml_end,
const double y_pml_start,
const double y_pml_end,
TimeStepper time_stepper_pt = &Mesh::Default_TimeStepper 
)
inline

Constructor: Create the underlying rectangular quad mesh.

197  : RectangularQuadMesh<ELEMENT>(n_pml_x,
198  n_pml_y,
199  x_pml_start,
200  x_pml_end,
201  y_pml_start,
202  y_pml_end,
203  time_stepper_pt)
204  {
205  }

Member Function Documentation

◆ pml_locate_zeta()

template<class ELEMENT >
void oomph::PMLQuadMeshBase< ELEMENT >::pml_locate_zeta ( const Vector< double > &  x,
FiniteElement *&  coarse_mesh_el_pt 
)
inlinevirtual

Overloaded function to allow the user to locate an element in mesh given the (Eulerian) position of a point. Note, the result may be nonunique if the given point lies on the boundary of two elements in the mesh. Additionally, we assume that the PML mesh is axis aligned when deciding if the given point lies inside the mesh

Implements oomph::PMLMeshBase.

214  {
215  //------------------------------------------
216  // Make sure the point lies inside the mesh:
217  //------------------------------------------
218  // Get the lower x limit of the mesh
219  const double x_min = this->x_min();
220 
221  // Get the upper x limit of the mesh
222  const double x_max = this->x_max();
223 
224  // Get the lower y limit of the mesh
225  const double y_min = this->y_min();
226 
227  // Get the upper y limit of the mesh
228  const double y_max = this->y_max();
229 
230  // Due to roundoff error we will choose a small tolerance which will be
231  // used to determine whether or not the point lies in the mesh
232  double tol = 1.0e-12;
233 
234  // Assuming the mesh is axis-aligned we merely need to check if the
235  // x-coordinate of the input lies in the range [x_min,x_max] and the
236  // y-coordinate of the input lies in the range [y_min,y_max]
237  if ((x[0] < x_min - tol) || (x[0] > x_max + tol) ||
238  (x[1] < y_min - tol) || (x[1] > y_max + tol))
239  {
240  // Open an output string stream
241  std::ostringstream error_message_stream;
242 
243  // Create an error message
244  error_message_stream << "Point does not lie in the mesh." << std::endl;
245 
246  // Throw the error message
247  throw OomphLibError(error_message_stream.str(),
250  }
251 
252  //-----------------------------------------
253  // Collect some information about the mesh:
254  //-----------------------------------------
255  // Find out how many elements there are in the x-direction
256  const unsigned nx = this->nx();
257 
258  // Find out how many elements there are in the y-direction
259  const unsigned ny = this->ny();
260 
261  // Find out how many nodes there are in one direction in each element
262  unsigned nnode_1d = this->finite_element_pt(0)->nnode_1d();
263 
264  // Find out how many nodes there are in each element
265  unsigned nnode = this->finite_element_pt(0)->nnode();
266 
267  // Create a pointer to store the element pointer as a FiniteElement
268  FiniteElement* el_pt = 0;
269 
270  // Vector to store the coordinate of each element corner node which
271  // lies on the bottom boundary of the mesh and varies, i.e. if we're
272  // on the bottom boundary of the PML mesh then the y-coordinate will
273  // be fixed therefore we need only store the x-coordinates of the
274  // corner nodes of each element attached to this boundary
275  Vector<double> bottom_boundary_x_coordinates(nx + 1, 0.0);
276 
277  // Vector to store the coordinate of each element corner node which lies
278  // on the right boundary of the mesh and varies
279  Vector<double> right_boundary_y_coordinates(ny + 1, 0.0);
280 
281  // Recall, we already know the start and end coordinates of the mesh
282  // in the x-direction so we can assign these entries straight away.
283  // Note, these values are exactly the same as those had we instead
284  // considered the upper boundary so it is only necessary to store
285  // the information of the one of these boundaries. A similar argument
286  // holds for the left and right boundaries.
287 
288  // The first entry of bottom_boundary_x_coordinates is:
289  bottom_boundary_x_coordinates[0] = x_min;
290 
291  // The last entry is:
292  bottom_boundary_x_coordinates[nx] = x_max;
293 
294  // The first entry of right_boundary_y_coordinates is:
295  right_boundary_y_coordinates[0] = y_min;
296 
297  // The last entry is:
298  right_boundary_y_coordinates[ny] = y_max;
299 
300  // To collect the remaining entries we need to loop over all of the
301  // elements on the bottom boundary and collect the x-coordinate of
302  // the bottom right node in the element. To avoid assigning the
303  // last entry twice we ignore the last element.
304 
305  // Store the lower boundary ID
306  unsigned lower_boundary_id = 0;
307 
308  // Store the right boundary ID
309  unsigned right_boundary_id = 1;
310 
311  // Loop over the elements on the bottom boundary
312  for (unsigned i = 0; i < nx; i++)
313  {
314  // Assign the element pointer
315  el_pt = this->boundary_element_pt(lower_boundary_id, i);
316 
317  // Store the x-coordinate of the bottom right node in this element
318  bottom_boundary_x_coordinates[i + 1] =
319  el_pt->node_pt(nnode_1d - 1)->x(0);
320  }
321 
322  // To collect the remaining entries we need to loop over all of the
323  // elements on the right boundary and collect the y-coordinate of
324  // the top right node in the element. To avoid assigning the
325  // last entry twice we ignore the last element.
326 
327  // Loop over the elements on the bottom boundary
328  for (unsigned i = 0; i < ny; i++)
329  {
330  // Assign the element pointer
331  el_pt = this->boundary_element_pt(right_boundary_id, i);
332 
333  // Store the y-coordinate of the top right node in this element
334  right_boundary_y_coordinates[i + 1] = el_pt->node_pt(nnode - 1)->x(1);
335  }
336 
337  //---------------------------
338  // Find the matching element:
339  //---------------------------
340  // Boolean variable to indicate that the ID of the element in the
341  // x-direction has been found. Note, the value of this ID must lie
342  // in the range [0,nx]
343  bool element_x_id_has_been_found = false;
344 
345  // Boolean variable to indicate that the ID of the element in the
346  // y-direction has been found. Note, the value of this ID must lie
347  // in the range [0,ny]
348  bool element_y_id_has_been_found = false;
349 
350  // Initialise the ID of the element in the x-direction
351  unsigned element_x_id = 0;
352 
353  // Initialise the ID of the element in the y-direction
354  unsigned element_y_id = 0;
355 
356  // Loop over all of the entries in bottom_boundary_x_coordinates
357  for (unsigned i = 0; i < nx; i++)
358  {
359  // Check if the x-coordinate of the input lies in this interval
360  if ((x[0] >= bottom_boundary_x_coordinates[i]) &&
361  (x[0] <= bottom_boundary_x_coordinates[i + 1]))
362  {
363  // The point lies in the i-th interval so the element ID is i
364  element_x_id = i;
365 
366  // Indicate that the element ID has been found
367  element_x_id_has_been_found = true;
368  }
369  } // for (unsigned i=0;i<nx;i++)
370 
371  // If the element ID hasn't been found
372  if (!element_x_id_has_been_found)
373  {
374  // Open an output string stream
375  std::ostringstream error_message_stream;
376 
377  // Create an error message
378  error_message_stream << "The ID of the element in the x-direction "
379  << "has not been found." << std::endl;
380 
381  // Throw the error message
382  throw OomphLibError(error_message_stream.str(),
385  }
386 
387  // Loop over all of the entries in right_boundary_y_coordinates
388  for (unsigned i = 0; i < ny; i++)
389  {
390  // Check if the y-coordinate of the input lies in this interval
391  if ((x[1] >= right_boundary_y_coordinates[i]) &&
392  (x[1] <= right_boundary_y_coordinates[i + 1]))
393  {
394  // The point lies in the i-th interval so the element ID is i
395  element_y_id = i;
396 
397  // Indicate that the element ID has been found
398  element_y_id_has_been_found = true;
399  }
400  } // for (unsigned i=0;i<ny;i++)
401 
402  // If the element ID hasn't been found
403  if (!element_y_id_has_been_found)
404  {
405  // Open an output string stream
406  std::ostringstream error_message_stream;
407 
408  // Create an error message
409  error_message_stream << "The ID of the element in the y-direction "
410  << "has not been found." << std::endl;
411 
412  // Throw the error message
413  throw OomphLibError(error_message_stream.str(),
416  }
417 
418  // Calculate the element number in the mesh
419  unsigned el_id = element_y_id * nx + element_x_id;
420 
421  // Set the pointer to this element
422  coarse_mesh_el_pt = dynamic_cast<FiniteElement*>(this->element_pt(el_id));
423  } // End of pml_locate_zeta
int i
Definition: BiCGSTAB_step_by_step.cpp:9
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
virtual unsigned nnode_1d() const
Definition: elements.h:2218
FiniteElement * finite_element_pt(const unsigned &e) const
Definition: mesh.h:473
FiniteElement * boundary_element_pt(const unsigned &b, const unsigned &e) const
Return pointer to e-th finite element on boundary b.
Definition: mesh.h:840
unsigned long nnode() const
Return number of nodes in the mesh.
Definition: mesh.h:596
const Vector< GeneralisedElement * > & element_pt() const
Return reference to the Vector of elements.
Definition: mesh.h:460
const double y_min() const
Return the minimum value of y coordinate.
Definition: rectangular_quadmesh.template.h:252
const double x_max() const
Return the maximum value of x coordinate.
Definition: rectangular_quadmesh.template.h:245
const double y_max() const
Return the maximum value of y coordinate.
Definition: rectangular_quadmesh.template.h:259
const unsigned & ny() const
Return number of elements in y direction.
Definition: rectangular_quadmesh.template.h:231
const double x_min() const
Return the minimum value of x coordinate.
Definition: rectangular_quadmesh.template.h:238
const unsigned & nx() const
Return number of elements in x direction.
Definition: rectangular_quadmesh.template.h:224
list x
Definition: plotDoE.py:28
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References oomph::Mesh::boundary_element_pt(), oomph::Mesh::element_pt(), oomph::Mesh::finite_element_pt(), i, oomph::FiniteElement::nnode(), oomph::Mesh::nnode(), oomph::FiniteElement::nnode_1d(), oomph::FiniteElement::node_pt(), oomph::RectangularQuadMesh< ELEMENT >::nx(), oomph::RectangularQuadMesh< ELEMENT >::ny(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, plotDoE::x, oomph::Node::x(), oomph::RectangularQuadMesh< ELEMENT >::x_max(), oomph::RectangularQuadMesh< ELEMENT >::x_min(), oomph::RectangularQuadMesh< ELEMENT >::y_max(), and oomph::RectangularQuadMesh< ELEMENT >::y_min().


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