Skip to content
Snippets Groups Projects
Commit a2012935 authored by Dr. Felix Tobias Schindler's avatar Dr. Felix Tobias Schindler
Browse files

[localoperator.codim1] modernize

* move traits to internal namespace
* do not guard dune includes
* update evaluation handling
* add safe integer conversion
* mark ctors explicit
parent 88489979
No related branches found
No related tags found
No related merge requests found
...@@ -13,13 +13,9 @@ ...@@ -13,13 +13,9 @@
#include <boost/numeric/conversion/cast.hpp> #include <boost/numeric/conversion/cast.hpp>
#include <dune/stuff/common/disable_warnings.hh>
#include <dune/common/densematrix.hh> #include <dune/common/densematrix.hh>
#include <dune/stuff/common/reenable_warnings.hh>
#include <dune/stuff/common/disable_warnings.hh>
#include <dune/geometry/quadraturerules.hh> #include <dune/geometry/quadraturerules.hh>
#include <dune/stuff/common/reenable_warnings.hh>
#include <dune/stuff/functions/interfaces.hh> #include <dune/stuff/functions/interfaces.hh>
...@@ -31,10 +27,16 @@ namespace GDT { ...@@ -31,10 +27,16 @@ namespace GDT {
namespace LocalOperator { namespace LocalOperator {
// forward, to be used in the traits // forwards
template <class QuaternaryEvaluationImp> template <class QuaternaryEvaluationImp>
class Codim1CouplingIntegral; class Codim1CouplingIntegral;
template <class BinaryEvaluationImp>
class Codim1BoundaryIntegral;
namespace internal {
template <class QuaternaryEvaluationImp> template <class QuaternaryEvaluationImp>
class Codim1CouplingIntegralTraits class Codim1CouplingIntegralTraits
...@@ -45,38 +47,52 @@ class Codim1CouplingIntegralTraits ...@@ -45,38 +47,52 @@ class Codim1CouplingIntegralTraits
public: public:
typedef Codim1CouplingIntegral<QuaternaryEvaluationImp> derived_type; typedef Codim1CouplingIntegral<QuaternaryEvaluationImp> derived_type;
typedef LocalEvaluation::Codim1Interface<typename QuaternaryEvaluationImp::Traits, 4> QuaternaryEvaluationType;
}; };
template <class BinaryEvaluationImp>
class Codim1BoundaryIntegralTraits
{
static_assert(std::is_base_of<LocalEvaluation::Codim1Interface<typename BinaryEvaluationImp::Traits, 2>,
BinaryEvaluationImp>::value,
"BinaryEvaluationImp has to be derived from LocalEvaluation::Codim1Interface< ..., 2 >!");
public:
typedef Codim1BoundaryIntegral<BinaryEvaluationImp> derived_type;
};
} // namespace internal
template <class QuaternaryEvaluationImp> template <class QuaternaryEvaluationImp>
class Codim1CouplingIntegral class Codim1CouplingIntegral
: public LocalOperator::Codim1CouplingInterface<Codim1CouplingIntegralTraits<QuaternaryEvaluationImp>> : public LocalOperator::Codim1CouplingInterface<internal::Codim1CouplingIntegralTraits<QuaternaryEvaluationImp>>
{ {
public: public:
typedef Codim1CouplingIntegralTraits<QuaternaryEvaluationImp> Traits; typedef internal::Codim1CouplingIntegralTraits<QuaternaryEvaluationImp> Traits;
typedef typename Traits::QuaternaryEvaluationType QuaternaryEvaluationType; typedef QuaternaryEvaluationImp QuaternaryEvaluationType;
private: private:
static const size_t numTmpObjectsRequired_ = 4; static const size_t numTmpObjectsRequired_ = 4;
public: public:
template <class... Args> template <class... Args>
Codim1CouplingIntegral(Args&&... args) explicit Codim1CouplingIntegral(Args&&... args)
: evaluation_(std::forward<Args>(args)...) : evaluation_(std::forward<Args>(args)...)
, over_integrate_(0) , over_integrate_(0)
{ {
} }
template <class... Args> template <class... Args>
Codim1CouplingIntegral(const size_t over_integrate, Args&&... args) explicit Codim1CouplingIntegral(const size_t over_integrate, Args&&... args)
: evaluation_(std::forward<Args>(args)...) : evaluation_(std::forward<Args>(args)...)
, over_integrate_(over_integrate) , over_integrate_(over_integrate)
{ {
} }
template <class... Args> template <class... Args>
Codim1CouplingIntegral(const int over_integrate, Args&&... args) explicit Codim1CouplingIntegral(const int over_integrate, Args&&... args)
: evaluation_(std::forward<Args>(args)...) : evaluation_(std::forward<Args>(args)...)
, over_integrate_(boost::numeric_cast<size_t>(over_integrate)) , over_integrate_(boost::numeric_cast<size_t>(over_integrate))
{ {
...@@ -103,11 +119,11 @@ public: ...@@ -103,11 +119,11 @@ public:
const auto localFunctionsNe = evaluation_.localFunctions(neighbor); const auto localFunctionsNe = evaluation_.localFunctions(neighbor);
// quadrature // quadrature
const size_t integrand_order = const size_t integrand_order =
evaluation().order( evaluation_.order(
localFunctionsEn, localFunctionsNe, entityTestBase, entityAnsatzBase, neighborTestBase, neighborAnsatzBase) localFunctionsEn, localFunctionsNe, entityTestBase, entityAnsatzBase, neighborTestBase, neighborAnsatzBase)
+ over_integrate_; + over_integrate_;
assert(integrand_order < std::numeric_limits<int>::max()); const auto& faceQuadrature =
const auto& faceQuadrature = QuadratureRules<D, d - 1>::rule(intersection.type(), int(integrand_order)); QuadratureRules<D, d - 1>::rule(intersection.type(), boost::numeric_cast<int>(integrand_order));
// check matrices // check matrices
entityEntityRet *= 0.0; entityEntityRet *= 0.0;
neighborNeighborRet *= 0.0; neighborNeighborRet *= 0.0;
...@@ -126,28 +142,28 @@ public: ...@@ -126,28 +142,28 @@ public:
assert(neighborEntityRet.rows() >= rowsEn); assert(neighborEntityRet.rows() >= rowsEn);
assert(neighborEntityRet.cols() >= colsEn); assert(neighborEntityRet.cols() >= colsEn);
assert(tmpLocalMatrices.size() >= numTmpObjectsRequired_); assert(tmpLocalMatrices.size() >= numTmpObjectsRequired_);
Dune::DynamicMatrix<R>& entityEntityVals = tmpLocalMatrices[0]; auto& entityEntityVals = tmpLocalMatrices[0];
Dune::DynamicMatrix<R>& neighborNeighborVals = tmpLocalMatrices[1]; auto& neighborNeighborVals = tmpLocalMatrices[1];
Dune::DynamicMatrix<R>& entityNeighborVals = tmpLocalMatrices[2]; auto& entityNeighborVals = tmpLocalMatrices[2];
Dune::DynamicMatrix<R>& neighborEntityVals = tmpLocalMatrices[3]; auto& neighborEntityVals = tmpLocalMatrices[3];
// loop over all quadrature points // loop over all quadrature points
for (auto quadPoint = faceQuadrature.begin(); quadPoint != faceQuadrature.end(); ++quadPoint) { for (auto quadPoint = faceQuadrature.begin(); quadPoint != faceQuadrature.end(); ++quadPoint) {
const Dune::FieldVector<D, d - 1> localPoint = quadPoint->position(); const Dune::FieldVector<D, d - 1> localPoint = quadPoint->position();
const R integrationFactor = intersection.geometry().integrationElement(localPoint); const auto integrationFactor = intersection.geometry().integrationElement(localPoint);
const R quadratureWeight = quadPoint->weight(); const auto quadratureWeight = quadPoint->weight();
// evaluate local // evaluate local
evaluation().evaluate(localFunctionsEn, evaluation_.evaluate(localFunctionsEn,
localFunctionsNe, localFunctionsNe,
entityTestBase, entityTestBase,
entityAnsatzBase, entityAnsatzBase,
neighborTestBase, neighborTestBase,
neighborAnsatzBase, neighborAnsatzBase,
intersection, intersection,
localPoint, localPoint,
entityEntityVals, entityEntityVals,
neighborNeighborVals, neighborNeighborVals,
entityNeighborVals, entityNeighborVals,
neighborEntityVals); neighborEntityVals);
// compute integral // compute integral
assert(entityEntityVals.rows() >= rowsEn); assert(entityEntityVals.rows() >= rowsEn);
assert(entityEntityVals.cols() >= colsEn); assert(entityEntityVals.cols() >= colsEn);
...@@ -191,41 +207,18 @@ public: ...@@ -191,41 +207,18 @@ public:
} // void apply(...) const } // void apply(...) const
private: private:
const QuaternaryEvaluationType& evaluation() const const QuaternaryEvaluationType evaluation_;
{
return evaluation_;
}
const QuaternaryEvaluationImp evaluation_;
const size_t over_integrate_; const size_t over_integrate_;
}; // class Codim1CouplingIntegral }; // class Codim1CouplingIntegral
// forward, to be used in the traits
template <class BinaryEvaluationImp>
class Codim1BoundaryIntegral;
template <class BinaryEvaluationImp>
class Codim1BoundaryIntegralTraits
{
static_assert(std::is_base_of<LocalEvaluation::Codim1Interface<typename BinaryEvaluationImp::Traits, 2>,
BinaryEvaluationImp>::value,
"BinaryEvaluationImp has to be derived from LocalEvaluation::Codim1Interface< ..., 2 >!");
public:
typedef Codim1BoundaryIntegral<BinaryEvaluationImp> derived_type;
typedef LocalEvaluation::Codim1Interface<typename BinaryEvaluationImp::Traits, 2> BinaryEvaluationType;
};
template <class BinaryEvaluationImp> template <class BinaryEvaluationImp>
class Codim1BoundaryIntegral class Codim1BoundaryIntegral
: public LocalOperator::Codim1BoundaryInterface<Codim1BoundaryIntegralTraits<BinaryEvaluationImp>> : public LocalOperator::Codim1BoundaryInterface<internal::Codim1BoundaryIntegralTraits<BinaryEvaluationImp>>
{ {
public: public:
typedef Codim1BoundaryIntegralTraits<BinaryEvaluationImp> Traits; typedef internal::Codim1BoundaryIntegralTraits<BinaryEvaluationImp> Traits;
typedef typename Traits::BinaryEvaluationType BinaryEvaluationType; typedef BinaryEvaluationImp BinaryEvaluationType;
private: private:
static const size_t numTmpObjectsRequired_ = 1; static const size_t numTmpObjectsRequired_ = 1;
...@@ -269,9 +262,9 @@ public: ...@@ -269,9 +262,9 @@ public:
// quadrature // quadrature
typedef Dune::QuadratureRules<D, d - 1> FaceQuadratureRules; typedef Dune::QuadratureRules<D, d - 1> FaceQuadratureRules;
typedef Dune::QuadratureRule<D, d - 1> FaceQuadratureType; typedef Dune::QuadratureRule<D, d - 1> FaceQuadratureType;
const size_t integrand_order = evaluation().order(localFunctions, testBase, ansatzBase) + over_integrate_; const auto integrand_order = evaluation_.order(localFunctions, testBase, ansatzBase) + over_integrate_;
assert(integrand_order < std::numeric_limits<int>::max()); const FaceQuadratureType& faceQuadrature =
const FaceQuadratureType& faceQuadrature = FaceQuadratureRules::rule(intersection.type(), int(integrand_order)); FaceQuadratureRules::rule(intersection.type(), boost::numeric_cast<int>(integrand_order));
// check matrix and tmp storage // check matrix and tmp storage
ret *= 0.0; ret *= 0.0;
const size_t rows = testBase.size(); const size_t rows = testBase.size();
...@@ -286,7 +279,7 @@ public: ...@@ -286,7 +279,7 @@ public:
const R integrationFactor = intersection.geometry().integrationElement(localPoint); const R integrationFactor = intersection.geometry().integrationElement(localPoint);
const R quadratureWeight = quadPoint->weight(); const R quadratureWeight = quadPoint->weight();
// evaluate local // evaluate local
evaluation().evaluate(localFunctions, testBase, ansatzBase, intersection, localPoint, localMatrix); evaluation_.evaluate(localFunctions, testBase, ansatzBase, intersection, localPoint, localMatrix);
// compute integral // compute integral
assert(localMatrix.rows() >= rows); assert(localMatrix.rows() >= rows);
assert(localMatrix.cols() >= cols); assert(localMatrix.cols() >= cols);
...@@ -303,12 +296,7 @@ public: ...@@ -303,12 +296,7 @@ public:
} // void apply(...) const } // void apply(...) const
private: private:
const BinaryEvaluationType& evaluation() const const BinaryEvaluationType evaluation_;
{
return evaluation_;
}
const BinaryEvaluationImp evaluation_;
const size_t over_integrate_; const size_t over_integrate_;
}; // class Codim1BoundaryIntegral }; // class Codim1BoundaryIntegral
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment