CMS 3D CMS Logo

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

#include <Node.h>

Public Member Functions

 Node (const edm::ParameterSet &iConfig)
 
std::vector< l1t::TrackerMuonprocessEvent (const std::vector< edm::Ptr< l1t::TrackerMuon::L1TTTrackType > > &tracks, const l1t::ObjectRefBxCollection< l1t::RegionalMuonCand > &muonTracks, const l1t::MuonStubRefVector &stubs)
 
 ~Node ()
 

Private Member Functions

l1t::MuonStubRefVector associateStubsWithNonant (const l1t::MuonStubRefVector &allStubs, uint processor)
 
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< ROITempAssociatorroi_assoc_
 
std::unique_ptr< Tauto3Mutauto3mu_
 
std::unique_ptr< TrackMuonMatchAlgorithmtrack_mu_match_
 
std::unique_ptr< TrackConvertertt_track_converter_
 
int verbose_
 

Detailed Description

Definition at line 13 of file Node.h.

Constructor & Destructor Documentation

◆ Node()

Phase2L1GMT::Node::Node ( const edm::ParameterSet iConfig)
inline

Definition at line 15 of file Node.h.

16  : verbose_(iConfig.getParameter<int>("verbose")),
17  tt_track_converter_(new TrackConverter(iConfig.getParameter<edm::ParameterSet>("trackConverter"))),
18  roi_assoc_(new ROITempAssociator(iConfig.getParameter<edm::ParameterSet>("roiTrackAssociator"))),
19  track_mu_match_(new TrackMuonMatchAlgorithm(iConfig.getParameter<edm::ParameterSet>("trackMatching"))),
20  isolation_(new Isolation(iConfig.getParameter<edm::ParameterSet>("isolation"))),
21  tauto3mu_(new Tauto3Mu(iConfig.getParameter<edm::ParameterSet>("tauto3mu"))) {}
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
std::unique_ptr< Isolation > isolation_
Definition: Node.h:140
std::unique_ptr< Tauto3Mu > tauto3mu_
Definition: Node.h:141
std::unique_ptr< TrackMuonMatchAlgorithm > track_mu_match_
Definition: Node.h:139
std::unique_ptr< ROITempAssociator > roi_assoc_
Definition: Node.h:138
std::unique_ptr< TrackConverter > tt_track_converter_
Definition: Node.h:137

◆ ~Node()

Phase2L1GMT::Node::~Node ( )
inline

Definition at line 23 of file Node.h.

23 {}

Member Function Documentation

◆ associateStubsWithNonant()

l1t::MuonStubRefVector Phase2L1GMT::Node::associateStubsWithNonant ( const l1t::MuonStubRefVector allStubs,
uint  processor 
)
inlineprivate

Definition at line 154 of file Node.h.

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

Referenced by processEvent().

154  {
156 
157  ap_int<BITSSTUBCOORD> center = ap_int<BITSSTUBCOORD>((processor * 910) / 32);
158 
159  for (const auto& s : allStubs) {
160  ap_int<BITSSTUBCOORD> phi = 0;
161  if (s->quality() & 0x1)
162  phi = s->coord1();
163  else
164  phi = s->coord2();
165 
166  ap_int<BITSSTUBCOORD> deltaPhi = phi - center;
167  ap_uint<BITSSTUBCOORD - 1> absDeltaPhi =
168  (deltaPhi < 0) ? ap_uint<BITSSTUBCOORD - 1>(-deltaPhi) : ap_uint<BITSSTUBCOORD - 1>(deltaPhi);
169  if (absDeltaPhi < 42)
170  out.push_back(s);
171 
172  /* if (processor==0 && phi>=-3000/32 && phi<=3000/32 ) */
173  /* out.push_back(s); */
174  /* else if (processor==1 && (phi>=-1000/32 && phi<=5000/32) ) */
175  /* out.push_back(s); */
176  /* else if (processor==2 && (phi>=500/32 && phi<=6500/32) ) */
177  /* out.push_back(s); */
178  /* else if (processor==3 && (phi>=2000/32 || phi<=-8000/32) ) */
179  /* out.push_back(s); */
180  /* else if (processor==4 && (phi>=4500/32 || phi<=-6000/32) ) */
181  /* out.push_back(s); */
182  /* else if (processor==5 && (phi>=6000/32 || phi<=-4500/32) ) */
183  /* out.push_back(s); */
184  /* else if (processor==6 && (phi>=8000/32 || phi<=-2000/32) ) */
185  /* out.push_back(s); */
186  /* else if (processor==7 && (phi>=-7000/32 && phi<=0) ) */
187  /* out.push_back(s); */
188  /* else if (processor==8 && (phi>=-4500/32 && phi<=1000/32) ) */
189  /* out.push_back(s); */
190  }
191  return out;
192  }
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> > Phase2L1GMT::Node::associateTracksWithNonant ( const std::vector< edm::Ptr< l1t::TrackerMuon::L1TTTrackType > > &  tracks,
uint  processor 
)
inlineprivate

