oomph::PerturbedSpineMesh Class Referenceabstract

#include <perturbed_spines.h>

+ Inheritance diagram for oomph::PerturbedSpineMesh:

Public Member Functions

virtual ~PerturbedSpineMesh ()
 Destructor to clean up the memory allocated to the perturbed spines. More...
 
PerturbedSpine *const & perturbed_spine_pt (const unsigned long &i) const
 Return the i-th perturbed spine in the mesh. More...
 
unsigned long nspine () const
 Return the number of perturbed spines in the mesh. More...
 
void add_perturbed_spine_pt (PerturbedSpine *const &spine_pt)
 Add a perturbed spine to the mesh. More...
 
PerturbedSpineNodenode_pt (const unsigned long &n)
 Return a pointer to the n-th global PerturbedSpineNode. More...
 
PerturbedSpineNodeelement_node_pt (const unsigned long &e, const unsigned &n)
 
unsigned long assign_global_eqn_numbers (Vector< double * > &Dof_pt)
 Assign equation numbers for perturbed spines. More...
 
void node_update (const bool &update_all_solid_nodes=false)
 
virtual void perturbed_spine_node_update (PerturbedSpineNode *spine_node_pt)=0
 
void dump (std::ofstream &dump_file, const bool &use_old_ordering=true) 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...
 
- 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 setup_boundary_element_info ()
 
virtual void setup_boundary_element_info (std::ostream &outfile)
 
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
 
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)
 
virtual void set_mesh_level_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
void set_consistent_pinned_values_for_continuation (ContinuationStorageScheme *const &continuation_stepper_pt)
 Set consistent values for pinned data in continuation. More...
 
bool does_pointer_correspond_to_mesh_data (double *const &parameter_pt)
 Does the double pointer correspond to any mesh data. More...
 
void set_nodal_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 Set the timestepper associated with the nodal data in the mesh. More...
 
void set_elemental_internal_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
virtual void compute_norm (double &norm)
 
