Skip to content
Snippets Groups Projects
Commit b8182f00 authored by Max Kahnt's avatar Max Kahnt
Browse files

Use the field type from the base type.

Norm derives the field_type from the FieldTraits struct
and overriding methods should be consistent with that.
parent 39e4fb15
Branches
No related tags found
1 merge request!17Cleanup/norms
Pipeline #
......@@ -10,10 +10,16 @@
#include "norm.hh"
//! A norm for blocks of vectors
template <class VectorType>
class BlockNorm: public Norm<VectorType>
template <class V>
class BlockNorm: public Norm<V>
{
public:
using VectorType = V;
using Base = Norm<V>;
/** \brief The type used for the result */
using typename Base::field_type;
// typedef typename std::vector<const Norm<typename VectorType::block_type>* > NormPointerVector;
typedef std::vector< std::shared_ptr<const Norm<typename VectorType::block_type> > > NormPointerVector;
......@@ -22,24 +28,24 @@ class BlockNorm: public Norm<VectorType>
{}
//! Compute the norm of the given vector
double operator()(const VectorType &v) const
field_type operator()(const VectorType &v) const override
{
double r = 0.0;
field_type r = 0.0;
for (int i=0; i<norms_.size(); ++i)
{
double ri = (*norms_[i])(v[i]);
auto ri = (*norms_[i])(v[i]);
r += ri*ri;
}
return std::sqrt(r);
}
//! Compute the norm of the difference of two vectors
double diff(const VectorType &v1, const VectorType &v2) const
field_type diff(const VectorType &v1, const VectorType &v2) const override
{
double r = 0.0;
field_type r = 0.0;
for (int i=0; i<norms_.size(); ++i)
{
double ri = (*norms_[i]).diff(v1[i], v2[i]);
auto ri = (*norms_[i]).diff(v1[i], v2[i]);
r += ri*ri;
}
return std::sqrt(r);
......
......@@ -21,7 +21,10 @@ class DiagNorm:
{
public:
typedef V VectorType;
typedef typename VectorType::field_type field_type;
using Base = Norm<V>;
/** \brief The type used for the result */
using typename Base::field_type;
private:
typedef typename VectorType::size_type SizeType;
......@@ -32,13 +35,13 @@ class DiagNorm:
{}
//! Compute the norm of the given vector
field_type operator()(const VectorType &f) const
field_type operator()(const VectorType &f) const override
{
return std::sqrt(normSquared(f));
}
//! Compute the square of the norm of the given vector
virtual field_type normSquared(const VectorType& f) const
virtual field_type normSquared(const VectorType& f) const override
{
field_type r = 0.0;
......@@ -53,7 +56,7 @@ class DiagNorm:
}
//! Compute the norm of the difference of two vectors
field_type diff(const VectorType &f1, const VectorType &f2) const
field_type diff(const VectorType &f1, const VectorType &f2) const override
{
field_type r = 0.0;
......@@ -79,8 +82,11 @@ class DiagNorm<Dune::BlockVector<Dune::FieldVector <double,1> >, Dune::BlockVect
public Norm<Dune::BlockVector<Dune::FieldVector <double,1> > >
{
public:
typedef double field_type;
typedef Dune::BlockVector<Dune::FieldVector <field_type,1> > VectorType;
typedef Dune::BlockVector<Dune::FieldVector <double,1> > VectorType;
using Base = Norm<VectorType>;
/** \brief The type used for the result */
using typename Base::field_type;
private:
typedef VectorType::size_type SizeType;
......@@ -91,13 +97,13 @@ class DiagNorm<Dune::BlockVector<Dune::FieldVector <double,1> >, Dune::BlockVect
{}
//! Compute the norm of the given vector
field_type operator()(const VectorType &f) const
field_type operator()(const VectorType &f) const override
{
return std::sqrt(normSquared(f));
}
//! Compute the square of the norm of the given vector
virtual field_type normSquared(const VectorType& f) const
virtual field_type normSquared(const VectorType& f) const override
{
field_type r = 0.0;
......@@ -108,7 +114,7 @@ class DiagNorm<Dune::BlockVector<Dune::FieldVector <double,1> >, Dune::BlockVect
}
//! Compute the norm of the difference of two vectors
field_type diff(const VectorType &f1, const VectorType &f2) const
field_type diff(const VectorType &f1, const VectorType &f2) const override
{
field_type r = 0.0;
......
......@@ -34,17 +34,18 @@ namespace Solvers {
{
public:
typedef V VectorType;
using Base = Norm<V>;
/** \brief The type used for the result */
typedef typename VectorType::field_type field_type;
using typename Base::field_type;
EnergyNorm(const double tol=1e-10 ) : iterationStep_(NULL), matrix_(NULL), tol_(tol) {}
EnergyNorm(const field_type tol=1e-10 ) : iterationStep_(NULL), matrix_(NULL), tol_(tol) {}
EnergyNorm(LinearIterationStep<MatrixType, VectorType>& it, const double tol=1e-10)
EnergyNorm(LinearIterationStep<MatrixType, VectorType>& it, const field_type tol=1e-10)
: iterationStep_(&it), matrix_(NULL), tol_(tol)
{}
EnergyNorm(const MatrixType& matrix, const double tol=1e-10)
EnergyNorm(const MatrixType& matrix, const field_type tol=1e-10)
: iterationStep_(NULL), matrix_(&matrix), tol_(tol)
{}
......@@ -57,7 +58,7 @@ namespace Solvers {
}
//! Compute the norm of the difference of two vectors
field_type diff(const VectorType& f1, const VectorType& f2) const {
field_type diff(const VectorType& f1, const VectorType& f2) const override {
if (iterationStep_ == NULL && matrix_ == NULL)
DUNE_THROW(Dune::Exception, "You have supplied neither a matrix nor an IterationStep to the EnergyNorm!");
......@@ -67,13 +68,13 @@ namespace Solvers {
}
//! Compute the norm of the given vector
field_type operator()(const VectorType& f) const
field_type operator()(const VectorType& f) const override
{
return std::sqrt(normSquared(f));
}
// \brief Compute the square of the norm of the given vector
virtual field_type normSquared(const VectorType& f) const
virtual field_type normSquared(const VectorType& f) const override
{
if (iterationStep_ == NULL && matrix_ == NULL)
DUNE_THROW(Dune::Exception, "You have supplied neither a matrix nor an IterationStep to the EnergyNorm!");
......@@ -110,7 +111,7 @@ namespace Solvers {
return ret;
};
}
protected:
......@@ -118,7 +119,7 @@ namespace Solvers {
const MatrixType* matrix_;
const double tol_;
const field_type tol_;
};
......
......@@ -23,7 +23,10 @@ class FullNorm: public Norm<V>
{
public:
typedef V VectorType;
typedef typename VectorType::field_type field_type;
using Base = Norm<V>;
/** \brief The type used for the result */
using typename Base::field_type;
FullNorm(const field_type alpha, const LowRankFactor &lowRankFactor) :
lowRankFactor_(lowRankFactor),
......@@ -31,13 +34,13 @@ class FullNorm: public Norm<V>
{}
//! Compute the norm of the given vector
field_type operator()(const VectorType &f) const
field_type operator()(const VectorType &f) const override
{
return std::sqrt(normSquared(f));
}
//! Compute the square of the norm of the given vector
virtual field_type normSquared(const VectorType& f) const
virtual field_type normSquared(const VectorType& f) const override
{
VectorType r(lowRankFactor_.N());
r = 0.0;
......@@ -48,7 +51,7 @@ class FullNorm: public Norm<V>
}
//! Compute the norm of the difference of two vectors
field_type diff(const VectorType &f1, const VectorType &f2) const
field_type diff(const VectorType &f1, const VectorType &f2) const override
{
VectorType r(lowRankFactor_.N());
r = 0.0;
......@@ -72,8 +75,11 @@ class FullNorm<Dune::BlockVector<Dune::FieldVector<double,1> >, Dune::BlockVecto
public Norm<Dune::BlockVector<Dune::FieldVector<double,1> > >
{
public:
typedef double field_type;
typedef Dune::BlockVector<Dune::FieldVector<field_type,1> > VectorType;
typedef Dune::BlockVector<Dune::FieldVector<double,1> > VectorType;
using Base = Norm<VectorType>;
/** \brief The type used for the result */
using typename Base::field_type;
private:
typedef VectorType::size_type SizeType;
public:
......@@ -84,13 +90,13 @@ class FullNorm<Dune::BlockVector<Dune::FieldVector<double,1> >, Dune::BlockVecto
{}
//! Compute the norm of the given vector
field_type operator()(const VectorType &f) const
field_type operator()(const VectorType &f) const override
{
return std::sqrt(normSquared(f));
}
//! Compute the square of the norm of the given vector
virtual field_type normSquared(const VectorType& f) const
virtual field_type normSquared(const VectorType& f) const override
{
field_type r = 0.0;
......@@ -101,7 +107,7 @@ class FullNorm<Dune::BlockVector<Dune::FieldVector<double,1> >, Dune::BlockVecto
}
//! Compute the norm of the difference of two vectors
field_type diff(const VectorType &f1, const VectorType &f2) const
field_type diff(const VectorType &f1, const VectorType &f2) const override
{
field_type r = 0.0;
......
......@@ -16,6 +16,13 @@ class H1SemiNorm : public Norm<V>
{
public:
typedef V VectorType;
using Base = Norm<V>;
using MatrixType = M;
/** \brief The type used for the result */
using typename Base::field_type;
H1SemiNorm() : matrix_(NULL) {}
H1SemiNorm(const MatrixType& matrix)
......@@ -25,12 +32,12 @@ public:
}
//! Compute the norm of the difference of two vectors
double diff(const VectorType& u1, const VectorType& u2) const
field_type diff(const VectorType& u1, const VectorType& u2) const override
{
assert(u1.size()==u2.size());
assert(u1.size()==matrix_->N());
double sum = 0;
field_type sum = 0;
for (size_t i=0; i<matrix_->N(); i++)
{
......@@ -46,7 +53,7 @@ public:
}
//! Compute the norm of the given vector
double operator()(const VectorType& u) const
field_type operator()(const VectorType& u) const override
{
if (matrix_ == NULL)
DUNE_THROW(Dune::Exception, "You have not supplied neither a matrix nor an IterationStep to the EnergyNorm routine!");
......@@ -55,7 +62,7 @@ public:
// we compute sqrt{uAu^t}. We have implemented this by hand because the matrix is
// always scalar but the vectors may not be
double sum = 0;
field_type sum = 0;
for (size_t i=0; i<matrix_->N(); i++) {
......
......@@ -14,13 +14,18 @@ template <class V>
class PNorm: public Norm<V>
{
public:
typedef V Vector;
PNorm(int p=2, double alpha=1.0):
typedef V VectorType;
using Base = Norm<V>;
/** \brief The type used for the result */
using typename Base::field_type;
PNorm(int p=2, field_type alpha=1.0):
p(p),
alpha(alpha)
{}
double operator()(const Vector &v) const
field_type operator()(const VectorType &v) const override
{
double r = 0.0;
......@@ -54,15 +59,15 @@ class PNorm: public Norm<V>
return alpha*r;
}
double diff(const Vector &v1, const Vector &v2) const
field_type diff(const VectorType &v1, const VectorType &v2) const override
{
double r = 0.0;
field_type r = 0.0;
if (p<1)
{
for(int row = 0; row < v1.size(); ++row)
{
double z = std::abs(v1[row]-v2[row]);
field_type z = std::abs(v1[row]-v2[row]);
if (z>r)
r = z;
}
......@@ -89,7 +94,7 @@ class PNorm: public Norm<V>
}
private:
double alpha;
field_type alpha;
int p;
};
......
......@@ -13,10 +13,16 @@
#include "norm.hh"
//! A norm for blocks of interlaced vectors
template <class VectorType, class ReorderedLocalVector>
class ReorderedBlockNorm: public Norm<VectorType>
template <class V, class ReorderedLocalVector>
class ReorderedBlockNorm: public Norm<V>
{
public:
typedef V VectorType;
using Base = Norm<V>;
/** \brief The type used for the result */
using typename Base::field_type;
typedef Dune::BlockVector<ReorderedLocalVector> ReorderedVector;
typedef std::vector< std::shared_ptr<const Norm<ReorderedLocalVector> > > NormPointerVector;
......@@ -26,13 +32,13 @@ class ReorderedBlockNorm: public Norm<VectorType>
{}
//! Compute the norm of the given vector
double operator()(const VectorType &v) const
field_type operator()(const VectorType &v) const override
{
return std::sqrt(normSquared(v));
}
//! Compute the square of the norm of the given vector
double normSquared(const VectorType& v) const
field_type normSquared(const VectorType& v) const override
{
double r = 0.0;
......@@ -48,7 +54,7 @@ class ReorderedBlockNorm: public Norm<VectorType>
}
//! Compute the norm of the difference of two vectors
double diff(const VectorType &v1, const VectorType &v2) const
field_type diff(const VectorType &v1, const VectorType &v2) const override
{
double r = 0.0;
......
......@@ -12,7 +12,10 @@ class SumNorm: public Norm<V>
{
public:
typedef V VectorType;
typedef typename VectorType::field_type field_type;
using Base = Norm<V>;
/** \brief The type used for the result */
using typename Base::field_type;
SumNorm(field_type _alpha1, const Norm<VectorType> &_norm1, field_type _alpha2, const Norm<VectorType> &_norm2) :
alpha1(_alpha1),
......@@ -22,13 +25,13 @@ class SumNorm: public Norm<V>
{}
//! Compute the norm of the given vector
field_type operator()(const VectorType &f) const
field_type operator()(const VectorType &f) const override
{
return std::sqrt(normSquared(f));
}
//! Compute the square of the norm of the given vector
virtual field_type normSquared(const VectorType& f) const
virtual field_type normSquared(const VectorType& f) const override
{
field_type r1 = norm1.normSquared(f);
field_type r2 = norm2.normSquared(f);
......@@ -37,7 +40,7 @@ class SumNorm: public Norm<V>
}
//! Compute the norm of the difference of two vectors
field_type diff(const VectorType &f1, const VectorType &f2) const
field_type diff(const VectorType &f1, const VectorType &f2) const override
{
field_type r1 = norm1.diff(f1,f2);
field_type r2 = norm2.diff(f1,f2);
......
......@@ -11,23 +11,27 @@
#include "norm.hh"
//! Wrapper around the two_norm() method of a vector class
template <class VectorType>
class TwoNorm : public Norm<VectorType>
template <class V>
class TwoNorm : public Norm<V>
{
public:
typedef V VectorType;
using Base = Norm<V>;
/** \brief The type used for the result */
using typename Base::field_type;
/** \brief Destructor, doing nothing */
virtual ~TwoNorm() {}
//! Compute the norm of the given vector
virtual double operator()(const VectorType& f) const
virtual field_type operator()(const VectorType& f) const override
{
return f.two_norm();
}
//! Compute the square of the norm of the given vector
virtual double normSquared(const VectorType& f) const
virtual field_type normSquared(const VectorType& f) const override
{
return f.two_norm2();
}
......@@ -36,7 +40,7 @@ class TwoNorm : public Norm<VectorType>
// assumed (above) anyway. Then also the specialization for
// FieldVector (below) can go.
//! Compute the norm of the difference of two vectors
virtual double diff(const VectorType& f1, const VectorType& f2) const
virtual field_type diff(const VectorType& f1, const VectorType& f2) const override
{
assert(f1.size() == f2.size());
double r = 0.0;
......@@ -66,19 +70,19 @@ class TwoNorm<Dune::FieldVector<T, dimension> >
virtual ~TwoNorm() {};
//! Compute the norm of the given vector
virtual double operator()(const VectorType& f) const
virtual T operator()(const VectorType& f) const override
{
return f.two_norm();
}
//! Compute the square of the norm of the given vector
virtual double normSquared(const VectorType& f) const
virtual T normSquared(const VectorType& f) const override
{
return f.two_norm2();
}
//! Compute the norm of the difference of two vectors
virtual double diff(const VectorType& f1, const VectorType& f2) const
virtual T diff(const VectorType& f1, const VectorType& f2) const override
{
VectorType tmp = f1;
tmp -= f2;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment