oomph::GeompackQuadScaffoldMesh Class Reference

#include <geompack_scaffold_mesh.h>

+ Inheritance diagram for oomph::GeompackQuadScaffoldMesh:

Public Member Functions

 GeompackQuadScaffoldMesh ()
 Empty constructor. More...
 
 GeompackQuadScaffoldMesh (const std::string &mesh_file_name, const std::string &curve_file_name)
 Constructor: Pass the filename of the mesh files. More...
 
 GeompackQuadScaffoldMesh (const GeompackQuadScaffoldMesh &)=delete
 Broken copy constructor. More...
 
void operator= (const GeompackQuadScaffoldMesh &)=delete
 Broken assignment operator. More...
 
 ~GeompackQuadScaffoldMesh ()
 Empty destructor. 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 node_update (const bool &update_all_solid_nodes=false)
 
virtual void reorder_nodes (const bool &use_old_ordering=true)
 
virtual void get_node_reordering (Vector< Node * > &reordering, const bool &use_old_ordering=true) const
 
template<class BULK_ELEMENT , template< class > class FACE_ELEMENT>
void build_face_mesh (const unsigned &b, Mesh *const &face_mesh_pt)
 
unsigned self_test ()
 Self-test: Check elements and nodes. Return 0 for OK. More...
 
void max_and_min_element_size (double &max_size, double &min_size)
 
double total_size ()
 
void check_inverted_elements (bool &mesh_has_inverted_elements, std::ofstream &inverted_element_file)
 
void check_inverted_elements (bool &mesh_has_inverted_elements)
 
unsigned check_for_repeated_nodes (const double &epsilon=1.0e-12)
 
Vector< Node * > prune_dead_nodes ()
 
unsigned nboundary () const
 Return number of boundaries. More...
 
unsigned long nboundary_node (const unsigned &ibound) const
 Return number of nodes on a particular boundary. More...
 
FiniteElementboundary_element_pt (const unsigned &b, const unsigned &e) const
 Return pointer to e-th finite element on boundary b. More...
 
Nodeget_some_non_boundary_node () const
 
unsigned nboundary_element (const unsigned &b) const
 Return number of finite elements that are adjacent to boundary b. More...
 
int face_index_at_boundary (const unsigned &b, const unsigned &e) const
 
virtual void dump (std::ofstream &dump_file, const bool &use_old_ordering=true) const
 Dump the data in the mesh into a file for restart. More...
 
void dump (const std::string &dump_file_name, const bool &use_old_ordering=true) const
 Dump the data in the mesh into a file for restart. More...
 
virtual void read (std::ifstream &restart_file)
 Read solution from restart file. More...
 
void output_paraview (std::ofstream &file_out, const unsigned &nplot) const
 
void output_fct_paraview (std::ofstream &file_out, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const
 
void output_fct_paraview (std::ofstream &file_out, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) const
 
void output (std::ostream &outfile)
 Output for all elements. More...
 
void output (std::ostream &outfile, const unsigned &n_plot)
 Output at f(n_plot) points in each element. More...
 
void output (FILE *file_pt)
 Output for all elements (C-style output) More...
 
void output (FILE *file_pt, const unsigned &nplot)
 Output at f(n_plot) points in each element (C-style output) More...
 
void output (const std::string &output_filename)
 Output for all elements. More...
 
void output (const std::string &output_filename, const unsigned &n_plot)
 Output at f(n_plot) points in each element. More...
 
void output_fct (std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt)
 Output a given Vector function at f(n_plot) points in each element. More...
 
void output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt)
 
void output_boundaries (std::ostream &outfile)
 Output the nodes on the boundaries (into separate tecplot zones) More...
 
void output_boundaries (const std::string &output_filename)
 
void assign_initial_values_impulsive ()
 Assign initial values for an impulsive start. More...
 
void shift_time_values ()
 
void calculate_predictions ()
 
