oomph::RefineableTetMeshBase Class Reference

Base class for refineable tet meshes. More...

#include <refineable_mesh.h>

+ Inheritance diagram for oomph::RefineableTetMeshBase:

Public Member Functions

doublemax_element_size ()
 Max element size allowed during adaptation. More...
 
doublemin_element_size ()
 Min element size allowed during adaptation. More...
 
doublemax_permitted_edge_ratio ()
 Min edge ratio before remesh gets triggered. More...
 
void doc_adaptivity_targets (std::ostream &outfile)
 Doc the targets for mesh adaptation. More...
 
double compute_volume_target (const Vector< double > &elem_error, Vector< double > &target_volume)
 
- Public Member Functions inherited from oomph::RefineableMeshBase
bool adapt_flag ()
 
 RefineableMeshBase ()
 
 RefineableMeshBase (const RefineableMeshBase &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const RefineableMeshBase &)=delete
 Broken assignment operator. More...
 
virtual ~RefineableMeshBase ()
 Empty Destructor: More...
 
unsigned nrefined ()
 Access fct for number of elements that were refined. More...
 
unsigned nunrefined ()
 Access fct for number of elements that were unrefined. More...
 
unsignednrefinement_overruled ()
 
unsignedmax_keep_unrefined ()
 
ErrorEstimator *& spatial_error_estimator_pt ()
 Access to spatial error estimator. More...
 
