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 : 19 : 20 : #ifndef LIBMESH_JUMP_ERROR_ESTIMATOR_H 21 : #define LIBMESH_JUMP_ERROR_ESTIMATOR_H 22 : 23 : // Local Includes 24 : #include "libmesh/dense_vector.h" 25 : #include "libmesh/error_estimator.h" 26 : #include "libmesh/fem_context.h" 27 : 28 : // C++ includes 29 : #include <cstddef> 30 : #include <string> 31 : #include <vector> 32 : #include <memory> 33 : 34 : namespace libMesh 35 : { 36 : 37 : // Forward Declarations 38 : class Point; 39 : class Elem; 40 : 41 : /** 42 : * This abstract base class implements utility functions for error estimators 43 : * which are based on integrated jumps between elements. 44 : * 45 : * \author Roy H. Stogner 46 : * \date 2006 47 : */ 48 : class JumpErrorEstimator : public ErrorEstimator 49 : { 50 : public: 51 : 52 : /** 53 : * Constructor. 54 : */ 55 21558 : JumpErrorEstimator() 56 21558 : : ErrorEstimator(), 57 20150 : scale_by_n_flux_faces(false), 58 20150 : use_unweighted_quadrature_rules(false), 59 20150 : integrate_slits(false), 60 20150 : integrate_boundary_sides(false), 61 : fine_context(), 62 : coarse_context(), 63 20150 : fine_error(0), 64 21558 : coarse_error(0) {} 65 : 66 : /** 67 : * This class cannot be (default) copy constructed/assigned because 68 : * it has unique_ptr members. Explicitly deleting these functions is 69 : * the best way to document this fact. 70 : */ 71 : JumpErrorEstimator (const JumpErrorEstimator &) = delete; 72 : JumpErrorEstimator & operator= (const JumpErrorEstimator &) = delete; 73 : 74 : /** 75 : * Defaulted move ctor, move assignment operator, and destructor. 76 : */ 77 : JumpErrorEstimator (JumpErrorEstimator &&) = default; 78 : JumpErrorEstimator & operator= (JumpErrorEstimator &&) = default; 79 4 : virtual ~JumpErrorEstimator() = default; 80 : 81 : /** 82 : * This function uses the derived class's jump error 83 : * estimate formula to estimate the error on each cell. 84 : * The estimated error is output in the vector 85 : * \p error_per_cell 86 : */ 87 : virtual void estimate_error (const System & system, 88 : ErrorVector & error_per_cell, 89 : const NumericVector<Number> * solution_vector = nullptr, 90 : bool estimate_parent_error = false) override; 91 : 92 : /** 93 : * This boolean flag allows you to scale the error indicator 94 : * result for each element by the number of "flux faces" the element 95 : * actually has. This tends to weight more evenly cells which are 96 : * on the boundaries and thus have fewer contributions to their flux. 97 : * The value is initialized to false, simply set it to true if you 98 : * want to use the feature. 99 : */ 100 : bool scale_by_n_flux_faces; 101 : 102 : /** 103 : * This boolean flag allows you to use "unweighted" quadrature rules 104 : * (sized to exactly integrate unweighted shape functions in master 105 : * element space) rather than "default" quadrature rules (sized to 106 : * exactly integrate polynomials of one higher degree than mass 107 : * matrix terms). The results with the former, lower-order rules 108 : * will be somewhat less accurate in many cases but will be much 109 : * cheaper to compute. 110 : * 111 : * The value is initialized to false, simply set it to true if you 112 : * want to use the feature. 113 : */ 114 : bool use_unweighted_quadrature_rules; 115 : 116 : /** 117 : * A boolean flag, by default false, to be set to true if integrations 118 : * should be performed on "slits" where two elements' faces overlap 119 : * even if those elements are not connected by neighbor links. 120 : * 121 : * This may only be useful for flex-IGA meshes, where 122 : * highly-nonconforming meshes are given pseudo-conforming solutions 123 : * via nodal constraints. 124 : * 125 : * Note that, to safely use this option in parallel, it is also 126 : * necessary to expand the default algebraic ghosting requirements 127 : * to include elements on the opposite sides of slits from local 128 : * elements. 129 : */ 130 : bool integrate_slits; 131 : 132 : protected: 133 : /** 134 : * A utility function to reinit the finite element data on elements sharing a 135 : * side 136 : */ 137 : void reinit_sides(); 138 : 139 : /** 140 : * A utility function to correctly increase n_flux_faces for the coarse element 141 : */ 142 : float coarse_n_flux_faces_increment(); 143 : 144 : /** 145 : * An initialization function, to give derived classes a chance to 146 : * request specific data from the FE objects 147 : */ 148 : virtual void init_context(FEMContext & c); 149 : 150 : /** 151 : * The function, to be implemented by derived classes, which calculates an error 152 : * term on an internal side 153 : */ 154 : virtual void internal_side_integration() = 0; 155 : 156 : /** 157 : * The function, to be implemented by derived classes, which 158 : * calculates an error term on a boundary side. 159 : * 160 : * \returns \p true if the flux bc function is in fact defined on 161 : * the current side. 162 : */ 163 0 : virtual bool boundary_side_integration() { return false; } 164 : 165 : /** 166 : * A boolean flag, by default false, to be set to true if integrations 167 : * with boundary_side_integration() should be performed 168 : */ 169 : bool integrate_boundary_sides; 170 : 171 : /** 172 : * Context objects for integrating on the fine and coarse elements 173 : * sharing a face 174 : */ 175 : std::unique_ptr<FEMContext> fine_context, coarse_context; 176 : 177 : /** 178 : * The fine and coarse error values to be set by each side_integration(); 179 : */ 180 : Real fine_error, coarse_error; 181 : 182 : /** 183 : * The variable number currently being evaluated 184 : */ 185 : unsigned int var; 186 : }; 187 : 188 : 189 : } // namespace libMesh 190 : 191 : #endif // LIBMESH_JUMP_ERROR_ESTIMATOR_H