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