CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
Phase2L1GMT::TPS Class Reference

#include <TPS.h>

Public Member Functions

std::vector< l1t::TrackerMuonprocessEvent (const std::vector< edm::Ptr< l1t::TrackerMuon::L1TTTrackType > > &, const l1t::MuonStubRefVector &)
 
 TPS (const edm::ParameterSet &iConfig)
 
 ~TPS ()=default
 

Private Member Functions

l1t::SAMuonRefVector associateMuonsWithNonant (const l1t::SAMuonRefVector &, uint)
 
l1t::MuonStubRefVector associateStubsWithNonant (const l1t::MuonStubRefVector &, uint)
 
std::vector< edm::Ptr< l1t::TrackerMuon::L1TTTrackType > > associateTracksWithNonant (const std::vector< edm::Ptr< l1t::TrackerMuon::L1TTTrackType > > &tracks, uint processor)
 

Private Attributes

std::unique_ptr< Isolationisolation_
 
std::unique_ptr< TPSAlgorithmtps_
 
std::unique_ptr< TrackConvertertt_track_converter_
 
int verbose_
 

Detailed Description

Definition at line 11 of file TPS.h.

Constructor & Destructor Documentation

◆ TPS()

TPS::TPS ( const edm::ParameterSet iConfig)

Definition at line 7 of file TPS.cc.

8  : verbose_(iConfig.getParameter<int>("verbose")),
9  tt_track_converter_(new TrackConverter(iConfig.getParameter<edm::ParameterSet>("trackConverter"))),
10  tps_(new TPSAlgorithm(iConfig.getParameter<edm::ParameterSet>("trackMatching"))),
11  isolation_(new Isolation(iConfig.getParameter<edm::ParameterSet>("isolation"))) {}
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
int verbose_
Definition: TPS.h:19
std::unique_ptr< Isolation > isolation_
Definition: TPS.h:22
std::unique_ptr< TrackConverter > tt_track_converter_
Definition: TPS.h:20
std::unique_ptr< TPSAlgorithm > tps_
Definition: TPS.h:21

◆ ~TPS()

Phase2L1GMT::TPS::~TPS ( )
default

Member Function Documentation

◆ associateMuonsWithNonant()

l1t::SAMuonRefVector TPS::associateMuonsWithNonant ( const l1t::SAMuonRefVector muons,
uint  processor 
)
private

Definition at line 88 of file TPS.cc.

References Phase2L1GMT::BITSPHI, SiPixelRawToDigiRegional_cfi::deltaPhi, DiMuonV_cfg::muons, MillePedeFileConverter_cfg::out, and alignCSCRings::s.

88  {
90 
91  ap_int<BITSPHI> center = ap_int<BITSPHI>(processor * 910);
92 
93  for (const auto& s : muons) {
94  ap_int<BITSSTUBCOORD> deltaPhi = s->hwPhi() - center;
95  ap_uint<BITSPHI - 1> absDeltaPhi =
96  (deltaPhi < 0) ? ap_uint<BITSPHI - 1>(-deltaPhi) : ap_uint<BITSPHI - 1>(deltaPhi);
97  if (absDeltaPhi < 683)
98  out.push_back(s);
99  }
100  return out;
101 }
std::vector< edm::Ref< SAMuonCollection > > SAMuonRefVector
Definition: SAMuon.h:19
muons
the two sets of parameters below are mutually exclusive, depending if RECO or ALCARECO is used the us...
Definition: DiMuonV_cfg.py:214
const int BITSPHI
Definition: Constants.h:25

◆ associateStubsWithNonant()

l1t::MuonStubRefVector TPS::associateStubsWithNonant ( const l1t::MuonStubRefVector allStubs,
uint  processor 
)
private

Definition at line 103 of file TPS.cc.

References Phase2L1GMT::BITSSTUBCOORD, SiPixelRawToDigiRegional_cfi::deltaPhi, MillePedeFileConverter_cfg::out, and alignCSCRings::s.

Referenced by processEvent().

103  {
105 
106  ap_int<BITSSTUBCOORD> center = ap_int<BITSSTUBCOORD>((processor * 910) / 8); //was 32
107 
108  for (const auto& s : allStubs) {
109  ap_int<BITSSTUBCOORD> phi = 0;
110  if (s->quality() & 0x1)
111  phi = s->coord1();
112  else
113  phi = s->coord2();
114 
115  ap_int<BITSSTUBCOORD> deltaPhi = phi - center;
116  ap_uint<BITSSTUBCOORD - 1> absDeltaPhi =
117  (deltaPhi < 0) ? ap_uint<BITSSTUBCOORD - 1>(-deltaPhi) : ap_uint<BITSSTUBCOORD - 1>(deltaPhi);
118  if (absDeltaPhi < 168) //was 42
119  out.push_back(s);
120  }
121  return out;
122 }
const int BITSSTUBCOORD
Definition: Constants.h:31
std::vector< edm::Ref< MuonStubCollection > > MuonStubRefVector
Definition: MuonStub.h:44

◆ associateTracksWithNonant()

std::vector< edm::Ptr< l1t::TrackerMuon::L1TTTrackType > > TPS::associateTracksWithNonant ( const std::vector< edm::Ptr< l1t::TrackerMuon::L1TTTrackType > > &  tracks,
uint  processor 
)
private

Definition at line 77 of file TPS.cc.

References MillePedeFileConverter_cfg::out, HLT_2024v13_cff::track, and DiMuonV_cfg::tracks.

Referenced by processEvent().

78  {
79  std::vector<edm::Ptr<l1t::TrackerMuon::L1TTTrackType> > out;
80  for (const auto& track : tracks) {
81  if (track->phiSector() == processor) {
82  out.push_back(track);
83  }
84  }
85  return out;
86 }

