libMesh
Public Member Functions | Protected Attributes | Private Attributes | List of all members
libMesh::VariationalMeshSmoother Class Reference

This is an implementation of Larisa Branets' smoothing algorithms. More...

#include <mesh_smoother_vsmoother.h>

Inheritance diagram for libMesh::VariationalMeshSmoother:
[legend]

Public Member Functions

 VariationalMeshSmoother (UnstructuredMesh &mesh, Real dilation_weight=0.5, const bool preserve_subdomain_boundaries=true)
 Simple constructor to use for smoothing purposes. More...
 
virtual ~VariationalMeshSmoother ()=default
 Destructor. More...
 
virtual void smooth () override
 Redefinition of the smooth function from the base class. More...
 
void smooth (unsigned int n_iterations)
 The actual smoothing function, gets called whenever the user specifies an actual number of smoothing iterations. More...
 

Protected Attributes

UnstructuredMesh_mesh
 

Private Attributes

const Real _dilation_weight
 Smoother control variables. More...
 
const bool _preserve_subdomain_boundaries
 Whether subdomain boundaries are subject to change via smoothing. More...
 

Detailed Description

This is an implementation of Larisa Branets' smoothing algorithms.

The initial implementation was done by her, the adaptation to libmesh was completed by Derek Gaston. The code was heavily refactored into something more closely resembling C++ by John Peterson in 2014.

Here are the relevant publications: 1) L. Branets, G. Carey, "Extension of a mesh quality metric for elements with a curved boundary edge or surface", Journal of Computing and Information Science in Engineering, vol. 5(4), pp.302-308, 2005.

2) L. Branets, G. Carey, "A local cell quality metric and variational grid smoothing algorithm", Engineering with Computers, vol. 21, pp.19-28, 2005.

3) L. Branets, "A variational grid optimization algorithm based on a local cell quality metric", Ph.D. thesis, The University of Texas at Austin, 2005.

Author
Derek R. Gaston
Date
2006

Definition at line 65 of file mesh_smoother_vsmoother.h.

Constructor & Destructor Documentation

◆ VariationalMeshSmoother()

libMesh::VariationalMeshSmoother::VariationalMeshSmoother ( UnstructuredMesh mesh,
Real  dilation_weight = 0.5,
const bool  preserve_subdomain_boundaries = true 
)

Simple constructor to use for smoothing purposes.

Definition at line 52 of file mesh_smoother_vsmoother.C.

54  :
56  _dilation_weight(dilation_weight),
57  _preserve_subdomain_boundaries(preserve_subdomain_boundaries)
58 {}
const Real _dilation_weight
Smoother control variables.
MeshBase & mesh
MeshSmoother(UnstructuredMesh &mesh)
Constructor.
Definition: mesh_smoother.h:46
const bool _preserve_subdomain_boundaries
Whether subdomain boundaries are subject to change via smoothing.

◆ ~VariationalMeshSmoother()

virtual libMesh::VariationalMeshSmoother::~VariationalMeshSmoother ( )
virtualdefault

Destructor.

Member Function Documentation

◆ smooth() [1/2]

virtual void libMesh::VariationalMeshSmoother::smooth ( )
inlineoverridevirtual

Redefinition of the smooth function from the base class.

All this does is call the smooth function in this class which takes an int, using a default value of 1.

Implements libMesh::MeshSmoother.

Definition at line 87 of file mesh_smoother_vsmoother.h.

References smooth().

Referenced by main(), and smooth().

87 {this->smooth(1); }
virtual void smooth() override
Redefinition of the smooth function from the base class.

◆ smooth() [2/2]

void libMesh::VariationalMeshSmoother::smooth ( unsigned int  n_iterations)

The actual smoothing function, gets called whenever the user specifies an actual number of smoothing iterations.

Definition at line 61 of file mesh_smoother_vsmoother.C.

References _dilation_weight, libMesh::MeshSmoother::_mesh, _preserve_subdomain_boundaries, libMesh::EquationSystems::add_system(), libMesh::System::attach_constraint_object(), libMesh::ParallelObject::comm(), libMesh::MeshBase::elem_dimensions(), libMesh::VariationalSmootherSystem::get_dilation_weight(), libMesh::EquationSystems::init(), libMesh::make_range(), libMesh::MeshBase::mesh_dimension(), libMesh::MeshBase::node_ref(), libMesh::FEMSystem::solve(), libMesh::Parallel::sync_dofobject_data_by_id(), libMesh::DifferentiableSystem::time_solver, and libMesh::TOLERANCE.

