///////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////// More...
#include <refineable_mesh.h>
Classes | |
struct | HangHelperStruct |
Helper struct to collate data required during TreeBasedRefineableMeshBase::synchronise_hanging_nodes. More... | |
Public Member Functions | |
TreeBasedRefineableMeshBase () | |
Constructor. More... | |
TreeBasedRefineableMeshBase (const TreeBasedRefineableMeshBase &dummy)=delete | |
Broken copy constructor. More... | |
void | operator= (const TreeBasedRefineableMeshBase &)=delete |
Broken assignment operator. More... | |
virtual | ~TreeBasedRefineableMeshBase () |
Empty Destructor: More... | |
void | adapt (const Vector< double > &elemental_error) |
Adapt mesh: Refine elements whose error is lager than err_max and (try to) unrefine those whose error is smaller than err_min. More... | |
void | p_adapt (const Vector< double > &elemental_error) |
p-adapt mesh: Refine elements whose error is lager than err_max and (try to) unrefine those whose error is smaller than err_min More... | |
void | refine_uniformly (DocInfo &doc_info) |
Refine mesh uniformly and doc process. More... | |
void | refine_uniformly () |
Refine mesh uniformly. More... | |
void | p_refine_uniformly (DocInfo &doc_info) |
p-refine mesh uniformly and doc process More... | |
void | p_refine_uniformly () |
p-refine mesh uniformly More... | |
unsigned | unrefine_uniformly () |
Unrefine mesh uniformly: Return 0 for success, 1 for failure (if unrefinement has reached the coarsest permitted level) More... | |
void | p_unrefine_uniformly (DocInfo &doc_info) |
p-unrefine mesh uniformly More... | |
virtual void | setup_tree_forest ()=0 |
Set up the tree forest associated with the Mesh (if any) More... | |
TreeForest * | forest_pt () |
Return pointer to the Forest represenation of the mesh. More... | |
void | doc_adaptivity_targets (std::ostream &outfile) |
Doc the targets for mesh adaptation. More... | |
unsigned & | max_refinement_level () |
Access fct for max. permissible refinement level (relative to base mesh) More... | |
unsigned & | min_refinement_level () |
Access fct for min. permissible refinement level (relative to base mesh) More... | |
unsigned & | max_p_refinement_level () |
Access fct for max. permissible p-refinement level (relative to base mesh) More... | |
unsigned & | min_p_refinement_level () |
Access fct for min. permissible p-refinement level (relative to base mesh) More... | |
virtual void | adapt_mesh (DocInfo &doc_info) |
Perform the actual tree-based mesh adaptation, documenting the progress in the directory specified in DocInfo object. More... | |
virtual void | adapt_mesh () |
Perform the actual tree-based mesh adaptation. A simple wrapper to call the function without documentation. More... | |
void | p_adapt_mesh (DocInfo &doc_info) |
Perform the actual tree-based mesh p-adaptation, documenting the progress in the directory specified in DocInfo object. More... | |
void | p_adapt_mesh () |
Perform the actual tree-based mesh p-adaptation. A simple wrapper to call the function without documentation. More... | |
virtual void | refine_selected_elements (const Vector< unsigned > &elements_to_be_refined) |
Refine mesh by splitting the elements identified by their numbers. More... | |
virtual void | refine_selected_elements (const Vector< RefineableElement * > &elements_to_be_refined) |
Refine mesh by splitting the elements identified by their pointers. More... | |
void | p_refine_selected_elements (const Vector< unsigned > &elements_to_be_refined) |
p-refine mesh by refining the elements identified by their numbers. More... | |
void | p_refine_selected_elements (const Vector< PRefineableElement * > &elements_to_be_refined_pt) |
p-refine mesh by refining the elements identified by their pointers. More... | |
virtual void | refine_base_mesh_as_in_reference_mesh (TreeBasedRefineableMeshBase *const &ref_mesh_pt) |
Refine base mesh to same degree as reference mesh (relative to original unrefined mesh). More... | |
virtual bool | refine_base_mesh_as_in_reference_mesh_minus_one (TreeBasedRefineableMeshBase *const &ref_mesh_pt) |
Refine base mesh to same degree as reference mesh minus one level of refinement (relative to original unrefined mesh). Useful function for multigrid solvers; allows the easy copy of a mesh to the level of refinement just below the current one. More... | |
virtual void | refine_as_in_reference_mesh (TreeBasedRefineableMeshBase *const &ref_mesh_pt) |
Refine mesh once so that its topology etc becomes that of the (finer!) reference mesh – if possible! Useful for meshes in multigrid hierarchies. If the meshes are too different and the conversion cannot be performed, the code dies (provided PARANOID is enabled). More... | |
virtual void | get_refinement_levels (unsigned &min_refinement_level, unsigned &max_refinement_level) |
Get max/min refinement levels in mesh. More... | |
virtual void | get_elements_at_refinement_level (unsigned &refinement_level, Vector< RefineableElement * > &level_elements) |
Extract the elements at a particular refinement level in the refinement pattern - used in Mesh::redistribute or whatever it's going to be called (RefineableMeshBase::reduce_halo_layers or something) More... | |
virtual void | get_refinement_pattern (Vector< Vector< unsigned >> &to_be_refined) |
Extract refinement pattern: Consider the hypothetical mesh obtained by truncating the refinement of the current mesh to a given level (where level=0 is the un-refined base mesh). To advance to the next refinement level, we need to refine (split) the to_be_refined [level].size() elements identified by the element numbers contained in vector to_be_refined[level][...]. More... | |
void | refine_base_mesh (Vector< Vector< unsigned >> &to_be_refined) |
Refine base mesh according to specified refinement pattern. More... | |
virtual void | refine (std::ifstream &restart_file) |
Refine mesh according to refinement pattern in restart file. More... | |
virtual void | dump_refinement (std::ostream &outfile) |
Dump refinement pattern to allow for rebuild. More... | |
virtual void | read_refinement (std::ifstream &restart_file, Vector< Vector< unsigned >> &to_be_refined) |
Read refinement pattern to allow for rebuild. More... | |
unsigned | uniform_refinement_level_when_pruned () const |
Level to which the mesh was uniformly refined when it was pruned (const version) More... | |
unsigned & | uniform_refinement_level_when_pruned () |
Level to which the mesh was uniformly refined when it was pruned. More... | |
void | classify_halo_and_haloed_nodes (DocInfo &doc_info, const bool &report_stats) |
Classify all halo and haloed information in the mesh (overloaded version from Mesh base class. Calls that one first, then synchronises hanging nodes) More... | |
void | classify_halo_and_haloed_nodes (const bool &report_stats=false) |
Classify the halo and haloed nodes in the mesh. More... | |
Public Member Functions inherited from oomph::RefineableMeshBase | |
bool | adapt_flag () |
RefineableMeshBase () | |
Constructor sets default values for refinement targets etc. and initialises pointer to spatial error estimator to NULL. More... | |
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... | |
unsigned & | nrefinement_overruled () |
Number of elements that would have liked to be refined further but can't because they've reached the max. refinement level. More... | |
unsigned & | max_keep_unrefined () |
Max. number of elements that we allow to remain unrefined if no other mesh adaptation is required (to avoid mesh-adaptations that would only unrefine a few elements and then force a new solve – this can't be worth our while!) More... | |
ErrorEstimator *& | spatial_error_estimator_pt () |
Access to spatial error estimator. More... | |
ErrorEstimator * | spatial_error_estimator_pt () const |
Access to spatial error estimator (const version. More... | |
double & | min_permitted_error () |
Access fct for min. error (i.e. (try to) merge elements if their error is smaller) More... | |
double & | max_permitted_error () |
Access fct for max. error (i.e. split elements if their error is larger) More... | |
double & | min_error () |
Access fct for min. actual error in present solution (i.e. before re-solve on adapted mesh) More... | |
double & | max_error () |
Access fct for max. actual error in present solution (i.e. before re-solve on adapted mesh) More... | |
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... | |
void | p_unrefine_uniformly (DocInfo &doc_info) |
p-unrefine mesh uniformly More... | |
Public Member Functions inherited from oomph::Mesh | |
void | resize_halo_nodes () |
Helper function that resizes halo nodes to the same size as their non-halo counterparts if required. (A discrepancy can arise if a FaceElement that introduces additional unknowns are attached to a bulk element that shares a node with a haloed element. In that case the joint node between haloed and non-haloed element is resized on that processor but not on the one that holds the halo counterpart (because no FaceElement is attached to the halo element) More... | |
Mesh () | |
Default constructor. More... | |
Mesh (const Vector< Mesh * > &sub_mesh_pt) | |
Constructor builds combined mesh from the meshes specified. Note: This simply merges the meshes' elements and nodes (ignoring duplicates; no boundary information etc. is created). More... | |
void | merge_meshes (const Vector< Mesh * > &sub_mesh_pt) |
Merge meshes. Note: This simply merges the meshes' elements and nodes (ignoring duplicates; no boundary information etc. is created). More... | |
virtual void | setup_boundary_element_info () |
Interface for function that is used to setup the boundary information (Empty virtual function – implement this for specific Mesh classes) More... | |
virtual void | setup_boundary_element_info (std::ostream &outfile) |
Setup lookup schemes which establish whic elements are located next to mesh's boundaries. Doc in outfile (if it's open). (Empty virtual function – implement this for specific Mesh classes) More... | |
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) |
Output boundary coordinates on boundary b – template argument specifies the bulk element type (needed to create FaceElement of appropriate type on mesh boundary). More... | |
virtual void | scale_mesh (const double &factor) |
Scale all nodal coordinates by given factor. Virtual so it can be overloaded in SolidMesh class where it also re-assigns the Lagrangian coordinates. More... | |
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 () |
Flush storage for elements and nodes by emptying the vectors that store the pointers to them. This is useful if a particular mesh is only built to generate a small part of a bigger mesh. Once the elements and nodes have been created, they are typically copied into the new mesh and the auxiliary mesh can be deleted. However, if we simply call the destructor of the auxiliary mesh, it will also wipe out the nodes and elements, because it still "thinks" it's in charge of these... More... | |
void | flush_element_storage () |
Flush storage for elements (only) by emptying the vectors that store the pointers to them. This is useful if a particular mesh is only built to generate a small part of a bigger mesh. Once the elements and nodes have been created, they are typically copied into the new mesh and the auxiliary mesh can be deleted. However, if we simply call the destructor of the auxiliary mesh, it will also wipe out the nodes and elements, because it still "thinks" it's in charge of these... More... | |
void | flush_node_storage () |
Flush storage for nodes (only) by emptying the vectors that store the pointers to them. More... | |
Node *& | node_pt (const unsigned long &n) |
Return pointer to global node n. More... | |
Node * | node_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... | |
GeneralisedElement * | element_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... | |
FiniteElement * | finite_element_pt (const unsigned &e) const |
Upcast (downcast?) to FiniteElement (needed to access FiniteElement member functions). More... | |
Node *& | boundary_node_pt (const unsigned &b, const unsigned &n) |
Return pointer to node n on boundary b. More... | |
Node * | boundary_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) |
Remove all information about nodes stored on the b-th boundary of the mesh. More... | |
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 () |
Replace existing boundary node lookup schemes with new schemes created using the boundary data stored in the nodes. More... | |
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) |
Update nodal positions in response to changes in the domain shape. Uses the FiniteElement::get_x(...) function for FiniteElements and doesn't do anything for other element types. If a MacroElement pointer has been set for a FiniteElement, the MacroElement representation is used to update the nodal positions; if not get_x(...) uses the FE interpolation and thus leaves the nodal positions unchanged. Virtual, so it can be overloaded by specific meshes, such as AlgebraicMeshes or SpineMeshes. Generally, this function updates the position of all nodes in response to changes in the boundary position. However, we ignore all SolidNodes since their position is computed as part of the solution – unless the bool flag is set to true. Such calls are typically made when the initial mesh is created and/or after a mesh has been refined repeatedly before the start of the computation. More... | |
virtual void | reorder_nodes (const bool &use_old_ordering=true) |
Re-order nodes in the order in which they appear in elements – can be overloaded for more efficient re-ordering. More... | |
virtual void | get_node_reordering (Vector< Node * > &reordering, const bool &use_old_ordering=true) const |
Get a reordering of the nodes in the order in which they appear in elements – can be overloaded for more efficient re-ordering. More... | |
template<class BULK_ELEMENT , template< class > class FACE_ELEMENT> | |
void | build_face_mesh (const unsigned &b, Mesh *const &face_mesh_pt) |
Constuct a Mesh of FACE_ELEMENTs along the b-th boundary of the mesh (which contains elements of type BULK_ELEMENT) More... | |
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) |
Determine max and min area for all FiniteElements in the mesh (non-FiniteElements are ignored) More... | |
double | total_size () |
Determine the sum of all "sizes" of the FiniteElements in the mesh (non-FiniteElements are ignored). This gives the length/area/volume occupied by the mesh. More... | |
void | check_inverted_elements (bool &mesh_has_inverted_elements, std::ofstream &inverted_element_file) |
Check for inverted elements and report outcome in boolean variable. This visits all elements at their integration points and checks if the Jacobian of the mapping between local and global coordinates is positive – using the same test that would be carried out (but only in PARANOID mode) during the assembly of the elements' Jacobian matrices. Inverted elements are output in inverted_element_file (if the stream is open). More... | |
void | check_inverted_elements (bool &mesh_has_inverted_elements) |
Check for inverted elements and report outcome in boolean variable. This visits all elements at their integration points and checks if the Jacobian of the mapping between local and global coordinates is positive – using the same test that would be carried out (but only in PARANOID mode) during the assembly of the elements' Jacobian matrices. More... | |
unsigned | check_for_repeated_nodes (const double &epsilon=1.0e-12) |
Check for repeated nodes within a given spatial tolerance. Return (0/1) for (pass/fail). More... | |
Vector< Node * > | prune_dead_nodes () |
Prune nodes. Nodes that have been marked as obsolete are removed from the mesh (and its boundary-node scheme). Returns vector of pointers to deleted nodes. More... | |
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... | |
FiniteElement * | boundary_element_pt (const unsigned &b, const unsigned &e) const |
Return pointer to e-th finite element on boundary b. More... | |
Node * | get_some_non_boundary_node () const |
Find a node not on any boundary in mesh_pt (useful for pinning a single node in a purely Neumann problem so that it is fully determined). More... | |
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 |
For the e-th finite element on boundary b, return int to indicate the face_index of the face adjacent to the boundary. This is consistent with input required during the generation of FaceElements. More... | |
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 |
Output in paraview format into specified file. Breaks up each element into sub-elements for plotting purposes. We assume that all elements are of the same type (fct will break break (in paranoid mode) if paraview output fcts of the elements are inconsistent). More... | |
void | output_fct_paraview (std::ofstream &file_out, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const |
Output in paraview format into specified file. Breaks up each element into sub-elements for plotting purposes. We assume that all elements are of the same type (fct will break break (in paranoid mode) if paraview output fcts of the elements are inconsistent). More... | |
void | output_fct_paraview (std::ofstream &file_out, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) const |
Output in paraview format into specified file. Breaks up each element into sub-elements for plotting purposes. We assume that all elements are of the same type (fct will break break (in paranoid mode) if paraview output fcts of the elements are inconsistent). More... | |
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) |
Output a given time-dep. Vector function at f(n_plot) points in each element. More... | |
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) |
Output the nodes on the boundaries (into separate tecplot zones). Specify filename. More... | |
void | assign_initial_values_impulsive () |
Assign initial values for an impulsive start. More... | |
void | shift_time_values () |
Shift time-dependent data along for next timestep: Deal with nodal Data/positions and the element's internal Data. More... | |
void | calculate_predictions () |
Calculate predictions for all Data and positions associated with the mesh, usually used in adaptive time-stepping. More... | |
void | set_nodal_and_elemental_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data) |
Set the timestepper associated with all nodal and elemental data stored in the mesh. More... | |
virtual void | set_mesh_level_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data) |
Function that can be used to set any additional timestepper data stored at the Mesh (as opposed to nodal and elemental) levels. This is virtual so that it can be overloaded in the appropriate Meshes. Examples include the SpineMeshes and adaptive triangle and tet meshes. More... | |
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 ¶meter_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) |
Set the timestepper associated with the internal data stored within elements in the meah. More... | |
virtual void | compute_norm (double &norm) |
Compute norm of solution by summing contributions of compute_norm(...) for all constituent elements in the mesh. What that norm means depends on what's defined in the element's function; may need to take the square root afterwards if the elements compute the square of the L2 norm, say. More... | |
virtual void | compute_norm (Vector< double > &norm) |
Compute norm of solution by summing contributions of compute_norm(...) for all constituent elements in the mesh. What that norm means depends on what's defined in the element's function; may need to take the square root afterwards if the elements compute the square of the L2 norm, say. More... | |
virtual void | compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm) |
Plot error when compared against a given exact solution. Also returns the norm of the error and that of the exact solution. More... | |
virtual void | compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm) |
Plot error when compared against a given time-depdendent exact solution. Also returns the norm of the error and that of the exact solution. More... | |
virtual void | compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm) |
Plot error when compared against a given time-dependent exact solution. Also returns the norm of the error and that of the exact solution. More... | |
virtual void | compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm) |
Plot error when compared against a given time-dependent exact solution. Also returns the norm of the error and that of the exact solution. More... | |
virtual void | compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm) |
Plot error when compared against a given time-depdendent exact solution. Also returns the norm of the error and that of the exact solution. Version with vectors of norms and errors so that different variables' norms and errors can be returned individually. More... | |
virtual void | compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm) |
Plot error when compared against a given time-depdendent exact solution. Also returns the norm of the error and that of the exact solution. Version with vectors of norms and errors so that different variables' norms and errors can be returned individually. More... | |
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) |
Returns the norm of the error and that of the exact solution. Version with vectors of norms and errors so that different variables' norms and errors can be returned individually. More... | |
bool | is_mesh_distributed () const |
Boolean to indicate if Mesh has been distributed. More... | |
OomphCommunicator * | communicator_pt () const |
Read-only access fct to communicator (Null if mesh is not distributed, i.e. if we don't have mpi). More... | |
void | set_communicator_pt (OomphCommunicator *comm_pt) |
Function to set communicator (mesh is assumed to be distributed if the communicator pointer is non-null). Only defined if mpi is enabled becaus Comm_pt itself is only defined when mpi is enabled. More... | |
void | set_keep_all_elements_as_halos () |
Call this function to keep all the elements as halo elements. More... | |
void | unset_keep_all_elements_as_halos () |
Calll this function to unset the flag that keeps all elements in the mesh as halo elements. More... | |
virtual void | distribute (OomphCommunicator *comm_pt, const Vector< unsigned > &element_domain, Vector< GeneralisedElement * > &deleted_element_pt, DocInfo &doc_info, const bool &report_stats, const bool &overrule_keep_as_halo_element_status) |
Distribute the problem and doc; make this virtual to allow overloading for particular meshes where further work is required. Add to vector of pointers to deleted elements. More... | |
void | distribute (OomphCommunicator *comm_pt, const Vector< unsigned > &element_domain, Vector< GeneralisedElement * > &deleted_element_pt, const bool &report_stats=false) |
Distribute the problem Add to vector of pointers to deleted elements. More... | |
void | prune_halo_elements_and_nodes (Vector< GeneralisedElement * > &deleted_element_pt, const bool &report_stats=false) |
(Irreversibly) prune halo(ed) elements and nodes, usually after another round of refinement, to get rid of excessively wide halo layers. Note that the current mesh will be now regarded as the base mesh and no unrefinement relative to it will be possible once this function has been called. More... | |
void | prune_halo_elements_and_nodes (Vector< GeneralisedElement * > &deleted_element_pt, DocInfo &doc_info, const bool &report_stats) |
(Irreversibly) prune halo(ed) elements and nodes, usually after another round of refinement, to get rid of excessively wide halo layers. Note that the current mesh will be now regarded as the base mesh and no unrefinement relative to it will be possible once this function has been called. More... | |
void | get_efficiency_of_mesh_distribution (double &av_efficiency, double &max_efficiency, double &min_efficiency) |
Get efficiency of mesh distribution: In an ideal distribution without halo overhead, each processor would only hold its own elements. Efficieny per processor = (number of non-halo elements)/ (total number of elements). More... | |
void | doc_mesh_distribution (DocInfo &doc_info) |
Doc the mesh distribution, to be processed with tecplot macros. More... | |
void | check_halo_schemes (DocInfo &doc_info, double &max_permitted_error_for_halo_check) |
Check halo and shared schemes on the mesh. More... | |
void | synchronise_shared_nodes (const bool &report_stats) |
Synchronise shared node lookup schemes to cater for the the case where: (1) a certain node on the current processor is halo with proc p (i.e. its non-halo counterpart lives on processor p) (2) that node is also exists (also as a halo) on another processor (q, say) where its non-halo counter part is also known to be on processor p. However, without calling this function the current processor does not necessarily know that it shares a node with processor q. This information can be required, e.g. when synchronising hanging node schemes over all processors. More... | |
void | get_all_halo_data (std::map< unsigned, double * > &map_of_halo_data) |
Get all the halo data stored in the mesh and add pointers to the data to the map, indexed by global equation number. More... | |
Vector< GeneralisedElement * > | halo_element_pt (const unsigned &p) |
Return vector of halo elements in this Mesh whose non-halo counterpart is held on processor p. More... | |
Vector< GeneralisedElement * > | haloed_element_pt (const unsigned &p) |
Return vector of haloed elements in this Mesh whose haloing counterpart is held on processor p. More... | |
unsigned | nnon_halo_element () |
Total number of non-halo elements in this mesh (Costly call computes result on the fly) More... | |
unsigned | nroot_halo_element () |
Total number of root halo elements in this Mesh. More... | |
unsigned | nroot_halo_element (const unsigned &p) |
Number of root halo elements in this Mesh whose non-halo counterpart is held on processor p. More... | |
Vector< GeneralisedElement * > | root_halo_element_pt (const unsigned &p) |
Vector of pointers to root halo elements in this Mesh whose non-halo counterpart is held on processor p. More... | |
GeneralisedElement *& | root_halo_element_pt (const unsigned &p, const unsigned &e) |
Access fct to the e-th root halo element in this Mesh whose non-halo counterpart is held on processor p. More... | |
void | add_root_halo_element_pt (const unsigned &p, GeneralisedElement *&el_pt) |
Add root halo element whose non-halo counterpart is held on processor p to this Mesh. More... | |
unsigned | nhalo_node () |
Total number of halo nodes in this Mesh. More... | |
unsigned | nhalo_node (const unsigned &p) |
Number of halo nodes in this Mesh whose non-halo counterpart is held on processor p. More... | |
void | add_halo_node_pt (const unsigned &p, Node *&nod_pt) |
Add halo node whose non-halo counterpart is held on processor p to the storage scheme for halo nodes. More... | |
Node * | halo_node_pt (const unsigned &p, const unsigned &j) |
Access fct to the j-th halo node in this Mesh whose non-halo counterpart is held on processor p. More... | |
unsigned | nroot_haloed_element () |
Total number of root haloed elements in this Mesh. More... | |
unsigned | nroot_haloed_element (const unsigned &p) |
Number of root haloed elements in this Mesh whose non-halo counterpart is held on processor p. More... | |
Vector< GeneralisedElement * > | root_haloed_element_pt (const unsigned &p) |
Vector of pointers to root haloed elements in this Mesh whose non-halo counterpart is held on processor p. More... | |
GeneralisedElement *& | root_haloed_element_pt (const unsigned &p, const unsigned &e) |
Access fct to the e-th root haloed element in this Mesh whose non-halo counterpart is held on processor p. More... | |
void | add_root_haloed_element_pt (const unsigned &p, GeneralisedElement *&el_pt) |
Add root haloed element whose non-halo counterpart is held on processor p to the storage scheme for haloed elements. Note: This does not add the element to the storage scheme for elements as it's understood to naturally live on this processor anyway! More... | |
unsigned | nhaloed_node () |
Total number of haloed nodes in this Mesh. More... | |
unsigned | nhaloed_node (const unsigned &p) |
Number of haloed nodes in this Mesh whose haloed counterpart is held on processor p. More... | |
Node * | haloed_node_pt (const unsigned &p, const unsigned &j) |
Access fct to the j-th haloed node in this Mesh whose halo counterpart is held on processor p. More... | |
void | add_haloed_node_pt (const unsigned &p, Node *&nod_pt) |
Add haloed node whose halo counterpart is held on processor p to the storage scheme for haloed nodes. More... | |
void | disable_resizing_of_halo_nodes () |
Function to suppress resizing of halo nodes – optmisation but call it at your own risk! More... | |
void | enable_resizing_of_halo_nodes () |
Function to (re-)enable resizing of halo nodes – this returns things to the default behaviour. More... | |
void | disable_output_of_halo_elements () |
Function to disable halo element output. More... | |
void | enable_output_of_halo_elements () |
Function to enable halo element output. More... | |
unsigned | nshared_node () |
Total number of shared nodes in this Mesh. More... | |
void | doc_shared_nodes () |
Doc shared nodes. More... | |
unsigned | nshared_node (const unsigned &p) |
Number of shared nodes in this Mesh who have a counterpart on processor p. More... | |
Node * | shared_node_pt (const unsigned &p, const unsigned &j) |
Access fct to the j-th shared node in this Mesh who has a counterpart on processor p. More... | |
void | get_shared_node_pt (const unsigned &p, Vector< Node * > &shared_node_pt) |
Get vector of pointers to shared nodes with processor p. Required for faster search in Missing_masters_functions::add_external_haloed_node_helper() and Missing_masters_functions::add_external_haloed_master_node_helper() More... | |
void | add_shared_node_pt (const unsigned &p, Node *&nod_pt) |
Add shared node whose counterpart is held on processor p to the storage scheme for shared nodes. (NB: ensure that this routine is called twice, once for each process) More... | |
void | get_halo_node_stats (double &av_number, unsigned &max_number, unsigned &min_number) |
Get halo node stats for this distributed mesh: Average/max/min number of halo nodes over all processors. Careful: Involves MPI Broadcasts and must therefore be called on all processors! More... | |
void | get_haloed_node_stats (double &av_number, unsigned &max_number, unsigned &min_number) |
Get haloed node stats for this distributed mesh: Average/max/min number of haloed nodes over all processors. Careful: Involves MPI Broadcasts and must therefore be called on all processors! More... | |
void | output_external_halo_elements (std::ostream &outfile, const unsigned &n_plot=5) |
Output all external halo elements. More... | |
void | output_external_halo_elements (const unsigned &p, std::ostream &outfile, const unsigned &n_plot=5) |
Output all external halo elements with processor p. More... | |
void | output_external_haloed_elements (std::ostream &outfile, const unsigned &n_plot=5) |
Output all external haloed elements. More... | |
void | output_external_haloed_elements (const unsigned &p, std::ostream &outfile, const unsigned &n_plot=5) |
Output all external haloed elements with processor p. More... | |
unsigned | nexternal_halo_element () |
Total number of external halo elements in this Mesh. More... | |
unsigned | nexternal_halo_element (const unsigned &p) |
Number of external halo elements in this Mesh whose non-halo counterpart is held on processor p. More... | |
GeneralisedElement *& | external_halo_element_pt (const unsigned &p, const unsigned &e) |
Access fct to the e-th external halo element in this Mesh whose non-halo counterpart is held on processor p. More... | |
void | add_external_halo_element_pt (const unsigned &p, GeneralisedElement *&el_pt) |
Add external halo element whose non-halo counterpart is held on processor p to this Mesh. More... | |
unsigned | nexternal_haloed_element () |
Total number of external haloed elements in this Mesh. More... | |
unsigned | nexternal_haloed_element (const unsigned &p) |
Number of external haloed elements in this Mesh whose non-halo counterpart is held on processor p. More... | |
GeneralisedElement *& | external_haloed_element_pt (const unsigned &p, const unsigned &e) |
Access fct to the e-th external haloed element in this Mesh whose non-halo counterpart is held on processor p. More... | |
unsigned | add_external_haloed_element_pt (const unsigned &p, GeneralisedElement *&el_pt) |
Add external haloed element whose non-halo counterpart is held on processor p to the storage scheme for haloed elements. More... | |
unsigned | nexternal_halo_node () |
Total number of external halo nodes in this Mesh. More... | |
void | get_external_halo_node_pt (Vector< Node * > &external_halo_node_pt) |
Get vector of pointers to all external halo nodes. More... | |
unsigned | nexternal_halo_node (const unsigned &p) |
Number of external halo nodes in this Mesh whose non-halo (external) counterpart is held on processor p. More... | |
void | add_external_halo_node_pt (const unsigned &p, Node *&nod_pt) |
Add external halo node whose non-halo (external) counterpart is held on processor p to the storage scheme for halo nodes. More... | |
Node *& | external_halo_node_pt (const unsigned &p, const unsigned &j) |
Access fct to the j-th external halo node in this Mesh whose non-halo external counterpart is held on processor p. More... | |
Vector< Node * > | external_halo_node_pt (const unsigned &p) |
Access fct to vector of external halo node in this Mesh whose non-halo external counterpart is held on processor p. (read only) More... | |
void | set_external_halo_node_pt (const unsigned &p, const Vector< Node * > &external_halo_node_pt) |
Set vector of external halo node in this Mesh whose non-halo external counterpart is held on processor p. More... | |
void | null_external_halo_node (const unsigned &p, Node *nod_pt) |
Null out specified external halo node (used when deleting duplicates) More... | |
void | remove_null_pointers_from_external_halo_node_storage () |
Consolidate external halo node storage by removing nulled out pointes in external halo and haloed schemes. More... | |
unsigned | nexternal_haloed_node () |
Total number of external haloed nodes in this Mesh. More... | |
unsigned | nexternal_haloed_node (const unsigned &p) |
Number of external haloed nodes in this Mesh whose halo (external) counterpart is held on processor p. More... | |
Node *& | external_haloed_node_pt (const unsigned &p, const unsigned &j) |
Access fct to the j-th external haloed node in this Mesh whose halo external counterpart is held on processor p. More... | |
unsigned | add_external_haloed_node_pt (const unsigned &p, Node *&nod_pt) |
Add external haloed node whose halo (external) counterpart is held on processor p to the storage scheme for haloed nodes. More... | |
Vector< Node * > | external_haloed_node_pt (const unsigned &p) |
Access fct to vector of external haloed node in this Mesh whose halo external counterpart is held on processor p. (read only) More... | |
void | set_external_haloed_node_pt (const unsigned &p, const Vector< Node * > &external_haloed_node_pt) |
Set vector of external haloed node in this Mesh whose halo external counterpart is held on processor p. More... | |
std::set< int > | external_halo_proc () |
Return the set of processors that hold external halo nodes. This is required to avoid having to pass a communicator into the node_update functions for Algebraic-based and MacroElement-based Meshes. More... | |
virtual void | create_shared_boundaries (OomphCommunicator *comm_pt, const Vector< unsigned > &element_domain, const Vector< GeneralisedElement * > &backed_up_el_pt, const Vector< FiniteElement * > &backed_up_f_el_pt, std::map< Data *, std::set< unsigned >> &processors_associated_with_data, const bool &overrule_keep_as_halo_element_status) |
Creates the shared boundaries, only used in unstructured meshes In this case with the "TriangleMesh" class. More... | |
virtual unsigned | try_to_add_root_haloed_element_pt (const unsigned &p, GeneralisedElement *&el_pt) |
virtual unsigned | try_to_add_haloed_node_pt (const unsigned &p, Node *&nod_pt) |
void | delete_all_external_storage () |
Wipe the storage for all externally-based elements. More... | |
Protected Member Functions | |
void | synchronise_hanging_nodes (const unsigned &ncont_interpolated_values) |
Synchronise the hanging nodes if the mesh is distributed. More... | |
virtual void | additional_synchronise_hanging_nodes (const unsigned &ncont_interpolated_values)=0 |
Additional synchronisation of hanging nodes Required for reconcilliation of hanging nodes on the outer edge of the halo layer when using elements with nonuniformly spaced nodes Must be implemented in templated derived class because ELEMENT template parameter is required for the node-creation functions in the Missing_masters_functions namespace. More... | |
void | synchronise_nonhanging_nodes () |
Synchronise the positions of non-hanging nodes that depend on non-existent neighbours (e.g. h-refinement of neighbouring elements with different p-orders where the shared edge is on the outer edge of the halo layer) More... | |
virtual void | split_elements_if_required ()=0 |
Split all the elements in the mesh if required. This template free interface will be overloaded in RefineableMesh<ELEMENT> so that any new elements that are created will be of the correct type. More... | |
virtual void | p_refine_elements_if_required ()=0 |
p-refine all the elements in the mesh if required. This template free interface will be overloaded in RefineableMesh<ELEMENT> so that any temporary copies of the element that are created will be of the correct type. More... | |
void | complete_hanging_nodes (const int &ncont_interpolated_values) |
Complete the hanging node scheme recursively. More... | |
void | complete_hanging_nodes_recursively (Node *&nod_pt, Vector< Node * > &master_nodes, Vector< double > &hang_weights, const int &ival) |
Auxiliary routine for recursive hanging node completion. More... | |
Protected Member Functions inherited from oomph::Mesh | |
void | setup_shared_node_scheme () |
Setup shared node scheme. More... | |
unsigned long | assign_global_eqn_numbers (Vector< double * > &Dof_pt) |
Assign the global equation numbers in the Data stored at the nodes and also internal element Data. Also, build (via push_back) the Vector of pointers to the dofs (variables). More... | |
void | describe_dofs (std::ostream &out, const std::string ¤t_string) const |
Function to describe the dofs of the Mesh. The ostream specifies the output stream to which the description is written; the string stores the currently assembled output that is ultimately written to the output stream by Data::describe_dofs(...); it is typically built up incrementally as we descend through the call hierarchy of this function when called from Problem::describe_dofs(...) More... | |
void | describe_local_dofs (std::ostream &out, const std::string ¤t_string) const |
Function to describe the local dofs of the elements. The ostream specifies the output stream to which the description is written; the string stores the currently assembled output that is ultimately written to the output stream by Data::describe_dofs(...); it is typically built up incrementally as we descend through the call hierarchy of this function when called from Problem::describe_dofs(...) More... | |
void | assign_local_eqn_numbers (const bool &store_local_dof_pt) |
Assign the local equation numbers in all elements If the boolean argument is true then also store pointers to dofs. More... | |
void | convert_to_boundary_node (Node *&node_pt, const Vector< FiniteElement * > &finite_element_pt) |
A function that upgrades an ordinary node to a boundary node We shouldn't ever really use this, but it does make life that bit easier for the lazy mesh writer. The pointer to the node is replaced by a pointer to the new boundary node in all element look-up schemes and in the mesh's Node_pt vector. The new node is also addressed by node_pt on return from the function. More... | |
void | convert_to_boundary_node (Node *&node_pt) |
A function that upgrades an ordinary node to a boundary node. All pointers to the node from the mesh's elements are found. and replaced by pointers to the new boundary node. If the node is present in the mesh's list of nodes, that pointer is also replaced. Finally, the pointer argument node_pt addresses the new node on return from the function. We shouldn't ever really use this, but it does make life that bit easier for the lazy mesh writer. More... | |
Protected Attributes | |
unsigned | Uniform_refinement_level_when_pruned |
Level to which the mesh was uniformly refined when it was pruned. More... | |
unsigned | Max_refinement_level |
Max. permissible refinement level (relative to base mesh) More... | |
unsigned | Min_refinement_level |
Min. permissible refinement level (relative to base mesh) More... | |
unsigned | Max_p_refinement_level |
Max. permissible p-refinement level (relative to base mesh) More... | |
unsigned | Min_p_refinement_level |
Min. permissible p-refinement level (relative to base mesh) More... | |
TreeForest * | Forest_pt |
Forest representation of the mesh. More... | |
Protected Attributes inherited from oomph::RefineableMeshBase | |
ErrorEstimator * | Spatial_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... | |
DocInfo * | Doc_info_pt |
Pointer to DocInfo. More... | |
unsigned | Max_keep_unrefined |
Max. number of elements that can remain unrefined if no other mesh adaptation is required (to avoid mesh-adaptations that would only unrefine a few elements and then force a new solve – this can't be worth our while!) More... | |
unsigned | Nrefinement_overruled |
Number of elements that would like to be refined further but can't because they've reached the max. refinement level. More... | |
Protected Attributes inherited from oomph::Mesh | |
Vector< Vector< Node * > > | Boundary_node_pt |
Vector of Vector of pointers to nodes on the boundaries: Boundary_node_pt(b,n). Note that this is private to force the use of the add_boundary_node() function, which ensures that the reverse look-up schemes for the nodes are set up. More... | |
bool | Lookup_for_elements_next_boundary_is_setup |
Flag to indicate that the lookup schemes for elements that are adjacent to the boundaries has been set up. More... | |
Vector< Vector< FiniteElement * > > | Boundary_element_pt |
Vector of Vector of pointers to elements on the boundaries: Boundary_element_pt(b,e) More... | |
Vector< Vector< int > > | Face_index_at_boundary |
For the e-th finite element on boundary b, this is the index of the face that lies along that boundary. More... | |
std::map< unsigned, Vector< GeneralisedElement * > > | Root_halo_element_pt |
Map of vectors holding the pointers to the root halo elements. More... | |
std::map< unsigned, Vector< GeneralisedElement * > > | Root_haloed_element_pt |
Map of vectors holding the pointers to the root haloed elements. More... | |
std::map< unsigned, Vector< Node * > > | Halo_node_pt |
Map of vectors holding the pointers to the halo nodes. More... | |
std::map< unsigned, Vector< Node * > > | Haloed_node_pt |
Map of vectors holding the pointers to the haloed nodes. More... | |
std::map< unsigned, Vector< Node * > > | Shared_node_pt |
Map of vectors holding the pointers to the shared nodes. These are all the nodes that are on two "neighbouring" processes (the halo(ed) lookup scheme depends upon which processor is in charge. More... | |
OomphCommunicator * | Comm_pt |
Pointer to communicator – set to NULL if mesh is not distributed. More... | |
std::map< unsigned, Vector< GeneralisedElement * > > | External_halo_element_pt |
External halo(ed) elements are created as and when they are needed to act as source elements for the particular process's mesh. The storage is wiped and rebuilt every time the mesh is refined. More... | |
std::map< unsigned, Vector< GeneralisedElement * > > | External_haloed_element_pt |
Map of vectors holding the pointers to the external haloed elements. More... | |
std::map< unsigned, Vector< Node * > > | External_halo_node_pt |
Map of vectors holding the pointers to the external halo nodes. More... | |
std::map< unsigned, Vector< Node * > > | External_haloed_node_pt |
Map of vectors holding the pointers to the external haloed nodes. More... | |
bool | Keep_all_elements_as_halos |
bool to indicate whether to keep all elements in a mesh as halos or not More... | |
bool | Resize_halo_nodes_not_required |
Set this to true to suppress resizing of halo nodes (at your own risk!) More... | |
Vector< Node * > | Node_pt |
Vector of pointers to nodes. More... | |
Vector< GeneralisedElement * > | Element_pt |
Vector of pointers to generalised elements. More... | |
std::vector< bool > | Boundary_coordinate_exists |
Vector of boolean data that indicates whether the boundary coordinates have been set for the boundary. More... | |
Additional Inherited Members | |
Public Types inherited from oomph::Mesh | |
typedef void(FiniteElement::* | SteadyExactSolutionFctPt) (const Vector< double > &x, Vector< double > &soln) |
Typedef for function pointer to function that computes steady exact solution. More... | |
typedef void(FiniteElement::* | UnsteadyExactSolutionFctPt) (const double &time, const Vector< double > &x, Vector< double > &soln) |
Typedef for function pointer to function that computes unsteady exact solution. More... | |
Public Attributes inherited from oomph::Mesh | |
bool | Output_halo_elements |
Bool for output of halo elements. More... | |
Static Public Attributes inherited from oomph::Mesh | |
static Steady< 0 > | Default_TimeStepper |
Default Steady Timestepper, to be used in default arguments to Mesh constructors. More... | |
static bool | Suppress_warning_about_empty_mesh_level_time_stepper_function |
Boolean used to control warning about empty mesh level timestepper function. More... | |
///////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////
Base class for tree-based refineable meshes.
Definition at line 375 of file refineable_mesh.h.
|
inline |
Constructor.
Definition at line 379 of file refineable_mesh.h.
References oomph::RefineableMeshBase::Doc_info_pt, Forest_pt, Max_p_refinement_level, Max_refinement_level, Min_p_refinement_level, Min_refinement_level, oomph::RefineableMeshBase::Nrefined, oomph::RefineableMeshBase::Nunrefined, and Uniform_refinement_level_when_pruned.
|
delete |
Broken copy constructor.
|
inlinevirtual |
|
virtual |
Adapt mesh: Refine elements whose error is lager than err_max and (try to) unrefine those whose error is smaller than err_min.
Do adaptive refinement for mesh.
Implements oomph::RefineableMeshBase.
Reimplemented in oomph::RefineableQuadFromTriangleMesh< ELEMENT >.
Definition at line 307 of file refineable_mesh.cc.
References adapt_mesh(), classify_halo_and_haloed_nodes(), oomph::Mesh::Comm_pt, oomph::Mesh::delete_all_external_storage(), oomph::RefineableElement::deselect_for_refinement(), oomph::RefineableElement::deselect_sons_for_unrefinement(), oomph::DocInfo::disable_doc(), oomph::RefineableMeshBase::doc_info(), oomph::RefineableMeshBase::doc_info_pt(), e, oomph::Mesh::element_pt(), oomph::Tree::father_pt(), oomph::Mesh::halo_element_pt(), oomph::Mesh::haloed_element_pt(), oomph::DocInfo::is_doc_enabled(), oomph::Tree::is_leaf(), oomph::Mesh::is_mesh_distributed(), oomph::OcTreeNames::LDB, oomph::RefineableMeshBase::max_keep_unrefined(), oomph::RefineableMeshBase::max_permitted_error(), max_refinement_level(), oomph::RefineableMeshBase::min_permitted_error(), min_refinement_level(), oomph::Mesh::nelement(), oomph::RefineableMeshBase::Nrefined, oomph::RefineableMeshBase::nrefinement_overruled(), oomph::Tree::nsons(), oomph::RefineableMeshBase::Nunrefined, oomph::Tree::object_pt(), oomph::oomph_info, oomph::RefineableElement::refinement_is_enabled(), oomph::RefineableElement::refinement_level(), oomph::Mesh::reorder_nodes(), oomph::RefineableElement::select_for_refinement(), oomph::RefineableElement::select_sons_for_unrefinement(), oomph::Tree::son_pt(), oomph::Tree::son_type(), and oomph::RefineableElement::tree_pt().
Referenced by unrefine_uniformly().
|
inlinevirtual |
Perform the actual tree-based mesh adaptation. A simple wrapper to call the function without documentation.
Definition at line 517 of file refineable_mesh.h.
References oomph::DocInfo::directory(), oomph::DocInfo::disable_doc(), and oomph::RefineableMeshBase::doc_info().
Referenced by adapt(), refine_base_mesh(), refine_selected_elements(), and refine_uniformly().
|
virtual |
Perform the actual tree-based mesh adaptation, documenting the progress in the directory specified in DocInfo object.
Adapt mesh, which exists in two representations, namely as:
Refinement/derefinement process is documented (in tecplot-able form) if requested.
Procedure:
After adaptation, all nodes (whether new or old) have up-to-date current and previous values.
If refinement process is being documented, the following information is documented:
can be viewed with
can be viewed with
to check the hanging node status.
and can be viewed with
Update the boundary element info – this can be a costly procedure and for this reason the mesh writer might have decided not to set up this scheme. If so, we won't change this and suppress its creation...
Definition at line 926 of file refineable_mesh.cc.
References oomph::Mesh::boundary_element_pt(), oomph::TreeForest::check_all_neighbours(), classify_halo_and_haloed_nodes(), oomph::TreeForest::close_hanging_node_files(), complete_hanging_nodes(), oomph::Tree::deactivate_object(), oomph::Mesh::delete_all_external_storage(), oomph::FiniteElement::dim(), oomph::DocInfo::directory(), oomph::Global_timings::Doc_comprehensive_timings, oomph::RefineableMeshBase::doc_info(), e, oomph::Mesh::Element_pt, oomph::Mesh::finite_element_pt(), Forest_pt, oomph::Node::get_boundaries_pt(), oomph::FiniteElement::get_x(), oomph::SolidFiniteElement::get_x_and_xi(), oomph::Node::hanging_pt(), i, oomph::DocInfo::is_doc_enabled(), oomph::Node::is_hanging(), oomph::Mesh::is_mesh_distributed(), oomph::Node::is_on_boundary(), oomph::Data::is_pinned(), oomph::SolidNode::lagrangian_position(), oomph::FiniteElement::local_coordinate_of_node(), oomph::Mesh::Lookup_for_elements_next_boundary_is_setup, oomph::FiniteElement::macro_elem_pt(), oomph::HangInfo::master_node_pt(), oomph::RefineableMeshBase::max_error(), oomph::RefineableElement::max_integrity_tolerance(), oomph::Tree::merge_sons_if_required(), oomph::Mesh::nboundary(), oomph::Mesh::nboundary_element(), oomph::Node::ndim(), oomph::Mesh::nelement(), oomph::SolidNode::nlagrangian(), oomph::HangInfo::nmaster(), oomph::FiniteElement::nnode(), oomph::Mesh::nnode(), oomph::FiniteElement::node_pt(), oomph::Mesh::node_pt(), oomph::AlgebraicNode::node_update(), oomph::TreeForest::ntree(), oomph::Data::ntstorage(), oomph::DocInfo::number(), oomph::Data::nvalue(), oomph::oomph_info, oomph::TreeForest::open_hanging_node_files(), oomph::Node::position(), oomph::Mesh::prune_dead_nodes(), oomph::Node::raw_value(), oomph::Mesh::reorder_nodes(), s, oomph::Node::set_non_obsolete(), oomph::Node::set_nonhanging(), oomph::Node::set_obsolete(), oomph::Data::set_value(), oomph::Mesh::setup_boundary_element_info(), split_elements_if_required(), oomph::TreeForest::stick_all_tree_nodes_into_vector(), oomph::TreeForest::stick_leaves_into_vector(), t, oomph::TimingHelpers::timer(), oomph::Tree::traverse_all(), oomph::Tree::traverse_all_but_leaves(), oomph::TreeForest::tree_pt(), oomph::SolidFiniteElement::undeformed_macro_elem_pt(), oomph::Node::value(), oomph::Node::x(), and oomph::SolidNode::xi().
|
protectedpure virtual |
Additional synchronisation of hanging nodes Required for reconcilliation of hanging nodes on the outer edge of the halo layer when using elements with nonuniformly spaced nodes Must be implemented in templated derived class because ELEMENT template parameter is required for the node-creation functions in the Missing_masters_functions namespace.
Implemented in oomph::TreeBasedRefineableMesh< ELEMENT >.
Referenced by synchronise_hanging_nodes().
|
inlinevirtual |
Classify the halo and haloed nodes in the mesh.
Reimplemented from oomph::Mesh.
Definition at line 699 of file refineable_mesh.h.
References classify_halo_and_haloed_nodes(), oomph::DocInfo::disable_doc(), and oomph::RefineableMeshBase::doc_info().
|
inlinevirtual |
Classify all halo and haloed information in the mesh (overloaded version from Mesh base class. Calls that one first, then synchronises hanging nodes)
Reimplemented from oomph::Mesh.
Definition at line 635 of file refineable_mesh.h.
References oomph::Mesh::classify_halo_and_haloed_nodes(), oomph::Mesh::Comm_pt, oomph::Global_timings::Doc_comprehensive_timings, oomph::RefineableMeshBase::doc_info(), oomph::Mesh::element_pt(), oomph::Mesh::nelement(), oomph::oomph_info, oomph::Mesh::resize_halo_nodes(), oomph::Mesh::Resize_halo_nodes_not_required, synchronise_hanging_nodes(), synchronise_nonhanging_nodes(), and oomph::TimingHelpers::timer().
Referenced by adapt(), adapt_mesh(), classify_halo_and_haloed_nodes(), p_adapt(), and p_adapt_mesh().
|
protected |
Complete the hanging node scheme recursively.
Complete the hanging node scheme recursively. After the initial markup scheme, hanging nodes can depend on other hanging nodes —> AAAAAAAAARGH! Need to translate this into a scheme where all hanging nodes only depend on non-hanging nodes...
Definition at line 2271 of file refineable_mesh.cc.
References complete_hanging_nodes_recursively(), oomph::Node::hanging_pt(), i, oomph::Node::is_hanging(), oomph::HangInfo::master_weight(), oomph::HangInfo::nmaster(), oomph::Mesh::nnode(), oomph::Mesh::node_pt(), oomph::oomph_info, oomph::Node::set_hanging_pt(), and oomph::HangInfo::set_master_node_pt().
Referenced by adapt_mesh(), and p_adapt_mesh().
|
protected |
Auxiliary routine for recursive hanging node completion.
Given a node, return a vector of pointers to master nodes and a vector of the associated weights. This is done recursively, so if a node is not hanging, the node is regarded as its own master node which has weight 1.0.
Definition at line 2214 of file refineable_mesh.cc.
References oomph::Node::hanging_pt(), i, oomph::Node::is_hanging(), oomph::HangInfo::master_node_pt(), oomph::HangInfo::master_weight(), and oomph::HangInfo::nmaster().
Referenced by complete_hanging_nodes().
|
inlinevirtual |
Doc the targets for mesh adaptation.
Reimplemented from oomph::RefineableMeshBase.
Definition at line 467 of file refineable_mesh.h.
References oomph::RefineableMeshBase::Max_keep_unrefined, Max_p_refinement_level, oomph::RefineableMeshBase::Max_permitted_error, Max_refinement_level, Min_p_refinement_level, oomph::RefineableMeshBase::Min_permitted_error, and Min_refinement_level.
|
virtual |
Dump refinement pattern to allow for rebuild.
Definition at line 218 of file refineable_mesh.cc.
References get_refinement_pattern(), and i.
|
inline |
Return pointer to the Forest represenation of the mesh.
Definition at line 460 of file refineable_mesh.h.
References Forest_pt.
Referenced by get_elements_at_refinement_level(), get_refinement_pattern(), p_adapt_mesh(), and refine_as_in_reference_mesh().
|
virtual |
Extract the elements at a particular refinement level in the refinement pattern - used in Mesh::redistribute or whatever it's going to be called (RefineableMeshBase::reduce_halo_layers or something)
Extract the elements at a particular refinement level in the refinement pattern (used in Mesh::redistribute or whatever it's going to be called (RefineableMeshBase::reduce_halo_layers or something)
Definition at line 113 of file refineable_mesh.cc.
References e, forest_pt(), and oomph::TreeForest::stick_all_tree_nodes_into_vector().
Referenced by oomph::Mesh::prune_halo_elements_and_nodes().
|
virtual |
Get max/min refinement levels in mesh.
Get max/min refinement level.
Definition at line 830 of file refineable_mesh.cc.
References e, oomph::Mesh::element_pt(), max_refinement_level(), min_refinement_level(), and oomph::Mesh::nelement().
Referenced by oomph::Problem::load_balance(), oomph::Mesh::prune_halo_elements_and_nodes(), oomph::Problem::read(), refine_as_in_reference_mesh(), refine_base_mesh(), oomph::RefineableBrickMesh< ELEMENT >::setup_octree_forest(), and oomph::RefineableQuadMesh< ELEMENT >::setup_quadtree_forest().
|
virtual |
Extract refinement pattern: Consider the hypothetical mesh obtained by truncating the refinement of the current mesh to a given level (where level=0
is the un-refined base mesh). To advance to the next refinement level, we need to refine (split) the to_be_refined
[level].size() elements identified by the element numbers contained in vector
to_be_refined[level][...].
Get refinement pattern of mesh: Consider the hypothetical mesh obtained by truncating the refinement of the current mesh to a given level (where level=0
is the un-refined base mesh). To advance to the next refinement level, we need to refine (split) the to_be_refined
[level].size() elements identified by the element numbers contained in vector
to_be_refined[level][...].
Definition at line 50 of file refineable_mesh.cc.
References e, forest_pt(), and oomph::TreeForest::stick_all_tree_nodes_into_vector().
Referenced by dump_refinement(), oomph::Mesh::prune_halo_elements_and_nodes(), refine_base_mesh_as_in_reference_mesh(), and refine_base_mesh_as_in_reference_mesh_minus_one().
|
inline |
Access fct for max. permissible p-refinement level (relative to base mesh)
Definition at line 499 of file refineable_mesh.h.
References Max_p_refinement_level.
Referenced by p_adapt().
|
inline |
Access fct for max. permissible refinement level (relative to base mesh)
Definition at line 486 of file refineable_mesh.h.
References Max_refinement_level.
Referenced by adapt(), and get_refinement_levels().
|
inline |
Access fct for min. permissible p-refinement level (relative to base mesh)
Definition at line 506 of file refineable_mesh.h.
References Min_p_refinement_level.
|
inline |
Access fct for min. permissible refinement level (relative to base mesh)
Definition at line 492 of file refineable_mesh.h.
References Min_refinement_level.
Referenced by adapt(), and get_refinement_levels().
|
delete |
Broken assignment operator.
|
virtual |
p-adapt mesh: Refine elements whose error is lager than err_max and (try to) unrefine those whose error is smaller than err_min
///////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////
Do adaptive p-refinement for mesh.
/ ... otherwise mark it as having been over-ruled
Reimplemented from oomph::RefineableMeshBase.
Definition at line 3951 of file refineable_mesh.cc.
References classify_halo_and_haloed_nodes(), oomph::Mesh::Comm_pt, oomph::Mesh::delete_all_external_storage(), oomph::PRefineableElement::deselect_for_p_refinement(), oomph::PRefineableElement::deselect_for_p_unrefinement(), oomph::DocInfo::disable_doc(), oomph::RefineableMeshBase::doc_info(), oomph::RefineableMeshBase::doc_info_pt(), e, oomph::Mesh::element_pt(), oomph::Mesh::halo_element_pt(), oomph::Mesh::haloed_element_pt(), oomph::PRefineableElement::initial_p_order(), oomph::DocInfo::is_doc_enabled(), oomph::Mesh::is_mesh_distributed(), oomph::RefineableMeshBase::max_keep_unrefined(), max_p_refinement_level(), oomph::RefineableMeshBase::max_permitted_error(), oomph::RefineableMeshBase::min_permitted_error(), oomph::Mesh::nelement(), oomph::RefineableMeshBase::Nrefined, oomph::RefineableMeshBase::nrefinement_overruled(), oomph::RefineableMeshBase::Nunrefined, oomph::oomph_info, p_adapt_mesh(), oomph::PRefineableElement::p_order(), oomph::PRefineableElement::p_refinement_is_enabled(), oomph::Mesh::reorder_nodes(), oomph::PRefineableElement::select_for_p_refinement(), and oomph::PRefineableElement::select_for_p_unrefinement().
|
inline |
Perform the actual tree-based mesh p-adaptation. A simple wrapper to call the function without documentation.
Definition at line 533 of file refineable_mesh.h.
References oomph::DocInfo::directory(), oomph::DocInfo::disable_doc(), and oomph::RefineableMeshBase::doc_info().
Referenced by p_adapt(), p_refine_selected_elements(), p_refine_uniformly(), and p_unrefine_uniformly().
void oomph::TreeBasedRefineableMeshBase::p_adapt_mesh | ( | DocInfo & | doc_info | ) |
Perform the actual tree-based mesh p-adaptation, documenting the progress in the directory specified in DocInfo object.
p-adapt mesh, which exists in two representations, namely as:
p-refinement/derefinement process is documented (in tecplot-able form) if requested.
Procedure:
After adaptation, all nodes (whether new or old) have up-to-date current and previous values.
If refinement process is being documented, the following information is documented:
can be viewed with
can be viewed with
to check the hanging node status.
and can be viewed with
Update the boundary element info – this can be a costly procedure and for this reason the mesh writer might have decided not to set up this scheme. If so, we won't change this and suppress its creation...
/BENFLAG: Check that all the nodes belong to their update elements
Definition at line 4433 of file refineable_mesh.cc.
References oomph::Mesh::boundary_element_pt(), oomph::TreeForest::check_all_neighbours(), classify_halo_and_haloed_nodes(), oomph::TreeForest::close_hanging_node_files(), complete_hanging_nodes(), oomph::Tree::deactivate_object(), oomph::Mesh::delete_all_external_storage(), oomph::PRefineableElement::deselect_for_p_refinement(), oomph::PRefineableElement::deselect_for_p_unrefinement(), oomph::FiniteElement::dim(), oomph::DocInfo::directory(), oomph::Global_timings::Doc_comprehensive_timings, oomph::RefineableMeshBase::doc_info(), e, oomph::Mesh::element_pt(), oomph::Mesh::finite_element_pt(), forest_pt(), oomph::MacroElementNodeUpdateElementBase::geom_object_pt(), oomph::Node::get_boundaries_pt(), oomph::FiniteElement::get_x(), oomph::SolidFiniteElement::get_x_and_xi(), oomph::Node::hanging_pt(), i, oomph::DocInfo::is_doc_enabled(), oomph::Node::is_hanging(), oomph::Mesh::is_mesh_distributed(), oomph::Node::is_on_boundary(), oomph::Data::is_pinned(), oomph::SolidNode::lagrangian_position(), oomph::FiniteElement::local_coordinate_of_node(), oomph::Mesh::Lookup_for_elements_next_boundary_is_setup, oomph::FiniteElement::macro_elem_pt(), oomph::HangInfo::master_node_pt(), oomph::RefineableMeshBase::max_error(), oomph::RefineableElement::max_integrity_tolerance(), oomph::Mesh::nboundary(), oomph::Mesh::nboundary_element(), oomph::Node::ndim(), oomph::Mesh::nelement(), oomph::SolidNode::nlagrangian(), oomph::HangInfo::nmaster(), oomph::FiniteElement::nnode(), oomph::Mesh::nnode(), oomph::FiniteElement::node_pt(), oomph::Mesh::node_pt(), oomph::AlgebraicNode::node_update(), oomph::TreeForest::ntree(), oomph::Data::ntstorage(), oomph::DocInfo::number(), oomph::Data::nvalue(), oomph::oomph_info, oomph::TreeForest::open_hanging_node_files(), p_refine_elements_if_required(), oomph::Node::position(), oomph::Mesh::prune_dead_nodes(), oomph::Node::raw_value(), oomph::Mesh::reorder_nodes(), s, oomph::MacroElementNodeUpdateElementBase::set_node_update_info(), oomph::Node::set_non_obsolete(), oomph::Node::set_nonhanging(), oomph::Node::set_obsolete(), oomph::Data::set_value(), oomph::Mesh::setup_boundary_element_info(), oomph::TreeForest::stick_all_tree_nodes_into_vector(), oomph::TreeForest::stick_leaves_into_vector(), t, oomph::TimingHelpers::timer(), oomph::Tree::traverse_all(), oomph::TreeForest::tree_pt(), oomph::SolidFiniteElement::undeformed_macro_elem_pt(), oomph::Node::value(), oomph::Node::x(), and oomph::SolidNode::xi().
|
protectedpure virtual |
p-refine all the elements in the mesh if required. This template free interface will be overloaded in RefineableMesh<ELEMENT> so that any temporary copies of the element that are created will be of the correct type.
Implemented in oomph::TreeBasedRefineableMesh< ELEMENT >.
Referenced by p_adapt_mesh().
void oomph::TreeBasedRefineableMeshBase::p_refine_selected_elements | ( | const Vector< PRefineableElement * > & | elements_to_be_refined_pt | ) |
p-refine mesh by refining the elements identified by their pointers.
Definition at line 5304 of file refineable_mesh.cc.
References e, oomph::Mesh::is_mesh_distributed(), and p_adapt_mesh().
void oomph::TreeBasedRefineableMeshBase::p_refine_selected_elements | ( | const Vector< unsigned > & | elements_to_be_refined | ) |
p-refine mesh by refining the elements identified by their numbers.
Definition at line 5264 of file refineable_mesh.cc.
References e, oomph::Mesh::element_pt(), oomph::Mesh::is_mesh_distributed(), p_adapt_mesh(), and oomph::PRefineableElement::select_for_p_refinement().
|
inlinevirtual |
p-refine mesh uniformly
Reimplemented from oomph::RefineableMeshBase.
Definition at line 443 of file refineable_mesh.h.
References oomph::RefineableMeshBase::p_refine_uniformly().
|
virtual |
p-refine mesh uniformly and doc process
p-refine mesh uniformly
Reimplemented from oomph::RefineableMeshBase.
Definition at line 1790 of file refineable_mesh.cc.
References oomph::RefineableMeshBase::doc_info(), e, oomph::Mesh::element_pt(), oomph::Mesh::nelement(), p_adapt_mesh(), and oomph::PRefineableElement::select_for_p_refinement().
void oomph::TreeBasedRefineableMeshBase::p_unrefine_uniformly | ( | DocInfo & | doc_info | ) |
p-unrefine mesh uniformly
p-unrefine mesh uniformly Unlike in h-refinement, we can simply p-unrefine each element in the mesh
Definition at line 5240 of file refineable_mesh.cc.
References oomph::RefineableMeshBase::doc_info(), e, oomph::Mesh::element_pt(), oomph::Mesh::nelement(), p_adapt_mesh(), and oomph::PRefineableElement::select_for_p_unrefinement().
|
virtual |
Read refinement pattern to allow for rebuild.
Definition at line 251 of file refineable_mesh.cc.
References i, and oomph::Global_string_for_annotation::string().
Referenced by refine().
|
virtual |
Refine mesh according to refinement pattern in restart file.
Refine base mesh according to refinement pattern in restart file.
Definition at line 201 of file refineable_mesh.cc.
References read_refinement(), and refine_base_mesh().
|
virtual |
Refine mesh once so that its topology etc becomes that of the (finer!) reference mesh – if possible! Useful for meshes in multigrid hierarchies. If the meshes are too different and the conversion cannot be performed, the code dies (provided PARANOID is enabled).
Definition at line 1947 of file refineable_mesh.cc.
References e, oomph::Tree::father_pt(), oomph::Mesh::finite_element_pt(), forest_pt(), get_refinement_levels(), i, oomph::Tree::is_leaf(), oomph::Node::ndim(), oomph::Mesh::nelement(), oomph::FiniteElement::nnode(), oomph::FiniteElement::node_pt(), oomph::Mesh::node_pt(), oomph::Tree::nsons(), oomph::oomph_info, oomph::Mesh::output(), oomph::pause(), refine_selected_elements(), oomph::Tree::son_pt(), oomph::TreeForest::stick_leaves_into_vector(), and oomph::Node::x().
void oomph::TreeBasedRefineableMeshBase::refine_base_mesh | ( | Vector< Vector< unsigned >> & | to_be_refined | ) |
Refine base mesh according to specified refinement pattern.
Refine original, unrefined mesh according to specified refinement pattern (relative to original, unrefined mesh).
Definition at line 137 of file refineable_mesh.cc.
References adapt_mesh(), oomph::Mesh::Comm_pt, oomph::Mesh::element_pt(), get_refinement_levels(), i, oomph::Mesh::is_mesh_distributed(), and unrefine_uniformly().
Referenced by refine(), refine_base_mesh_as_in_reference_mesh(), refine_base_mesh_as_in_reference_mesh_minus_one(), and oomph::Problem::refine_distributed_base_mesh().
|
virtual |
Refine base mesh to same degree as reference mesh (relative to original unrefined mesh).
Refine to same degree as the reference mesh.
Definition at line 1885 of file refineable_mesh.cc.
References get_refinement_pattern(), and refine_base_mesh().
|
virtual |
Refine base mesh to same degree as reference mesh minus one level of refinement (relative to original unrefined mesh). Useful function for multigrid solvers; allows the easy copy of a mesh to the level of refinement just below the current one.
Refine to same degree as the reference mesh minus one. Useful function for multigrid solvers; allows the easy copy of a mesh to the level of refinement just below the current one. Returns a boolean variable which indicates if the reference mesh has not been refined at all.
Definition at line 1906 of file refineable_mesh.cc.
References get_refinement_pattern(), and refine_base_mesh().
Referenced by oomph::MGSolver< DIM >::setup_mg_hierarchy(), and oomph::HelmholtzMGPreconditioner< DIM >::setup_mg_hierarchy().
|
virtual |
Refine mesh by splitting the elements identified by their pointers.
Definition at line 1852 of file refineable_mesh.cc.
References adapt_mesh(), e, and oomph::Mesh::is_mesh_distributed().
|
virtual |
Refine mesh by splitting the elements identified by their numbers.
Definition at line 1816 of file refineable_mesh.cc.
References adapt_mesh(), e, oomph::Mesh::element_pt(), and oomph::Mesh::is_mesh_distributed().
Referenced by refine_as_in_reference_mesh().
|
inlinevirtual |
Refine mesh uniformly.
Reimplemented from oomph::RefineableMeshBase.
Reimplemented in oomph::RefineableQuadFromTriangleMesh< ELEMENT >.
Definition at line 434 of file refineable_mesh.h.
References oomph::RefineableMeshBase::refine_uniformly().
|
virtual |
Refine mesh uniformly and doc process.
Refine mesh uniformly.
Implements oomph::RefineableMeshBase.
Reimplemented in oomph::RefineableQuadFromTriangleMesh< ELEMENT >.
Definition at line 1772 of file refineable_mesh.cc.
References adapt_mesh(), oomph::RefineableMeshBase::doc_info(), e, oomph::Mesh::element_pt(), and oomph::Mesh::nelement().
Referenced by oomph::Problem::load_balance().
|
pure virtual |
Set up the tree forest associated with the Mesh (if any)
Implemented in oomph::RefineableQuadMesh< ELEMENT >, oomph::RefineableLineMesh< ELEMENT >, and oomph::RefineableBrickMesh< ELEMENT >.
Referenced by oomph::Mesh::distribute(), and oomph::Mesh::prune_halo_elements_and_nodes().
|
protectedpure virtual |
Split all the elements in the mesh if required. This template free interface will be overloaded in RefineableMesh<ELEMENT> so that any new elements that are created will be of the correct type.
Implemented in oomph::TreeBasedRefineableMesh< ELEMENT >.
Referenced by adapt_mesh().
|
protected |
Synchronise the hanging nodes if the mesh is distributed.
Deal with nodes that are hanging on one process but not another (i.e. the hanging status of the haloed and halo layers disagrees)
Definition at line 2375 of file refineable_mesh.cc.
References additional_synchronise_hanging_nodes(), oomph::Mesh::Comm_pt, oomph::Global_timings::Doc_comprehensive_timings, oomph::Mesh::halo_node_pt(), oomph::Mesh::haloed_node_pt(), oomph::TreeBasedRefineableMeshBase::HangHelperStruct::Hang_pt, oomph::Node::hanging_pt(), i, oomph::TreeBasedRefineableMeshBase::HangHelperStruct::icont, oomph::Node::is_hanging(), oomph::TreeBasedRefineableMeshBase::HangHelperStruct::Master_node_index, oomph::HangInfo::master_node_pt(), oomph::HangInfo::master_weight(), oomph::Mesh::nhalo_node(), oomph::Mesh::nhaloed_node(), oomph::HangInfo::nmaster(), oomph::TreeBasedRefineableMeshBase::HangHelperStruct::Node_pt, oomph::Data::non_halo_proc_ID(), oomph::oomph_info, oomph::TreeBasedRefineableMeshBase::HangHelperStruct::Sending_processor, oomph::Node::set_hanging_pt(), oomph::HangInfo::set_master_node_pt(), oomph::TreeBasedRefineableMeshBase::HangHelperStruct::Shared_node_id_on_sending_processor, oomph::TreeBasedRefineableMeshBase::HangHelperStruct::Shared_node_proc, oomph::Mesh::Shared_node_pt, oomph::Mesh::shared_node_pt(), oomph::TimingHelpers::timer(), and oomph::TreeBasedRefineableMeshBase::HangHelperStruct::Weight.
Referenced by classify_halo_and_haloed_nodes().
|
protected |
Synchronise the positions of non-hanging nodes that depend on non-existent neighbours (e.g. h-refinement of neighbouring elements with different p-orders where the shared edge is on the outer edge of the halo layer)
Definition at line 3650 of file refineable_mesh.cc.
References oomph::Mesh::Comm_pt, oomph::FiniteElement::dim(), oomph::Global_timings::Doc_comprehensive_timings, e, oomph::FiniteElement::get_x(), oomph::Mesh::halo_element_pt(), oomph::Mesh::haloed_element_pt(), oomph::Node::is_hanging(), oomph::FiniteElement::local_coordinate_of_node(), oomph::FiniteElement::nnode(), oomph::FiniteElement::node_pt(), oomph::Data::ntstorage(), oomph::oomph_info, s, t, oomph::TimingHelpers::timer(), and oomph::Node::x().
Referenced by classify_halo_and_haloed_nodes().
|
inline |
Level to which the mesh was uniformly refined when it was pruned.
Definition at line 625 of file refineable_mesh.h.
References Uniform_refinement_level_when_pruned.
|
inline |
Level to which the mesh was uniformly refined when it was pruned (const version)
Definition at line 618 of file refineable_mesh.h.
References Uniform_refinement_level_when_pruned.
Referenced by oomph::Problem::load_balance(), and oomph::Mesh::prune_halo_elements_and_nodes().
|
virtual |
Unrefine mesh uniformly: Return 0 for success, 1 for failure (if unrefinement has reached the coarsest permitted level)
Unrefine mesh uniformly. Return 0 for success, 1 for failure (if unrefinement has reached the coarsest permitted level)
Implements oomph::RefineableMeshBase.
Definition at line 2150 of file refineable_mesh.cc.
References adapt(), e, oomph::RefineableMeshBase::max_keep_unrefined(), oomph::RefineableMeshBase::max_permitted_error(), oomph::RefineableMeshBase::min_permitted_error(), and oomph::Mesh::nelement().
Referenced by refine_base_mesh().
|
protected |
Forest representation of the mesh.
Definition at line 768 of file refineable_mesh.h.
Referenced by adapt_mesh(), oomph::ElasticRefineableQuarterPipeMesh< ELEMENT >::ElasticRefineableQuarterPipeMesh(), forest_pt(), oomph::TreeBasedRefineableMesh< ELEMENT >::p_refine_elements_if_required(), oomph::RefineableEighthSphereMesh< ELEMENT >::RefineableEighthSphereMesh(), oomph::RefineableFullCircleMesh< ELEMENT >::RefineableFullCircleMesh(), oomph::RefineableQuarterPipeMesh< ELEMENT >::RefineableQuarterPipeMesh(), oomph::RefineableQuarterTubeMesh< ELEMENT >::RefineableQuarterTubeMesh(), oomph::RefineableTubeMesh< ELEMENT >::RefineableTubeMesh(), oomph::RefineableLineMesh< ELEMENT >::setup_binary_tree_forest(), oomph::RefineableBrickMesh< ELEMENT >::setup_octree_forest(), oomph::RefineableQuadMesh< ELEMENT >::setup_quadtree_forest(), oomph::TreeBasedRefineableMesh< ELEMENT >::split_elements_if_required(), TreeBasedRefineableMeshBase(), oomph::RefineableQuarterPipeMesh< ELEMENT >::~RefineableQuarterPipeMesh(), and ~TreeBasedRefineableMeshBase().
|
protected |
Max. permissible p-refinement level (relative to base mesh)
Definition at line 762 of file refineable_mesh.h.
Referenced by doc_adaptivity_targets(), max_p_refinement_level(), and TreeBasedRefineableMeshBase().
|
protected |
Max. permissible refinement level (relative to base mesh)
Definition at line 756 of file refineable_mesh.h.
Referenced by doc_adaptivity_targets(), max_refinement_level(), and TreeBasedRefineableMeshBase().
|
protected |
Min. permissible p-refinement level (relative to base mesh)
Definition at line 765 of file refineable_mesh.h.
Referenced by doc_adaptivity_targets(), min_p_refinement_level(), and TreeBasedRefineableMeshBase().
|
protected |
Min. permissible refinement level (relative to base mesh)
Definition at line 759 of file refineable_mesh.h.
Referenced by doc_adaptivity_targets(), min_refinement_level(), and TreeBasedRefineableMeshBase().
|
protected |
Level to which the mesh was uniformly refined when it was pruned.
Definition at line 753 of file refineable_mesh.h.
Referenced by TreeBasedRefineableMeshBase(), and uniform_refinement_level_when_pruned().