◆ processEvent()

std::vector< l1t::TrackerMuon > TPS::processEvent ( const std::vector< edm::Ptr< l1t::TrackerMuon::L1TTTrackType > > &  tracks,
const l1t::MuonStubRefVector muonStubs 
)

Definition at line 13 of file TPS.cc.

References associateStubsWithNonant(), associateTracksWithNonant(), l1tGTTFileWriter_cfi::convertedTracks, mps_fire::i, isolation_, SimL1Emulator_cff::stubs, tps_, DiMuonV_cfg::tracks, tt_track_converter_, and findQualityFiles::v.

14  {
15  //Split tracks to the links as they come
16  std::array<std::vector<edm::Ptr<l1t::TrackerMuon::L1TTTrackType> >, 9> loctracks;
17  for (unsigned i = 0; i < 9; ++i)
18  loctracks[i] = associateTracksWithNonant(tracks, i);
19 
20  //Convert TT tracks to our internal tracking format
21  std::array<std::vector<ConvertedTTTrack>, 9> convertedTracks;
22  for (unsigned i = 0; i < 9; ++i)
23  convertedTracks[i] = tt_track_converter_->convertTracks(loctracks.at(i));
24 
25  //Transition stubs to different nonants with overlap
26  std::array<l1t::MuonStubRefVector, 9> stubs;
27  for (int i = 0; i < 9; ++i)
28  stubs[i] = associateStubsWithNonant(muonStubs, i);
29 
30  //run track - muon matching per nonant
31  std::array<std::vector<PreTrackMatchedMuon>, 9> mus;
32  for (int i = 0; i < 9; ++i)
33  mus[i] = tps_->processNonant(convertedTracks.at(i), stubs.at(i));
34 
35  //clean neighboring nonants
36  std::vector<std::vector<PreTrackMatchedMuon> > muCleaneds;
37  muCleaneds.push_back(tps_->cleanNeighbor(mus[0], mus[8], mus[1], true));
38  muCleaneds.push_back(tps_->cleanNeighbor(mus[1], mus[0], mus[2], false));
39  muCleaneds.push_back(tps_->cleanNeighbor(mus[2], mus[1], mus[3], true));
40  muCleaneds.push_back(tps_->cleanNeighbor(mus[3], mus[2], mus[4], false));
41  muCleaneds.push_back(tps_->cleanNeighbor(mus[4], mus[3], mus[5], true));
42  muCleaneds.push_back(tps_->cleanNeighbor(mus[5], mus[4], mus[6], false));
43  muCleaneds.push_back(tps_->cleanNeighbor(mus[6], mus[5], mus[7], true));
44  muCleaneds.push_back(tps_->cleanNeighbor(mus[7], mus[6], mus[8], false));
45  //ARGH! 9 sectors - so some duplicates very rarely
46  muCleaneds.push_back(tps_->cleanNeighbor(mus[8], mus[7], mus[0], false));
47 
48  //merge all the collections
49  std::vector<PreTrackMatchedMuon> mergedCleaned;
50  for (auto&& v : muCleaneds) {
51  mergedCleaned.insert(mergedCleaned.end(), v.begin(), v.end());
52  }
53 
54  std::vector<l1t::TrackerMuon> trackMatchedMuonsNoIso = tps_->convert(mergedCleaned, 32);
55 
56  //sorter here:
57  std::vector<l1t::TrackerMuon> sortedTrackMuonsNoIso = tps_->sort(trackMatchedMuonsNoIso, 12);
58 
59  tps_->SetQualityBits(sortedTrackMuonsNoIso);
60 
61  //Isolation and tau3mu will read those muons and all 9 collections of convertedTracks*
62  std::vector<ConvertedTTTrack> mergedconvertedTracks;
63  for (auto&& v : convertedTracks) {
64  mergedconvertedTracks.insert(mergedconvertedTracks.end(), v.begin(), v.end());
65  }
66 
67  isolation_->isolation_allmu_alltrk(sortedTrackMuonsNoIso, mergedconvertedTracks);
68 
69  // To be enable when tau->3mu rate study is ready
70  //tauto3mu_->GetTau3Mu(sortedTrackMuonsNoIso, mergedconvertedTracks);
71 
72  tps_->outputGT(sortedTrackMuonsNoIso);
73 
74  return sortedTrackMuonsNoIso; //when we add more collections like tau3mu etc we change that
75 }
std::vector< edm::Ptr< l1t::TrackerMuon::L1TTTrackType > > associateTracksWithNonant(const std::vector< edm::Ptr< l1t::TrackerMuon::L1TTTrackType > > &tracks, uint processor)
Definition: TPS.cc:77
l1t::MuonStubRefVector associateStubsWithNonant(const l1t::MuonStubRefVector &, uint)
Definition: TPS.cc:103
std::unique_ptr< Isolation > isolation_
Definition: TPS.h:22
std::unique_ptr< TrackConverter > tt_track_converter_
Definition: TPS.h:20
std::unique_ptr< TPSAlgorithm > tps_
Definition: TPS.h:21

Member Data Documentation

◆ isolation_

std::unique_ptr<Isolation> Phase2L1GMT::TPS::isolation_
private

Definition at line 22 of file TPS.h.

Referenced by processEvent().

◆ tps_

std::unique_ptr<TPSAlgorithm> Phase2L1GMT::TPS::tps_
private

Definition at line 21 of file TPS.h.

Referenced by processEvent().

◆ tt_track_converter_

std::unique_ptr<TrackConverter> Phase2L1GMT::TPS::tt_track_converter_
private

Definition at line 20 of file TPS.h.

Referenced by processEvent().

◆ verbose_

int Phase2L1GMT::TPS::verbose_
private

Definition at line 19 of file TPS.h.