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_ERROR_ESTIMATOR_H 21 : #define LIBMESH_ERROR_ESTIMATOR_H 22 : 23 : // Local Includes 24 : #include "libmesh/libmesh_common.h" 25 : #include "libmesh/system_norm.h" 26 : 27 : // C++ includes 28 : #include <cstddef> 29 : #include <map> 30 : #include <string> 31 : #include <vector> 32 : #include <memory> 33 : 34 : namespace libMesh 35 : { 36 : 37 : // Forward Declarations 38 : class ErrorVector; 39 : class EquationSystems; 40 : class System; 41 : template <typename T> class NumericVector; 42 : 43 : namespace Parallel { 44 : class Communicator; 45 : } 46 : 47 : enum ErrorEstimatorType : int; 48 : 49 : /** 50 : * This class holds functions that will estimate the error 51 : * in a finite element solution on a given mesh. These error 52 : * estimates can be useful in their own right, or may be used 53 : * to guide adaptive mesh refinement. 54 : * 55 : * \note The computed errors are stored as floats rather 56 : * than doubles since the required precision is low. 57 : * 58 : * \author Benjamin S. Kirk 59 : * \date 2003 60 : */ 61 : class ErrorEstimator 62 : { 63 : public: 64 : 65 : /** 66 : * Constructor. Empty. Derived classes should reset error_norm as 67 : * appropriate. 68 : */ 69 35034 : ErrorEstimator() = default; 70 : 71 : /** 72 : * Copy/move ctor, copy/move assignment operator, and destructor are 73 : * all explicitly defaulted for this simple class. 74 : */ 75 : ErrorEstimator (const ErrorEstimator &) = default; 76 : ErrorEstimator (ErrorEstimator &&) = default; 77 : ErrorEstimator & operator= (const ErrorEstimator &) = default; 78 : ErrorEstimator & operator= (ErrorEstimator &&) = default; 79 5456 : virtual ~ErrorEstimator() = default; 80 : 81 : 82 : /** 83 : * This pure virtual function must be redefined 84 : * in derived classes to compute the error for each 85 : * active element and place it in the "error_per_cell" vector. 86 : * 87 : * If solution_vector is not nullptr, the estimator will 88 : * (if able) attempt to estimate an error in that field 89 : * instead of in system.solution. 90 : * 91 : * If estimate_parent_error is not false, the estimator will (if 92 : * able) attempt to give a consistent estimate of errors in parent 93 : * elements that would be generated by coarsening. 94 : */ 95 : virtual void estimate_error (const System & system, 96 : ErrorVector & error_per_cell, 97 : const NumericVector<Number> * solution_vector = nullptr, 98 : bool estimate_parent_error = false) = 0; 99 : 100 : /** 101 : * This virtual function can be redefined 102 : * in derived classes, but by default computes the sum of 103 : * the error_per_cell for each system in the equation_systems. 104 : * 105 : * Currently this function ignores the error_norm member variable, 106 : * and uses the function argument error_norms instead. 107 : * 108 : * This function is named estimate_errors instead of estimate_error 109 : * because otherwise C++ can get confused. 110 : */ 111 : virtual void estimate_errors (const EquationSystems & equation_systems, 112 : ErrorVector & error_per_cell, 113 : const std::map<const System *, SystemNorm> & error_norms, 114 : const std::map<const System *, const NumericVector<Number> *> * solution_vectors = nullptr, 115 : bool estimate_parent_error = false); 116 : 117 : /** 118 : * When calculating many error vectors at once, we need a data structure to 119 : * hold them all 120 : */ 121 : typedef std::map<std::pair<const System *, unsigned int>, std::unique_ptr<ErrorVector>> ErrorMap; 122 : 123 : /** 124 : * This virtual function can be redefined 125 : * in derived classes, but by default it calls estimate_error 126 : * repeatedly to calculate the requested error vectors. 127 : * 128 : * Currently this function ignores the error_norm.weight() values 129 : * because it calculates each variable's error individually, unscaled. 130 : * 131 : * The user selects which errors get computed by filling a map with error 132 : * vectors: If errors_per_cell[&system][v] exists, it will be filled with the 133 : * error values in variable \p v of \p system 134 : */ 135 : virtual void estimate_errors (const EquationSystems & equation_systems, 136 : ErrorMap & errors_per_cell, 137 : const std::map<const System *, const NumericVector<Number> *> * solution_vectors = nullptr, 138 : bool estimate_parent_error = false); 139 : 140 : /** 141 : * \returns The type for the ErrorEstimator subclass. 142 : */ 143 : virtual ErrorEstimatorType type() const = 0; 144 : 145 : /** 146 : * When estimating the error in a single system, the \p error_norm 147 : * is used to control the scaling and norm choice for each variable. 148 : * Not all estimators will support all norm choices. The default 149 : * scaling is for all variables to be weighted equally. The default 150 : * norm choice depends on the error estimator. 151 : * 152 : * Part of this functionality was supported via component_scale and 153 : * sobolev_order in older libMesh versions, and a small part was 154 : * supported via component_mask in even older versions. Hopefully 155 : * the encapsulation here will allow us to avoid changing this API 156 : * again. 157 : */ 158 : SystemNorm error_norm; 159 : 160 : protected: 161 : 162 : /** 163 : * This method takes the local error contributions in 164 : * \p error_per_cell from each processor and combines 165 : * them to get the global error vector. 166 : */ 167 : void reduce_error (std::vector<ErrorVectorReal> & error_per_cell, 168 : const Parallel::Communicator & comm) const; 169 : }; 170 : 171 : 172 : } // namespace libMesh 173 : 174 : #endif // LIBMESH_ERROR_ESTIMATOR_H