oomph::AnnularSpineMesh< ELEMENT > Class Template Reference
+ Inheritance diagram for oomph::AnnularSpineMesh< ELEMENT >:

Public Member Functions

 AnnularSpineMesh (const unsigned &n_r, const unsigned &n_y, const unsigned &n_theta, const double &r_min, const double &r_max, const double &l_y, const double &theta_min, const double &theta_max, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
virtual void spine_node_update (SpineNode *spine_node_pt)
 
- Public Member Functions inherited from oomph::SingleLayerCubicSpineMesh< ELEMENT >
 SingleLayerCubicSpineMesh (const unsigned &nx, const unsigned &ny, const unsigned &nz, const double &lx, const double &ly, const double &h, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
- Public Member Functions inherited from oomph::SimpleCubicMesh< ELEMENT >
 SimpleCubicMesh (const unsigned &nx, const unsigned &ny, const unsigned &nz, const double &lx, const double &ly, const double &lz, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
 SimpleCubicMesh (const unsigned &nx, const unsigned &ny, const unsigned &nz, const double &xmin, const double &xmax, const double &ymin, const double &ymax, const double &zmin, const double &zmax, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
const unsignednx () const
 Access function for number of elements in x directions. More...
 
const unsignedny () const
 Access function for number of elements in y directions. More...
 
const unsignednz () const
 Access function for number of elements in y directions. More...
 
- Public Member Functions inherited from oomph::BrickMeshBase
 BrickMeshBase ()
 Constructor (empty) More...
 
 BrickMeshBase (const BrickMeshBase &)=delete
 Broken copy constructor. More...
 
virtual ~BrickMeshBase ()
 Broken assignment operator. 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...
 

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::SingleLayerCubicSpineMesh< ELEMENT >
virtual void build_single_layer_mesh (TimeStepper *time_stepper_pt)
 
- Protected Member Functions inherited from oomph::SimpleCubicMesh< ELEMENT >
void build_mesh (TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 Generic mesh construction function: contains all the hard work. More...
 
- Protected Member Functions inherited from oomph::Mesh
unsigned long assign_global_eqn_numbers (Vector< double * > &Dof_pt)
 Assign (global) equation numbers to the nodes. More...
 
void describe_dofs (std::ostream &out, const std::string &current_string) const
 
void describe_local_dofs (std::ostream &out, const std::string &current_string) const
 
void assign_local_eqn_numbers (const bool &store_local_dof_pt)
 Assign local equation numbers in all elements. More...
 
void convert_to_boundary_node (Node *&node_pt, const Vector< FiniteElement * > &finite_element_pt)
 
void convert_to_boundary_node (Node *&node_pt)
 
- Protected Attributes inherited from oomph::SimpleCubicMesh< ELEMENT >
unsigned Nx
 Number of elements in x direction. More...
 
unsigned Ny
 Number of elements in y direction. More...
 
unsigned Nz
 Number of elements in y direction. 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
 Minimum value of y coordinate. More...
 
double Zmin
 Minimum value of z coordinate. More...
 
double Zmax
 Maximum value of z coordinate. 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
 
- Protected Attributes inherited from oomph::SpineMesh
Vector< Spine * > Spine_pt
 A Spine mesh contains a Vector of pointers to spines. More...
 

Detailed Description

template<class ELEMENT>
class oomph::AnnularSpineMesh< ELEMENT >

Deform the existing cubic spine mesh into a annular section with spines directed radially inwards from the wall

Constructor & Destructor Documentation

◆ AnnularSpineMesh()

template<class ELEMENT >
oomph::AnnularSpineMesh< ELEMENT >::AnnularSpineMesh ( const unsigned n_r,
const unsigned n_y,
const unsigned n_theta,
const double r_min,
const double r_max,
const double l_y,
const double theta_min,
const double theta_max,
TimeStepper time_stepper_pt = &Mesh::Default_TimeStepper 
)
inline
129  :
130  //This will make a cubic mesh with n_r in the x-direction
131  //n_y in the y-direction and n_theta in the z-direction
132  //The coordinates will run from 0 to r_max, 0 to l_y and 0 to theta_max
133  SingleLayerCubicSpineMesh<ELEMENT>(n_theta,n_y,n_r,r_max,l_y,theta_max,
134  time_stepper_pt)
135  {
136 
137  //Find out how many nodes there are
138  unsigned n_node = this->nnode();
139 
140  //Loop over all nodes
141  for (unsigned n=0;n<n_node;n++)
142  {
143  //pointer to node
144  Node* nod_pt=this->node_pt(n);
145  SpineNode* spine_node_pt=dynamic_cast<SpineNode*>(nod_pt);
146  //Get x/y/z coordinates
147  double x_old=nod_pt->x(0);
148  double y_old=nod_pt->x(1);
149  double z_old=nod_pt->x(2);
150 
151  //Mapping
152 
153  double r = r_min + (r_max-r_min)*z_old;
154  double theta = (theta_min + (theta_max-theta_min)*x_old);
155  double y = y_old;
156 
157 
158  if(spine_node_pt->spine_pt()->ngeom_parameter()==0)
159  {
160  spine_node_pt->h() =
162  spine_node_pt->spine_pt()->add_geom_parameter(theta);
163  }
164 
165  //cout << spine_node_pt->spine_pt()->ngeom_parameter() << std::endl;
166  //Set new nodal coordinates
167  nod_pt->x(0)=r*cos(theta);
168  nod_pt->x(2)=r*sin(theta);
169  nod_pt->x(1)=y;
170  //cout << "one" << theta << std::endl;
171  }
172 
173  }
AnnoyingScalar cos(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:136
AnnoyingScalar sin(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:137
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
unsigned long nnode() const
Return number of nodes in the mesh.
Definition: mesh.h:596
Definition: nodes.h:906
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition: nodes.h:1060
Definition: single_layer_cubic_spine_mesh.template.h:48
SpineNode * node_pt(const unsigned long &n)
Return a pointer to the n-th global SpineNode.
Definition: spines.h:648
Definition: spines.h:328
double & h()
Access function to spine height.
Definition: spines.h:397
Spine *& spine_pt()
Access function to spine.
Definition: spines.h:372
void add_geom_parameter(const double &geom_parameter)
Definition: spines.h:280
unsigned ngeom_parameter()
Definition: spines.h:265
Scalar * y
Definition: level1_cplx_impl.h:128
double r_min
Definition: two_d_biharmonic.cc:237
double r_max
Definition: two_d_biharmonic.cc:238
double theta
Definition: two_d_biharmonic.cc:236
double Film_Thickness
Definition: 3d_rayleigh_instability_surfactant.cc:54
r
Definition: UniformPSDSelfTest.py:20

References oomph::Spine::add_geom_parameter(), cos(), Global_Physical_Variables::Film_Thickness, oomph::SpineNode::h(), n, oomph::Spine::ngeom_parameter(), UniformPSDSelfTest::r, BiharmonicTestFunctions2::r_max, BiharmonicTestFunctions2::r_min, sin(), oomph::SpineNode::spine_pt(), BiharmonicTestFunctions2::theta, oomph::Node::x(), and y.

Member Function Documentation

◆ spine_node_update()

template<class ELEMENT >
virtual void oomph::AnnularSpineMesh< 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

Reimplemented from oomph::SingleLayerCubicSpineMesh< ELEMENT >.

177  {
178  //Get fraction along the spine
179  double W = spine_node_pt->fraction();
180  //Get theta along the spine
181  double theta = spine_node_pt->spine_pt()->geom_parameter(0);
182  //Get spine height
183  double H = spine_node_pt->h();
184  //Set the value of y
185  spine_node_pt->x(0) = (1.0-W*H)*cos(theta);
186  spine_node_pt->x(2) = (1.0-W*H)*sin(theta);
187  }
MatrixXf H
Definition: HessenbergDecomposition_matrixH.cpp:4
double & fraction()
Set reference to fraction along spine.
Definition: spines.h:378
double & geom_parameter(const unsigned &i)
Definition: spines.h:287
@ W
Definition: quadtree.h:63

References cos(), oomph::SpineNode::fraction(), oomph::Spine::geom_parameter(), H, oomph::SpineNode::h(), sin(), oomph::SpineNode::spine_pt(), BiharmonicTestFunctions2::theta, oomph::QuadTreeNames::W, and oomph::Node::x().


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