CMS 3D CMS Logo

AlignmentGlobalTrackSelector.cc
Go to the documentation of this file.
1 //Framework
6 
7 //DataFormats
12 
14 
15 #include <DataFormats/RecoCandidate/interface/RecoCandidate.h> //for the get<TrackRef>() Call
16 
18 
19 //STL
20 #include <cmath>
21 
23 
24 using namespace std;
25 using namespace edm;
26 
27 // constructor ----------------------------------------------------------------
29  : theGMFilterSwitch(cfg.getParameter<bool>("applyGlobalMuonFilter")),
30  theIsoFilterSwitch(cfg.getParameter<bool>("applyIsolationtest")),
31  theJetCountFilterSwitch(cfg.getParameter<bool>("applyJetCountFilter")) {
33  LogDebug("Alignment") << "> applying global Trackfilter ...";
34 
35  if (theGMFilterSwitch) {
36  edm::InputTag theMuonSource = cfg.getParameter<InputTag>("muonSource");
37  theMuonToken = iC.consumes<reco::MuonCollection>(theMuonSource);
38  theMaxTrackDeltaR = cfg.getParameter<double>("maxTrackDeltaR");
39  theMinGlobalMuonCount = cfg.getParameter<int>("minGlobalMuonCount");
40  LogDebug("Alignment") << "> GlobalMuonFilter : source, maxTrackDeltaR, min. Count : " << theMuonSource
41  << " , " << theMaxTrackDeltaR << " , " << theMinIsolatedCount;
42  }
43 
44  if (theIsoFilterSwitch) {
45  edm::InputTag theJetIsoSource = cfg.getParameter<InputTag>("jetIsoSource");
46  theJetIsoToken = iC.consumes<reco::CaloJetCollection>(theJetIsoSource);
47  theMaxJetPt = cfg.getParameter<double>("maxJetPt");
48  theMinJetDeltaR = cfg.getParameter<double>("minJetDeltaR");
49  theMinIsolatedCount = cfg.getParameter<int>("minIsolatedCount");
50  LogDebug("Alignment") << "> Isolationtest : source, maxJetPt, minJetDeltaR, min. Count: " << theJetIsoSource
51  << " , " << theMaxJetPt << " ," << theMinJetDeltaR << " ," << theMinGlobalMuonCount;
52  }
53 
55  edm::InputTag theJetCountSource = cfg.getParameter<InputTag>("jetCountSource");
56  theJetCountToken = iC.consumes<reco::CaloJetCollection>(theJetCountSource);
57  theMinJetPt = cfg.getParameter<double>("minJetPt");
58  theMaxJetCount = cfg.getParameter<int>("maxJetCount");
59  LogDebug("Alignment") << "> JetCountFilter : source, minJetPt, maxJetCount : " << theJetCountSource
60  << " , " << theMinJetPt << " ," << theMaxJetCount;
61  }
62 }
63 
64 // destructor -----------------------------------------------------------------
66 
70 }
71 
72 // do selection ---------------------------------------------------------------
74  const edm::Event& iEvent,
75  const edm::EventSetup& iSetup) {
77 
84  LogDebug("Alignment") << "> Global: tracks all, kept: " << tracks.size() << ", " << result.size();
85  // LogDebug("Alignment")<<"> o kept:";
86  // printTracks(result);
87 
88  return result;
89 }
90 
93  const edm::Event& iEvent) const {
94  Tracks result;
96 
97  //fill globalMuons with muons
99  iEvent.getByToken(theMuonToken, muons);
100 
101  if (muons.isValid()) {
102  for (reco::MuonCollection::const_iterator itMuon = muons->begin(); itMuon != muons->end(); ++itMuon) {
103  const reco::Track* muonTrack = (*itMuon).get<reco::TrackRef>().get();
104  if (!muonTrack) {
105  LogDebug("Alignment") << "@SUB=AlignmentGlobalTrackSelector::findMuons"
106  << "Found muon without track: Standalone Muon!";
107  } else {
108  globalMuons.push_back(muonTrack);
109  }
110  }
111  } else {
112  LogError("Alignment") << "@SUB=AlignmentGlobalTrackSelector::findMuons"
113  << "> could not optain mounCollection!";
114  }
115 
116  result = this->matchTracks(tracks, globalMuons);
117 
118  if (static_cast<int>(result.size()) < theMinGlobalMuonCount)
119  result.clear();
120 
121  return result;
122 }
123 
126  const edm::Event& iEvent) const {
127  Tracks result;
128  result.clear();
129 
131  iEvent.getByToken(theJetIsoToken, jets);
132 
133  if (jets.isValid()) {
134  for (Tracks::const_iterator it = cands.begin(); it != cands.end(); ++it) {
135  bool isolated = true;
136  for (reco::CaloJetCollection::const_iterator itJet = jets->begin(); itJet != jets->end(); ++itJet)
137  isolated &= !((*itJet).pt() > theMaxJetPt && deltaR(*(*it), (*itJet)) < theMinJetDeltaR);
138 
139  if (isolated)
140  result.push_back(*it);
141  }
142  // LogDebug("Alignment") << "D Found "<<result.size()<<" isolated of "<< cands.size()<<" Tracks!";
143 
144  } else
145  LogError("Alignment") << "@SUB=AlignmentGlobalTrackSelector::checkIsolation"
146  << "> could not optain jetCollection!";
147 
148  if (static_cast<int>(result.size()) < theMinIsolatedCount)
149  result.clear();
150 
151  return result;
152 }
153 
156  const edm::Event& iEvent) const {
157  Tracks result;
158  result.clear();
159 
161  iEvent.getByToken(theJetCountToken, jets);
162 
163  if (jets.isValid()) {
164  int jetCount = 0;
165  for (reco::CaloJetCollection::const_iterator itJet = jets->begin(); itJet != jets->end(); ++itJet) {
166  if ((*itJet).pt() > theMinJetPt)
167  jetCount++;
168  }
169 
170  if (jetCount <= theMaxJetCount)
171  result = tracks;
172 
173  LogDebug("Alignment") << "> found " << jetCount << " Jets";
174  } else
175  LogError("Alignment") << "@SUB=AlignmentGlobalTrackSelector::checkJetCount"
176  << "> could not optain jetCollection!";
177 
178  return result;
179 }
180 
181 //===================HELPERS===================
182 
185  const Tracks& comp) const {
186  Tracks result;
187  for (Tracks::const_iterator itComp = comp.begin(); itComp != comp.end(); ++itComp) {
188  int match = -1;
189  double min = theMaxTrackDeltaR;
190  for (unsigned int i = 0; i < src.size(); i++) {
191  // LogDebug("Alignment") << "> Trackmatch dist: "<<deltaR(src.at(i),*itComp);
192  if (min > deltaR(*(src.at(i)), *(*itComp))) {
193  min = deltaR(*(src.at(i)), *(*itComp));
194  match = static_cast<int>(i);
195  }
196  }
197  if (match > -1)
198  result.push_back(src.at(match));
199  }
200  return result;
201 }
202 
205  int count = 0;
206  LogDebug("Alignment") << ">......................................";
207  for (Tracks::const_iterator it = col.begin(); it < col.end(); ++it, ++count) {
208  LogDebug("Alignment") << "> Track No. " << count << ": p = (" << (*it)->px() << "," << (*it)->py() << ","
209  << (*it)->pz() << ")\n"
210  << "> pT = " << (*it)->pt() << " eta = " << (*it)->eta()
211  << " charge = " << (*it)->charge();
212  }
213  LogDebug("Alignment") << ">......................................";
214 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
AlignmentGlobalTrackSelector(const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
constructor
std::vector< const reco::Track * > Tracks
Log< level::Error, false > LogError
Tracks findMuons(const Tracks &tracks, const edm::Event &iEvent) const
filter for Tracks that match the Track of a global Muon
edm::EDGetTokenT< reco::CaloJetCollection > theJetCountToken
muons
the two sets of parameters below are mutually exclusive, depending if RECO or ALCARECO is used the us...
Definition: DiMuonV_cfg.py:212
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
int iEvent
Definition: GenABIO.cc:224
Tracks select(const Tracks &tracks, const edm::Event &iEvent, const edm::EventSetup &eSetup)
select tracks
edm::EDGetTokenT< reco::MuonCollection > theMuonToken
Tracks matchTracks(const Tracks &src, const Tracks &comp) const
matches [src] with [comp] returns collection with matching Tracks coming from [src] ...
Tracks checkIsolation(const Tracks &cands, const edm::Event &iEvent) const
returns only isolated tracks in [cands]
HLT enums.
col
Definition: cuy.py:1009
void printTracks(const Tracks &col) const
print Information on Track-Collection
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
Tracks checkJetCount(const Tracks &cands, const edm::Event &iEvent) const
returns [tracks] if there are less than theMaxCount Jets with theMinJetPt and an empty set if not ...
bool useThisFilter()
returns if any of the Filters is used.
edm::EDGetTokenT< reco::CaloJetCollection > theJetIsoToken
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects
#define LogDebug(id)