62 {
63  // Check for multiple dimensions
64  if (_mesh.elem_dimensions().size() > 1)
65  libmesh_not_implemented_msg("Meshes containing elements of differing dimension are not yet supported.");
66 
67  /*
68  Ideally we'd want to update _mesh directly even if it already has an
69  EquationSystems attached and doing work on it ... and that's thwarted by the
70  problems:
71 
72  1. We can't easily tell if there's already an EquationSystems attached.
73  2. We can't attach a second EquationSystems safely (if it already has a system 0)
74  because the DoF indexing will need to be overwritten by our system.
75  3. The destructor of es won't even clean up after itself (we generally expect
76  a mesh to go unused after its EquationSystems is destroyed), much less know
77  how to restore anything from a previous EquationSystems.
78 
79  To avoid these issues, we'll just construct a new DistributedMesh mesh_copy
80  from _mesh, then do the solve on mesh_copy, then copy its node locations back
81  to _mesh after the solve is done. That'll be slightly less memory-efficient
82  and somewhat more CPU-efficient in the case where _mesh is serial, though
83  it'll be significantly less memory-efficient when _mesh is already distributed,
84  but either way the robustness is probably worth it.
85  */
86 
87  // Create a new mesh, EquationSystems, and System
88  DistributedMesh mesh_copy(_mesh);
89  EquationSystems es(mesh_copy);
90  VariationalSmootherSystem & sys = es.add_system<VariationalSmootherSystem>("variational_smoother_system");
91 
92  // Set this to something > 0 to add more quadrature points than the default
93  // rule that integrates order 2 * fe_order + 1 polynomials exactly.
94  // Using higher quadrature orders has not had a significant effect on observed solutions.
95  //sys.extra_quadrature_order = 0;
96 
97  // Uncomment these to debug
98  //sys.print_element_solutions=true;
99  //sys.print_element_residuals=true;
100  //sys.print_element_jacobians=true;
101 
102  // Add boundary node and hanging node constraints
103  VariationalSmootherConstraint constraint(sys, _preserve_subdomain_boundaries);
104  sys.attach_constraint_object(constraint);
105 
106  // Set system parameters
107  sys.get_dilation_weight() = _dilation_weight;
108 
109  // Set up solver
110  sys.time_solver =
111  std::make_unique<SteadySolver>(sys);
112 
113  // Uncomment this line and use -snes_test_jacobian and -snes_test_jacobian_view
114  // flags to compare the hand-coded jacobian in VariationalSmootherSystem
115  // to finite difference jacobians.
116  //sys.time_solver->diff_solver() = std::make_unique<PetscDiffSolver>(sys);
117 
118  es.init();
119 
120  // More debugging options
121  //DiffSolver & solver = *(sys.time_solver->diff_solver().get());
122  //solver.quiet = false;
123  //solver.verbose = true;
124 
125  sys.time_solver->diff_solver()->relative_residual_tolerance = TOLERANCE*TOLERANCE;
126 
127  sys.solve();
128 
129  // Update _mesh from mesh_copy
130  for (auto * node_copy : mesh_copy.local_node_ptr_range())
131  {
132  auto & node = _mesh.node_ref(node_copy->id());
133  for (const auto d : make_range(mesh_copy.mesh_dimension()))
134  node(d) = (*node_copy)(d);
135  }
136 
137  SyncNodalPositions sync_object(_mesh);
138  Parallel::sync_dofobject_data_by_id (_mesh.comm(), _mesh.nodes_begin(), _mesh.nodes_end(), sync_object);
139 }
const Real _dilation_weight
Smoother control variables.
static constexpr Real TOLERANCE
const Parallel::Communicator & comm() const
UnstructuredMesh & _mesh
Definition: mesh_smoother.h:61
const std::set< unsigned char > & elem_dimensions() const
Definition: mesh_base.h:282
void sync_dofobject_data_by_id(const Communicator &comm, const Iterator &range_begin, const Iterator &range_end, SyncFunctor &sync)
Request data about a range of ghost dofobjects uniquely identified by their id.
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:140
virtual const Node & node_ref(const dof_id_type i) const
Definition: mesh_base.h:596
const bool _preserve_subdomain_boundaries
Whether subdomain boundaries are subject to change via smoothing.

Member Data Documentation

◆ _dilation_weight

const Real libMesh::VariationalMeshSmoother::_dilation_weight
private

Smoother control variables.

Definition at line 102 of file mesh_smoother_vsmoother.h.

Referenced by smooth().

◆ _mesh

UnstructuredMesh& libMesh::MeshSmoother::_mesh
protectedinherited

◆ _preserve_subdomain_boundaries

const bool libMesh::VariationalMeshSmoother::_preserve_subdomain_boundaries
private

Whether subdomain boundaries are subject to change via smoothing.

Definition at line 105 of file mesh_smoother_vsmoother.h.

Referenced by smooth().


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