virtual void compute_norm (Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 Returns the norm of the error and that of the exact solution. More...
 
virtual void compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
bool is_mesh_distributed () const
 Boolean to indicate if Mesh has been distributed. More...
 
OomphCommunicatorcommunicator_pt () const
 
void delete_all_external_storage ()
 Wipe the storage for all externally-based elements. More...
 

Protected Attributes

Vector< PerturbedSpine * > PerturbedSpine_pt
 A PerturbedSpine mesh contains a Vector of pointers to perturbed spines. More...
 
- Protected Attributes inherited from oomph::Mesh
Vector< Vector< Node * > > Boundary_node_pt
 
bool Lookup_for_elements_next_boundary_is_setup
 
Vector< Vector< FiniteElement * > > Boundary_element_pt
 
Vector< Vector< int > > Face_index_at_boundary
 
Vector< Node * > Node_pt
 Vector of pointers to nodes. More...
 
Vector< GeneralisedElement * > Element_pt
 Vector of pointers to generalised elements. More...
 
std::vector< boolBoundary_coordinate_exists
 

Additional Inherited Members

- Public Types inherited from oomph::Mesh
typedef void(FiniteElement::* SteadyExactSolutionFctPt) (const Vector< double > &x, Vector< double > &soln)
 
typedef void(FiniteElement::* UnsteadyExactSolutionFctPt) (const double &time, const Vector< double > &x, Vector< double > &soln)
 
- Static Public Attributes inherited from oomph::Mesh
static Steady< 0 > Default_TimeStepper
 The Steady Timestepper. More...
 
static bool Suppress_warning_about_empty_mesh_level_time_stepper_function
 Static boolean flag to control warning about mesh level timesteppers. More...
 
- Protected Member Functions inherited from oomph::Mesh
unsigned long assign_global_eqn_numbers (Vector< double * > &Dof_pt)
 Assign (global) equation numbers to the nodes. More...
 
void describe_dofs (std::ostream &out, const std::string &current_string) const
 
void describe_local_dofs (std::ostream &out, const std::string &current_string) const
 
void assign_local_eqn_numbers (const bool &store_local_dof_pt)
 Assign local equation numbers in all elements. More...
 
void convert_to_boundary_node (Node *&node_pt, const Vector< FiniteElement * > &finite_element_pt)
 
void convert_to_boundary_node (Node *&node_pt)
 

Detailed Description

/////////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////////// General PerturbedSpineMesh class.

Derived from Mesh with virtual so that perturbed spine meshes can be derived from general meshes, without multiple copies of Mesh objects.

Constructor & Destructor Documentation

◆ ~PerturbedSpineMesh()

virtual oomph::PerturbedSpineMesh::~PerturbedSpineMesh ( )
inlinevirtual

Destructor to clean up the memory allocated to the perturbed spines.

563  {
564  // Set the range of PerturbedSpine_pt
565  unsigned long PerturbedSpine_pt_range = PerturbedSpine_pt.size();
566 
567  // Loop over the entries in reverse and free memory
568  for(unsigned long i=PerturbedSpine_pt_range;i>0;i--)
569  {
570  delete PerturbedSpine_pt[i-1]; PerturbedSpine_pt[i-1] = 0;
571  }
572  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Vector< PerturbedSpine * > PerturbedSpine_pt
A PerturbedSpine mesh contains a Vector of pointers to perturbed spines.
Definition: perturbed_spines.h:557

References i, and PerturbedSpine_pt.

Member Function Documentation

◆ add_perturbed_spine_pt()

void oomph::PerturbedSpineMesh::add_perturbed_spine_pt ( PerturbedSpine *const &  spine_pt)
inline

Add a perturbed spine to the mesh.

585  {
586  PerturbedSpine_pt.push_back(spine_pt);
587  }

References PerturbedSpine_pt.

◆ assign_global_eqn_numbers()

unsigned long oomph::PerturbedSpineMesh::assign_global_eqn_numbers ( Vector< double * > &  Dof_pt)

Assign equation numbers for perturbed spines.

Assign (global) equation numbers to spines, nodes and elements.

193 {
194  // First call the general mesh assign_eqn_numbers
195  unsigned long equation_number = Mesh::assign_global_eqn_numbers(Dof_pt);
196 
197  // Loop over spines and set global equation numbers for the spine heights
198  // (they are the only Data items whose global eqn numbers are assigned
199  // here)
200  unsigned long PerturbedSpine_pt_range = PerturbedSpine_pt.size();
201  for(unsigned long i=0;i<PerturbedSpine_pt_range;i++)
202  {
203  PerturbedSpine_pt[i]->height_pt()->assign_eqn_numbers(equation_number,
204  Dof_pt);
205  }
206 
207  // Return total number of equations
208  return(equation_number);
209 }
unsigned long assign_global_eqn_numbers(Vector< double * > &Dof_pt)
Assign (global) equation numbers to the nodes.
Definition: mesh.cc:677

References oomph::Mesh::assign_global_eqn_numbers(), i, and PerturbedSpine_pt.

◆ dump()

void oomph::PerturbedSpineMesh::dump ( std::ofstream &  dump_file,
const bool use_old_ordering = true 
) const
virtual

Overload the dump function so that the spine data is dumped.

Overload the dump function so that the spine data is also dumped.

Reimplemented from oomph::Mesh.

218 {
219  // Call the standard mesh dump function
220  Mesh::dump(dump_file,use_old_ordering);
221 
222  // Now loop over the spine data and dump the spine height data
223  // The ASSUMPTION is that the geometric data is stored elsewhere and will
224  // be dumped elsewhere
225 
226  // Find the number of perturbed spines
227  const unsigned long n_spine = nspine();
228 
229  // Doc number of spines
230  dump_file << n_spine << " # number of spines " << std::endl;
231 
232  // Loop over the perturbed spines
233  for(unsigned long s=0;s<n_spine;s++)
234  {
235  perturbed_spine_pt(s)->height_pt()->dump(dump_file);
236  }
237 }
void dump(std::ostream &dump_file) const
Dump the data object to a file.
Definition: nodes.cc:645
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.
Definition: mesh.cc:1088
unsigned long nspine() const
Return the number of perturbed spines in the mesh.
Definition: perturbed_spines.h:581
PerturbedSpine *const & perturbed_spine_pt(const unsigned long &i) const
Return the i-th perturbed spine in the mesh.
Definition: perturbed_spines.h:575
Data *& height_pt()
Access function to Data object that stores the spine "heights".
Definition: perturbed_spines.h:114
RealScalar s
Definition: level1_cplx_impl.h:130

References oomph::Mesh::dump(), oomph::Data::dump(), oomph::PerturbedSpine::height_pt(), nspine(), perturbed_spine_pt(), and s.

◆ element_node_pt()

PerturbedSpineNode* oomph::PerturbedSpineMesh::element_node_pt ( const unsigned long &  e,
const unsigned n 
)
inline

Return the n-th local PerturbedSpineNode in element e. This is required to cast the nodes in a spine mesh to be PerturbedSpineNodes and therefore allow access to the extra PerturbedSpineNode data.

615  {
616 #ifdef PARANOID
617  // Try to cast to FiniteElement
618  FiniteElement* el_pt=dynamic_cast<FiniteElement*>(Element_pt[e]);
619  if (el_pt==0)
620  {
621  throw OomphLibError(
622  "Can't execute element_node_pt(...) for non FiniteElements",
623  "PerturbedSpineMesh::element_node_pt()",
625  }
626  if(!dynamic_cast<PerturbedSpineNode*>(el_pt->node_pt(n)))
627  {
628  std::ostringstream error_message;
629  error_message << "Node " << n << "is a "
630  << typeid(el_pt->node_pt(n)).name()
631  << ", not a PerturbedSpineNode" << std::endl;
632 
633  throw OomphLibError(error_message.str(),
634  "PerturbedSpineMesh::node_pt()",
636  }
637 #endif
638  // Return a cast to the node pointer
639  return(dynamic_cast<PerturbedSpineNode*>(
640  dynamic_cast<FiniteElement*>(Element_pt[e])->node_pt(n)));
641  }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Array< double, 1, 3 > e(1./3., 0.5, 2.)
Vector< GeneralisedElement * > Element_pt
Vector of pointers to generalised elements.
Definition: mesh.h:186
PerturbedSpineNode * node_pt(const unsigned long &n)
Return a pointer to the n-th global PerturbedSpineNode.
Definition: perturbed_spines.h:590
string name
Definition: plotDoE.py:33
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61

References e(), oomph::Mesh::Element_pt, n, plotDoE::name, oomph::FiniteElement::node_pt(), node_pt(), and OOMPH_EXCEPTION_LOCATION.

◆ node_pt()

PerturbedSpineNode* oomph::PerturbedSpineMesh::node_pt ( const unsigned long &  n)
inline

Return a pointer to the n-th global PerturbedSpineNode.

591  {
592 #ifdef PARANOID
593  if(!dynamic_cast<PerturbedSpineNode*>(Node_pt[n]))
594  {
595  std::ostringstream error_message;
596  error_message << "Node " << n << "is a "
597  << typeid(Node_pt[n]).name()
598  << ", not a PerturbedSpineNode" << std::endl;
599 
600  throw OomphLibError(error_message.str(),
601  "PerturbedSpineMesh::node_pt()",
603  }
604 #endif
605 
606  // Return a cast to the pointer to the node
607  return (dynamic_cast<PerturbedSpineNode*>(Node_pt[n]));
608  }
Vector< Node * > Node_pt
Vector of pointers to nodes.
Definition: mesh.h:183

References n, plotDoE::name, oomph::Mesh::Node_pt, and OOMPH_EXCEPTION_LOCATION.

Referenced by element_node_pt(), and PerturbedStateProblem< BASE_ELEMENT, PERTURBED_ELEMENT >::PerturbedStateProblem().

◆ node_update()

void oomph::PerturbedSpineMesh::node_update ( const bool update_all_solid_nodes = false)
virtual

Update function to update all nodes of mesh [Doesn't make sense to use this mesh with SolidElements anyway, so we buffer the case if update_all_solid_nodes is set to true.]

////////////////////////////////////////////////////////////////////// ////////////////////////////////////////////////////////////////////// ////////////////////////////////////////////////////////////////////// ////////////////////////////////////////////////////////////////////// Update function to update all nodes of mesh. [Doesn't make sense to use this mesh with SolidElements anyway, so we buffer the case if update_all_solid_nodes (which defaults to false) is set to true.]

Reimplemented from oomph::Mesh.

148 {
149 #ifdef PARANOID
150  if (update_all_solid_nodes)
151  {
152  std::string error_message =
153  "Doesn't make sense to use an PerturbedSpineMesh with\n";
154  error_message +=
155  "SolidElements so specifying update_all_solid_nodes=true\n";
156  error_message += "doesn't make sense either\n";
157 
158  throw OomphLibError(error_message,
159  "PerturbedSpineMesh::node_update()",
161  }
162 #endif
163 
164  // Loop over all the nodes
165  unsigned long Node_pt_range = Node_pt.size();
166  for(unsigned long l=0;l<Node_pt_range;l++)
167  {
168 #ifdef PARANOID
169  if(!dynamic_cast<PerturbedSpineNode*>(Node_pt[l]))
170  {
171  std::ostringstream error_stream;
172  error_stream
173  << "Error: Node " << l << "is a " << typeid(Node_pt[l]).name()
174  << ", not a PerturbedSpineNode" << std::endl;
175  throw OomphLibError(error_stream.str(),
176  "PerturbedSpineMesh::node_update()",
178  }
179 #endif
180 
181  // Need to cast to spine node to get to update function
182  dynamic_cast<PerturbedSpineNode*>(Node_pt[l])->node_update();
183  }
184 }
void node_update(const bool &update_all_solid_nodes=false)
Definition: perturbed_spines.cc:147
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286

References plotDoE::name, oomph::Mesh::Node_pt, OOMPH_EXCEPTION_LOCATION, and oomph::Global_string_for_annotation::string().

◆ nspine()

unsigned long oomph::PerturbedSpineMesh::nspine ( ) const
inline

Return the number of perturbed spines in the mesh.

581 { return PerturbedSpine_pt.size(); }

References PerturbedSpine_pt.

Referenced by dump(), PerturbedStateProblem< BASE_ELEMENT, PERTURBED_ELEMENT >::PerturbedStateProblem(), and read().

◆ perturbed_spine_node_update()

virtual void oomph::PerturbedSpineMesh::perturbed_spine_node_update ( PerturbedSpineNode spine_node_pt)
pure virtual

Update function for given spine node – this must be implemented by all specific PerturbedSpineMeshes.

Implemented in oomph::TwoLayerPerturbedSpineMesh< ELEMENT >, and oomph::TwoLayerPerturbedSpineMesh< PERTURBED_ELEMENT >.

Referenced by oomph::PerturbedSpineNode::node_update().

◆ perturbed_spine_pt()

PerturbedSpine* const& oomph::PerturbedSpineMesh::perturbed_spine_pt ( const unsigned long &  i) const
inline

Return the i-th perturbed spine in the mesh.

576  {
577  return PerturbedSpine_pt[i];
578  }

References i, and PerturbedSpine_pt.

Referenced by dump(), PerturbedStateProblem< BASE_ELEMENT, PERTURBED_ELEMENT >::PerturbedStateProblem(), and read().

◆ read()

void oomph::PerturbedSpineMesh::read ( std::ifstream &  restart_file)
virtual

Overload the read function so that the spine data is also read.

Overload the read function so that the spine data is read from the restart file

Reimplemented from oomph::Mesh.

245 {
246  // Call the standard mesh read function
247  Mesh::read(restart_file);
248 
249  // Now loop over the spine data and dump the spine height data
250  // The ASSUMPTION is that the geometric data is stored elsewhere and will
251  // be dumped elsewhere
252 
253  // Get the number of perturbed spines
254  const unsigned long n_spine = nspine();
255 
256  std::string input_string;
257  // Read line up to termination sign
258  getline(restart_file,input_string,'#');
259  // Ignore the restr of the line
260  restart_file.ignore(80,'\n');
261 
262  // Check the number of spines
263  const unsigned long check_n_spine = atoi(input_string.c_str());
264 
265  if(check_n_spine != n_spine)
266  {
267  std::ostringstream error_stream;
268  error_stream
269  << "Number of spines in the restart file, " << check_n_spine
270  << std::endl << "does not equal the number of spines in the mesh "
271  << n_spine << std::endl;
272 
273  throw OomphLibError(error_stream.str(),
274  "PerturbedSpineMesh::read()",
276  }
277 
278  // Loop over the spines and read the data
279  for(unsigned long s=0;s<n_spine;s++)
280  {
281  perturbed_spine_pt(s)->height_pt()->read(restart_file);
282  }
283 }
void read(std::ifstream &restart_file)
Read data object from a file.
Definition: nodes.cc:672
virtual void read(std::ifstream &restart_file)
Read solution from restart file.
Definition: mesh.cc:1130

References oomph::PerturbedSpine::height_pt(), nspine(), OOMPH_EXCEPTION_LOCATION, perturbed_spine_pt(), oomph::Mesh::read(), oomph::Data::read(), s, and oomph::Global_string_for_annotation::string().

Member Data Documentation

◆ PerturbedSpine_pt

Vector<PerturbedSpine*> oomph::PerturbedSpineMesh::PerturbedSpine_pt
protected

A PerturbedSpine mesh contains a Vector of pointers to perturbed spines.

Referenced by add_perturbed_spine_pt(), assign_global_eqn_numbers(), nspine(), perturbed_spine_pt(), and ~PerturbedSpineMesh().


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