ErrorEstimatorspatial_error_estimator_pt () const
 Access to spatial error estimator (const version. More...
 
doublemin_permitted_error ()
 
doublemax_permitted_error ()
 
doublemin_error ()
 
doublemax_error ()
 
DocInfo *& doc_info_pt ()
 Access fct for pointer to DocInfo. More...
 
void enable_adaptation ()
 Enable adaptation. More...
 
void disable_adaptation ()
 Disable adaptation. More...
 
void enable_p_adaptation ()
 Enable adaptation. More...
 
void disable_p_adaptation ()
 Disable adaptation. More...
 
void enable_additional_synchronisation_of_hanging_nodes ()
 Enable additional synchronisation of hanging nodes. More...
 
void disable_additional_synchronisation_of_hanging_nodes ()
 Disable additional synchronisation of hanging nodes. More...
 
bool is_adaptation_enabled () const
 Return whether the mesh is to be adapted. More...
 
bool is_p_adaptation_enabled () const
 Return whether the mesh is to be adapted. More...
 
bool is_additional_synchronisation_of_hanging_nodes_disabled () const
 Return whether additional synchronisation is enabled. More...
 
DocInfo doc_info ()
 Access fct for DocInfo. More...
 
virtual void adapt (const Vector< double > &elemental_error)=0
 
virtual void p_adapt (const Vector< double > &elemental_error)
 
virtual void refine_uniformly (DocInfo &doc_info)=0
 Refine mesh uniformly and doc process. More...
 
virtual void refine_uniformly ()
 Refine mesh uniformly. More...
 
virtual void p_refine_uniformly (DocInfo &doc_info)
 p-refine mesh uniformly and doc process More...
 
virtual void p_refine_uniformly ()
 p-refine mesh uniformly More...
 
virtual unsigned unrefine_uniformly ()=0
 
void p_unrefine_uniformly (DocInfo &doc_info)
 p-unrefine mesh uniformly 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...
 

Public Attributes

double Max_element_size
 Max permitted element size. More...
 
double Min_element_size
 Min permitted element size. More...
 
double Max_permitted_edge_ratio
 Max edge ratio before remesh gets triggered. 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::RefineableMeshBase
ErrorEstimatorSpatial_error_estimator_pt
 Pointer to spatial error estimator. More...
 
double Max_permitted_error
 Max. error (i.e. split elements if their error is larger) More...
 
double Min_permitted_error
 Min. error (i.e. (try to) merge elements if their error is smaller) More...
 
double Min_error
 Min.actual error. More...
 
double Max_error
 Max. actual error. More...
 
unsigned Nrefined
 Stats: Number of elements that were refined. More...
 
unsigned Nunrefined
 Stats: Number of elements that were unrefined. More...
 
bool Adapt_flag
 Flag that requests adaptation. More...
 
bool P_adapt_flag
 Flag that requests p-adaptation. More...
 
bool Additional_synchronisation_of_hanging_nodes_not_required
 Flag that disables additional synchronisation of hanging nodes. More...
 
DocInfoDoc_info_pt
 Pointer to DocInfo. More...
 
unsigned Max_keep_unrefined
 
unsigned Nrefinement_overruled
 
- 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

Base class for refineable tet meshes.

///////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////

Member Function Documentation

◆ compute_volume_target()

double oomph::RefineableTetMeshBase::compute_volume_target ( const Vector< double > &  elem_error,
Vector< double > &  target_volume 
)
inline

Compute target volume based on the elements' error and the error target; return max edge ratio

924  {
925  double max_edge_ratio = 0.0;
926  unsigned count_unrefined = 0;
927  unsigned count_refined = 0;
928  this->Nrefinement_overruled = 0;
929 
930  unsigned nel = this->nelement();
931  for (unsigned e = 0; e < nel; e++)
932  {
933  // Get element
934  FiniteElement* el_pt = this->finite_element_pt(e);
935 
936  // Calculate the volume of the element
937  double volume = el_pt->size();
938 
939  // Find the vertex coordinates
940  // (vertices are enumerated first)
941  double vertex[4][3];
942  for (unsigned n = 0; n < 4; ++n)
943  {
944  for (unsigned i = 0; i < 3; ++i)
945  {
946  vertex[n][i] = el_pt->node_pt(n)->x(i);
947  }
948  }
949 
950  // Compute the radius of the circumsphere of the tetrahedron
951  // Algorithm stolen from tetgen for consistency
952  DenseDoubleMatrix A(3);
953  for (unsigned i = 0; i < 3; ++i)
954  {
955  A(0, i) = vertex[1][i] - vertex[0][i];
956  A(1, i) = vertex[2][i] - vertex[0][i];
957  A(2, i) = vertex[3][i] - vertex[0][i];
958  }
959 
960  Vector<double> rhs(3);
961  // Compute the right hand side vector b (3x1).
962  for (unsigned i = 0; i < 3; ++i)
963  {
964  rhs[i] = 0.0;
965  for (unsigned k = 0; k < 3; ++k)
966  {
967  rhs[i] += A(i, k) * A(i, k);
968  }
969  rhs[i] *= 0.5;
970  }
971 
972  // Solve the linear system, in which the rhs is over-written with
973  // the solution
974  A.solve(rhs);
975 
976  // Calculate the circum-radius
977  double circum_radius =
978  sqrt(rhs[0] * rhs[0] + rhs[1] * rhs[1] + rhs[2] * rhs[2]);
979 
980  // Now find the shortest edge length
981  Vector<double> edge(3);
982  double min_length = DBL_MAX;
983  for (unsigned start = 0; start < 4; ++start)
984  {
985  for (unsigned end = start + 1; end < 4; ++end)
986  {
987  for (unsigned i = 0; i < 3; ++i)
988  {
989  edge[i] = vertex[start][i] - vertex[end][i];
990  }
991  double length =
992  sqrt(edge[0] * edge[0] + edge[1] * edge[1] + edge[2] * edge[2]);
993  if (length < min_length)
994  {
995  min_length = length;
996  }
997  }
998  }
999 
1000  // Now calculate the minimum edge ratio for this element
1001  double local_max_edge_ratio = circum_radius / min_length;
1002  if (local_max_edge_ratio > max_edge_ratio)
1003  {
1004  max_edge_ratio = local_max_edge_ratio;
1005  }
1006 
1007  // Mimick refinement in tree-based procedure: Target volumes
1008  // for elements that exceed permitted error is 1/4 of their
1009  // current volume, corresponding to a uniform sub-division.
1010  if (elem_error[e] > this->max_permitted_error())
1011  {
1012  target_volume[e] = std::max(volume / 4.0, Min_element_size);
1013  if (target_volume[e] != Min_element_size)
1014  {
1015  count_refined++;
1016  }
1017  else
1018  {
1019  this->Nrefinement_overruled++;
1020  }
1021  }
1022  else if (elem_error[e] < this->min_permitted_error())
1023  {
1024  target_volume[e] = std::min(4.0 * volume, Max_element_size);
1025  if (target_volume[e] != Max_element_size)
1026  {
1027  count_unrefined++;
1028  }
1029  }
1030  else
1031  {
1032  // Leave it alone
1033  target_volume[e] = std::max(volume, Min_element_size);
1034  }
1035 
1036  } // End of loop over elements
1037 
1038  // Tell everybody
1039  this->Nrefined = count_refined;
1040  this->Nunrefined = count_unrefined;
1041 
1042  return max_edge_ratio;
1043  }
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
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.)
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
FiniteElement * finite_element_pt(const unsigned &e) const
Definition: mesh.h:473
unsigned long nelement() const
Return number of elements in the mesh.
Definition: mesh.h:590
unsigned Nrefined
Stats: Number of elements that were refined.
Definition: refineable_mesh.h:338
unsigned Nunrefined
Stats: Number of elements that were unrefined.
Definition: refineable_mesh.h:341
double & min_permitted_error()
Definition: refineable_mesh.h:156
unsigned Nrefinement_overruled
Definition: refineable_mesh.h:363
double & max_permitted_error()
Definition: refineable_mesh.h:163
double Max_element_size
Max permitted element size.
Definition: refineable_mesh.h:1046
double Min_element_size
Min permitted element size.
Definition: refineable_mesh.h:1049
#define min(a, b)
Definition: datatypes.h:22
#define max(a, b)
Definition: datatypes.h:23
static constexpr lastp1_t end
Definition: IndexedViewHelper.h:79
char char char int int * k
Definition: level2_impl.h:374
void start(const unsigned &i)
(Re-)start i-th timer
Definition: oomph_utilities.cc:243

