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/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 52 852 : VariationalMeshSmoother::VariationalMeshSmoother(UnstructuredMesh & mesh, 53 : Real dilation_weight, 54 852 : const bool preserve_subdomain_boundaries) : 55 : MeshSmoother(mesh), 56 804 : _dilation_weight(dilation_weight), 57 852 : _preserve_subdomain_boundaries(preserve_subdomain_boundaries) 58 852 : {} 59 : 60 : 61 852 : void VariationalMeshSmoother::smooth(unsigned int) 62 : { 63 : // Check for multiple dimensions 64 876 : if (_mesh.elem_dimensions().size() > 1) 65 0 : 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 900 : DistributedMesh mesh_copy(_mesh); 89 900 : EquationSystems es(mesh_copy); 90 852 : 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 852 : VariationalSmootherConstraint constraint(sys, _preserve_subdomain_boundaries); 104 852 : sys.attach_constraint_object(constraint); 105 : 106 : // Set system parameters 107 852 : sys.get_dilation_weight() = _dilation_weight; 108 : 109 : // Set up solver 110 : sys.time_solver = 111 852 : 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 852 : 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 852 : sys.time_solver->diff_solver()->relative_residual_tolerance = TOLERANCE*TOLERANCE; 126 : 127 852 : sys.solve(); 128 : 129 : // Update _mesh from mesh_copy 130 88074 : for (auto * node_copy : mesh_copy.local_node_ptr_range()) 131 : { 132 47124 : auto & node = _mesh.node_ref(node_copy->id()); 133 223914 : for (const auto d : make_range(mesh_copy.mesh_dimension())) 134 137520 : node(d) = (*node_copy)(d); 135 804 : } 136 : 137 852 : SyncNodalPositions sync_object(_mesh); 138 1680 : Parallel::sync_dofobject_data_by_id (_mesh.comm(), _mesh.nodes_begin(), _mesh.nodes_end(), sync_object); 139 852 : } 140 : 141 : } // namespace libMesh 142 : 143 : #endif // defined(LIBMESH_ENABLE_VSMOOTHER)