oomph::ChannelSpineMesh< ELEMENT > Class Template Reference

#include <channel_spine_mesh.template.h>

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

Public Member Functions

 ChannelSpineMesh (const unsigned &nx0, const unsigned &nx1, const unsigned &nx2, const unsigned &ny, const double &lx0, const double &lx1, const double &lx2, const double &h, GeomObject *wall_pt, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
 ChannelSpineMesh (const unsigned &nx0, const unsigned &nx1, const unsigned &nx2, const unsigned &ny, const double &lx0, const double &lx1, const double &lx2, const double &h, GeomObject *wall_pt, const bool &periodic_in_x, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
FiniteElement *& left_element_pt (const unsigned long &i)
 
FiniteElement *& centre_element_pt (const unsigned long &i)
 
FiniteElement *& right_element_pt (const unsigned long &i)
 
unsigned long nleft () const
 Number of elements in left region. More...
 
unsigned long ncentre () const
 Number of elements in centre region. More...
 
unsigned long nright () const
 Number of elements in right region. More...
 
unsigned long nbulk () const
 Number of elements in bulk. More...
 
void element_reorder ()
 
virtual void spine_node_update (SpineNode *spine_node_pt)
 
virtual double x_spacing_function (unsigned xelement, unsigned xnode, unsigned yelement, unsigned ynode)
 
Spine *& left_spine_pt (const unsigned long &i)
 Access function for spines in left region. More...
 
Spine *& centre_spine_pt (const unsigned long &i)
 Access function for spines in centre region. More...
 
Spine *& right_spine_pt (const unsigned long &i)
 Access function for spines in right region. More...
 
unsigned nleft_spine ()
 Access function for the number of spines in the left region. More...
 
unsigned ncentre_spine ()
 Access function for the number of spines in the centre region. More...
 
unsigned nright_spine ()
 Access function for the number of spines in the right region. More...
 
GeomObjectwall_pt ()
 Access function to the GeomObject for upper wall. More...
 
GeomObjectstraight_wall_pt ()
 Access function to the GeomObject for the straight upper wall. More...
 
- 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 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 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...
 
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)
 
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...
 
- Public Member Functions inherited from oomph::SpineMesh
virtual ~SpineMesh ()
 Destructor to clean up the memory allocated to the spines. More...
 
Spine *& spine_pt (const unsigned long &i)
 Return the i-th spine in the mesh. More...
 
const Spinespine_pt (const unsigned long &i) const
 Return the i-th spine in the mesh (const version) More...
 
unsigned long nspine () const
 Return the number of spines in the mesh. More...
 
void add_spine_pt (Spine *const &spine_pt)
 Add a spine to the mesh. More...
 
SpineNodenode_pt (const unsigned long &n)
 Return a pointer to the n-th global SpineNode. More...
 
SpineNodeelement_node_pt (const unsigned long &e, const unsigned &n)
 
unsigned long assign_global_spine_eqn_numbers (Vector< double * > &Dof_pt)
 Assign spines to Spine_pt vector of element. More...
 
void describe_spine_dofs (std::ostream &out, const std::string &current_string) const
 
void set_mesh_level_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
void set_spine_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 Assign time stepper to spines data. More...
 
void set_consistent_pinned_spine_values_for_continuation (ContinuationStorageScheme *const &continuation_stepper_pt)
 
bool does_pointer_correspond_to_spine_data (double *const &parameter_pt)
 
void node_update (const bool &update_all_solid_nodes=false)
 
void dump (std::ofstream &dump_file) const
 Overload the dump function so that the spine data is dumped. More...
 
void read (std::ifstream &restart_file)
 Overload the read function so that the spine data is also read. More...
 

Protected Member Functions

virtual void build_channel_spine_mesh (TimeStepper *time_stepper_pt)
 
- 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

Vector< FiniteElement * > Left_element_pt
 Vector of pointers to element in the left region. More...
 
Vector< FiniteElement * > Centre_element_pt
 Vector of pointers to element in the centre region. More...
 
Vector< FiniteElement * > Right_element_pt
 Vector of pointers to element in the right region. More...
 
unsigned Nx0
 Number of elements in the left region. More...
 
unsigned Nx1
 Number of elements in the centre region. More...
 
unsigned Nx2
 Number of elements in the right region. More...
 
double Lx0
 Length of left region. More...
 
double Lx1
 Length of centre region. More...
 
double Lx2
 Length of right region. More...
 
unsigned Nleft_spine
 Number of spines in left region. More...
 
unsigned Ncentre_spine
 Number of spines in centre region. More...
 
unsigned Nright_spine
 Number of spines in right region. More...
 
GeomObjectWall_pt
 GeomObject for upper wall. More...
 
GeomObjectStraight_wall_pt
 GeomObject for the straight upper wall. More...
 
- 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
 
- Protected Attributes inherited from oomph::SpineMesh
Vector< Spine * > Spine_pt
 A Spine mesh contains a Vector of pointers to spines. 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...
 

Detailed Description

template<class ELEMENT>
class oomph::ChannelSpineMesh< ELEMENT >

Spine mesh class derived from standard 2D mesh. The mesh contains a StraightLine GeomObject which defines the height of the left and right regions (0,2) and another GeomObject is passed to the constructor to define the height in the central region.

Constructor & Destructor Documentation

◆ ChannelSpineMesh() [1/2]

template<class ELEMENT >
oomph::ChannelSpineMesh< ELEMENT >::ChannelSpineMesh ( const unsigned nx0,
const unsigned nx1,
const unsigned nx2,
const unsigned ny,
const double lx0,
const double lx1,
const double lx2,
const double h,
GeomObject wall_pt,
TimeStepper time_stepper_pt = &Mesh::Default_TimeStepper 
)

Constructor: Pass number of elements in x-direction in regions 0,1 and 2, number of elements in y-direction, length in x direction in regions 0,1 and 2, height mesh, pointer to the GeomObject defining the heightof the central region and pointer to timestepper (defaults to Steady timestepper)

Constructor for spine 2D mesh: Pass number of elements in x-direction, number of elements in y-direction, axial length and height of layer, and pointer to timestepper (defaults to Static timestepper).

The mesh contains a layer of spinified fluid elements (of type ELEMENT; e.g SpineElement<QCrouzeixRaviartElement<2>) and a surface layer of corresponding Spine interface elements of type INTERFACE_ELEMENT, e.g. SpineLineFluidInterfaceElement<ELEMENT> for 2D planar problems.

58  : RectangularQuadMesh<ELEMENT>(nx0 + nx1 + nx2,
59  ny,
60  0.0,
61  lx0 + lx1 + lx2,
62  0.0,
63  h,
64  false,
65  false,
66  time_stepper_pt),
67  Nx0(nx0),
68  Nx1(nx1),
69  Nx2(nx2),
70  Lx0(lx0),
71  Lx1(lx1),
72  Lx2(lx2),
74  {
75  // Mesh can only be built with 2D Qelements.
76  MeshChecker::assert_geometric_element<QElementGeometricBase, ELEMENT>(2);
77 
78  // Mesh can only be built with spine elements
79  MeshChecker::assert_geometric_element<SpineFiniteElement, ELEMENT>(2);
80 
81  // We've called the "generic" constructor for the RectangularQuadMesh
82  // which doesn't do much...
83 
84  // Build the straight line object
85  Straight_wall_pt = new StraightLine(h);
86 
87  // Now build the mesh:
88  build_channel_spine_mesh(time_stepper_pt);
89  }
virtual void build_channel_spine_mesh(TimeStepper *time_stepper_pt)
Definition: channel_spine_mesh.template.cc:153
unsigned Nx0
Number of elements in the left region.
Definition: channel_spine_mesh.template.h:288
GeomObject * Straight_wall_pt
GeomObject for the straight upper wall.
Definition: channel_spine_mesh.template.h:318
double Lx1
Length of centre region.
Definition: channel_spine_mesh.template.h:300
double Lx0
Length of left region.
Definition: channel_spine_mesh.template.h:297
unsigned Nx2
Number of elements in the right region.
Definition: channel_spine_mesh.template.h:294
double Lx2
Length of right region.
Definition: channel_spine_mesh.template.h:303
unsigned Nx1
Number of elements in the centre region.
Definition: channel_spine_mesh.template.h:291
GeomObject * wall_pt()
Access function to the GeomObject for upper wall.
Definition: channel_spine_mesh.template.h:262
GeomObject * Wall_pt
GeomObject for upper wall.
Definition: channel_spine_mesh.template.h:315
const unsigned & ny() const
Return number of elements in y direction.
Definition: rectangular_quadmesh.template.h:231

References oomph::ChannelSpineMesh< ELEMENT >::build_channel_spine_mesh(), and oomph::ChannelSpineMesh< ELEMENT >::Straight_wall_pt.

◆ ChannelSpineMesh() [2/2]

template<class ELEMENT >
oomph::ChannelSpineMesh< ELEMENT >::ChannelSpineMesh ( const unsigned nx0,
const unsigned nx1,
const unsigned nx2,
const unsigned ny,
const double lx0,
const double lx1,
const double lx2,
const double h,
GeomObject wall_pt,
const bool periodic_in_x,
TimeStepper time_stepper_pt = &Mesh::Default_TimeStepper 
)

Constructor: Pass number of elements in x-direction in regions 0,1 and 2, number of elements in y-direction, length in x direction in regions 0,1 and 2, height mesh, pointer to the GeomObject defining the heightof the central region, a boolean flag to indicate whether or not the mesh is periodic and pointer to timestepper (defaults to Steady timestepper)

Constuctor for spine 2D mesh: Pass number of elements in x-direction, number of elements in y-direction, axial length and height of layer, a boolean flag to make the mesh periodic in the x-direction, and pointer to timestepper (defaults to Static timestepper).

The mesh contains a layer of elements (of type ELEMENT) and a surface layer of corresponding Spine interface elements (of type SpineLineFluidInterfaceElement<ELEMENT>).

114  : RectangularQuadMesh<ELEMENT>(nx0 + nx1 + nx2,
115  ny,
116  0.0,
117  lx0 + lx1 + lx2,
118  0.0,
119  h,
120  periodic_in_x,
121  false,
122  time_stepper_pt),
123  Nx0(nx0),
124  Nx1(nx1),
125  Nx2(nx2),
126  Lx0(lx0),
127  Lx1(lx1),
128  Lx2(lx2),
130  {
131  // Mesh can only be built with 2D Qelements.
132  MeshChecker::assert_geometric_element<QElementGeometricBase, ELEMENT>(2);
133 
134  // Mesh can only be built with spine elements
135  MeshChecker::assert_geometric_element<SpineFiniteElement, ELEMENT>(2);
136 
137  // We've called the "generic" constructor for the RectangularQuadMesh
138  // which doesn't do much...
139 
140  // Build the straight line object
141  Straight_wall_pt = new StraightLine(h);
142 
143  // Now build the mesh:
144  build_channel_spine_mesh(time_stepper_pt);
145  }

References oomph::ChannelSpineMesh< ELEMENT >::build_channel_spine_mesh(), and oomph::ChannelSpineMesh< ELEMENT >::Straight_wall_pt.

Member Function Documentation

◆ build_channel_spine_mesh()

template<class ELEMENT >
void oomph::ChannelSpineMesh< ELEMENT >::build_channel_spine_mesh ( TimeStepper time_stepper_pt)
protectedvirtual

Helper function to actually build the channel-spine mesh (called from various constructors)

Helper function that actually builds the channel-spine mesh based on the parameters set in the various constructors

155  {
156  // Build the underlying quad mesh:
158 
159  // Read out the number of elements in the x-direction and y-direction
160  // and in each of the left, centre and right regions
161  unsigned n_x = this->Nx;
162  unsigned n_y = this->Ny;
163  unsigned n_x0 = this->Nx0;
164  unsigned n_x1 = this->Nx1;
165  unsigned n_x2 = this->Nx2;
166 
167  // Set up the pointers to elements in the left region
168  unsigned nleft = n_x0 * n_y;
169  ;
170  Left_element_pt.reserve(nleft);
171  unsigned ncentre = n_x1 * n_y;
172  ;
173  Centre_element_pt.reserve(ncentre);
174  unsigned nright = n_x2 * n_y;
175  ;
176  Right_element_pt.reserve(nright);
177  for (unsigned irow = 0; irow < n_y; irow++)
178  {
179  for (unsigned e = 0; e < n_x0; e++)
180  {
181  Left_element_pt.push_back(this->finite_element_pt(irow * n_x + e));
182  }
183  for (unsigned e = 0; e < n_x1; e++)
184  {
185  Centre_element_pt.push_back(
186  this->finite_element_pt(irow * n_x + (n_x0 + e)));
187  }
188  for (unsigned e = 0; e < n_x2; e++)
189  {
190  Right_element_pt.push_back(
191  this->finite_element_pt(irow * n_x + (n_x0 + n_x1 + e)));
192  }
193  }
194 
195 #ifdef PARANOID
196  // Check that we have the correct number of elements
197  if (nelement() != nleft + ncentre + nright)
198  {
199  throw OomphLibError("Incorrect number of element pointers!",
202  }
203 #endif
204 
205  // Allocate memory for the spines and fractions along spines
206  //---------------------------------------------------------
207 
208  // Read out number of linear points in the element
209  unsigned n_p = dynamic_cast<ELEMENT*>(finite_element_pt(0))->nnode_1d();
210 
211  unsigned nspine;
212  // Allocate store for the spines:
213  if (this->Xperiodic)
214  {
215  nspine = (n_p - 1) * n_x;
216  Spine_pt.reserve(nspine);
217  // Number of spines in each region
218  // NOTE that boundary spines are in both regions
219  Nleft_spine = (n_p - 1) * n_x0 + 1;
220  Ncentre_spine = (n_p - 1) * n_x1 + 1;
221  Nright_spine = (n_p - 1) * n_x2;
222  }
223  else
224  {
225  nspine = (n_p - 1) * n_x + 1;
226  Spine_pt.reserve(nspine);
227  // Number of spines in each region
228  // NOTE that boundary spines are in both regions
229  Nleft_spine = (n_p - 1) * n_x0 + 1;
230  Ncentre_spine = (n_p - 1) * n_x1 + 1;
231  Nright_spine = (n_p - 1) * n_x2 + 1;
232  }
233 
234  // end Allocating memory
235 
236 
237  // set up the vectors of geometric data & objects for building spines
238  Vector<double> r_wall(2), zeta(1), s_wall(1);
239  GeomObject* geometric_object_pt = 0;
240 
241  // LEFT REGION
242  // ===========
243 
244  // SPINES IN LEFT REGION
245  // ---------------------
246 
247  // Set up zeta increments
248  double zeta_lo = 0.0;
249  double dzeta = Lx0 / n_x0;
250 
251  // Initialise number of elements in previous regions:
252  unsigned n_prev_elements = 0;
253 
254 
255  // FIRST SPINE
256  //-----------
257 
258  // Element 0
259  // Node 0
260  // Assign the new spine with unit length
261  Spine* new_spine_pt = new Spine(1.0);
262  new_spine_pt->spine_height_pt()->pin(0);
263  Spine_pt.push_back(new_spine_pt);
264 
265  // Get pointer to node
266  SpineNode* nod_pt = element_node_pt(0, 0);
267  // Set the pointer to the spine
268  nod_pt->spine_pt() = new_spine_pt;
269  // Set the fraction
270  nod_pt->fraction() = 0.0;
271  // Pointer to the mesh that implements the update fct
272  nod_pt->spine_mesh_pt() = this;
273  // Set update fct id
274  nod_pt->node_update_fct_id() = 0;
275 
276  // Provide spine with additional storage for wall coordinate
277  // and wall geom object:
278 
279  {
280  // Get the Lagrangian coordinate in the Lower Wall
281  zeta[0] = 0.0;
282  // Get the geometric object and local coordinate
283  Straight_wall_pt->locate_zeta(zeta, geometric_object_pt, s_wall);
284 
285  // The local coordinate is a geometric parameter
286  // This needs to be set (rather than added) because the
287  // same spine may be visited more than once
288  Vector<double> parameters(1, s_wall[0]);
289  nod_pt->spine_pt()->set_geom_parameter(parameters);
290 
291  // Get position of wall
292  Straight_wall_pt->position(s_wall, r_wall);
293 
294  // Adjust spine height
295  nod_pt->spine_pt()->height() = r_wall[1];
296 
297  // The sub geom object is one (and only) geom object
298  // for spine:
299  Vector<GeomObject*> geom_object_pt(1);
300  geom_object_pt[0] = geometric_object_pt;
301 
302  // Pass geom object(s) to spine
303  nod_pt->spine_pt()->set_geom_object_pt(geom_object_pt);
304  }
305 
306  // Loop vertically along the spine
307  // Loop over the elements
308  for (unsigned long i = 0; i < n_y; i++)
309  {
310  // Loop over the vertical nodes, apart from the first
311  for (unsigned l1 = 1; l1 < n_p; l1++)
312  {
313  // Get pointer to node
314  SpineNode* nod_pt = element_node_pt(i * n_x, l1 * n_p);
315  // Set the pointer to the spine
316  nod_pt->spine_pt() = new_spine_pt;
317  // Set the fraction
318  nod_pt->fraction() =
319  (double(i) + double(l1) / double(n_p - 1)) / double(n_y);
320  // Pointer to the mesh that implements the update fct
321  nod_pt->spine_mesh_pt() = this;
322  // Set update fct id
323  nod_pt->node_update_fct_id() = 0;
324  }
325  } // end loop over elements
326 
327 
328  // LOOP OVER OTHER SPINES
329  //----------------------
330 
331  // Now loop over the elements horizontally
332  for (unsigned long j = 0; j < n_x0; j++)
333  {
334  // Loop over the nodes in the elements horizontally, ignoring
335  // the first column
336  unsigned n_pmax = n_p;
337  for (unsigned l2 = 1; l2 < n_pmax; l2++)
338  {
339  // Assign the new spine with unit height
340  new_spine_pt = new Spine(1.0);
341  new_spine_pt->spine_height_pt()->pin(0);
342  Spine_pt.push_back(new_spine_pt);
343 
344  // Get the node
345  SpineNode* nod_pt = element_node_pt(j, l2);
346  // Set the pointer to spine
347  nod_pt->spine_pt() = new_spine_pt;
348  // Set the fraction
349  nod_pt->fraction() = 0.0;
350  // Pointer to the mesh that implements the update fct
351  nod_pt->spine_mesh_pt() = this;
352  // Set update fct id
353  nod_pt->node_update_fct_id() = 0;
354 
355  {
356  // Provide spine with additional storage for wall coordinate
357  // and wall geom object:
358 
359  // Get the Lagrangian coordinate in the Lower Wall
360  zeta[0] =
361  zeta_lo + double(j) * dzeta + double(l2) * dzeta / (n_p - 1.0);
362  // Get the geometric object and local coordinate
363  Straight_wall_pt->locate_zeta(zeta, geometric_object_pt, s_wall);
364 
365  // The local coordinate is a geometric parameter
366  // This needs to be set (rather than added) because the
367  // same spine may be visited more than once
368  Vector<double> parameters(1, s_wall[0]);
369  nod_pt->spine_pt()->set_geom_parameter(parameters);
370 
371  // Get position of wall
372  Straight_wall_pt->position(s_wall, r_wall);
373 
374  // Adjust spine height
375  nod_pt->spine_pt()->height() = r_wall[1];
376 
377  // The sub geom object is one (and only) geom object
378  // for spine:
379  Vector<GeomObject*> geom_object_pt(1);
380  geom_object_pt[0] = geometric_object_pt;
381 
382  // Pass geom object(s) to spine
383  nod_pt->spine_pt()->set_geom_object_pt(geom_object_pt);
384  }
385 
386  // Loop vertically along the spine
387  // Loop over the elements
388  for (unsigned long i = 0; i < n_y; i++)
389  {
390  // Loop over the vertical nodes, apart from the first
391  for (unsigned l1 = 1; l1 < n_p; l1++)
392  {
393  // Get the node
394  SpineNode* nod_pt = element_node_pt(i * n_x + j, l1 * n_p + l2);
395  // Set the pointer to the spine
396  nod_pt->spine_pt() = new_spine_pt;
397  // Set the fraction
398  nod_pt->fraction() =
399  (double(i) + double(l1) / double(n_p - 1)) / double(n_y);
400  // Pointer to the mesh that implements the update fct
401  nod_pt->spine_mesh_pt() = this;
402  // Set update fct id
403  nod_pt->node_update_fct_id() = 0;
404  }
405  }
406  }
407  }
408 
409  // Increment number of previous elements
410  n_prev_elements += n_x0 * n_y;
411 
412 
413  // CENTRE REGION
414  // ===========
415 
416  zeta_lo = Lx0;
417  dzeta = Lx1 / n_x1;
418 
419  // SPINES IN LEFT REGION
420  // ---------------------
421 
422  // LOOP OVER OTHER SPINES
423  //----------------------
424 
425  // Now loop over the elements horizontally
426  for (unsigned long j = n_x0; j < n_x0 + n_x1; j++)
427  {
428  // Loop over the nodes in the elements horizontally, ignoring
429  // the first column
430  unsigned n_pmax = n_p;
431  for (unsigned l2 = 1; l2 < n_pmax; l2++)
432  {
433  // Assign the new spine with unit height
434  new_spine_pt = new Spine(1.0);
435  new_spine_pt->spine_height_pt()->pin(0);
436  Spine_pt.push_back(new_spine_pt);
437 
438  // Get the node
439  SpineNode* nod_pt = element_node_pt(j, l2);
440  // Set the pointer to spine
441  nod_pt->spine_pt() = new_spine_pt;
442  // Set the fraction
443  nod_pt->fraction() = 0.0;
444  // Pointer to the mesh that implements the update fct
445  nod_pt->spine_mesh_pt() = this;
446  // Set update fct id
447  nod_pt->node_update_fct_id() = 1;
448 
449  {
450  // Provide spine with additional storage for wall coordinate
451  // and wall geom object:
452 
453  // Get the Lagrangian coordinate in the Lower Wall
454  zeta[0] = zeta_lo + double(j - n_x0) * dzeta +
455  double(l2) * dzeta / (n_p - 1.0);
456  // Get the geometric object and local coordinate
457  Wall_pt->locate_zeta(zeta, geometric_object_pt, s_wall);
458 
459  // The local coordinate is a geometric parameter
460  // This needs to be set (rather than added) because the
461  // same spine may be visited more than once
462  Vector<double> parameters(1, s_wall[0]);
463  nod_pt->spine_pt()->set_geom_parameter(parameters);
464 
465  // Get position of wall
466  Wall_pt->position(s_wall, r_wall);
467 
468  // Adjust spine height
469  nod_pt->spine_pt()->height() = r_wall[1];
470 
471  // The sub geom object is one (and only) geom object
472  // for spine:
473  Vector<GeomObject*> geom_object_pt(1);
474  geom_object_pt[0] = geometric_object_pt;
475 
476  // Pass geom object(s) to spine
477  nod_pt->spine_pt()->set_geom_object_pt(geom_object_pt);
478  }
479 
480  // Loop vertically along the spine
481  // Loop over the elements
482  for (unsigned long i = 0; i < n_y; i++)
483  {
484  // Loop over the vertical nodes, apart from the first
485  for (unsigned l1 = 1; l1 < n_p; l1++)
486  {
487  // Get the node
488  SpineNode* nod_pt = element_node_pt(i * n_x + j, l1 * n_p + l2);
489  // Set the pointer to the spine
490  nod_pt->spine_pt() = new_spine_pt;
491  // Set the fraction
492  nod_pt->fraction() =
493  (double(i) + double(l1) / double(n_p - 1)) / double(n_y);
494  // Pointer to the mesh that implements the update fct
495  nod_pt->spine_mesh_pt() = this;
496  // Set update fct id
497  nod_pt->node_update_fct_id() = 1;
498  }
499  }
500  }
501  }
502 
503  // Increment number of previous elements
504  n_prev_elements += n_x1 * n_y;
505 
506 
507  // RIGHT REGION
508  // ===========
509 
510  // SPINES IN RIGHT REGION
511  // ----------------------
512 
513  // Set up zeta increments
514  zeta_lo = Lx0 + Lx1;
515  dzeta = Lx2 / n_x2;
516 
517  // LOOP OVER OTHER SPINES
518  //----------------------
519 
520  // Now loop over the elements horizontally
521  for (unsigned long j = n_x0 + n_x1; j < n_x0 + n_x1 + n_x2; j++)
522  {
523  // Need to copy last spine in previous element over to first spine
524  // in next elements, for all elements other than the first
525  if (j > 0)
526  {
527  // For first spine, re-assign the geometric objects, since
528  // we treat it as part of the right region.
529  if (j == n_x0 + n_x1)
530  {
531  SpineNode* nod_pt = element_node_pt(j, 0);
532  // Set update fct id
533  nod_pt->node_update_fct_id() = 2;
534  {
535  // Provide spine with additional storage for wall coordinate
536  // and wall geom object:
537 
538  // Get the Lagrangian coordinate in the Lower Wall
539  zeta[0] = zeta_lo + double(j - n_x0 - n_x1) * dzeta;
540  // Get the geometric object and local coordinate
541  Straight_wall_pt->locate_zeta(zeta, geometric_object_pt, s_wall);
542 
543  // The local coordinate is a geometric parameter
544  // This needs to be set (rather than added) because the
545  // same spine may be visited more than once
546  Vector<double> parameters(1, s_wall[0]);
547  nod_pt->spine_pt()->set_geom_parameter(parameters);
548 
549  // Get position of wall
550  Straight_wall_pt->position(s_wall, r_wall);
551 
552  // Adjust spine height
553  nod_pt->spine_pt()->height() = r_wall[1];
554 
555  // The sub geom object is one (and only) geom object
556  // for spine:
557  Vector<GeomObject*> geom_object_pt(1);
558  geom_object_pt[0] = geometric_object_pt;
559 
560  // Pass geom object(s) to spine
561  nod_pt->spine_pt()->set_geom_object_pt(geom_object_pt);
562  }
563  }
564  }
565  // Loop over the nodes in the elements horizontally, ignoring
566  // the first column
567 
568  // Last spine needs special treatment in x-periodic meshes:
569  unsigned n_pmax = n_p;
570  if ((this->Xperiodic) && (j == n_x - 1)) n_pmax = n_p - 1;
571 
572  for (unsigned l2 = 1; l2 < n_pmax; l2++)
573  {
574  // Assign the new spine with unit height
575  new_spine_pt = new Spine(1.0);
576  new_spine_pt->spine_height_pt()->pin(0);
577  Spine_pt.push_back(new_spine_pt);
578 
579  // Get the node
580  SpineNode* nod_pt = element_node_pt(j, l2);
581  // Set the pointer to spine
582  nod_pt->spine_pt() = new_spine_pt;
583  // Set the fraction
584  nod_pt->fraction() = 0.0;
585  // Pointer to the mesh that implements the update fct
586  nod_pt->spine_mesh_pt() = this;
587  // Set update fct id
588  nod_pt->node_update_fct_id() = 2;
589 
590  {
591  // Provide spine with additional storage for wall coordinate
592  // and wall geom object:
593 
594  // Get the Lagrangian coordinate in the Lower Wall
595  zeta[0] = zeta_lo + double(j - n_x0 - n_x1) * dzeta +
596  double(l2) * dzeta / (n_p - 1.0);
597  // Get the geometric object and local coordinate
598  Straight_wall_pt->locate_zeta(zeta, geometric_object_pt, s_wall);
599 
600  // The local coordinate is a geometric parameter
601  // This needs to be set (rather than added) because the
602  // same spine may be visited more than once
603  Vector<double> parameters(1, s_wall[0]);
604  nod_pt->spine_pt()->set_geom_parameter(parameters);
605 
606  // Get position of wall
607  Straight_wall_pt->position(s_wall, r_wall);
608 
609  // Adjust spine height
610  nod_pt->spine_pt()->height() = r_wall[1];
611 
612  // The sub geom object is one (and only) geom object
613  // for spine:
614  Vector<GeomObject*> geom_object_pt(1);
615  geom_object_pt[0] = geometric_object_pt;
616 
617  // Pass geom object(s) to spine
618  nod_pt->spine_pt()->set_geom_object_pt(geom_object_pt);
619  }
620 
621  // Loop vertically along the spine
622  // Loop over the elements
623  for (unsigned long i = 0; i < n_y; i++)
624  {
625  // Loop over the vertical nodes, apart from the first
626  for (unsigned l1 = 1; l1 < n_p; l1++)
627  {
628  // Get the node
629  SpineNode* nod_pt = element_node_pt(i * n_x + j, l1 * n_p + l2);
630  // Set the pointer to the spine
631  nod_pt->spine_pt() = new_spine_pt;
632  // Set the fraction
633  nod_pt->fraction() =
634  (double(i) + double(l1) / double(n_p - 1)) / double(n_y);
635  // Pointer to the mesh that implements the update fct
636  nod_pt->spine_mesh_pt() = this;
637  // Set update fct id
638  nod_pt->node_update_fct_id() = 2;
639  }
640  }
641  }
642  }
643 
644  // Increment number of previous elements
645  n_prev_elements += n_x2 * n_y;
646 
647  // Last spine needs special treatment for periodic meshes
648  // because it's the same as the first one...
649  if (this->Xperiodic)
650  {
651  // Last spine is the same as first one...
652  Spine* final_spine_pt = Spine_pt[0];
653 
654  // Get the node
655  SpineNode* nod_pt = element_node_pt((n_x - 1), (n_p - 1));
656 
657  // Set the pointer for the first node
658  nod_pt->spine_pt() = final_spine_pt;
659  // Set the fraction to be the same as for the nodes on the first row
660  nod_pt->fraction() = element_node_pt(0, 0)->fraction();
661  // Pointer to the mesh that implements the update fct
662  nod_pt->spine_mesh_pt() = element_node_pt(0, 0)->spine_mesh_pt();
663 
664  // Now loop vertically along the spine
665  for (unsigned i = 0; i < n_y; i++)
666  {
667  // Loop over the vertical nodes, apart from the first
668  for (unsigned l1 = 1; l1 < n_p; l1++)
669  {
670  // Get the node
671  SpineNode* nod_pt =
672  element_node_pt(i * n_x + (n_x - 1), l1 * n_p + (n_p - 1));
673 
674  // Set the pointer to the spine
675  nod_pt->spine_pt() = final_spine_pt;
676  // Set the fraction to be the same as in first row
677  nod_pt->fraction() = element_node_pt(i * n_x, l1 * n_p)->fraction();
678  // Pointer to the mesh that implements the update fct
679  nod_pt->spine_mesh_pt() =
680  element_node_pt(i * n_x, l1 * n_p)->spine_mesh_pt();
681  }
682  }
683  }
684 
685 
686  } // end of build_channel_spine_mesh
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Array< double, 1, 3 > e(1./3., 0.5, 2.)
Vector< FiniteElement * > Centre_element_pt
Vector of pointers to element in the centre region.
Definition: channel_spine_mesh.template.h:278
unsigned long ncentre() const
Number of elements in centre region.
Definition: channel_spine_mesh.template.h:107
unsigned long nright() const
Number of elements in right region.
Definition: channel_spine_mesh.template.h:113
unsigned Nright_spine
Number of spines in right region.
Definition: channel_spine_mesh.template.h:312
unsigned Nleft_spine
Number of spines in left region.
Definition: channel_spine_mesh.template.h:306
unsigned long nleft() const
Number of elements in left region.
Definition: channel_spine_mesh.template.h:101
Vector< FiniteElement * > Right_element_pt
Vector of pointers to element in the right region.
Definition: channel_spine_mesh.template.h:281
Vector< FiniteElement * > Left_element_pt
Vector of pointers to element in the left region.
Definition: channel_spine_mesh.template.h:275
unsigned Ncentre_spine
Number of spines in centre region.
Definition: channel_spine_mesh.template.h:309
virtual void position(const Vector< double > &zeta, Vector< double > &r) const =0
Parametrised position on object at current time: r(zeta).
virtual void locate_zeta(const Vector< double > &zeta, GeomObject *&sub_geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
Definition: geom_objects.h:381
FiniteElement * finite_element_pt(const unsigned &e) const
Definition: mesh.h:473
unsigned long nelement() const
Return number of elements in the mesh.
Definition: mesh.h:590
unsigned Nx
Nx: number of elements in x-direction.
Definition: rectangular_quadmesh.template.h:63
unsigned Ny
Ny: number of elements in y-direction.
Definition: rectangular_quadmesh.template.h:65
bool Xperiodic
Definition: rectangular_quadmesh.template.h:81
void build_mesh(TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Generic mesh construction function: contains all the hard work.
Definition: rectangular_quadmesh.template.cc:43
Vector< Spine * > Spine_pt
A Spine mesh contains a Vector of pointers to spines.
Definition: spines.h:616
unsigned long nspine() const
Return the number of spines in the mesh.
Definition: spines.h:635
SpineNode * element_node_pt(const unsigned long &e, const unsigned &n)
Definition: spines.h:669
SpineMesh *& spine_mesh_pt()
Definition: spines.h:391
double & fraction()
Set reference to fraction along spine.
Definition: spines.h:378
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
#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::RectangularQuadMesh< ELEMENT >::build_mesh(), e(), oomph::SpineNode::fraction(), oomph::Spine::height(), i, j, oomph::SpineNode::node_update_fct_id(), GlobalParameters::Nx, GlobalParameters::Ny, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::Data::pin(), oomph::Spine::set_geom_object_pt(), oomph::Spine::set_geom_parameter(), oomph::Spine::spine_height_pt(), oomph::SpineNode::spine_mesh_pt(), oomph::SpineNode::spine_pt(), and Eigen::zeta().

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

◆ centre_element_pt()

template<class ELEMENT >
FiniteElement*& oomph::ChannelSpineMesh< ELEMENT >::centre_element_pt ( const unsigned long &  i)
inline

Access functions for pointers to the \( i \) -th element in the centre region.

89  {
90  return Centre_element_pt[i];
91  }

References oomph::ChannelSpineMesh< ELEMENT >::Centre_element_pt, and i.

◆ centre_spine_pt()

template<class ELEMENT >
Spine*& oomph::ChannelSpineMesh< ELEMENT >::centre_spine_pt ( const unsigned long &  i)
inline

Access function for spines in centre region.

215  {
216  if (i > Ncentre_spine)
217  {
218  throw OomphLibError("Arguemnt out of range",
221  }
222  else
223  {
224  return Spine_pt[Nleft_spine + i];
225  }
226  }

References i, oomph::ChannelSpineMesh< ELEMENT >::Ncentre_spine, oomph::ChannelSpineMesh< ELEMENT >::Nleft_spine, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::SpineMesh::Spine_pt.

◆ element_reorder()

template<class ELEMENT >
void oomph::ChannelSpineMesh< ELEMENT >::element_reorder
virtual

Reorder the elements so we loop over them vertically first (advantageous in "wide" domains if a frontal solver is used).

Reorder the elements, so we loop over them vertically first (advantageous in "wide" domains if a frontal solver is used).

Reimplemented from oomph::RectangularQuadMesh< ELEMENT >.

695  {
696  unsigned n_x = this->Nx;
697  unsigned n_y = this->Ny;
698  // Find out how many elements there are
699  unsigned long Nelement = nelement();
700  // Find out how many fluid elements there are
701  unsigned long Nfluid = n_x * n_y;
702  // Create a dummy array of elements
703  Vector<FiniteElement*> dummy;
704 
705  // Loop over the elements in horizontal order
706  for (unsigned long j = 0; j < n_x; j++)
707  {
708  // Loop over the elements in lower layer vertically
709  for (unsigned long i = 0; i < n_y; i++)
710  {
711  // Push back onto the new stack
712  dummy.push_back(finite_element_pt(n_x * i + j));
713  }
714 
715  // Push back the line element onto the stack
716  dummy.push_back(finite_element_pt(Nfluid + j));
717  }
718 
719  // Now copy the reordered elements into the element_pt
720  for (unsigned long e = 0; e < Nelement; e++)
721  {
722  Element_pt[e] = dummy[e];
723  }
724 
725  } // end of element_reorder
Vector< GeneralisedElement * > Element_pt
Vector of pointers to generalised elements.
Definition: mesh.h:186

References e(), i, j, GlobalParameters::Nx, and GlobalParameters::Ny.

◆ left_element_pt()

template<class ELEMENT >
FiniteElement*& oomph::ChannelSpineMesh< ELEMENT >::left_element_pt ( const unsigned long &  i)
inline

Access functions for pointers to the \( i \) -th element in the left region.

82  {
83  return Left_element_pt[i];
84  }

References i, and oomph::ChannelSpineMesh< ELEMENT >::Left_element_pt.

◆ left_spine_pt()

template<class ELEMENT >
Spine*& oomph::ChannelSpineMesh< ELEMENT >::left_spine_pt ( const unsigned long &  i)
inline

Access function for spines in left region.

201  {
202 #ifdef PARNOID
203  if (i > Nleft_spine)
204  {
205  throw OomphLibError("Arguemnt out of range",
208  }
209 #endif
210  return Spine_pt[i];
211  }

References i, oomph::ChannelSpineMesh< ELEMENT >::Nleft_spine, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::SpineMesh::Spine_pt.

◆ nbulk()

template<class ELEMENT >
unsigned long oomph::ChannelSpineMesh< ELEMENT >::nbulk ( ) const
inline

Number of elements in bulk.

120  {
121  unsigned long Nbulk = Left_element_pt.size() + Centre_element_pt.size() +
122  Right_element_pt.size();
123  return Nbulk;
124  }

References oomph::ChannelSpineMesh< ELEMENT >::Centre_element_pt, oomph::ChannelSpineMesh< ELEMENT >::Left_element_pt, and oomph::ChannelSpineMesh< ELEMENT >::Right_element_pt.

◆ ncentre()

template<class ELEMENT >
unsigned long oomph::ChannelSpineMesh< ELEMENT >::ncentre ( ) const
inline

Number of elements in centre region.

108  {
109  return Centre_element_pt.size();
110  }

References oomph::ChannelSpineMesh< ELEMENT >::Centre_element_pt.

◆ ncentre_spine()

template<class ELEMENT >
unsigned oomph::ChannelSpineMesh< ELEMENT >::ncentre_spine ( )
inline

Access function for the number of spines in the centre region.

251  {
252  return Ncentre_spine;
253  }

References oomph::ChannelSpineMesh< ELEMENT >::Ncentre_spine.

◆ nleft()

template<class ELEMENT >
unsigned long oomph::ChannelSpineMesh< ELEMENT >::nleft ( ) const
inline

Number of elements in left region.

102  {
103  return Left_element_pt.size();
104  }

References oomph::ChannelSpineMesh< ELEMENT >::Left_element_pt.

◆ nleft_spine()

template<class ELEMENT >
unsigned oomph::ChannelSpineMesh< ELEMENT >::nleft_spine ( )
inline

Access function for the number of spines in the left region.

245  {
246  return Nleft_spine;
247  }

References oomph::ChannelSpineMesh< ELEMENT >::Nleft_spine.

◆ nright()

template<class ELEMENT >
unsigned long oomph::ChannelSpineMesh< ELEMENT >::nright ( ) const
inline

Number of elements in right region.

114  {
115  return Right_element_pt.size();
116  }

References oomph::ChannelSpineMesh< ELEMENT >::Right_element_pt.

◆ nright_spine()

template<class ELEMENT >
unsigned oomph::ChannelSpineMesh< ELEMENT >::nright_spine ( )
inline

Access function for the number of spines in the right region.

257  {
258  return Nright_spine;
259  }

References oomph::ChannelSpineMesh< ELEMENT >::Nright_spine.

◆ right_element_pt()

template<class ELEMENT >
FiniteElement*& oomph::ChannelSpineMesh< ELEMENT >::right_element_pt ( const unsigned long &  i)
inline

Access functions for pointers to the \( i \) -th element in the right region.

96  {
97  return Right_element_pt[i];
98  }

References i, and oomph::ChannelSpineMesh< ELEMENT >::Right_element_pt.

◆ right_spine_pt()

template<class ELEMENT >
Spine*& oomph::ChannelSpineMesh< ELEMENT >::right_spine_pt ( const unsigned long &  i)
inline

Access function for spines in right region.

230  {
231  if (i > Nright_spine)
232  {
233  throw OomphLibError("Arguemnt out of range",
236  }
237  else
238  {
239  return Spine_pt[Nleft_spine + Ncentre_spine - 1 + i];
240  }
241  }

References i, oomph::ChannelSpineMesh< ELEMENT >::Ncentre_spine, oomph::ChannelSpineMesh< ELEMENT >::Nleft_spine, oomph::ChannelSpineMesh< ELEMENT >::Nright_spine, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::SpineMesh::Spine_pt.

◆ spine_node_update()

template<class ELEMENT >
virtual void oomph::ChannelSpineMesh< ELEMENT >::spine_node_update ( SpineNode spine_node_pt)
inlinevirtual

General node update function implements pure virtual function defined in SpineMesh base class and performs specific node update actions: along vertical spines

Implements oomph::SpineMesh.

134  {
135  // Get spine node's fraction along the spine
136  double W = spine_node_pt->fraction();
137 
138  // Get local coordinates
139  Vector<double> s_wall(1);
140  s_wall[0] = spine_node_pt->spine_pt()->geom_parameter(0);
141 
142  // Get position vector to wall
143  Vector<double> position(2);
144  spine_node_pt->spine_pt()->geom_object_pt(0)->position(s_wall, position);
145 
146  // Set the value of y
147  spine_node_pt->x(1) = this->Ymin + W * position[1];
148  }
double Ymin
Minimum value of y coordinate.
Definition: rectangular_quadmesh.template.h:75
@ W
Definition: quadtree.h:63

References oomph::SpineNode::fraction(), oomph::Spine::geom_object_pt(), oomph::Spine::geom_parameter(), oomph::GeomObject::position(), oomph::SpineNode::spine_pt(), oomph::QuadTreeNames::W, oomph::Node::x(), and oomph::RectangularQuadMesh< ELEMENT >::Ymin.

◆ straight_wall_pt()

template<class ELEMENT >
GeomObject* oomph::ChannelSpineMesh< ELEMENT >::straight_wall_pt ( )
inline

Access function to the GeomObject for the straight upper wall.

269  {
270  return Straight_wall_pt;
271  }

References oomph::ChannelSpineMesh< ELEMENT >::Straight_wall_pt.

◆ wall_pt()

template<class ELEMENT >
GeomObject* oomph::ChannelSpineMesh< ELEMENT >::wall_pt ( )
inline

Access function to the GeomObject for upper wall.

263  {
264  return Wall_pt;
265  }

References oomph::ChannelSpineMesh< ELEMENT >::Wall_pt.

Referenced by doc_sparse_node_update().

◆ x_spacing_function()

template<class ELEMENT >
virtual double oomph::ChannelSpineMesh< ELEMENT >::x_spacing_function ( unsigned  xelement,
unsigned  xnode,
unsigned  yelement,
unsigned  ynode 
)
inlinevirtual

Return the value of the x-coordinate at the node given by the local node number (xnode, ynode) in the element (xelement,yelement). The description is in a "psudeo" two-dimensional coordinate system, so the range of xelement is [0,Nx-1], yelement is [0,Ny-1], and that of xnode and ynode is [0,Np-1]. The default is to return nodes that are equally spaced in the x coodinate.

Reimplemented from oomph::RectangularQuadMesh< ELEMENT >.

160  {
161  // Calculate the values of equal increments in nodal values in left region
162  double xstep1 = Lx0 / ((this->Np - 1) * Nx0);
163  // Calculate the values of equal increments in nodal values in centre
164  // region
165  double xstep2 = Lx1 / ((this->Np - 1) * Nx1);
166  // Calculate the values of equal increments in nodal values in right
167  // region
168  double xstep3 = Lx2 / ((this->Np - 1) * Nx2);
169 
170  // left region
171  if (xelement < Nx0)
172  {
173  // Return the appropriate value
174  return (this->Xmin + xstep1 * ((this->Np - 1) * xelement + xnode));
175  }
176  // centre region
177  else if (xelement < Nx0 + Nx1)
178  {
179  // Return the appropriate value
180  return (Lx0 + xstep2 * ((this->Np - 1) * (xelement - Nx0) + xnode));
181  }
182  // right region
183  else if (xelement < Nx0 + Nx1 + Nx2)
184  {
185  // Return the appropriate value
186  return (Lx0 + Lx1 +
187  xstep3 * ((this->Np - 1) * (xelement - Nx0 - Nx1) + xnode));
188  }
189  else
190  {
191  throw OomphLibError("Should not have got here",
194  }
195  // Dummy return to keep compiler from barking
196  return 0.0;
197  }
unsigned Np
Np: number of (linear) points in the element.
Definition: rectangular_quadmesh.template.h:67
double Xmin
Minimum value of x coordinate.
Definition: rectangular_quadmesh.template.h:70

References oomph::ChannelSpineMesh< ELEMENT >::Lx0, oomph::ChannelSpineMesh< ELEMENT >::Lx1, oomph::ChannelSpineMesh< ELEMENT >::Lx2, oomph::RectangularQuadMesh< ELEMENT >::Np, oomph::ChannelSpineMesh< ELEMENT >::Nx0, oomph::ChannelSpineMesh< ELEMENT >::Nx1, oomph::ChannelSpineMesh< ELEMENT >::Nx2, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::RectangularQuadMesh< ELEMENT >::Xmin.

Member Data Documentation

◆ Centre_element_pt

template<class ELEMENT >
Vector<FiniteElement*> oomph::ChannelSpineMesh< ELEMENT >::Centre_element_pt
protected

◆ Left_element_pt

template<class ELEMENT >
Vector<FiniteElement*> oomph::ChannelSpineMesh< ELEMENT >::Left_element_pt
protected

◆ Lx0

template<class ELEMENT >
double oomph::ChannelSpineMesh< ELEMENT >::Lx0
protected

Length of left region.

Referenced by oomph::ChannelSpineMesh< ELEMENT >::x_spacing_function().

◆ Lx1

template<class ELEMENT >
double oomph::ChannelSpineMesh< ELEMENT >::Lx1
protected

Length of centre region.

Referenced by oomph::ChannelSpineMesh< ELEMENT >::x_spacing_function().

◆ Lx2

template<class ELEMENT >
double oomph::ChannelSpineMesh< ELEMENT >::Lx2
protected

Length of right region.

Referenced by oomph::ChannelSpineMesh< ELEMENT >::x_spacing_function().

◆ Ncentre_spine

◆ Nleft_spine

◆ Nright_spine

template<class ELEMENT >
unsigned oomph::ChannelSpineMesh< ELEMENT >::Nright_spine
protected

◆ Nx0

template<class ELEMENT >
unsigned oomph::ChannelSpineMesh< ELEMENT >::Nx0
protected

Number of elements in the left region.

Referenced by oomph::ChannelSpineMesh< ELEMENT >::x_spacing_function().

◆ Nx1

template<class ELEMENT >
unsigned oomph::ChannelSpineMesh< ELEMENT >::Nx1
protected

Number of elements in the centre region.

Referenced by oomph::ChannelSpineMesh< ELEMENT >::x_spacing_function().

◆ Nx2

template<class ELEMENT >
unsigned oomph::ChannelSpineMesh< ELEMENT >::Nx2
protected

Number of elements in the right region.

Referenced by oomph::ChannelSpineMesh< ELEMENT >::x_spacing_function().

◆ Right_element_pt

template<class ELEMENT >
Vector<FiniteElement*> oomph::ChannelSpineMesh< ELEMENT >::Right_element_pt
protected

◆ Straight_wall_pt

template<class ELEMENT >
GeomObject* oomph::ChannelSpineMesh< ELEMENT >::Straight_wall_pt
protected

◆ Wall_pt

template<class ELEMENT >
GeomObject* oomph::ChannelSpineMesh< ELEMENT >::Wall_pt
protected

GeomObject for upper wall.

Referenced by oomph::ChannelSpineMesh< ELEMENT >::wall_pt().


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