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_EXACT_ERROR_ESTIMATOR_H 21 : #define LIBMESH_EXACT_ERROR_ESTIMATOR_H 22 : 23 : // Local Includes 24 : #include "libmesh/error_estimator.h" 25 : #include "libmesh/function_base.h" 26 : #include "libmesh/vector_value.h" 27 : #include "libmesh/tensor_value.h" 28 : 29 : // C++ includes 30 : #include <cstddef> 31 : #include <string> 32 : #include <vector> 33 : 34 : namespace libMesh 35 : { 36 : 37 : // Forward Declarations 38 : class Elem; 39 : template <typename T> class FEGenericBase; 40 : typedef FEGenericBase<Real> FEBase; 41 : class MeshFunction; 42 : class Point; 43 : class Parameters; 44 : template <typename T> class DenseVector; 45 : 46 : /** 47 : * This class implements an "error estimator" 48 : * based on the difference between the approximate 49 : * and exact solution. In theory the quadrature error 50 : * in this estimate should be much lower than the 51 : * approximation error in other estimates, so this 52 : * estimator can be used to calculate effectivity. 53 : * 54 : * \author Roy Stogner 55 : * \date 2006 56 : */ 57 : class ExactErrorEstimator : public ErrorEstimator 58 : { 59 : public: 60 : 61 : /** 62 : * Constructor. Responsible for initializing the _bc_function function 63 : * pointer to nullptr, and defaulting the norm type to H1. 64 : */ 65 : ExactErrorEstimator(); 66 : 67 : /** 68 : * This class cannot be (default) copy constructed/assigned because 69 : * it has containers of unique_ptrs. Explicitly deleting these functions is 70 : * the best way to document this fact. 71 : */ 72 : ExactErrorEstimator (const ExactErrorEstimator &) = delete; 73 : ExactErrorEstimator & operator= (const ExactErrorEstimator &) = delete; 74 : 75 : /** 76 : * Defaulted move ctor, move assignment operator, and destructor. 77 : */ 78 : ExactErrorEstimator (ExactErrorEstimator &&) = default; 79 : ExactErrorEstimator & operator= (ExactErrorEstimator &&) = default; 80 0 : virtual ~ExactErrorEstimator() = default; 81 : 82 : /** 83 : * Clone and attach arbitrary functors which compute the exact 84 : * values of the EquationSystems' solutions at any point. 85 : */ 86 : void attach_exact_values (std::vector<FunctionBase<Number> *> f); 87 : 88 : /** 89 : * Clone and attach an arbitrary functor which computes the exact 90 : * value of the system \p sys_num solution at any point. 91 : */ 92 : void attach_exact_value (unsigned int sys_num, 93 : FunctionBase<Number> * f); 94 : 95 : /** 96 : * Attach an arbitrary function which computes the exact value of 97 : * the solution at any point. 98 : */ 99 : typedef Number (*ValueFunctionPointer)(const Point & p, 100 : const Parameters & Parameters, 101 : const std::string & sys_name, 102 : const std::string & unknown_name); 103 : void attach_exact_value (ValueFunctionPointer fptr); 104 : 105 : /** 106 : * Clone and attach arbitrary functors which compute the exact 107 : * gradients of the EquationSystems' solutions at any point. 108 : */ 109 : void attach_exact_derivs (std::vector<FunctionBase<Gradient> *> g); 110 : 111 : /** 112 : * Clone and attach an arbitrary functor which computes the exact 113 : * gradient of the system \p sys_num solution at any point. 114 : */ 115 : void attach_exact_deriv (unsigned int sys_num, 116 : FunctionBase<Gradient> * g); 117 : 118 : /** 119 : * Attach an arbitrary function which computes the exact gradient of 120 : * the solution at any point. 121 : */ 122 : typedef Gradient (*GradientFunctionPointer)(const Point & p, 123 : const Parameters & parameters, 124 : const std::string & sys_name, 125 : const std::string & unknown_name); 126 : void attach_exact_deriv (GradientFunctionPointer gptr); 127 : 128 : /** 129 : * Clone and attach arbitrary functors which compute the exact 130 : * second derivatives of the EquationSystems' solutions at any point. 131 : */ 132 : void attach_exact_hessians (std::vector<FunctionBase<Tensor> *> h); 133 : 134 : /** 135 : * Clone and attach an arbitrary functor which computes the exact 136 : * second derivatives of the system \p sys_num solution at any point. 137 : */ 138 : void attach_exact_hessian (unsigned int sys_num, 139 : FunctionBase<Tensor> * h); 140 : 141 : /** 142 : * Attach an arbitrary function which computes the exact second 143 : * derivatives of the solution at any point. 144 : */ 145 : typedef Tensor (*HessianFunctionPointer)(const Point & p, 146 : const Parameters & parameters, 147 : const std::string & sys_name, 148 : const std::string & unknown_name); 149 : void attach_exact_hessian (HessianFunctionPointer hptr); 150 : 151 : /** 152 : * Attach function similar to system.h which 153 : * allows the user to attach a second EquationSystems 154 : * object with a reference fine grid solution. 155 : */ 156 : void attach_reference_solution (EquationSystems * es_fine); 157 : 158 : /** 159 : * Increases or decreases the order of the quadrature rule used for numerical 160 : * integration. The default \p extraorder is 1, because properly 161 : * integrating L2 error requires integrating the squares of terms 162 : * with order p+1, and 2p+2 is 1 higher than what we default to 163 : * using for reasonable mass matrix integration. 164 : */ 165 : void extra_quadrature_order (const int extraorder) 166 : { _extra_order = extraorder; } 167 : 168 : // Bring the base class functionality into the name lookup 169 : // procedure. This allows for alternative calling formats 170 : // defined in the base class. Thanks Wolfgang. 171 : // GCC 2.95.3 cannot compile such code. Since it was not really 172 : // essential to the functioning of this class, it's been removed. 173 : // using ErrorEstimator::estimate_error; 174 : 175 : /** 176 : * This function uses the exact solution function 177 : * to estimate the error on each cell. 178 : * The estimated error is output in the vector 179 : * \p error_per_cell 180 : */ 181 : virtual void estimate_error (const System & system, 182 : ErrorVector & error_per_cell, 183 : const NumericVector<Number> * solution_vector = nullptr, 184 : bool estimate_parent_error = false) override; 185 : 186 : virtual ErrorEstimatorType type() const override; 187 : 188 : private: 189 : 190 : /** 191 : * Function pointer to user-provided function which 192 : * computes the exact value of the solution. 193 : */ 194 : ValueFunctionPointer _exact_value; 195 : 196 : /** 197 : * Function pointer to user-provided function which 198 : * computes the exact derivative of the solution. 199 : */ 200 : GradientFunctionPointer _exact_deriv; 201 : 202 : /** 203 : * Function pointer to user-provided function which 204 : * computes the exact hessian of the solution. 205 : */ 206 : HessianFunctionPointer _exact_hessian; 207 : 208 : /** 209 : * User-provided functors which compute the exact value of the 210 : * solution for each system. 211 : */ 212 : std::vector<std::unique_ptr<FunctionBase<Number>>> _exact_values; 213 : 214 : /** 215 : * User-provided functors which compute the exact derivative of the 216 : * solution for each system. 217 : */ 218 : std::vector<std::unique_ptr<FunctionBase<Gradient>>> _exact_derivs; 219 : 220 : /** 221 : * User-provided functors which compute the exact hessians of the 222 : * solution for each system. 223 : */ 224 : std::vector<std::unique_ptr<FunctionBase<Tensor>>> _exact_hessians; 225 : 226 : /** 227 : * Constant pointer to the \p EquationSystems object 228 : * containing a fine grid solution. 229 : */ 230 : EquationSystems * _equation_systems_fine; 231 : 232 : /** 233 : * Helper method for calculating on each element 234 : */ 235 : Real find_squared_element_error (const System & system, 236 : const std::string & var_name, 237 : const Elem * elem, 238 : const DenseVector<Number> & Uelem, 239 : FEBase * fe, 240 : MeshFunction * fine_values) const; 241 : 242 : /** 243 : * Helper method for cleanup 244 : */ 245 : void clear_functors (); 246 : 247 : /** 248 : * Extra order to use for quadrature rule 249 : */ 250 : int _extra_order; 251 : }; 252 : 253 : 254 : } // namespace libMesh 255 : 256 : #endif // LIBMESH_EXACT_ERROR_ESTIMATOR_H