Skip to content

Commit

Permalink
Intrepid2: fix C++20 deprecated code and potential type mismatch issue (
Browse files Browse the repository at this point in the history
#13394)

* Intrepid2: fix Inclusion checks for non-double point types

* Intrepid2: changes to avoid warnings with c++20 standard (fixes #12786)
  • Loading branch information
mperego authored Aug 27, 2024
1 parent 4c55a25 commit 3132245
Show file tree
Hide file tree
Showing 9 changed files with 92 additions and 74 deletions.
28 changes: 14 additions & 14 deletions packages/intrepid2/src/Cell/Intrepid2_CellData.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -337,21 +337,21 @@ template<unsigned CellTopologyKey>
*/
template<>
struct PointInclusion<shards::Line<>::key> {
template<typename PointViewType>
template<typename PointViewType, typename ScalarType>
KOKKOS_INLINE_FUNCTION
static bool
check(const PointViewType &point, const double threshold);
check(const PointViewType &point, const ScalarType threshold);
};

/**
\brief Triangle topology
*/
template<>
struct PointInclusion<shards::Triangle<>::key> {
template<typename PointViewType>
template<typename PointViewType, typename ScalarType>
KOKKOS_INLINE_FUNCTION
static bool
check(const PointViewType &point, const double threshold);
check(const PointViewType &point, const ScalarType threshold);
};

/**
Expand All @@ -360,54 +360,54 @@ template<unsigned CellTopologyKey>
template<>
struct PointInclusion<shards::Quadrilateral<>::key> {

template<typename PointViewType>
template<typename PointViewType, typename ScalarType>
KOKKOS_INLINE_FUNCTION
static bool
check(const PointViewType &point, const double threshold);
check(const PointViewType &point, const ScalarType threshold);
};

/**
\brief Tetrahedron topology
*/
template<>
struct PointInclusion<shards::Tetrahedron<>::key> {
template<typename PointViewType>
template<typename PointViewType, typename ScalarType>
KOKKOS_INLINE_FUNCTION
static bool
check(const PointViewType &point, const double threshold);
check(const PointViewType &point, const ScalarType threshold);
};

/**
\brief Hexahedron topology
*/
template<>
struct PointInclusion<shards::Hexahedron<>::key> {
template<typename PointViewType>
template<typename PointViewType, typename ScalarType>
KOKKOS_INLINE_FUNCTION
static bool
check(const PointViewType &point, const double threshold);
check(const PointViewType &point, const ScalarType threshold);
};

/**
\brief Pyramid topology
*/
template<>
struct PointInclusion<shards::Pyramid<>::key> {
template<typename PointViewType>
template<typename PointViewType, typename ScalarType>
KOKKOS_INLINE_FUNCTION
static bool
check(const PointViewType &point, const double threshold);
check(const PointViewType &point, const ScalarType threshold);
};

/**
\brief Wedge topology
*/
template<>
struct PointInclusion<shards::Wedge<>::key> {
template<typename PointViewType>
template<typename PointViewType, typename ScalarType>
KOKKOS_INLINE_FUNCTION
static bool
check(const PointViewType &point, const double threshold);
check(const PointViewType &point, const ScalarType threshold);
};

}
Expand Down
48 changes: 24 additions & 24 deletions packages/intrepid2/src/Cell/Intrepid2_CellDataDef.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -826,76 +826,76 @@ refCenterDataStatic_ = {
// Point Inclusion


template<typename PointViewType>
template<typename PointViewType, typename ScalarType>
KOKKOS_INLINE_FUNCTION
bool
PointInclusion<shards::Line<>::key>::
check(const PointViewType &point, const double threshold) {
const double minus_one = -1.0 - threshold, plus_one = 1.0 + threshold;
check(const PointViewType &point, const ScalarType threshold) {
const ScalarType minus_one = -1.0 - threshold, plus_one = 1.0 + threshold;
return (minus_one <= point(0) && point(0) <= plus_one);
}

template<typename PointViewType>
template<typename PointViewType, typename ScalarType>
KOKKOS_INLINE_FUNCTION
bool
PointInclusion<shards::Triangle<>::key>::
check(const PointViewType &point, const double threshold) {
const double distance = max( max( -point(0), -point(1) ), point(0) + point(1) - 1.0 );
check(const PointViewType &point, const ScalarType threshold) {
const ScalarType distance = max( max( -point(0), -point(1) ), point(0) + point(1) - 1.0 );
return distance < threshold;
}

template<typename PointViewType>
template<typename PointViewType, typename ScalarType>
KOKKOS_INLINE_FUNCTION
bool
PointInclusion<shards::Quadrilateral<>::key>::
check(const PointViewType &point,
const double threshold) {
const double minus_one = -1.0 - threshold, plus_one = 1.0 + threshold;
const ScalarType threshold) {
const ScalarType minus_one = -1.0 - threshold, plus_one = 1.0 + threshold;
return ((minus_one <= point(0) && point(0) <= plus_one) &&
(minus_one <= point(1) && point(1) <= plus_one));
}

template<typename PointViewType>
template<typename PointViewType, typename ScalarType>
KOKKOS_INLINE_FUNCTION
bool
PointInclusion<shards::Tetrahedron<>::key>::
check(const PointViewType &point, const double threshold) {
const double distance = max( max(-point(0),-point(1)),
check(const PointViewType &point, const ScalarType threshold) {
const ScalarType distance = max( max(-point(0),-point(1)),
max(-point(2), point(0) + point(1) + point(2) - 1) );
return distance < threshold;
}

template<typename PointViewType>
template<typename PointViewType, typename ScalarType>
KOKKOS_INLINE_FUNCTION
bool
PointInclusion<shards::Hexahedron<>::key>::
check(const PointViewType &point, const double threshold) {
const double minus_one = -1.0 - threshold, plus_one = 1.0 + threshold;
check(const PointViewType &point, const ScalarType threshold) {
const ScalarType minus_one = -1.0 - threshold, plus_one = 1.0 + threshold;
return ((minus_one <= point(0) && point(0) <= plus_one) &&
(minus_one <= point(1) && point(1) <= plus_one) &&
(minus_one <= point(2) && point(2) <= plus_one));
}

template<typename PointViewType>
template<typename PointViewType, typename ScalarType>
KOKKOS_INLINE_FUNCTION
bool
PointInclusion<shards::Pyramid<>::key>::
check(const PointViewType &point, const double threshold) {
const double minus_one = -1.0 - threshold, plus_one = 1.0 + threshold, minus_zero = -threshold;
const double left = minus_one + point(2);
const double right = plus_one - point(2);
check(const PointViewType &point, const ScalarType threshold) {
const ScalarType minus_one = -1.0 - threshold, plus_one = 1.0 + threshold, minus_zero = -threshold;
const ScalarType left = minus_one + point(2);
const ScalarType right = plus_one - point(2);
return ((left <= point(0) && point(0) <= right) &&
(left <= point(1) && point(1) <= right) &&
(minus_zero <= point(2) && point(2) <= plus_one));
}

template<typename PointViewType>
template<typename PointViewType, typename ScalarType>
KOKKOS_INLINE_FUNCTION
bool
PointInclusion<shards::Wedge<>::key>::
check(const PointViewType &point, const double threshold) {
const double minus_one = -1.0 - threshold, plus_one = 1.0 + threshold;
const double distance = max( max( -point(0), -point(1) ), point(0) + point(1) - 1 );
check(const PointViewType &point, const ScalarType threshold) {
const ScalarType minus_one = -1.0 - threshold, plus_one = 1.0 + threshold;
const ScalarType distance = max( max( -point(0), -point(1) ), point(0) + point(1) - 1 );
return (distance < threshold && (minus_one <= point(2) && point(2) <= plus_one));
}

Expand Down
17 changes: 11 additions & 6 deletions packages/intrepid2/src/Cell/Intrepid2_CellTools.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -1396,11 +1396,13 @@ namespace Intrepid2 {
\param threshold [in] - "tightness" of the inclusion test
\return true if the point is in the closure of the specified reference cell and false otherwise.
*/
template<typename pointViewType>
template<typename PointViewType>
static bool
checkPointInclusion( const pointViewType point,
checkPointInclusion( const PointViewType point,
const shards::CellTopology cellTopo,
const double thres = threshold() );
const typename ScalarTraits<typename PointViewType::value_type>::scalar_type thres =
threshold<typename ScalarTraits<typename PointViewType::value_type>::scalar_type>() );



/** \brief Checks every point for inclusion in the reference cell of a given topology.
Expand All @@ -1417,7 +1419,8 @@ namespace Intrepid2 {
typename InputViewType>
static void checkPointwiseInclusion( OutputViewType inCell,
const InputViewType points,
const double thresh = threshold());
const typename ScalarTraits<typename InputViewType::value_type>::scalar_type thresh =
threshold<typename ScalarTraits<typename InputViewType::value_type>::scalar_type>());



Expand All @@ -1434,7 +1437,8 @@ namespace Intrepid2 {
static void checkPointwiseInclusion( InCellViewType inCell,
const PointViewType points,
const shards::CellTopology cellTopo,
const double thres = threshold() );
const typename ScalarTraits<typename PointViewType::value_type>::scalar_type thres =
threshold<typename ScalarTraits<typename PointViewType::value_type>::scalar_type>() );

/** \brief Checks every points for inclusion in physical cells from a cell workset.
The points can belong to a global set and stored in a rank-2 (P,D) view,
Expand All @@ -1454,7 +1458,8 @@ namespace Intrepid2 {
const Kokkos::DynRankView<pointValueType,pointProperties...> points,
const Kokkos::DynRankView<cellWorksetValueType,cellWorksetProperties...> cellWorkset,
const shards::CellTopology cellTopo,
const double thres = threshold() );
const typename ScalarTraits<pointValueType>::scalar_type thres =
threshold<typename ScalarTraits<pointValueType>::scalar_type>() );


// //============================================================================================//
Expand Down
29 changes: 15 additions & 14 deletions packages/intrepid2/src/Cell/Intrepid2_CellToolsDefInclusion.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -34,9 +34,9 @@ namespace Intrepid2 {
template<typename PointViewType>
bool
CellTools<DeviceType>::
checkPointInclusion( const PointViewType point,
const shards::CellTopology cellTopo,
const double threshold) {
checkPointInclusion( const PointViewType point,
const shards::CellTopology cellTopo,
const typename ScalarTraits<typename PointViewType::value_type>::scalar_type threshold) {
#ifdef HAVE_INTREPID2_DEBUG
INTREPID2_TEST_FOR_EXCEPTION( point.rank() != 1, std::invalid_argument,
">>> ERROR (Intrepid2::CellTools::checkPointInclusion): Point must have rank 1. ");
Expand Down Expand Up @@ -94,12 +94,13 @@ namespace Intrepid2 {
struct checkPointInclusionFunctor {
OutputViewType output_;
InputViewType input_;
double threshold_;
using ScalarType = typename ScalarTraits<typename InputViewType::value_type>::scalar_type;
ScalarType threshold_;

KOKKOS_INLINE_FUNCTION
checkPointInclusionFunctor( OutputViewType output,
const InputViewType input,
const double threshold)
checkPointInclusionFunctor( OutputViewType output,
const InputViewType input,
const ScalarType threshold)
: output_(output),
input_(input),
threshold_(threshold) {}
Expand Down Expand Up @@ -129,7 +130,7 @@ namespace Intrepid2 {
void CellTools<DeviceType>::
checkPointwiseInclusion( OutputViewType inCell,
const InputViewType points,
const double threshold) {
const typename ScalarTraits<typename InputViewType::value_type>::scalar_type threshold) {

using FunctorType = checkPointInclusionFunctor<cellTopologyKey,decltype(inCell),decltype(points)>;
if (points.rank() == 2) { // inCell.rank() == 1
Expand All @@ -144,13 +145,13 @@ namespace Intrepid2 {

template<typename DeviceType>
template<typename InCellViewType,
typename PointViewType>
typename InputViewType>
void
CellTools<DeviceType>::
checkPointwiseInclusion( InCellViewType inCell,
const PointViewType points,
const shards::CellTopology cellTopo,
const double threshold ) {
checkPointwiseInclusion( InCellViewType inCell,
const InputViewType points,
const shards::CellTopology cellTopo,
const typename ScalarTraits<typename InputViewType::value_type>::scalar_type threshold ) {
#ifdef HAVE_INTREPID2_DEBUG
{
INTREPID2_TEST_FOR_EXCEPTION( (inCell.rank() != 1) && (inCell.rank() != 2), std::invalid_argument,
Expand Down Expand Up @@ -218,7 +219,7 @@ namespace Intrepid2 {
const Kokkos::DynRankView<pointValueType,pointProperties...> points,
const Kokkos::DynRankView<cellWorksetValueType,cellWorksetProperties...> cellWorkset,
const shards::CellTopology cellTopo,
const double threshold ) {
const typename ScalarTraits<pointValueType>::scalar_type threshold ) {
#ifdef HAVE_INTREPID2_DEBUG
{
const auto key = cellTopo.getBaseKey();
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -123,7 +123,7 @@ namespace Intrepid2 {
// prepare for allocation of temporary storage
// note: tempStorage goes "backward", starting from the final component, which needs just one entry

const bool allocateFadStorage = !std::is_pod<Scalar>::value;
const bool allocateFadStorage = !(std::is_standard_layout<Scalar>::value && std::is_trivial<Scalar>::value);
if (allocateFadStorage)
{
fad_size_output_ = dimension_scalar(integralView_);
Expand Down Expand Up @@ -1063,7 +1063,7 @@ namespace Intrepid2 {
// prepare for allocation of temporary storage
// note: tempStorage goes "backward", starting from the final component, which needs just one entry

const bool allocateFadStorage = !std::is_pod<Scalar>::value;
const bool allocateFadStorage = !(std::is_standard_layout<Scalar>::value && std::is_trivial<Scalar>::value);
if (allocateFadStorage)
{
fad_size_output_ = dimension_scalar(integralView_);
Expand Down Expand Up @@ -2127,7 +2127,7 @@ void IntegrationTools<DeviceType>::integrate(Data<Scalar,DeviceType> integrals,
{
ScalarView<Scalar,DeviceType> componentIntegralView;

const bool allocateFadStorage = !std::is_pod<Scalar>::value;
const bool allocateFadStorage = !(std::is_standard_layout<Scalar>::value && std::is_trivial<Scalar>::value);
if (allocateFadStorage)
{
auto fad_size_output = dimension_scalar(integrals.getUnderlyingView());
Expand Down
4 changes: 2 additions & 2 deletions packages/intrepid2/src/Shared/Intrepid2_TestUtils.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -203,7 +203,7 @@ namespace Intrepid2
template<typename ValueType, typename DeviceType, class ... DimArgs>
inline ViewType<ValueType,DeviceType> getView(const std::string &label, DimArgs... dims)
{
const bool allocateFadStorage = !std::is_pod<ValueType>::value;
const bool allocateFadStorage = !(std::is_standard_layout<ValueType>::value && std::is_trivial<ValueType>::value);
if (!allocateFadStorage)
{
return ViewType<ValueType,DeviceType>(label,dims...);
Expand All @@ -218,7 +218,7 @@ namespace Intrepid2
template<typename ValueType, class ... DimArgs>
inline FixedRankViewType< typename RankExpander<ValueType, sizeof...(DimArgs) >::value_type, DefaultTestDeviceType > getFixedRankView(const std::string &label, DimArgs... dims)
{
const bool allocateFadStorage = !std::is_pod<ValueType>::value;
const bool allocateFadStorage = !(std::is_standard_layout<ValueType>::value && std::is_trivial<ValueType>::value);
using value_type = typename RankExpander<ValueType, sizeof...(dims) >::value_type;
if (!allocateFadStorage)
{
Expand Down
16 changes: 14 additions & 2 deletions packages/intrepid2/src/Shared/Intrepid2_Types.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -65,14 +65,26 @@ namespace Intrepid2 {
return epsilon<double>();
}

template<typename ValueType>
KOKKOS_FORCEINLINE_FUNCTION
ValueType tolerence() {
return 100.0*epsilon<ValueType>();
}

KOKKOS_FORCEINLINE_FUNCTION
double tolerence() {
return 100.0*epsilon();
return tolerence<double>();
}

template<typename ValueType>
KOKKOS_FORCEINLINE_FUNCTION
ValueType threshold() {
return 10.0*epsilon<ValueType>();
}

KOKKOS_FORCEINLINE_FUNCTION
double threshold() {
return 10.0*epsilon();
return threshold<double>();
}

/// Define constants
Expand Down
Loading

0 comments on commit 3132245

Please sign in to comment.