void set_nodal_and_elemental_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
virtual void set_mesh_level_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
void set_consistent_pinned_values_for_continuation (ContinuationStorageScheme *const &continuation_stepper_pt)
 Set consistent values for pinned data in continuation. More...
 
bool does_pointer_correspond_to_mesh_data (double *const &parameter_pt)
 Does the double pointer correspond to any mesh data. More...
 
void set_nodal_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 Set the timestepper associated with the nodal data in the mesh. More...
 
void set_elemental_internal_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
virtual void compute_norm (double &norm)
 
virtual void compute_norm (Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 Returns the norm of the error and that of the exact solution. More...
 
virtual void compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
bool is_mesh_distributed () const
 Boolean to indicate if Mesh has been distributed. More...
 
OomphCommunicatorcommunicator_pt () const
 
void delete_all_external_storage ()
 Wipe the storage for all externally-based elements. More...
 

Additional Inherited Members

- Public Types inherited from oomph::Mesh
typedef void(FiniteElement::* SteadyExactSolutionFctPt) (const Vector< double > &x, Vector< double > &soln)
 
typedef void(FiniteElement::* UnsteadyExactSolutionFctPt) (const double &time, const Vector< double > &x, Vector< double > &soln)
 
- Static Public Attributes inherited from oomph::Mesh
static Steady< 0 > Default_TimeStepper
 The Steady Timestepper. More...
 
static bool Suppress_warning_about_empty_mesh_level_time_stepper_function
 Static boolean flag to control warning about mesh level timesteppers. More...
 
- Protected Member Functions inherited from oomph::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::Mesh
Vector< Vector< Node * > > Boundary_node_pt
 
bool Lookup_for_elements_next_boundary_is_setup
 
Vector< Vector< FiniteElement * > > Boundary_element_pt
 
Vector< Vector< int > > Face_index_at_boundary
 
Vector< Node * > Node_pt
 Vector of pointers to nodes. More...
 
Vector< GeneralisedElement * > Element_pt
 Vector of pointers to generalised elements. More...
 
std::vector< boolBoundary_coordinate_exists
 

Detailed Description

Mesh that is based on input files generated by the quadrilateral mesh generator Geompack.

Constructor & Destructor Documentation

◆ GeompackQuadScaffoldMesh() [1/3]

oomph::GeompackQuadScaffoldMesh::GeompackQuadScaffoldMesh ( )
inline

Empty constructor.

42 {}

◆ GeompackQuadScaffoldMesh() [2/3]

oomph::GeompackQuadScaffoldMesh::GeompackQuadScaffoldMesh ( const std::string &  mesh_file_name,
const std::string &  curve_file_name 
)

Constructor: Pass the filename of the mesh files.

35  {
36  // Read the output mesh file to find informations about the nodes
37  // and elements of the mesh
38 
39  // Process mesh file
40  //---------------------
41  std::ifstream mesh_file(mesh_file_name.c_str(), std::ios_base::in);
42 
43  // Number of nodes
44  unsigned n_node;
45  mesh_file >> n_node;
46 
47  // Coordinates of nodes and extra information index
48  Vector<double> x_node(n_node);
49  Vector<double> y_node(n_node);
50  Vector<int> vertinfo(n_node);
51  for (unsigned i = 0; i < n_node; i++)
52  {
53  mesh_file >> x_node[i];
54  mesh_file >> y_node[i];
55  mesh_file >> vertinfo[i];
56  }
57 
58  // Extra information (nodes lying on a curve)
59  unsigned n_vx;
60  mesh_file >> n_vx;
61 
62  // Dummy storage for the node code
63  int dummy_node_code;
64 
65  // Storage for the curve indice
66  Vector<int> icurv(n_vx); // it's negative if not used
67 
68  // Dummy storage for the location of the point on the curve
69  double dummy_ucurv;
70 
71  for (unsigned i = 0; i < n_vx; i++)
72  {
73  mesh_file >> dummy_node_code;
74  mesh_file >> icurv[i];
75  mesh_file >> dummy_ucurv;
76  }
77 
78  // Number of local nodes per element
79  unsigned n_local_node;
80  mesh_file >> n_local_node;
81 
82  // Number of elements
83  unsigned n_element;
84  mesh_file >> n_element;
85 
86  // Storage for global node numbers for all elements
87  Vector<unsigned> global_node(n_local_node * n_element);
88 
89  // Storage for edge information
90  // (needed for a possible construction of midside node
91  // in the following build from scaffold function)
92  Vector<int> edgeinfo(n_local_node * n_element);
93 
94  // Initialize counter
95  unsigned k = 0;
96 
97  // Read global node numbers for all elements
98  for (unsigned i = 0; i < n_element; i++)
99  {
100  for (unsigned j = 0; j < n_local_node; j++)
101  {
102  mesh_file >> global_node[k];
103  k++;
104  }
105  }
106 
107  // Initialize counter
108  unsigned l = 0;
109 
110  // Read the edge information
111  for (unsigned i = 0; i < n_element; i++)
112  {
113  for (unsigned j = 0; j < n_local_node; j++)
114  {
115  mesh_file >> edgeinfo[l];
116  l++;
117  }
118  }
119 
120  mesh_file.close();
121 
122  // Create a vector of boolean so as not to create the same node twice
123  std::vector<bool> done(n_node);
124  for (unsigned i = 0; i < n_node; i++)
125  {
126  done[i] = false;
127  }
128 
129  // Resize the Node vector
130  Node_pt.resize(n_node, 0);
131 
132  // Resize the Element vector
133  Element_pt.resize(n_element);
134 
135 
136  // Process curve file to extract information about boundaries
137  // ----------------------------------------------------------
138 
139  // Important: the input file must NOT have NURBS curve
140  std::ifstream curve_file(curve_file_name.c_str(), std::ios_base::in);
141 
142  // Number of curves
143  unsigned n_curv;
144  curve_file >> n_curv;
145 
146  // Storage of several information for each curve
147  Vector<Vector<int>> curv;
148 
149  // Resize to n_curv rows
150  curv.resize(n_curv);
151 
152  // Curve type
153  unsigned type;
154 
155  // Loop over the curves to extract information
156  for (unsigned i = 0; i < n_curv; i++)
157  {
158  curve_file >> type;
159  if (type == 1)
160  {
161  curv[i].resize(4);
162  curv[i][0] = type;
163  for (unsigned j = 1; j < 4; j++)
164  {
165  curve_file >> curv[i][j];
166  }
167  }
168  else if (type == 2)
169  {
170  curv[i].resize(5);
171  curv[i][0] = type;
172  for (unsigned j = 1; j < 5; j++)
173  {
174  curve_file >> curv[i][j];
175  }
176  }
177  else
178  {
179  std::ostringstream error_stream;
180  error_stream << "Current we can only process curves of\n"
181  << "type 1 (straight lines) and 2 (circular arcs\n"
182  << "You've specified: type " << type << std::endl;
183 
184  throw OomphLibError(
185  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
186  }
187  }
188 
189  curve_file.close();
190 
191  // Searching the number of boundaries
192  int d = 0;
193  for (unsigned i = 0; i < n_curv; i++)
194  {
195  if (d < curv[i][1])
196  {
197  d = curv[i][1]; // the boundary code is the 2nd value of each curve
198  }
199  }
200  oomph_info << "The number of boundaries is " << d << std::endl;
201 
202  // Set number of boundaries
203  if (d > 0)
204  {
205  set_nboundary(d);
206  }
207 
208  // Search the boundary number of node located on a boundary
209  // If after this, boundary_of_node[j][0] is -1 then node j
210  // is not located on any boundary.
211  // If boundary_of_node[j][0] is positive, the node is located
212  // on the boundary indicated by that number.
213  // If boundary_of_node[j][1] is also positive, the node is also
214  // located on that boundary. Note: We're ignoring the (remote)
215  // possibility that node is located on 3 or more boundaries
216  // as one of them would have to be an internal boundary which
217  // would be odd...
218  Vector<Vector<int>> boundary_of_node;
219  boundary_of_node.resize(n_node);
220  unsigned n;
221  for (unsigned i = 0; i < n_node; i++)
222  {
223  n = 0;
224  boundary_of_node[i].resize(2);
225  boundary_of_node[i][0] = -1;
226  boundary_of_node[i][1] = -1;
227  if (vertinfo[i] == 2) // it's an endpoint
228  {
229  for (unsigned j = 0; j < n_curv; j++)
230  {
231  for (unsigned m = 2; m < curv[j].size(); m++)
232  {
233  if (curv[j][m] ==
234  static_cast<int>(i + 1)) // node number begins at 1
235  { // in the mesh file !!!
236  boundary_of_node[i][n] = curv[j][1];
237  n++;
238  }
239  }
240  }
241  }
242  if (vertinfo[i] > 20)
243  {
244  int a = 0;
245  a = (vertinfo[i]) / 20;
246  int b;
247  b = icurv[a - 1]; // 1st value of vector at [0] !!
248  boundary_of_node[i][0] =
249  curv[b - 1][1]; // 1st value of vector at [0] !!
250  }
251  }
252 
253 
254  // Create the elements
255  //--------------------
256 
257  unsigned count = 0;
258  unsigned c;
259  for (unsigned e = 0; e < n_element; e++)
260  {
261  // Build simplex four node quad in the scaffold mesh
262  Element_pt[e] = new QElement<2, 2>;
263 
264  // Construction of the two first nodes of the element
265  for (unsigned j = 0; j < 2; j++)
266  {
267  c = global_node[count];
268  if (done[c - 1] == false) // c-1 because node number begins
269  // at 1 in the mesh file
270  {
271  // If the node is located on a boundary construct a boundary node
272  if ((d > 0) && ((boundary_of_node[c - 1][0] > 0) ||
273  (boundary_of_node[c - 1][1] > 0)))
274  {
275  // Construct a boundary node
277  // Add to the look=up schemes
278  if (boundary_of_node[c - 1][0] > 0)
279  {
280  add_boundary_node(boundary_of_node[c - 1][0] - 1, Node_pt[c - 1]);
281  }
282  if (boundary_of_node[c - 1][1] > 0)
283  {
284  add_boundary_node(boundary_of_node[c - 1][1] - 1, Node_pt[c - 1]);
285  }
286  }
287  // Otherwise construct a normal node
288  else
289  {
291  }
292  done[c - 1] = true;
293  Node_pt[c - 1]->x(0) = x_node[c - 1];
294  Node_pt[c - 1]->x(1) = y_node[c - 1];
295  }
296  else
297  {
298  finite_element_pt(e)->node_pt(j) = Node_pt[c - 1];
299  }
300  count++;
301  }
302 
303  // Construction of the third node not in anticlockwise order
304  c = global_node[count + 1];
305  if (done[c - 1] ==
306  false) // c-1 because node number begins at 1 in the mesh file
307  {
308  // If the node is on a boundary, construct a boundary node
309  if ((d > 0) && ((boundary_of_node[c - 1][0] > 0) ||
310  (boundary_of_node[c - 1][1] > 0)))
311  {
312  // Construct the node
314  // Add to the look-up schemes
315  if (boundary_of_node[c - 1][0] > 0)
316  {
317  add_boundary_node(boundary_of_node[c - 1][0] - 1, Node_pt[c - 1]);
318  }
319  if (boundary_of_node[c - 1][1] > 0)
320  {
321  add_boundary_node(boundary_of_node[c - 1][1] - 1, Node_pt[c - 1]);
322  }
323  }
324  // otherwise construct a normal node
325  else
326  {
328  }
329  done[c - 1] = true;
330  Node_pt[c - 1]->x(0) = x_node[c - 1];
331  Node_pt[c - 1]->x(1) = y_node[c - 1];
332  }
333  else
334  {
335  finite_element_pt(e)->node_pt(2) = Node_pt[c - 1];
336  }
337 
338  count++;
339 
340  // Construction of the fourth node
341  c = global_node[count - 1];
342  if (done[c - 1] ==
343  false) // c-1 because node number begins at 1 in the mesh file
344  {
345  // If the node is on a boundary, constuct a boundary node
346  if ((d > 0) && ((boundary_of_node[c - 1][0] > 0) ||
347  (boundary_of_node[c - 1][1] > 0)))
348  {
349  // Construct the boundary node
351  // Add to the look-up schemes
352  if (boundary_of_node[c - 1][0] > 0)
353  {
354  add_boundary_node(boundary_of_node[c - 1][0] - 1, Node_pt[c - 1]);
355  }
356  if (boundary_of_node[c - 1][1] > 0)
357  {
358  add_boundary_node(boundary_of_node[c - 1][1] - 1, Node_pt[c - 1]);
359  }
360  }
361  // Otherwise construct a normal node
362  else
363  {
365  }
366  done[c - 1] = true;
367  Node_pt[c - 1]->x(0) = x_node[c - 1];
368  Node_pt[c - 1]->x(1) = y_node[c - 1];
369  }
370  else
371  {
372  finite_element_pt(e)->node_pt(3) = Node_pt[c - 1];
373  }
374 
375  count++;
376  }
377  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Array< double, 1, 3 > e(1./3., 0.5, 2.)
Scalar * b
Definition: benchVecAdd.cpp:17
virtual Node * construct_node(const unsigned &n)
Definition: elements.h:2509
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
virtual Node * construct_boundary_node(const unsigned &n)
Definition: elements.h:2538
void add_boundary_node(const unsigned &b, Node *const &node_pt)
Add a (pointer to) a node to the b-th boundary.
Definition: mesh.cc:243
Vector< Node * > Node_pt
Vector of pointers to nodes.
Definition: mesh.h:183
FiniteElement * finite_element_pt(const unsigned &e) const
Definition: mesh.h:473
void set_nboundary(const unsigned &nbound)
Set the number of boundaries in the mesh.
Definition: mesh.h:505
Vector< GeneralisedElement * > Element_pt
Vector of pointers to generalised elements.
Definition: mesh.h:186
const Scalar * a
Definition: level2_cplx_impl.h:32
int * m
Definition: level2_cplx_impl.h:294
char char char int int * k
Definition: level2_impl.h:374
int c
Definition: calibrate.py:100
type
Definition: compute_granudrum_aor.py:141
OomphInfo oomph_info
Definition: oomph_definitions.cc:319
#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 a, oomph::Mesh::add_boundary_node(), b, calibrate::c, oomph::FiniteElement::construct_boundary_node(), oomph::FiniteElement::construct_node(), e(), oomph::Mesh::Element_pt, oomph::Mesh::finite_element_pt(), i, j, k, m, n, oomph::FiniteElement::node_pt(), oomph::Mesh::Node_pt, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::oomph_info, oomph::Mesh::set_nboundary(), and compute_granudrum_aor::type.

◆ GeompackQuadScaffoldMesh() [3/3]

oomph::GeompackQuadScaffoldMesh::GeompackQuadScaffoldMesh ( const GeompackQuadScaffoldMesh )
delete

Broken copy constructor.

◆ ~GeompackQuadScaffoldMesh()

oomph::GeompackQuadScaffoldMesh::~GeompackQuadScaffoldMesh ( )
inline

Empty destructor.

55 {}

Member Function Documentation

◆ operator=()

void oomph::GeompackQuadScaffoldMesh::operator= ( const GeompackQuadScaffoldMesh )
delete

Broken assignment operator.


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