References e(), Eigen::placeholders::end, oomph::Mesh::finite_element_pt(), i, k, max, Max_element_size, oomph::RefineableMeshBase::max_permitted_error(), min, Min_element_size, oomph::RefineableMeshBase::min_permitted_error(), n, oomph::Mesh::nelement(), oomph::FiniteElement::node_pt(), oomph::RefineableMeshBase::Nrefined, oomph::RefineableMeshBase::Nrefinement_overruled, oomph::RefineableMeshBase::Nunrefined, oomph::FiniteElement::size(), sqrt(), oomph::CumulativeTimings::start(), and oomph::Node::x().

◆ doc_adaptivity_targets()

void oomph::RefineableTetMeshBase::doc_adaptivity_targets ( std::ostream &  outfile)
inlinevirtual

Doc the targets for mesh adaptation.

Reimplemented from oomph::RefineableMeshBase.

904  {
905  outfile << std::endl;
906  outfile << "Targets for mesh adaptation: " << std::endl;
907  outfile << "---------------------------- " << std::endl;
908  outfile << "Target for max. error: " << Max_permitted_error << std::endl;
909  outfile << "Target for min. error: " << Min_permitted_error << std::endl;
910  outfile << "Target max edge ratio: " << Max_permitted_edge_ratio
911  << std::endl;
912  outfile << "Min. allowed element size: " << Min_element_size << std::endl;
913  outfile << "Max. allowed element size: " << Max_element_size << std::endl;
914  outfile << "Don't unrefine if less than " << Max_keep_unrefined
915  << " elements need unrefinement." << std::endl;
916  outfile << std::endl;
917  }
unsigned Max_keep_unrefined
Definition: refineable_mesh.h:359
double Min_permitted_error
Min. error (i.e. (try to) merge elements if their error is smaller)
Definition: refineable_mesh.h:329
double Max_permitted_error
Max. error (i.e. split elements if their error is larger)
Definition: refineable_mesh.h:326
double Max_permitted_edge_ratio
Max edge ratio before remesh gets triggered.
Definition: refineable_mesh.h:1052

References Max_element_size, oomph::RefineableMeshBase::Max_keep_unrefined, Max_permitted_edge_ratio, oomph::RefineableMeshBase::Max_permitted_error, Min_element_size, and oomph::RefineableMeshBase::Min_permitted_error.

◆ max_element_size()

double& oomph::RefineableTetMeshBase::max_element_size ( )
inline

Max element size allowed during adaptation.

885  {
886  return Max_element_size;
887  }

References Max_element_size.

◆ max_permitted_edge_ratio()

double& oomph::RefineableTetMeshBase::max_permitted_edge_ratio ( )
inline

Min edge ratio before remesh gets triggered.

897  {
899  }

References Max_permitted_edge_ratio.

◆ min_element_size()

double& oomph::RefineableTetMeshBase::min_element_size ( )
inline

Min element size allowed during adaptation.

891  {
892  return Min_element_size;
893  }

References Min_element_size.

Member Data Documentation

◆ Max_element_size

double oomph::RefineableTetMeshBase::Max_element_size

◆ Max_permitted_edge_ratio

double oomph::RefineableTetMeshBase::Max_permitted_edge_ratio

◆ Min_element_size

double oomph::RefineableTetMeshBase::Min_element_size

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