Definition at line 143 of file Node.h.

References MillePedeFileConverter_cfg::out, HLT_2023v12_cff::track, and pwdgSkimBPark_cfi::tracks.

Referenced by processEvent().

144  {
145  std::vector<edm::Ptr<l1t::TrackerMuon::L1TTTrackType> > out;
146  for (const auto& track : tracks) {
147  if (track->phiSector() == processor) {
148  out.push_back(track);
149  }
150  }
151  return out;
152  }

◆ processEvent()

std::vector<l1t::TrackerMuon> Phase2L1GMT::Node::processEvent ( const std::vector< edm::Ptr< l1t::TrackerMuon::L1TTTrackType > > &  tracks,
const l1t::ObjectRefBxCollection< l1t::RegionalMuonCand > &  muonTracks,
const l1t::MuonStubRefVector stubs 
)
inline

Definition at line 25 of file Node.h.

References associateStubsWithNonant(), associateTracksWithNonant(), l1tGTTFileWriter_cfi::convertedTracks, filterCSVwithJSON::copy, isolation_, TrackingDataMCValidation_Standalone_cff::muonTracks, roi_assoc_, track_mu_match_, pwdgSkimBPark_cfi::tracks, tt_track_converter_, and verbose_.

27  {
28  //Split tracks to the links as they come
29  std::vector<edm::Ptr<l1t::TrackerMuon::L1TTTrackType> > tracks0 = associateTracksWithNonant(tracks, 0);
30  std::vector<edm::Ptr<l1t::TrackerMuon::L1TTTrackType> > tracks1 = associateTracksWithNonant(tracks, 1);
31  std::vector<edm::Ptr<l1t::TrackerMuon::L1TTTrackType> > tracks2 = associateTracksWithNonant(tracks, 2);
32  std::vector<edm::Ptr<l1t::TrackerMuon::L1TTTrackType> > tracks3 = associateTracksWithNonant(tracks, 3);
33  std::vector<edm::Ptr<l1t::TrackerMuon::L1TTTrackType> > tracks4 = associateTracksWithNonant(tracks, 4);
34  std::vector<edm::Ptr<l1t::TrackerMuon::L1TTTrackType> > tracks5 = associateTracksWithNonant(tracks, 5);
35  std::vector<edm::Ptr<l1t::TrackerMuon::L1TTTrackType> > tracks6 = associateTracksWithNonant(tracks, 6);
36  std::vector<edm::Ptr<l1t::TrackerMuon::L1TTTrackType> > tracks7 = associateTracksWithNonant(tracks, 7);
37  std::vector<edm::Ptr<l1t::TrackerMuon::L1TTTrackType> > tracks8 = associateTracksWithNonant(tracks, 8);
38 
39  //Transition stubs to different nonants with overlap
49 
50  //Convert TT tracks to our internal tracking format
51  std::vector<ConvertedTTTrack> convertedTracks0 = tt_track_converter_->convertTracks(tracks0);
52  std::vector<ConvertedTTTrack> convertedTracks1 = tt_track_converter_->convertTracks(tracks1);
53  std::vector<ConvertedTTTrack> convertedTracks2 = tt_track_converter_->convertTracks(tracks2);
54  std::vector<ConvertedTTTrack> convertedTracks3 = tt_track_converter_->convertTracks(tracks3);
55  std::vector<ConvertedTTTrack> convertedTracks4 = tt_track_converter_->convertTracks(tracks4);
56  std::vector<ConvertedTTTrack> convertedTracks5 = tt_track_converter_->convertTracks(tracks5);
57  std::vector<ConvertedTTTrack> convertedTracks6 = tt_track_converter_->convertTracks(tracks6);
58  std::vector<ConvertedTTTrack> convertedTracks7 = tt_track_converter_->convertTracks(tracks7);
59  std::vector<ConvertedTTTrack> convertedTracks8 = tt_track_converter_->convertTracks(tracks8);
60 
61  //Build ROIs per nonant
62  std::vector<MuonROI> rois0 = roi_assoc_->associate(0, muonTracks, stubs0);
63  std::vector<MuonROI> rois1 = roi_assoc_->associate(0, muonTracks, stubs1);
64  std::vector<MuonROI> rois2 = roi_assoc_->associate(0, muonTracks, stubs2);
65  std::vector<MuonROI> rois3 = roi_assoc_->associate(0, muonTracks, stubs3);
66  std::vector<MuonROI> rois4 = roi_assoc_->associate(0, muonTracks, stubs4);
67  std::vector<MuonROI> rois5 = roi_assoc_->associate(0, muonTracks, stubs5);
68  std::vector<MuonROI> rois6 = roi_assoc_->associate(0, muonTracks, stubs6);
69  std::vector<MuonROI> rois7 = roi_assoc_->associate(0, muonTracks, stubs7);
70  std::vector<MuonROI> rois8 = roi_assoc_->associate(0, muonTracks, stubs8);
71 
72  //run track - muon matching per nonant
73  std::vector<PreTrackMatchedMuon> mu0 = track_mu_match_->processNonant(convertedTracks0, rois0);
74  std::vector<PreTrackMatchedMuon> mu1 = track_mu_match_->processNonant(convertedTracks1, rois1);
75  std::vector<PreTrackMatchedMuon> mu2 = track_mu_match_->processNonant(convertedTracks2, rois2);
76  std::vector<PreTrackMatchedMuon> mu3 = track_mu_match_->processNonant(convertedTracks3, rois3);
77  std::vector<PreTrackMatchedMuon> mu4 = track_mu_match_->processNonant(convertedTracks4, rois4);
78  std::vector<PreTrackMatchedMuon> mu5 = track_mu_match_->processNonant(convertedTracks5, rois5);
79  std::vector<PreTrackMatchedMuon> mu6 = track_mu_match_->processNonant(convertedTracks6, rois6);
80  std::vector<PreTrackMatchedMuon> mu7 = track_mu_match_->processNonant(convertedTracks7, rois7);
81  std::vector<PreTrackMatchedMuon> mu8 = track_mu_match_->processNonant(convertedTracks8, rois8);
82  if (verbose_)
83  printf("Matching Nonant 5 with %zu tracks and %zu rois and %zu stubs\n",
84  convertedTracks5.size(),
85  rois5.size(),
86  stubs5.size());
87 
88  //clean neighboring nonants
89  std::vector<PreTrackMatchedMuon> muCleaned = track_mu_match_->cleanNeighbor(mu0, mu8, mu1, true);
90  std::vector<PreTrackMatchedMuon> muCleaned1 = track_mu_match_->cleanNeighbor(mu1, mu0, mu2, false);
91  std::vector<PreTrackMatchedMuon> muCleaned2 = track_mu_match_->cleanNeighbor(mu2, mu1, mu3, true);
92  std::vector<PreTrackMatchedMuon> muCleaned3 = track_mu_match_->cleanNeighbor(mu3, mu2, mu4, false);
93  std::vector<PreTrackMatchedMuon> muCleaned4 = track_mu_match_->cleanNeighbor(mu4, mu3, mu5, true);
94  std::vector<PreTrackMatchedMuon> muCleaned5 = track_mu_match_->cleanNeighbor(mu5, mu4, mu6, false);
95  std::vector<PreTrackMatchedMuon> muCleaned6 = track_mu_match_->cleanNeighbor(mu6, mu5, mu7, true);
96  std::vector<PreTrackMatchedMuon> muCleaned7 = track_mu_match_->cleanNeighbor(mu7, mu6, mu8, false);
97  std::vector<PreTrackMatchedMuon> muCleaned8 =
98  track_mu_match_->cleanNeighbor(mu8, mu7, mu0, false); //ARGH! 9 sectors - so some duplicates very rarely
99 
100  //merge all the collections
101  std::copy(muCleaned1.begin(), muCleaned1.end(), std::back_inserter(muCleaned));
102  std::copy(muCleaned2.begin(), muCleaned2.end(), std::back_inserter(muCleaned));
103  std::copy(muCleaned3.begin(), muCleaned3.end(), std::back_inserter(muCleaned));
104  std::copy(muCleaned4.begin(), muCleaned4.end(), std::back_inserter(muCleaned));
105  std::copy(muCleaned5.begin(), muCleaned5.end(), std::back_inserter(muCleaned));
106  std::copy(muCleaned6.begin(), muCleaned6.end(), std::back_inserter(muCleaned));
107  std::copy(muCleaned7.begin(), muCleaned7.end(), std::back_inserter(muCleaned));
108  std::copy(muCleaned8.begin(), muCleaned8.end(), std::back_inserter(muCleaned));
109 
110  std::vector<l1t::TrackerMuon> trackMatchedMuonsNoIso = track_mu_match_->convert(muCleaned, 32);
111 
112  //Isolation and tau3mu will read those muons and all 9 collections of convertedTracks*
113  std::vector<ConvertedTTTrack> convertedTracks = convertedTracks0;
114  std::copy(convertedTracks1.begin(), convertedTracks1.end(), std::back_inserter(convertedTracks));
115  std::copy(convertedTracks2.begin(), convertedTracks2.end(), std::back_inserter(convertedTracks));
116  std::copy(convertedTracks3.begin(), convertedTracks3.end(), std::back_inserter(convertedTracks));
117  std::copy(convertedTracks4.begin(), convertedTracks4.end(), std::back_inserter(convertedTracks));
118  std::copy(convertedTracks5.begin(), convertedTracks5.end(), std::back_inserter(convertedTracks));
119  std::copy(convertedTracks6.begin(), convertedTracks6.end(), std::back_inserter(convertedTracks));
120  std::copy(convertedTracks7.begin(), convertedTracks7.end(), std::back_inserter(convertedTracks));
121  std::copy(convertedTracks8.begin(), convertedTracks8.end(), std::back_inserter(convertedTracks));
122 
123  //sorter here:
124  std::vector<l1t::TrackerMuon> sortedTrackMuonsNoIso = track_mu_match_->sort(trackMatchedMuonsNoIso, 12);
125 
126  isolation_->isolation_allmu_alltrk(sortedTrackMuonsNoIso, convertedTracks);
127 
128  //tauto3mu_->GetTau3Mu(sortedTrackMuonsNoIso, convertedTracks);
129 
130  track_mu_match_->outputGT(sortedTrackMuonsNoIso);
131 
132  return sortedTrackMuonsNoIso; //when we add more collections like tau3mu etc we change that
133  }
std::unique_ptr< Isolation > isolation_
Definition: Node.h:140
std::vector< edm::Ref< MuonStubCollection > > MuonStubRefVector
Definition: MuonStub.h:44
l1t::MuonStubRefVector associateStubsWithNonant(const l1t::MuonStubRefVector &allStubs, uint processor)
Definition: Node.h:154
std::vector< edm::Ptr< l1t::TrackerMuon::L1TTTrackType > > associateTracksWithNonant(const std::vector< edm::Ptr< l1t::TrackerMuon::L1TTTrackType > > &tracks, uint processor)
Definition: Node.h:143
std::unique_ptr< TrackMuonMatchAlgorithm > track_mu_match_
Definition: Node.h:139
std::unique_ptr< ROITempAssociator > roi_assoc_
Definition: Node.h:138
std::unique_ptr< TrackConverter > tt_track_converter_
Definition: Node.h:137

Member Data Documentation

◆ isolation_

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

Definition at line 140 of file Node.h.

Referenced by processEvent().

◆ roi_assoc_

std::unique_ptr<ROITempAssociator> Phase2L1GMT::Node::roi_assoc_
private

Definition at line 138 of file Node.h.

Referenced by processEvent().

◆ tauto3mu_

std::unique_ptr<Tauto3Mu> Phase2L1GMT::Node::tauto3mu_
private

Definition at line 141 of file Node.h.

◆ track_mu_match_

std::unique_ptr<TrackMuonMatchAlgorithm> Phase2L1GMT::Node::track_mu_match_
private

Definition at line 139 of file Node.h.

Referenced by processEvent().

◆ tt_track_converter_

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

Definition at line 137 of file Node.h.

Referenced by processEvent().

◆ verbose_

int Phase2L1GMT::Node::verbose_
private

Definition at line 136 of file Node.h.

Referenced by processEvent().