Skip to content

Commit

Permalink
Remove template parameter from FPArrayComparator
Browse files Browse the repository at this point in the history
  • Loading branch information
peterrum committed Sep 20, 2022
1 parent 51d619e commit a22baf0
Show file tree
Hide file tree
Showing 5 changed files with 104 additions and 67 deletions.
4 changes: 2 additions & 2 deletions include/deal.II/matrix_free/dof_info.h
Original file line number Diff line number Diff line change
Expand Up @@ -44,7 +44,7 @@ namespace internal
template <int dim>
class HangingNodes;

template <typename, typename>
template <typename>
struct FPArrayComparator;

struct TaskInfo;
Expand Down Expand Up @@ -111,7 +111,7 @@ namespace internal
std::pair<std::vector<Number>, types::global_dof_index> next_constraint;
std::map<std::vector<Number>,
types::global_dof_index,
FPArrayComparator<Number, VectorizedArray<Number>>>
FPArrayComparator<VectorizedArray<Number>>>
constraints;
};

Expand Down
2 changes: 1 addition & 1 deletion include/deal.II/matrix_free/dof_info.templates.h
Original file line number Diff line number Diff line change
Expand Up @@ -41,7 +41,7 @@ namespace internal
{
template <typename Number>
ConstraintValues<Number>::ConstraintValues()
: constraints(FPArrayComparator<Number>(1.))
: constraints(FPArrayComparator<VectorizedArray<Number>>(1.))
{}


Expand Down
27 changes: 12 additions & 15 deletions include/deal.II/matrix_free/mapping_info.h
Original file line number Diff line number Diff line change
Expand Up @@ -310,10 +310,14 @@ namespace internal
* satisfied). This is not a problem in the use cases for this class, but be
* careful when using it in other contexts.
*/
template <typename Number,
typename VectorizedArrayType = VectorizedArray<Number>>
template <typename VectorizedArrayType>
struct FPArrayComparator
{
using Number = typename internal::VectorizedArrayTrait<
VectorizedArrayType>::value_type;
static constexpr std::size_t width =
internal::VectorizedArrayTrait<VectorizedArrayType>::width;

FPArrayComparator(const Number scaling);

/**
Expand All @@ -328,33 +332,26 @@ namespace internal
* issues).
*/
bool
operator()(
const Tensor<1, VectorizedArrayType::size(), Number> &t1,
const Tensor<1, VectorizedArrayType::size(), Number> &t2) const;
operator()(const Tensor<1, width, Number> &t1,
const Tensor<1, width, Number> &t2) const;

/**
* Compare two rank-1 tensors of vectorized arrays (stored as tensors to
* avoid alignment issues).
*/
template <int dim>
bool
operator()(
const Tensor<1, dim, Tensor<1, VectorizedArrayType::size(), Number>>
&t1,
const Tensor<1, dim, Tensor<1, VectorizedArrayType::size(), Number>>
&t2) const;
operator()(const Tensor<1, dim, Tensor<1, width, Number>> &t1,
const Tensor<1, dim, Tensor<1, width, Number>> &t2) const;

/**
* Compare two rank-2 tensors of vectorized arrays (stored as tensors to
* avoid alignment issues).
*/
template <int dim>
bool
operator()(
const Tensor<2, dim, Tensor<1, VectorizedArrayType::size(), Number>>
&t1,
const Tensor<2, dim, Tensor<1, VectorizedArrayType::size(), Number>>
&t2) const;
operator()(const Tensor<2, dim, Tensor<1, width, Number>> &t1,
const Tensor<2, dim, Tensor<1, width, Number>> &t2) const;

/**
* Compare two arrays of tensors.
Expand Down
103 changes: 70 additions & 33 deletions include/deal.II/matrix_free/mapping_info.templates.h
Original file line number Diff line number Diff line change
Expand Up @@ -347,12 +347,12 @@ namespace internal
struct CompressedCellData
{
CompressedCellData(const double expected_size)
: data(FPArrayComparator<Number, VectorizedArrayType>(expected_size))
: data(FPArrayComparator<VectorizedArrayType>(expected_size))
{}

std::map<Tensor<2, dim, Tensor<1, VectorizedArrayType::size(), Number>>,
unsigned int,
FPArrayComparator<Number, VectorizedArrayType>>
FPArrayComparator<VectorizedArrayType>>
data;
};

Expand Down Expand Up @@ -1077,10 +1077,11 @@ namespace internal
// first quadrature point on the cell - we use a relatively coarse
// tolerance to account for some inaccuracies in the manifold
// evaluation
const FPArrayComparator<double> comparator(1e4 * jacobian_size);
const FPArrayComparator<VectorizedArray<double>> comparator(
1e4 * jacobian_size);
std::map<std::array<Tensor<2, dim>, dim + 1>,
unsigned int,
FPArrayComparator<double>>
FPArrayComparator<VectorizedArray<double>>>
compressed_jacobians(comparator);

unsigned int n_data_buckets = 0;
Expand Down Expand Up @@ -1503,8 +1504,7 @@ namespace internal
// CompressedCellData) and add another factor of 512 to account for
// some roundoff effects.
CompressedFaceData(const Number jacobian_size)
: data(FPArrayComparator<Number, VectorizedArrayType>(512. /
jacobian_size))
: data(FPArrayComparator<VectorizedArrayType>(512. / jacobian_size))
, jacobian_size(jacobian_size)
{}

Expand All @@ -1518,7 +1518,7 @@ namespace internal
2 * dim * dim + dim + 1,
Tensor<1, VectorizedArrayType::size(), Number>>,
unsigned int,
FPArrayComparator<Number, VectorizedArrayType>>
FPArrayComparator<VectorizedArrayType>>
data;

// Store the scaling factor
Expand Down Expand Up @@ -3351,19 +3351,21 @@ namespace internal

/* ------------------------------------------------------------------ */

template <typename Number, typename VectorizedArrayType>
FPArrayComparator<Number, VectorizedArrayType>::FPArrayComparator(
const Number scaling)
template <typename VectorizedArrayType>
FPArrayComparator<VectorizedArrayType>::FPArrayComparator(
const typename FPArrayComparator<VectorizedArrayType>::Number scaling)
: tolerance(scaling * std::numeric_limits<double>::epsilon() * 1024.)
{}



template <typename Number, typename VectorizedArrayType>
template <typename VectorizedArrayType>
bool
FPArrayComparator<Number, VectorizedArrayType>::operator()(
const std::vector<Number> &v1,
const std::vector<Number> &v2) const
FPArrayComparator<VectorizedArrayType>::operator()(
const std::vector<typename FPArrayComparator<VectorizedArrayType>::Number>
&v1,
const std::vector<typename FPArrayComparator<VectorizedArrayType>::Number>
&v2) const
{
const unsigned int s1 = v1.size(), s2 = v2.size();
if (s1 < s2)
Expand All @@ -3381,13 +3383,20 @@ namespace internal



template <typename Number, typename VectorizedArrayType>
template <typename VectorizedArrayType>
bool
FPArrayComparator<Number, VectorizedArrayType>::operator()(
const Tensor<1, VectorizedArrayType::size(), Number> &t1,
const Tensor<1, VectorizedArrayType::size(), Number> &t2) const
FPArrayComparator<VectorizedArrayType>::operator()(
const Tensor<1,
FPArrayComparator<VectorizedArrayType>::width,
typename FPArrayComparator<VectorizedArrayType>::Number> &t1,
const Tensor<1,
FPArrayComparator<VectorizedArrayType>::width,
typename FPArrayComparator<VectorizedArrayType>::Number> &t2)
const
{
for (unsigned int k = 0; k < VectorizedArrayType::size(); ++k)
for (unsigned int k = 0;
k < FPArrayComparator<VectorizedArrayType>::width;
++k)
if (t1[k] < t2[k] - tolerance)
return true;
else if (t1[k] > t2[k] + tolerance)
Expand All @@ -3397,16 +3406,28 @@ namespace internal



template <typename Number, typename VectorizedArrayType>
template <typename VectorizedArrayType>
template <int dim>
bool
FPArrayComparator<Number, VectorizedArrayType>::operator()(
const Tensor<1, dim, Tensor<1, VectorizedArrayType::size(), Number>> &t1,
const Tensor<1, dim, Tensor<1, VectorizedArrayType::size(), Number>> &t2)
FPArrayComparator<VectorizedArrayType>::operator()(
const Tensor<
1,
dim,
Tensor<1,
FPArrayComparator<VectorizedArrayType>::width,
typename FPArrayComparator<VectorizedArrayType>::Number>> &t1,
const Tensor<
1,
dim,
Tensor<1,
FPArrayComparator<VectorizedArrayType>::width,
typename FPArrayComparator<VectorizedArrayType>::Number>> &t2)
const
{
for (unsigned int d = 0; d < dim; ++d)
for (unsigned int k = 0; k < VectorizedArrayType::size(); ++k)
for (unsigned int k = 0;
k < FPArrayComparator<VectorizedArrayType>::width;
++k)
if (t1[d][k] < t2[d][k] - tolerance)
return true;
else if (t1[d][k] > t2[d][k] + tolerance)
Expand All @@ -3416,17 +3437,29 @@ namespace internal



template <typename Number, typename VectorizedArrayType>
template <typename VectorizedArrayType>
template <int dim>
bool
FPArrayComparator<Number, VectorizedArrayType>::operator()(
const Tensor<2, dim, Tensor<1, VectorizedArrayType::size(), Number>> &t1,
const Tensor<2, dim, Tensor<1, VectorizedArrayType::size(), Number>> &t2)
FPArrayComparator<VectorizedArrayType>::operator()(
const Tensor<
2,
dim,
Tensor<1,
FPArrayComparator<VectorizedArrayType>::width,
typename FPArrayComparator<VectorizedArrayType>::Number>> &t1,
const Tensor<
2,
dim,
Tensor<1,
FPArrayComparator<VectorizedArrayType>::width,
typename FPArrayComparator<VectorizedArrayType>::Number>> &t2)
const
{
for (unsigned int d = 0; d < dim; ++d)
for (unsigned int e = 0; e < dim; ++e)
for (unsigned int k = 0; k < VectorizedArrayType::size(); ++k)
for (unsigned int k = 0;
k < FPArrayComparator<VectorizedArrayType>::width;
++k)
if (t1[d][e][k] < t2[d][e][k] - tolerance)
return true;
else if (t1[d][e][k] > t2[d][e][k] + tolerance)
Expand All @@ -3436,12 +3469,16 @@ namespace internal



template <typename Number, typename VectorizedArrayType>
template <typename VectorizedArrayType>
template <int dim>
bool
FPArrayComparator<Number, VectorizedArrayType>::operator()(
const std::array<Tensor<2, dim, Number>, dim + 1> &t1,
const std::array<Tensor<2, dim, Number>, dim + 1> &t2) const
FPArrayComparator<VectorizedArrayType>::operator()(
const std::array<
Tensor<2, dim, typename FPArrayComparator<VectorizedArrayType>::Number>,
dim + 1> &t1,
const std::array<
Tensor<2, dim, typename FPArrayComparator<VectorizedArrayType>::Number>,
dim + 1> &t2) const
{
for (unsigned int i = 0; i < t1.size(); ++i)
for (unsigned int d = 0; d < dim; ++d)
Expand Down
35 changes: 19 additions & 16 deletions source/matrix_free/mapping_info.cc
Original file line number Diff line number Diff line change
Expand Up @@ -31,31 +31,34 @@ DEAL_II_NAMESPACE_OPEN

#if SPLIT_INSTANTIATIONS_INDEX == 0

template struct internal::MatrixFreeFunctions::
FPArrayComparator<double, VectorizedArray<double, 1>>;
template struct internal::MatrixFreeFunctions::
FPArrayComparator<float, VectorizedArray<float, 1>>;
template struct internal::MatrixFreeFunctions::FPArrayComparator<double>;
template struct internal::MatrixFreeFunctions::FPArrayComparator<float>;

template struct internal::MatrixFreeFunctions::FPArrayComparator<
VectorizedArray<double, 1>>;
template struct internal::MatrixFreeFunctions::FPArrayComparator<
VectorizedArray<float, 1>>;

# if (DEAL_II_VECTORIZATION_WIDTH_IN_BITS >= 128 && defined(__SSE2__)) || \
(DEAL_II_VECTORIZATION_WIDTH_IN_BITS >= 128 && defined(__ALTIVEC__))
template struct internal::MatrixFreeFunctions::
FPArrayComparator<double, VectorizedArray<double, 2>>;
template struct internal::MatrixFreeFunctions::
FPArrayComparator<float, VectorizedArray<float, 4>>;
template struct internal::MatrixFreeFunctions::FPArrayComparator<
VectorizedArray<double, 2>>;
template struct internal::MatrixFreeFunctions::FPArrayComparator<
VectorizedArray<float, 4>>;
# endif

# if DEAL_II_VECTORIZATION_WIDTH_IN_BITS >= 256 && defined(__AVX__)
template struct internal::MatrixFreeFunctions::
FPArrayComparator<double, VectorizedArray<double, 4>>;
template struct internal::MatrixFreeFunctions::
FPArrayComparator<float, VectorizedArray<float, 8>>;
template struct internal::MatrixFreeFunctions::FPArrayComparator<
VectorizedArray<double, 4>>;
template struct internal::MatrixFreeFunctions::FPArrayComparator<
VectorizedArray<float, 8>>;
# endif

# if DEAL_II_VECTORIZATION_WIDTH_IN_BITS >= 512 && defined(__AVX512F__)
template struct internal::MatrixFreeFunctions::
FPArrayComparator<double, VectorizedArray<double, 8>>;
template struct internal::MatrixFreeFunctions::
FPArrayComparator<float, VectorizedArray<float, 16>>;
template struct internal::MatrixFreeFunctions::FPArrayComparator<
VectorizedArray<double, 8>>;
template struct internal::MatrixFreeFunctions::FPArrayComparator<
VectorizedArray<float, 16>>;
# endif

#endif
Expand Down

0 comments on commit a22baf0

Please sign in to comment.