Skip to content

Commit

Permalink
Merge Pull Request #13428 from mperego/Trilinos/intrepid2_fixinclusion
Browse files Browse the repository at this point in the history
Automatically Merged using Trilinos Pull Request AutoTester
PR Title: b'Intrepid2: Fix inclusion test'
PR Author: mperego
  • Loading branch information
trilinos-autotester authored Sep 6, 2024
2 parents 446b5b0 + 13092ec commit 869d1df
Show file tree
Hide file tree
Showing 15 changed files with 164 additions and 156 deletions.
31 changes: 22 additions & 9 deletions packages/intrepid2/src/Cell/Intrepid2_CellDataDef.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -831,6 +831,7 @@ refCenterDataStatic_ = {
bool
PointInclusion<shards::Line<>::key>::
check(const PointViewType &point, const ScalarType threshold) {
//this implementation should work when PointType is a Sacado Fad<ScalarType>
const ScalarType minus_one = -1.0 - threshold, plus_one = 1.0 + threshold;
return (minus_one <= point(0) && point(0) <= plus_one);
}
Expand All @@ -840,16 +841,19 @@ refCenterDataStatic_ = {
bool
PointInclusion<shards::Triangle<>::key>::
check(const PointViewType &point, const ScalarType threshold) {
const ScalarType distance = max( max( -point(0), -point(1) ), ScalarType( point(0) + point(1) - 1.0 ) );
//this implementation should work when PointType is a Sacado Fad<ScalarType>
using PointType = typename PointViewType::value_type;
const PointType one = 1.0;
const PointType distance = max( max( -point(0), -point(1) ), point(0) + point(1) - one );
return distance < threshold;
}

template<typename PointViewType, typename ScalarType>
KOKKOS_INLINE_FUNCTION
bool
PointInclusion<shards::Quadrilateral<>::key>::
check(const PointViewType &point,
const ScalarType threshold) {
check(const PointViewType &point, const ScalarType threshold) {
//this implementation should work when PointType is a Sacado Fad<ScalarType>
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));
Expand All @@ -860,8 +864,11 @@ refCenterDataStatic_ = {
bool
PointInclusion<shards::Tetrahedron<>::key>::
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) );
//this implementation should work when PointType is a Sacado Fad<ScalarType>
using PointType = typename PointViewType::value_type;
const PointType one = 1.0;
const PointType distance = max( max(-point(0),-point(1)),
max(-point(2), point(0) + point(1) + point(2) - one) );
return distance < threshold;
}

Expand All @@ -870,6 +877,7 @@ refCenterDataStatic_ = {
bool
PointInclusion<shards::Hexahedron<>::key>::
check(const PointViewType &point, const ScalarType threshold) {
//this implementation should work when PointType is a Sacado Fad<ScalarType>
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) &&
Expand All @@ -881,9 +889,11 @@ refCenterDataStatic_ = {
bool
PointInclusion<shards::Pyramid<>::key>::
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);
//this implementation should work when PointType is a Sacado Fad<ScalarType>
using PointType = typename PointViewType::value_type;
const PointType minus_one = -1.0 - threshold, plus_one = 1.0 + threshold, minus_zero = -threshold;
const PointType left = minus_one + point(2);
const PointType 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));
Expand All @@ -894,8 +904,11 @@ refCenterDataStatic_ = {
bool
PointInclusion<shards::Wedge<>::key>::
check(const PointViewType &point, const ScalarType threshold) {
//this implementation should work when PointType is a Sacado Fad<ScalarType>
using PointType = typename PointViewType::value_type;
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 );
const PointType one = 1.0;
const PointType distance = max( max( -point(0), -point(1) ), point(0) + point(1) - one );
return (distance < threshold && (minus_one <= point(2) && point(2) <= plus_one));
}

