oomph::EighthSphereMesh< ELEMENT > Class Template Reference

#include <eighth_sphere_mesh.template.h>

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

Public Member Functions

 EighthSphereMesh (const double &radius, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
 
 ~EighthSphereMesh ()
 Destructor. 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 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...
 

Protected Attributes

DomainDomain_pt
 Pointer to the domain. More...
 
double Radius
 Radius of the sphere. 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

template<class ELEMENT>
class oomph::EighthSphereMesh< ELEMENT >

Eight of a sphere brick mesh, based on the EightSphereDomain Non-refineable version with four brick elements. The eighth-sphere is located in the positive octant, centred at the origin. The mesh boundaries are numbered as follows:

  • Boundary 0: Plane x=0
  • Boundary 1: Plane y=0
  • Boundary 2: Plane z=0
  • Boundary 3: The surface of the sphere.

Constructor & Destructor Documentation

◆ EighthSphereMesh()

template<class ELEMENT >
oomph::EighthSphereMesh< ELEMENT >::EighthSphereMesh ( const double radius,
TimeStepper time_stepper_pt = &Mesh::Default_TimeStepper 
)

Constructor: Pass radius and timestepper; defaults to static default timestepper

Constructor for the eighth of a sphere mesh: Pass timestepper; defaults to static default timestepper.

41  : Radius(radius)
42  {
43  // Mesh can only be built with 3D Qelements.
44  MeshChecker::assert_geometric_element<QElementGeometricBase, ELEMENT>(3);
45 
46  // Set the number of boundaries
47  this->set_nboundary(4);
48 
49  // Provide storage for the four elements
50  this->Element_pt.resize(4);
51 
52  // Set the domain pointer: Pass the radius of the sphere
53  Domain_pt = new EighthSphereDomain(Radius);
54 
55 
56  Vector<double> s(3), s_fraction(3);
57  Vector<double> r(3);
58 
59  // Create first element
60  //--------------------
61  this->Element_pt[0] = new ELEMENT;
62 
63  // Give it nodes
64 
65  // How many 1D nodes are there?
66  unsigned nnode1d = this->finite_element_pt(0)->nnode_1d();
67 
68  // Loop over nodes
69  for (unsigned i = 0; i < nnode1d; i++)
70  {
71  for (unsigned j = 0; j < nnode1d; j++)
72  {
73  for (unsigned k = 0; k < nnode1d; k++)
74  {
75  unsigned jnod = k * nnode1d * nnode1d + j * nnode1d + i;
76 
77  // If we're on a boundary, make a boundary node
78  if ((i == 0) || (j == 0) || (k == 0))
79  {
80  // Allocate the node according to the element's wishes
81  this->Node_pt.push_back(
82  this->finite_element_pt(0)->construct_boundary_node(
83  jnod, time_stepper_pt));
84  }
85  // Otherwise it's a normal node
86  else
87  {
88  // Allocate the node according to the element's wishes
89  this->Node_pt.push_back(this->finite_element_pt(0)->construct_node(
90  jnod, time_stepper_pt));
91  }
92 
93  // Work out the node's coordinates in the finite element's local
94  // coordinate system
95  this->finite_element_pt(0)->local_fraction_of_node(jnod, s_fraction);
96 
97 
98  // Get the position of the node from macro element mapping
99  s[0] = -1.0 + 2.0 * s_fraction[0];
100  s[1] = -1.0 + 2.0 * s_fraction[1];
101  s[2] = -1.0 + 2.0 * s_fraction[2];
103 
104  // Assign coordinates
105  this->finite_element_pt(0)->node_pt(jnod)->x(0) = r[0];
106  this->finite_element_pt(0)->node_pt(jnod)->x(1) = r[1];
107  this->finite_element_pt(0)->node_pt(jnod)->x(2) = r[2];
108 
109  // Add the node to the appropriate boundary
110  if (i == 0)
111  add_boundary_node(0, this->finite_element_pt(0)->node_pt(jnod));
112  if (j == 0)
113  add_boundary_node(1, this->finite_element_pt(0)->node_pt(jnod));
114  if (k == 0)
115  add_boundary_node(2, this->finite_element_pt(0)->node_pt(jnod));
116  }
117  }
118  }
119 
120 
121  // Create a second element
122  //------------------------
123  this->Element_pt[1] = new ELEMENT;
124  ;
125 
126  // Give it nodes
127 
128  // Loop over nodes
129  for (unsigned i = 0; i < nnode1d; i++)
130  {
131  for (unsigned j = 0; j < nnode1d; j++)
132  {
133  for (unsigned k = 0; k < nnode1d; k++)
134  {
135  unsigned jnod = k * nnode1d * nnode1d + j * nnode1d + i;
136 
137  // If node has not been yet created, create it
138  if (i > 0)
139  {
140  // If the node is on a boundary, make a boundary node
141  if ((i == nnode1d - 1) || (j == 0) || (k == 0))
142  {
143  this->Node_pt.push_back(
144  this->finite_element_pt(1)->construct_boundary_node(
145  jnod, time_stepper_pt));
146  }
147  // Otherwise make a normal node
148  else
149  {
150  this->Node_pt.push_back(
151  this->finite_element_pt(1)->construct_node(jnod,
152  time_stepper_pt));
153  }
154 
155  // Work out the node's coordinates in the finite element's local
156  // coordinate system
158  s_fraction);
159 
160  // Get the position of the node from macro element mapping
161  s[0] = -1.0 + 2.0 * s_fraction[0];
162  s[1] = -1.0 + 2.0 * s_fraction[1];
163  s[2] = -1.0 + 2.0 * s_fraction[2];
165 
166  // Assign coordinate
167  this->finite_element_pt(1)->node_pt(jnod)->x(0) = r[0];
168  this->finite_element_pt(1)->node_pt(jnod)->x(1) = r[1];
169  this->finite_element_pt(1)->node_pt(jnod)->x(2) = r[2];
170 
171  // Add the node to the approriate boundary
172  if (j == 0)
173  add_boundary_node(1, this->finite_element_pt(1)->node_pt(jnod));
174  if (k == 0)
175  add_boundary_node(2, this->finite_element_pt(1)->node_pt(jnod));
176  if (i == nnode1d - 1)
177  add_boundary_node(3, this->finite_element_pt(1)->node_pt(jnod));
178  }
179 
180  // ...else use the node already created
181  else
182  {
183  this->finite_element_pt(1)->node_pt(jnod) =
184  this->finite_element_pt(0)->node_pt(jnod + nnode1d - 1);
185  }
186  }
187  }
188  }
189 
190 
191  // Create a third element
192  //------------------------
193  this->Element_pt[2] = new ELEMENT;
194 
195  // Give it nodes
196 
197  // Loop over nodes
198  for (unsigned i = 0; i < nnode1d; i++)
199  {
200  for (unsigned j = 0; j < nnode1d; j++)
201  {
202  for (unsigned k = 0; k < nnode1d; k++)
203  {
204  unsigned jnod = k * nnode1d * nnode1d + j * nnode1d + i;
205 
206  // If the node has not been yet created, create it
207  if ((i < nnode1d - 1) && (j > 0))
208  {
209  // If it's on a boundary, make a boundary node
210  if ((i == 0) || (j == nnode1d - 1) || (k == 0))
211  {
212  // Allocate the node according to the element's wishes
213  this->Node_pt.push_back(
214  this->finite_element_pt(2)->construct_boundary_node(
215  jnod, time_stepper_pt));
216  }
217  // Otherwise allocate a normal node
218  else
219  {
220  // Allocate the node according to the element's wishes
221  this->Node_pt.push_back(
222  this->finite_element_pt(2)->construct_node(jnod,
223  time_stepper_pt));
224  }
225 
226  // Work out the node's coordinates in the finite element's local
227  // coordinate system
229  s_fraction);
230 
231  // Get the position of the node from macro element mapping
232  s[0] = -1.0 + 2.0 * s_fraction[0];
233  s[1] = -1.0 + 2.0 * s_fraction[1];
234  s[2] = -1.0 + 2.0 * s_fraction[2];
236 
237  // Assign coordinates
238  this->finite_element_pt(2)->node_pt(jnod)->x(0) = r[0];
239  this->finite_element_pt(2)->node_pt(jnod)->x(1) = r[1];
240  this->finite_element_pt(2)->node_pt(jnod)->x(2) = r[2];
241 
242  // Add the node to the appropriate boundary
243  if (i == 0)
244  add_boundary_node(0, this->finite_element_pt(2)->node_pt(jnod));
245  if (k == 0)
246  add_boundary_node(2, this->finite_element_pt(2)->node_pt(jnod));
247  if (j == nnode1d - 1)
248  add_boundary_node(3, this->finite_element_pt(2)->node_pt(jnod));
249  }
250 
251  // ...else use the nodes already created
252  else
253  {
254  // If the node belongs to the element 0
255  if (j == 0)
256  this->finite_element_pt(2)->node_pt(jnod) =
257  this->finite_element_pt(0)->node_pt(jnod +
258  nnode1d * (nnode1d - 1));
259 
260  // ...else it belongs to the element 1
261  else if (i == nnode1d - 1)
262  this->finite_element_pt(2)->node_pt(jnod) =
263  this->finite_element_pt(1)->node_pt(nnode1d * nnode1d * k + j +
264  i * nnode1d);
265  }
266  }
267  }
268  }
269 
270 
271  // Create the fourth element
272  //-------------------------
273  this->Element_pt[3] = new ELEMENT;
274 
275  // Give it nodes
276 
277  // Loop over nodes
278  for (unsigned i = 0; i < nnode1d; i++)
279  {
280  for (unsigned j = 0; j < nnode1d; j++)
281  {
282  for (unsigned k = 0; k < nnode1d; k++)
283  {
284  unsigned jnod = k * nnode1d * nnode1d + j * nnode1d + i;
285 
286  // If the node has not been yet created, create it
287  if ((k > 0) && (i < nnode1d - 1) && (j < nnode1d - 1))
288  {
289  // If it's on a boundary, allocate a boundary node
290  if ((i == 0) || (j == 0) || (k == nnode1d - 1))
291  {
292  // Allocate the node according to the element's wishes
293  this->Node_pt.push_back(
294  this->finite_element_pt(3)->construct_boundary_node(
295  jnod, time_stepper_pt));
296  }
297  // Otherwise allocate a normal node
298  else
299  {
300  // Allocate the node according to the element's wishes
301  this->Node_pt.push_back(
302  this->finite_element_pt(3)->construct_node(jnod,
303  time_stepper_pt));
304  }
305 
306  // Work out the node's coordinates in the finite element's local
307  // coordinate system
309  s_fraction);
310 
311  // Get the position of the node from macro element mapping
312  s[0] = -1.0 + 2.0 * s_fraction[0];
313  s[1] = -1.0 + 2.0 * s_fraction[1];
314  s[2] = -1.0 + 2.0 * s_fraction[2];
316 
317  // Assign coordinates
318  this->finite_element_pt(3)->node_pt(jnod)->x(0) = r[0];
319  this->finite_element_pt(3)->node_pt(jnod)->x(1) = r[1];
320  this->finite_element_pt(3)->node_pt(jnod)->x(2) = r[2];
321 
322  // Add the node to the appropriate boundary
323  if (i == 0)
324  add_boundary_node(0, this->finite_element_pt(3)->node_pt(jnod));
325  if (j == 0)
326  add_boundary_node(1, this->finite_element_pt(3)->node_pt(jnod));
327  if (k == nnode1d - 1)
328  add_boundary_node(3, this->finite_element_pt(3)->node_pt(jnod));
329  }
330 
331  // ...otherwise the node was already created: use it.
332  else
333  {
334  // if k=0 then the node belongs to element 0
335  if (k == 0)
336  {
337  this->finite_element_pt(3)->node_pt(jnod) =
338  this->finite_element_pt(0)->node_pt(jnod + nnode1d * nnode1d *
339  (nnode1d - 1));
340  }
341  else
342  {
343  // else if i==nnode1d-1 the node already exists in element 1
344  if (i == nnode1d - 1)
345  {
346  this->finite_element_pt(3)->node_pt(jnod) =
347  this->finite_element_pt(1)->node_pt(
348  k + i * nnode1d * nnode1d + j * nnode1d);
349  }
350  else
351  // else, the node exists in element 2
352  {
353  this->finite_element_pt(3)->node_pt(jnod) =
354  this->finite_element_pt(2)->node_pt(i + k * nnode1d +
355  j * nnode1d * nnode1d);
356  }
357  }
358  }
359  }
360  }
361  }
362 
363  // Setup boundary element lookup schemes
365  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
void setup_boundary_element_info()
Definition: brick_mesh.h:195
MacroElement * macro_element_pt(const unsigned &i)
Access to i-th macro element.
Definition: domain.h:116
Domain * Domain_pt
Pointer to the domain.
Definition: eighth_sphere_mesh.template.h:71
double Radius
Radius of the sphere.
Definition: eighth_sphere_mesh.template.h:74
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
virtual unsigned nnode_1d() const
Definition: elements.h:2218
virtual void local_fraction_of_node(const unsigned &j, Vector< double > &s_fraction)
Definition: elements.cc:3191
void macro_map(const Vector< double > &s, Vector< double > &r)
The mapping from local to global coordinates at the current time : r(s)
Definition: macro_element.h:126
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
Node *& node_pt(const unsigned long &n)
Return pointer to global node n.
Definition: mesh.h:436
Vector< GeneralisedElement * > Element_pt
Vector of pointers to generalised elements.
Definition: mesh.h:186
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition: nodes.h:1060
RealScalar s
Definition: level1_cplx_impl.h:130
char char char int int * k
Definition: level2_impl.h:374
r
Definition: UniformPSDSelfTest.py:20
radius
Definition: UniformPSDSelfTest.py:15
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References oomph::Mesh::add_boundary_node(), oomph::EighthSphereMesh< ELEMENT >::Domain_pt, oomph::Mesh::Element_pt, oomph::Mesh::finite_element_pt(), i, j, k, oomph::FiniteElement::local_fraction_of_node(), oomph::Domain::macro_element_pt(), oomph::MacroElement::macro_map(), oomph::FiniteElement::nnode_1d(), oomph::FiniteElement::node_pt(), oomph::Mesh::Node_pt, oomph::Mesh::node_pt(), UniformPSDSelfTest::r, oomph::EighthSphereMesh< ELEMENT >::Radius, s, oomph::Mesh::set_nboundary(), oomph::BrickMeshBase::setup_boundary_element_info(), and oomph::Node::x().

◆ ~EighthSphereMesh()

template<class ELEMENT >
oomph::EighthSphereMesh< ELEMENT >::~EighthSphereMesh ( )
inline

Destructor.

64  {
65  delete Domain_pt;
66  Domain_pt = 0;
67  }

References oomph::EighthSphereMesh< ELEMENT >::Domain_pt.

Member Data Documentation

◆ Domain_pt

◆ Radius

template<class ELEMENT >
double oomph::EighthSphereMesh< ELEMENT >::Radius
protected

Radius of the sphere.

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


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