3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
numericepsilon.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_ASSEMBLY_NUMERIC_EPSILON_HH
25#define DUMUX_ASSEMBLY_NUMERIC_EPSILON_HH
26
27#include <dune/common/fvector.hh>
30
31namespace Dumux {
32
39template<class Scalar, int numEq>
41{
42 using NumEqVector = Dune::FieldVector<Scalar, numEq>;
43
44public:
45 explicit NumericEpsilon(const std::string& paramGroup = "")
46 {
47 // get epsilons from input file with invalid default
48 baseEps_ = getParamFromGroup<Scalar>(paramGroup, "Assembly.NumericDifference.BaseEpsilon", 1e-10);
49 magnitude_ = getParamFromGroup<NumEqVector>(paramGroup, "Assembly.NumericDifference.PriVarMagnitude", NumEqVector(-1));
50 }
51
57 Scalar operator() (Scalar priVar, int priVarIdx) const noexcept
58 {
59 return magnitude_[priVarIdx] > 0.0 ? baseEps_*magnitude_[priVarIdx]
60 : NumericDifferentiation::epsilon(priVar, baseEps_);
61 }
62
63private:
64 Scalar baseEps_;
65 NumEqVector magnitude_;
66};
67
68} // end namespace Dumux
69
70#endif
A class for numeric differentiation.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
A helper class for local assemblers using numeric differentiation to determine the epsilon.
Definition: numericepsilon.hh:41
NumericEpsilon(const std::string &paramGroup="")
Definition: numericepsilon.hh:45
Scalar operator()(Scalar priVar, int priVarIdx) const noexcept
get the epsilon
Definition: numericepsilon.hh:57
static Scalar epsilon(const Scalar value, const Scalar baseEps=1e-10)
Computes the epsilon used for numeric differentiation.
Definition: numericdifferentiation.hh:47