Expand Down
12 changes: 5 additions & 7 deletions packages/intrepid2/src/Cell/Intrepid2_CellToolsDefInclusion.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -184,6 +184,10 @@ namespace Intrepid2 {
checkPointwiseInclusion<shards::Quadrilateral<>::key,decltype(inCell),decltype(points)>(inCell, points, threshold);
break;

case shards::Tetrahedron<>::key :
checkPointwiseInclusion<shards::Tetrahedron<>::key,decltype(inCell),decltype(points)>(inCell, points, threshold);
break;

case shards::Hexahedron<>::key :
checkPointwiseInclusion<shards::Hexahedron<>::key,decltype(inCell),decltype(points)>(inCell, points, threshold);
break;
Expand All @@ -197,13 +201,7 @@ namespace Intrepid2 {
break;

default:
INTREPID2_TEST_FOR_EXCEPTION( !( (key == shards::Line<>::key ) ||
(key == shards::Triangle<>::key) ||
(key == shards::Quadrilateral<>::key) ||
(key == shards::Tetrahedron<>::key) ||
(key == shards::Hexahedron<>::key) ||
(key == shards::Wedge<>::key) ||
(key == shards::Pyramid<>::key) ),
INTREPID2_TEST_FOR_EXCEPTION( false,
std::invalid_argument,
">>> ERROR (Intrepid2::CellTools::checkPointInclusion): Invalid cell topology. ");
}
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -18,7 +18,7 @@ namespace Intrepid2 {
template <typename DT, typename PT, typename WT>
CubatureDirectLineGaussJacobi20<DT,PT,WT>::
CubatureDirectLineGaussJacobi20(const ordinal_type degree)
: CubatureDirect<DT>(degree, 1) {
: CubatureDirect<DT,PT,WT>(degree, 1) {

INTREPID2_TEST_FOR_EXCEPTION( degree < 0 ||
degree > static_cast<ordinal_type>(Parameters::MaxCubatureDegreePyr), std::out_of_range,
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -18,7 +18,7 @@ namespace Intrepid2 {
template <typename DT, typename PT, typename WT>
CubatureDirectTetDefault<DT,PT,WT>::
CubatureDirectTetDefault(const ordinal_type degree)
: CubatureDirect<DT>(degree, 3) {
: CubatureDirect<DT,PT,WT>(degree, 3) {

INTREPID2_TEST_FOR_EXCEPTION( degree < 0 ||
degree > static_cast<ordinal_type>(Parameters::MaxCubatureDegreeTet), std::out_of_range,
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -18,7 +18,7 @@ namespace Intrepid2 {
template <typename DT, typename PT, typename WT>
CubatureDirectTetSymmetric<DT,PT,WT>::
CubatureDirectTetSymmetric(const ordinal_type degree)
: CubatureDirect<DT>(degree, 3) {
: CubatureDirect<DT,PT,WT>(degree, 3) {

INTREPID2_TEST_FOR_EXCEPTION( degree < 0 ||
degree > static_cast<ordinal_type>(Parameters::MaxCubatureDegreeTet), std::out_of_range,
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -18,7 +18,7 @@ namespace Intrepid2 {
template <typename DT, typename PT, typename WT>
CubatureDirectTriDefault<DT,PT,WT>::
CubatureDirectTriDefault(const ordinal_type degree)
: CubatureDirect<DT>(degree, 2) {
: CubatureDirect<DT,PT,WT>(degree, 2) {

INTREPID2_TEST_FOR_EXCEPTION( degree < 0 ||
degree >= cubatureDataStaticSize, std::out_of_range,
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -18,7 +18,7 @@ namespace Intrepid2 {
template <typename DT, typename PT, typename WT>
CubatureDirectTriSymmetric<DT,PT,WT>::
CubatureDirectTriSymmetric(const ordinal_type degree)
: CubatureDirect<DT>(degree, 2) {
: CubatureDirect<DT,PT,WT>(degree, 2) {

INTREPID2_TEST_FOR_EXCEPTION( degree < 0 ||
degree >= cubatureDataStaticSize, std::out_of_range,
Expand Down
2 changes: 1 addition & 1 deletion packages/intrepid2/src/Shared/Intrepid2_Types.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -56,7 +56,7 @@ namespace Intrepid2 {

flt_32 s;
s.f32 = 1;
s.f32++;
s.i32++;
return (s.i32 < 0 ? 1 - s.f32 : s.f32 - 1);
}

Expand Down
4 changes: 4 additions & 0 deletions packages/intrepid2/unit-test/Cell/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -8,6 +8,10 @@ LIST(APPEND Intrepid2_TEST_ETI_VALUETYPE_NAME "DOUBLE")
LIST(APPEND Intrepid2_TEST_ETI_VALUETYPE "double")
LIST(APPEND Intrepid2_TEST_ETI_SACADO "0")

LIST(APPEND Intrepid2_TEST_ETI_VALUETYPE_NAME "FLOAT")
LIST(APPEND Intrepid2_TEST_ETI_VALUETYPE "float")
LIST(APPEND Intrepid2_TEST_ETI_SACADO "0")

IF (HAVE_INTREPID2_SACADO)
# LIST(APPEND Intrepid2_TEST_ETI_VALUETYPE_NAME "DFAD_DOUBLE")
# LIST(APPEND Intrepid2_TEST_ETI_VALUETYPE "Sacado::Fad::DFad<double>")
Expand Down
2 changes: 1 addition & 1 deletion packages/intrepid2/unit-test/Cell/test_01.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -200,7 +200,7 @@ namespace Intrepid2 {
typedef CellTools<DeviceType> ct;
typedef Kokkos::DynRankView<value_type,DeviceType> DynRankView;

const value_type tol = tolerence()*100.0;
const value_type tol = tolerence<value_type>()*100.0;

int errorFlag = 0;

Expand Down
4 changes: 2 additions & 2 deletions packages/intrepid2/unit-test/Cell/test_02.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -81,7 +81,7 @@ namespace Intrepid2 {
using ct = CellTools<DeviceType>;
using DynRankView = Kokkos::DynRankView<ValueType,DeviceType>;

const ValueType tol = tolerence()*100.0;
const ValueType tol = tolerence<ValueType>()*100.0;

int errorFlag = 0;

Expand Down Expand Up @@ -112,7 +112,7 @@ namespace Intrepid2 {
{
// Define cubature on the edge parametrization domain:
const auto testAccuracy = 6;
CubatureDirectLineGauss<DeviceType> edgeCubature(testAccuracy);
CubatureDirectLineGauss<DeviceType,ValueType,ValueType> edgeCubature(testAccuracy);

const auto cubDim = edgeCubature.getDimension();
const auto numCubPoints = edgeCubature.getNumPoints();
Expand Down
2 changes: 1 addition & 1 deletion packages/intrepid2/unit-test/Cell/test_03.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -85,7 +85,7 @@ namespace Intrepid2 {
using DynRankView = Kokkos::DynRankView<ValueType,DeviceType>;
using DynRankViewHost = Kokkos::DynRankView<ValueType,Kokkos::HostSpace>;

const ValueType tol = tolerence()*100.0;
const ValueType tol = tolerence<ValueType>()*100.0;

int errorFlag = 0;

Expand Down
4 changes: 2 additions & 2 deletions packages/intrepid2/unit-test/Cell/test_04.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -79,7 +79,7 @@ namespace Intrepid2 {
using ct = CellTools<DeviceType>;
using DynRankView = Kokkos::DynRankView<ValueType,DeviceType>;

const ValueType tol = tolerence()*100.0;
const ValueType tol = tolerence<ValueType>()*100.0;

int errorFlag = 0;

Expand Down Expand Up @@ -132,7 +132,7 @@ namespace Intrepid2 {
continue;

// 1. Define a single reference point set using cubature factory with order 4 cubature
const auto cellCubature = cubFactory.create<DeviceType>(cell, testAccuracy);
const auto cellCubature = cubFactory.create<DeviceType,ValueType,ValueType>(cell, testAccuracy);
const auto cubDim = cellCubature->getDimension();
const auto cubNumPoints = cellCubature->getNumPoints();

Expand Down
2 changes: 1 addition & 1 deletion packages/intrepid2/unit-test/Cell/test_06.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -130,7 +130,7 @@ namespace Intrepid2 {
<< "| |\n"
<< "===============================================================================\n";

const ValueType tol = tolerence()*100.0;
const ValueType tol = tolerence<ValueType>()*100.0;

int errorFlag = 0;

Expand Down
Loading

0 comments on commit 869d1df

Please sign in to comment.