libMesh
mesh_smoother_vsmoother.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 #include "libmesh/libmesh_config.h"
19 #if defined(LIBMESH_ENABLE_VSMOOTHER)
20 
21 // Local includes
22 #include "libmesh/mesh_smoother_vsmoother.h"
23 #include "libmesh/mesh_tools.h"
24 #include "libmesh/elem.h"
25 #include "libmesh/unstructured_mesh.h"
26 #include "libmesh/utility.h"
27 #include "libmesh/boundary_info.h"
28 #include "libmesh/equation_systems.h"
29 #include "libmesh/distributed_mesh.h"
30 #include "libmesh/steady_solver.h"
31 #include "libmesh/diff_solver.h"
32 #include "libmesh/variational_smoother_constraint.h"
33 #include "libmesh/parallel_ghost_sync.h"
34 
35 // C++ includes
36 #include <time.h> // for clock_t, clock()
37 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC
38 #include <cmath>
39 #include <iomanip>
40 #include <limits>
41 
42 namespace libMesh
43 {
44 
45 // Optimization at -O2 or greater seem to break Intel's icc. So if we are
46 // being compiled with icc let's dumb-down the optimizations for this file
47 #ifdef __INTEL_COMPILER
48 # pragma optimize ( "", off )
49 #endif
50 
51 // Member functions for the Variational Smoother
53  Real dilation_weight,
54  const bool preserve_subdomain_boundaries) :
56  _dilation_weight(dilation_weight),
57  _preserve_subdomain_boundaries(preserve_subdomain_boundaries)
58 {}
59 
60 
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
104  sys.attach_constraint_object(constraint);
105 
106  // Set system parameters
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  sys.time_solver->diff_solver()->relative_residual_tolerance = TOLERANCE*TOLERANCE;
125 
126  sys.solve();
127 
128  // Update _mesh from mesh_copy
129  for (auto * node_copy : mesh_copy.local_node_ptr_range())
130  {
131  auto & node = _mesh.node_ref(node_copy->id());
132  for (const auto d : make_range(mesh_copy.mesh_dimension()))
133  node(d) = (*node_copy)(d);
134  }
135 
136  SyncNodalPositions sync_object(_mesh);
137  Parallel::sync_dofobject_data_by_id (_mesh.comm(), _mesh.nodes_begin(), _mesh.nodes_end(), sync_object);
138 }
139 
140 } // namespace libMesh
141 
142 #endif // defined(LIBMESH_ENABLE_VSMOOTHER)
This is the EquationSystems class.
const Real _dilation_weight
Smoother control variables.
static constexpr Real TOLERANCE
std::unique_ptr< TimeSolver > time_solver
A pointer to the solver object we&#39;re going to use.
Definition: diff_system.h:253
MeshBase & mesh
const Parallel::Communicator & comm() const
UnstructuredMesh & _mesh
Definition: mesh_smoother.h:61
The libMesh namespace provides an interface to certain functionality in the library.
void attach_constraint_object(Constraint &constrain)
Register a user object for imposing constraints.
Definition: system.C:2209
This is a generic class that defines a solver to handle ImplicitSystem classes, including NonlinearIm...
Definition: diff_solver.h:68
The UnstructuredMesh class is derived from the MeshBase class.
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.
The DistributedMesh class is derived from the MeshBase class, and is intended to provide identical fu...
bool verbose
The DiffSolver may print a lot more to libMesh::out if verbose is set to true; default is false...
Definition: diff_solver.h:166
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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
unsigned int mesh_dimension() const
Definition: mesh_base.C:372
virtual void solve() override
Invokes the solver associated with the system.
Definition: fem_system.C:1076
VariationalMeshSmoother(UnstructuredMesh &mesh, Real dilation_weight=0.5, const bool preserve_subdomain_boundaries=true)
Simple constructor to use for smoothing purposes.
virtual const Node & node_ref(const dof_id_type i) const
Definition: mesh_base.h:596
This class provides the necessary interface for mesh smoothing.
Definition: mesh_smoother.h:38
virtual void init()
Initialize all the systems.
virtual System & add_system(std::string_view system_type, std::string_view name)
Add the system of type system_type named name to the systems array.
const bool _preserve_subdomain_boundaries
Whether subdomain boundaries are subject to change via smoothing.
virtual void smooth() override
Redefinition of the smooth function from the base class.