Line data Source code
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/parallel_ghost_sync.h" 33 : 34 : // C++ includes 35 : #include <time.h> // for clock_t, clock() 36 : #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC 37 : #include <cmath> 38 : #include <iomanip> 39 : #include <limits> 40 : 41 : namespace libMesh 42 : { 43 : 44 : // Optimization at -O2 or greater seem to break Intel's icc. So if we are 45 : // being compiled with icc let's dumb-down the optimizations for this file 46 : #ifdef __INTEL_COMPILER 47 : # pragma optimize ( "", off ) 48 : #endif 49 : 50 : // Member functions for the Variational Smoother 51 2414 : VariationalMeshSmoother::VariationalMeshSmoother( 52 : UnstructuredMesh &mesh, Real dilation_weight, 53 : const bool preserve_subdomain_boundaries, 54 : const Real relative_residual_tolerance, 55 2414 : const Real absolute_residual_tolerance, const unsigned int verbosity) 56 2278 : : MeshSmoother(mesh), _verbosity(verbosity), 57 2278 : _dilation_weight(dilation_weight), 58 2278 : _preserve_subdomain_boundaries(preserve_subdomain_boundaries), 59 2278 : _setup_called(false), 60 2278 : _relative_residual_tolerance(relative_residual_tolerance), 61 2482 : _absolute_residual_tolerance(absolute_residual_tolerance) {} 62 : 63 2414 : void VariationalMeshSmoother::setup() 64 : { 65 : // Check for multiple dimensions 66 2482 : if (_mesh.elem_dimensions().size() > 1) 67 0 : libmesh_not_implemented_msg("Meshes containing elements of differing dimension are not yet supported."); 68 : 69 : /* 70 : Ideally we'd want to update _mesh directly even if it already has an 71 : EquationSystems attached and doing work on it ... and that's thwarted by the 72 : problems: 73 : 74 : 1. We can't easily tell if there's already an EquationSystems attached. 75 : 2. We can't attach a second EquationSystems safely (if it already has a system 0) 76 : because the DoF indexing will need to be overwritten by our system. 77 : 3. The destructor of es won't even clean up after itself (we generally expect 78 : a mesh to go unused after its EquationSystems is destroyed), much less know 79 : how to restore anything from a previous EquationSystems. 80 : 81 : To avoid these issues, we'll just construct a new DistributedMesh _mesh_copy 82 : from _mesh, then do the solve on _mesh_copy, then copy its node locations back 83 : to _mesh after the solve is done. That'll be slightly less memory-efficient 84 : and somewhat more CPU-efficient in the case where _mesh is serial, though 85 : it'll be significantly less memory-efficient when _mesh is already distributed, 86 : but either way the robustness is probably worth it. 87 : */ 88 : 89 : // Create a new mesh, EquationSystems, and System 90 4760 : _mesh_copy = std::make_unique<DistributedMesh>(_mesh); 91 4760 : _equation_systems = std::make_unique<EquationSystems>(*_mesh_copy); 92 2414 : _system = 93 2414 : &(_equation_systems->add_system<VariationalSmootherSystem>("variational_smoother_system")); 94 : 95 : // Set this to something > 0 to add more quadrature points than the default 96 : // rule that integrates order 2 * fe_order + 1 polynomials exactly. 97 : // Using higher quadrature orders has not had a significant effect on observed solutions. 98 : //system()->extra_quadrature_order = 0; 99 : 100 : // Uncomment these to debug 101 : //system()->print_element_solutions=true; 102 : //system()->print_element_residuals=true; 103 : //system()->print_element_jacobians=true; 104 : 105 : // Add boundary node and hanging node constraints 106 2482 : _constraint = std::make_unique<VariationalSmootherConstraint>( 107 2414 : *_system, _preserve_subdomain_boundaries, _verbosity); 108 2482 : system()->attach_constraint_object(*_constraint); 109 : 110 : // Set system parameters 111 2414 : system()->set_verbosity(_verbosity); 112 2414 : system()->get_dilation_weight() = _dilation_weight; 113 : 114 : // Set up solver 115 2482 : system()->time_solver = std::make_unique<SteadySolver>(*_system); 116 : 117 : // Uncomment this line and use -snes_test_jacobian and -snes_test_jacobian_view 118 : // flags to compare the hand-coded jacobian in VariationalSmootherSystem 119 : // to finite difference jacobians. 120 : //system()->time_solver->diff_solver() = std::make_unique<PetscDiffSolver>(*_system); 121 : 122 2414 : _equation_systems->init(); 123 : 124 : // Solver verbosity 125 2414 : if (_verbosity > 15) { 126 146 : DiffSolver &solver = *(system()->time_solver->diff_solver().get()); 127 142 : solver.quiet = false; 128 142 : solver.verbose = true; 129 : } 130 : 131 : // Solver convergence tolerances 132 2482 : system()->time_solver->diff_solver()->relative_residual_tolerance = _relative_residual_tolerance; 133 2482 : system()->time_solver->diff_solver()->absolute_residual_tolerance = _absolute_residual_tolerance; 134 : 135 2414 : _setup_called = true; 136 2414 : } 137 : 138 2414 : void VariationalMeshSmoother::smooth() 139 : { 140 2414 : if (!_setup_called) 141 2272 : setup(); 142 : 143 : libmesh_try 144 : { 145 2414 : system()->solve(); 146 : } 147 : 148 0 : libmesh_catch (ConvergenceFailure& e) 149 : { 150 : #ifdef LIBMESH_ENABLE_EXCEPTIONS 151 0 : throw ConvergenceFailure("The VariationalMeshSmoother solve failed to converge."); 152 : #endif 153 0 : } 154 : 155 : // Update _mesh from _mesh_copy 156 75952 : for (auto * node_copy : _mesh_copy->local_node_ptr_range()) 157 : { 158 38832 : auto & node = _mesh.node_ref(node_copy->id()); 159 151812 : for (const auto d : make_range(_mesh_copy->mesh_dimension())) 160 112980 : node(d) = (*node_copy)(d); 161 2278 : } 162 : 163 2414 : SyncNodalPositions sync_object(_mesh); 164 4760 : Parallel::sync_dofobject_data_by_id (_mesh.comm(), _mesh.nodes_begin(), _mesh.nodes_end(), sync_object); 165 : 166 : // Release memory occupied by _mesh_copy 167 : // Destruct this before _mesh_copy because it references _mesh_copy 168 68 : _equation_systems.reset(); 169 68 : _mesh_copy.reset(); 170 : // We'll need to call setup again since we'll have to reconstruct _mesh_copy 171 2414 : _setup_called = false; 172 2414 : } 173 : 174 142 : const MeshQualityInfo & VariationalMeshSmoother::get_mesh_info() const 175 : { 176 142 : libmesh_error_msg_if(!_setup_called, "Need to first call the setup() method of " 177 : << "this VariationalMeshSmoother object, then call get_mesh_info() prior " 178 : << "to calling smooth()."); 179 142 : return _system->get_mesh_info(); 180 : } 181 : 182 : } // namespace libMesh 183 : 184 : #endif // defined(LIBMESH_ENABLE_VSMOOTHER)