File indexing completed on 2025-07-12 07:53:36
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include <boost/test/unit_test.hpp>
0010
0011 #include <Acts/Plugins/FastJet/TrackJets.hpp>
0012
0013 class Track {
0014 public:
0015 static constexpr float mass = 139.57061 * Acts::UnitConstants::MeV;
0016
0017 Track(float pt, float eta, float phi) {
0018 Acts::Vector3 p3 = Acts::Vector3::Zero();
0019 p3[0] = pt * std::cos(phi);
0020 p3[1] = pt * std::sin(phi);
0021 p3[2] = pt * std::sinh(eta);
0022 float e = std::sqrt(mass * mass + p3.squaredNorm());
0023 m_fourMom[0] = p3[0];
0024 m_fourMom[1] = p3[1];
0025 m_fourMom[2] = p3[2];
0026 m_fourMom[3] = e;
0027 }
0028
0029 Acts::Vector4 fourMomentum() const { return m_fourMom; }
0030
0031 private:
0032 Acts::Vector4 m_fourMom{};
0033 };
0034
0035 bool operator==(Track const& lhs, Track const& rhs) {
0036 return lhs.fourMomentum() == rhs.fourMomentum();
0037 }
0038
0039 class TrackContainer {
0040 public:
0041 using TrackProxy = Track;
0042
0043 TrackContainer() = default;
0044 void insert(Track track) { m_vec.push_back(std::move(track)); }
0045 std::size_t size() { return m_vec.size(); }
0046
0047 using ConstTrackProxy = const Track&;
0048 ConstTrackProxy getTrack(std::size_t i) {
0049 if (i < size()) {
0050 return m_vec[i];
0051 }
0052 throw std::runtime_error("Too few tracks");
0053 }
0054
0055 private:
0056 std::vector<Track> m_vec{};
0057 };
0058
0059 BOOST_AUTO_TEST_CASE(SingleTrack) {
0060 TrackContainer tracks;
0061 tracks.insert(Track(100, 0, 0));
0062
0063 Acts::FastJet::InputTracks inputTracks(tracks);
0064
0065 Acts::FastJet::TrackJetSequence jetSeq =
0066 Acts::FastJet::TrackJetSequence::create(inputTracks.fourMomenta());
0067 std::vector<fastjet::PseudoJet> jets = jetSeq.jets();
0068
0069 BOOST_CHECK_EQUAL(jets.size(), 1);
0070 BOOST_CHECK_EQUAL(jets[0].constituents().size(), 1);
0071 BOOST_CHECK_EQUAL(jets[0].constituents()[0].user_index(), 0);
0072 BOOST_CHECK_CLOSE(jets[0].pt(), 100, 1e-3);
0073 BOOST_CHECK_CLOSE(jets[0].eta(), 0, 1e-3);
0074 BOOST_CHECK_CLOSE(jets[0].phi(), 0, 1e-3);
0075 BOOST_CHECK_CLOSE(jets[0].m(), Track::mass, 1);
0076 }
0077
0078 BOOST_AUTO_TEST_CASE(TwoTracksTwoJets) {
0079 TrackContainer tracks;
0080 tracks.insert(Track(100, 0, 0.0));
0081 tracks.insert(Track(100, 0, std::numbers::pi));
0082
0083 Acts::FastJet::InputTracks inputTracks(tracks);
0084
0085 Acts::FastJet::TrackJetSequence jetSeq =
0086 Acts::FastJet::TrackJetSequence::create(inputTracks.fourMomenta());
0087 std::vector<fastjet::PseudoJet> jets = jetSeq.jets();
0088
0089 BOOST_CHECK_EQUAL(jets.size(), 2);
0090
0091 std::vector<Track> trks_0 = inputTracks.tracksInJet(jets[0]);
0092 BOOST_CHECK_EQUAL(trks_0.size(), 1);
0093 BOOST_CHECK(trks_0[0] == tracks.getTrack(0) ||
0094 trks_0[0] == tracks.getTrack(1));
0095
0096 std::vector<Track> trks_1 = inputTracks.tracksInJet(jets[1]);
0097 BOOST_CHECK_EQUAL(trks_1.size(), 1);
0098 BOOST_CHECK(trks_1[0] == tracks.getTrack(0) ||
0099 trks_1[0] == tracks.getTrack(1));
0100 BOOST_CHECK(trks_0[0] != trks_1[0]);
0101 }
0102
0103 BOOST_AUTO_TEST_CASE(TwoTracksOneJet) {
0104 TrackContainer tracks;
0105 tracks.insert(Track(100, 0, 0.0));
0106 tracks.insert(Track(100, 0, 0.2));
0107
0108 Acts::FastJet::InputTracks inputTracks(tracks);
0109
0110 Acts::FastJet::TrackJetSequence jetSeq =
0111 Acts::FastJet::TrackJetSequence::create(inputTracks.fourMomenta());
0112 std::vector<fastjet::PseudoJet> jets = jetSeq.jets();
0113
0114 BOOST_CHECK_EQUAL(jets.size(), 1);
0115
0116 std::vector<Track> trks_0 = inputTracks.tracksInJet(jets[0]);
0117 BOOST_CHECK_EQUAL(trks_0.size(), 2);
0118 BOOST_CHECK(trks_0[0] == tracks.getTrack(0) ||
0119 trks_0[0] == tracks.getTrack(1));
0120 BOOST_CHECK(trks_0[1] == tracks.getTrack(0) ||
0121 trks_0[1] == tracks.getTrack(1));
0122 BOOST_CHECK(trks_0[0] != trks_0[1]);
0123 }
0124
0125 BOOST_AUTO_TEST_CASE(TracksInJetCore) {
0126 TrackContainer tracks;
0127 tracks.insert(Track(100, 0, 0));
0128 tracks.insert(Track(10, 0.05, 0));
0129 tracks.insert(Track(10, -0.05, 0));
0130 tracks.insert(Track(10, 0.2, 0));
0131 tracks.insert(Track(10, -0.2, 0));
0132
0133 Acts::FastJet::InputTracks inputTracks(tracks);
0134
0135 Acts::FastJet::TrackJetSequence jetSeq =
0136 Acts::FastJet::TrackJetSequence::create(inputTracks.fourMomenta());
0137 std::vector<fastjet::PseudoJet> jets = jetSeq.jets();
0138
0139 BOOST_REQUIRE_EQUAL(jets.size(), 1);
0140
0141 std::vector<Track> trks = inputTracks.tracksInJet(jets[0], 0.1);
0142 BOOST_CHECK_EQUAL(trks.size(), 3);
0143
0144 BOOST_CHECK(std::find(trks.begin(), trks.end(), tracks.getTrack(0)) !=
0145 trks.end());
0146 BOOST_CHECK(std::find(trks.begin(), trks.end(), tracks.getTrack(1)) !=
0147 trks.end());
0148 BOOST_CHECK(std::find(trks.begin(), trks.end(), tracks.getTrack(2)) !=
0149 trks.end());
0150 BOOST_CHECK(std::find(trks.begin(), trks.end(), tracks.getTrack(3)) ==
0151 trks.end());
0152 BOOST_CHECK(std::find(trks.begin(), trks.end(), tracks.getTrack(4)) ==
0153 trks.end());
0154 }
0155
0156 BOOST_AUTO_TEST_CASE(EmptyTrackContainer) {
0157 Acts::FastJet::TrackJetSequence jetSeq =
0158 Acts::FastJet::TrackJetSequence::create(
0159 std::vector<fastjet::PseudoJet>());
0160 BOOST_CHECK_EQUAL(jetSeq.jets().size(), 0);
0161 }
0162
0163 BOOST_AUTO_TEST_CASE(InvalidCoreRadius) {
0164 TrackContainer tracks;
0165 tracks.insert(Track(100, 0, 0));
0166 Acts::FastJet::InputTracks inputTracks(tracks);
0167 Acts::FastJet::TrackJetSequence jetSeq =
0168 Acts::FastJet::TrackJetSequence::create(inputTracks.fourMomenta());
0169 BOOST_CHECK_THROW(inputTracks.tracksInJet(jetSeq.jets()[0], -1.0),
0170 std::invalid_argument);
0171 }