forked from cms-sw/cmssw
-
Notifications
You must be signed in to change notification settings - Fork 1
Commit
This commit does not belong to any branch on this repository, and may belong to a fork outside of the repository.
Port alpaka phi functions to DataFormats
- Loading branch information
Showing
1 changed file
with
45 additions
and
0 deletions.
There are no files selected for viewing
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,45 @@ | ||
#ifndef DataFormats_PortableMath_deltaPhi_h | ||
#define DataFormats_PortableMath_deltaPhi_h | ||
|
||
/* functions to compute heterogeneous deltaPhi | ||
* Ported from the LST algorithm: | ||
* https://github.com/cms-sw/cmssw/pull/45117#discussion_r1749143439 | ||
*/ | ||
|
||
namespace cms::alpakatools { | ||
|
||
// reduce to [-pi,pi] | ||
template <typename TAcc, typename T> | ||
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr T reduceRange(TAcc const& acc, T x) { | ||
constexpr T o2pi = 1. / (2. * M_PI); | ||
if (alpaka::math::abs(acc, x) <= T(M_PI)) | ||
return x; | ||
T n = alpaka::math::round(acc, x * o2pi); | ||
return x - n * T(2. * M_PI); | ||
} | ||
|
||
template <typename TAcc, typename T> | ||
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr T phi(TAcc const& acc, T x, T y) { | ||
return reduceRange(acc, M_PI + alpaka::math::atan2(acc, -y, -x)); | ||
} | ||
|
||
template <typename TAcc, typename T> | ||
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr T deltaPhi(TAcc const& acc, T x1, T y1, T x2, T y2) { | ||
T phi1 = phi(acc, x1, y1); | ||
T phi2 = phi(acc, x2, y2); | ||
return reduceRange(acc, phi2 - phi1); | ||
} | ||
|
||
template <typename TAcc, typename T> | ||
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr T deltaPhi(T phi1, T phi2) { | ||
return reduceRange(phi1 - phi2); | ||
} | ||
|
||
template <typename TAcc, typename T> | ||
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr T deltaPhiChange(TAcc const& acc, T x1, T y1, T x2, T y2) { | ||
return deltaPhi(acc, x1, y1, x2 - x1, y2 - y1); | ||
} | ||
|
||
} // namespace cms::alpakatools | ||
|
||
#endif |