CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AlignmentGlobalTrackSelector.cc
Go to the documentation of this file.
1 //Framework
5 
6 //DataFormats
13 
14 #include <DataFormats/RecoCandidate/interface/RecoCandidate.h> //for the get<TrackRef>() Call
15 
17 
18 //STL
19 #include <math.h>
20 
22 
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")),
32  theMuonSource("muons"),
33  theJetIsoSource("fastjet6CaloJets"),
34  theJetCountSource("fastjet6CaloJets")
35 {
37  LogDebug("Alignment") << "> applying global Trackfilter ...";
38 
39  if (theGMFilterSwitch) {
40  theMuonSource = cfg.getParameter<InputTag>("muonSource");
41  theMaxTrackDeltaR =cfg.getParameter<double>("maxTrackDeltaR");
42  theMinGlobalMuonCount = cfg.getParameter<int>("minGlobalMuonCount");
43  LogDebug("Alignment") << "> GlobalMuonFilter : source, maxTrackDeltaR, min. Count : "
44  << theMuonSource << " , "
45  << theMaxTrackDeltaR << " , "
47  }
48 
49  if (theIsoFilterSwitch) {
50  theJetIsoSource = cfg.getParameter<InputTag>("jetIsoSource");
51  theMaxJetPt = cfg.getParameter<double>("maxJetPt");
52  theMinJetDeltaR = cfg.getParameter<double>("minJetDeltaR");
53  theMinIsolatedCount = cfg.getParameter<int>("minIsolatedCount");
54  LogDebug("Alignment") << "> Isolationtest : source, maxJetPt, minJetDeltaR, min. Count: "
55  << theJetIsoSource << " , "
56  << theMaxJetPt << " ,"
57  << theMinJetDeltaR << " ,"
59  }
60 
62  theJetCountSource = cfg.getParameter<InputTag>("jetCountSource");
63  theMinJetPt = cfg.getParameter<double>("minJetPt");
64  theMaxJetCount = cfg.getParameter<int>("maxJetCount");
65  LogDebug("Alignment") << "> JetCountFilter : source, minJetPt, maxJetCount : "
66  << theJetCountSource << " , "
67  << theMinJetPt << " ,"
68  << theMaxJetCount;
69  }
70 }
71 
72 // destructor -----------------------------------------------------------------
74 {}
75 
78 {
80 }
81 
82 // do selection ---------------------------------------------------------------
85 {
87 
88  if (theGMFilterSwitch) result = findMuons(result, iEvent);
89  if (theIsoFilterSwitch) result = checkIsolation(result, iEvent);
90  if (theJetCountFilterSwitch) result = checkJetCount(result, iEvent);
91  LogDebug("Alignment") << "> Global: tracks all, kept: " << tracks.size() << ", " << result.size();
92 // LogDebug("Alignment")<<"> o kept:";
93 // printTracks(result);
94 
95  return result;
96 }
97 
101 {
102  Tracks result;
104 
105  //fill globalMuons with muons
107  iEvent.getByLabel(theMuonSource, muons);
108 
109  if (muons.isValid()) {
110  for (reco::MuonCollection::const_iterator itMuon = muons->begin();
111  itMuon != muons->end();
112  ++itMuon) {
113  const reco::Track* muonTrack = (*itMuon).get<reco::TrackRef>().get();
114  if (!muonTrack) {
115  LogDebug("Alignment") << "@SUB=AlignmentGlobalTrackSelector::findMuons"
116  << "Found muon without track: Standalone Muon!";
117  } else {
118  globalMuons.push_back(muonTrack);
119  }
120  }
121  } else {
122  LogError("Alignment") << "@SUB=AlignmentGlobalTrackSelector::findMuons"
123  <<"> could not optain mounCollection!";
124  }
125 
126  result = this->matchTracks(tracks, globalMuons);
127 
128  if (static_cast<int>(result.size()) < theMinGlobalMuonCount) result.clear();
129 
130  return result;
131 }
132 
136 {
137  Tracks result; result.clear();
138 
140  iEvent.getByLabel(theJetIsoSource, jets);
141 
142  if (jets.isValid()) {
143  for (Tracks::const_iterator it = cands.begin();
144  it!=cands.end();
145  ++it) {
146  bool isolated = true;
147  for (reco::CaloJetCollection::const_iterator itJet = jets->begin();
148  itJet!=jets->end();
149  ++itJet)
150  isolated &= !((*itJet).pt() > theMaxJetPt && deltaR(*(*it),(*itJet)) < theMinJetDeltaR);
151 
152  if (isolated)
153  result.push_back(*it);
154  }
155  // LogDebug("Alignment") << "D Found "<<result.size()<<" isolated of "<< cands.size()<<" Tracks!";
156 
157  } else
158  LogError("Alignment")<< "@SUB=AlignmentGlobalTrackSelector::checkIsolation"
159  << "> could not optain jetCollection!";
160 
161  if (static_cast<int>(result.size()) < theMinIsolatedCount) result.clear();
162 
163  return result;
164 }
165 
169 {
170  Tracks result; result.clear();
171 
173  iEvent.getByLabel(theJetCountSource, jets);
174 
175  if (jets.isValid()) {
176  int jetCount = 0;
177  for (reco::CaloJetCollection::const_iterator itJet = jets->begin();
178  itJet!=jets->end();
179  ++itJet) {
180  if ((*itJet).pt() > theMinJetPt)
181  jetCount++;
182  }
183 
184  if (jetCount <= theMaxJetCount)
185  result = tracks;
186 
187  LogDebug("Alignment") << "> found " << jetCount << " Jets";
188  } else
189  LogError("Alignment") << "@SUB=AlignmentGlobalTrackSelector::checkJetCount"
190  << "> could not optain jetCollection!";
191 
192  return result;
193 }
194 
195 //===================HELPERS===================
196 
200 {
201  Tracks result;
202  for (Tracks::const_iterator itComp = comp.begin();
203  itComp!=comp.end();
204  ++itComp) {
205  int match = -1;
206  double min = theMaxTrackDeltaR;
207  for (unsigned int i=0;i<src.size();i++) {
208  // LogDebug("Alignment") << "> Trackmatch dist: "<<deltaR(src.at(i),*itComp);
209  if(min > deltaR(*(src.at(i)),*(*itComp))){
210  min = deltaR(*(src.at(i)),*(*itComp));
211  match = static_cast<int>(i);
212  }
213  }
214  if (match > -1)
215  result.push_back(src.at(match));
216  }
217  return result;
218 }
219 
222 {
223  int count = 0;
224  LogDebug("Alignment") << ">......................................";
225  for (Tracks::const_iterator it = col.begin();
226  it < col.end();
227  ++it,++count) {
228  LogDebug("Alignment")
229  << "> Track No. " << count << ": p = ("
230  << (*it)->px() << ","
231  << (*it)->py() << ","
232  << (*it)->pz() << ")\n"
233  << "> pT = "
234  << (*it)->pt() << " eta = "
235  << (*it)->eta() << " charge = "
236  << (*it)->charge();
237  }
238  LogDebug("Alignment") << ">......................................";
239 }
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
Tracks findMuons(const Tracks &tracks, const edm::Event &iEvent) const
filter for Tracks that match the Track of a global Muon
#define min(a, b)
Definition: mlp_lapack.h:161
std::vector< const reco::Track * > Tracks
int iEvent
Definition: GenABIO.cc:243
vector< PseudoJet > jets
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 ...
tuple result
Definition: query.py:137
void printTracks(const Tracks &col) const
print Information on Track-Collection
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
Tracks checkIsolation(const Tracks &cands, const edm::Event &iEvent) const
returns only isolated tracks in [cands]
tuple tracks
Definition: testEve_cfg.py:39
AlignmentGlobalTrackSelector(const edm::ParameterSet &cfg)
constructor
tuple muons
Definition: patZpeak.py:38
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:6
Tracks select(const Tracks &tracks, const edm::Event &iEvent)
select tracks
bool useThisFilter()
returns if any of the Filters is used.
Tracks matchTracks(const Tracks &src, const Tracks &comp) const
matches [src] with [comp] returns collection with matching Tracks coming from [src] ...