CMS 3D CMS Logo

Tau3MuReco.cc
Go to the documentation of this file.
1 #include <TMath.h>
2 
4 
6 
8 
9 Tau3MuReco::Tau3MuReco(const edm::ParameterSet& iConfig, edm::ConsumesCollector&& iC):m_kMatchingDeltaR(iConfig.getParameter<double>("RecoAnalysisMatchingDeltaR")),
10  m_kMatchingPt(iConfig.getParameter<double>("RecoAnalysisMatchingPt")),
11  m_kTauMassCut(iConfig.getParameter<double>("RecoAnalysisTauMassCut")),
12  m_kTauMass(iConfig.getParameter<double>("RecoAnalysisTauMass")),
13  m_kMuonMass(iConfig.getParameter<double>("RecoAnalysisMuonMass")),
14  m_kMuonSourceToken(iC.consumes<reco::MuonCollection>(iConfig.getParameter<edm::InputTag>("MuonSourceTag"))),
15  m_kTrackSourceToken(iC.consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("TrackSourceTag")))
16 {
17 }
18 
20 {
21 
22 }
23 
25 {
28 
31 
32  reco::MuonCollection::const_iterator muon;
33 
34  iEvent.getByToken(m_kMuonSourceToken, muons);
35  iEvent.getByToken(m_kTrackSourceToken, tracks);
36 
37  for( muon = muons->begin(); muon != muons->end(); ++muon )
38  {
39  m_TrackCollection->push_back(*(muon->track().get()));
40  m_MuonCollection->push_back(*muon);
41  }
42 
43  if( m_TrackCollection->size() > 3 )
44  {
45  //find the right three ones coming from tau
46  if(findCorrectPairing()==false)
47  { //maybe implement something like in ==3 (throw away ....)
48  LogDebug("Tau3MuReco") << "Could not find correct combination!" << std::endl;
49  return false;
50  }
51 
52  return true;
53  }
54 
55  if( m_TrackCollection->size() == 3 )
56  {
58  return true;
59  else//throw away one muon which don't match
61 
62  }
63 
64  if( m_TrackCollection->size() == 2)
65  {
66  //search the third track
67 
68  //get 3rd muon canidate from tracks
69  if(find3rdTrack(iEvent, iSetup, *(tracks.product()))==false)
70  {
71  LogDebug("Tau3MuReco") << "A 3rd Track can not be assigned!" << std::endl;
72  return false;
73  }
74  else
75  return true;
76  }
77 
78  // cannot use this event, because less than 2 muons have been found
79 
80  LogDebug("Tau3MuReco") << "Not enough (" << m_TrackCollection->size() << ") muons found! Event skipped!" << std::endl;
81 
82  return false;
83 }
84 
85 //private
87 {
88  reco::TrackCollection::const_iterator iter;
89 
90  for(iter = m_TrackCollection->begin(); iter!=m_TrackCollection->end(); iter++)
91  {
92  //check if the track has the right charge
93  //and check if dR is smaller than fMatchingDeltaR
94  if((*iter).charge() == track.charge()
95  && getDeltaR(*iter,track)<m_kMatchingDeltaR
96  && fabs(((*iter).pt())-(track.pt()))<=(track.pt()*m_kMatchingPt))
97  {
98  LogDebug("Tau3MuReco") << "Found muon track in Tracks with DeltaR: " << getDeltaR(*iter,track)
99  << " Pt: " << track.pt()
100  << " Muons Pt is: " << (*iter).pt() << std::endl;
101  return true;
102  }
103  }
104  return false;
105 }
106 
107 
109 {
110  //check size of TrackVector (should be two!)
111  if(m_TrackCollection->size()!=2)
112  return false;
113 
114  //more then two tracks should be in the event
115  if(Tracks.size()<=2)
116  return false;
117 
118  double SumDeltaR = 0;
119 
120  double MuonDeltaR = getDeltaR(m_TrackCollection->at(0),m_TrackCollection->at(1));
121 
122  //Loop overall tracks
123 
124  LogDebug("Tau3MuReco") << "Number of tracks found: " << Tracks.size() << std::endl;
125 
126  std::multimap<double, const reco::Track> TrackMultiMap;
127 
128  unsigned short muonCounter = 0;
129 
130  reco::TrackCollection::const_iterator track;
131 
132  for(track=Tracks.begin(); track!=Tracks.end(); track++)
133  {
134  if(check4MuonTrack(*track))
135  {
136  muonCounter++;
137  continue;
138  }
139 
140  SumDeltaR = MuonDeltaR;
141 
142  SumDeltaR += getDeltaR(m_TrackCollection->at(1), *track);
143  SumDeltaR += getDeltaR(*track,m_TrackCollection->at(0));
144 
145  std::pair<double, const reco::Track> actTrack(SumDeltaR,*track);
146 
147  TrackMultiMap.insert(actTrack);
148  }
149 
150  //two tracks should be clearly identified as muons by check4MuonTrack
151  //else event is not useable
152  if(muonCounter<2)
153  {
154  LogDebug("Tau3MuReco") << "Not enough muons (" << muonCounter << ") assigned to a track! Event skipped!" << std::endl;
155  return false;
156  }
157 
158  std::multimap<double, const reco::Track>::iterator it = TrackMultiMap.begin();
159 
160  if(it==TrackMultiMap.end())
161  {
162  LogDebug("Tau3MuReco") << "Not enough tracks (0) left! Event skipped!" << std::endl;
163  return false;
164  }
165 
166  //get 2mu+track with minimal DeltaR Sum (MultiMaps are sorted)
167  m_TrackCollection->push_back((*it).second);
168 
169  //and check charge of this track
170  //and check invariant mass of this combination
171  //and make a vertex fit
172 
173  char Charge = m_TrackCollection->at(0).charge() * m_TrackCollection->at(1).charge();
174 
175  unsigned int count = 0;
176 
177  //Charge > 0 means the two muons have same charge, so the third track has to have the opposit charge
178  while((Charge > 0 && ((*it).second).charge()==(m_TrackCollection->at(0)).charge())
180  )
181  {
182 
183  count++;
184 
185  LogDebug("Tau3MuReco") << "Track canidate: " << count << std::endl;
186 
187  if(Charge > 0 && ((*it).second).charge()!=(m_TrackCollection->at(0)).charge())
188  LogDebug("Tau3MuReco") << "\tWrong charge!" << std::endl;
189  LogDebug("Tau3MuReco") << "\tInvariant Mass deviation! " << fabs(getInvariantMass(m_TrackCollection)-m_kTauMass) << std::endl;
190  LogDebug("Tau3MuReco") << "\tTrack Pt: " << (*it).second.pt() << std::endl;
191  LogDebug("Tau3MuReco") << "\tDelta R: " << (*it).first << std::endl;
192  LogDebug("Tau3MuReco") << "\tChi2: " << ((*it).second).normalizedChi2() << std::endl;
193 
194  ++it;
195 
196  //was not the best canidate
197  m_TrackCollection->pop_back();
198 
199  if(it==TrackMultiMap.end())
200  return false;
201 
202  //get next to minimal (Delta R Sum) track
203  m_TrackCollection->push_back((*it).second);
204  }
205 
206  LogDebug("Tau3MuReco") << "Found corresponding 3rd track: " << std::endl;
207  LogDebug("Tau3MuReco") << "Track canidate: " << count << std::endl;
208  LogDebug("Tau3MuReco") << "\tInvariant Mass deviation! " << fabs(getInvariantMass(m_TrackCollection)-m_kTauMass) << std::endl;
209  LogDebug("Tau3MuReco") << "\tDelta R: " << (*it).first << std::endl;
210  LogDebug("Tau3MuReco") << "\tNormChi2: " << ((*it).second).normalizedChi2() << std::endl;
211 
212  //choose this track, because it is the best canidate
213  return true;
214 }
215 
217 {
218  Combinatorics myCombinatorics(m_TrackCollection->size(), 3);
219 
220  std::vector < std::vector<UInt_t> > CombinationVec = myCombinatorics.GetCombinations();
221 
222  std::vector< std::vector<UInt_t> >::iterator it = CombinationVec.begin();
223 
224  char Charge = 0;
225 
226  reco::TrackCollection tempTrackCollection;
227  reco::MuonCollection tempMuonCollection;
228 
229  do
230  {
231  if(it==CombinationVec.end())
232  return false;
233 
234  Charge = 0;
235 
236  tempMuonCollection.clear();
237  tempTrackCollection.clear();
238 
239  for(UInt_t i=0; i< (*it).size(); i++)
240  {
241  Charge += m_TrackCollection->at((*it).at(i)).charge();
242  tempTrackCollection.push_back(m_TrackCollection->at((*it).at(i)));
243  tempMuonCollection.push_back(m_MuonCollection->at((*it).at(i)));
244  }
245 
246  LogDebug("Tau3MuReco") << "Charge is: " << (int)Charge << " Have to be -1 or 1!!!" << std::endl;
247  LogDebug("Tau3MuReco") << "Invariant mass is: " << fabs(getInvariantMass(&tempTrackCollection)-m_kTauMass) << " Have to be smaller than " << m_kTauMassCut << std::endl;
248 
249  it++;
250  }
251  while(abs(Charge)!=1 || fabs(getInvariantMass(&tempTrackCollection)-m_kTauMass)>m_kTauMassCut);
252 
253  *m_MuonCollection = tempMuonCollection;
254  *m_TrackCollection = tempTrackCollection;
255 
256  return true;
257 }
258 
260 {
261 
262  double deltaR12 = getDeltaR(m_TrackCollection->at(0),m_TrackCollection->at(1));
263  double deltaR23 = getDeltaR(m_TrackCollection->at(1),m_TrackCollection->at(2));
264  double deltaR31 = getDeltaR(m_TrackCollection->at(2),m_TrackCollection->at(0));
265 
266  //if DeltaR12 is the smallest, than the 3rd one seems to be wrong
267  //if DeltaR23 is the smallest, than the 2nd one seems to be wrong
268  //if DeltaR31 is the smallest, than the 1st one seems to be wrong
269 
270  unsigned char temp;
271  double junk;
272 
273  deltaR12 < deltaR23 ? temp=3 : temp=1;
274  deltaR12 < deltaR23 ? junk=deltaR12 : junk=deltaR23;
275 
276  if(deltaR31 < junk)
277  temp=2;
278 
279  m_TrackCollection->erase(m_TrackCollection->begin()+temp-1);
280 
281  return true;
282 }
283 
285 {
286  unsigned int numOfParticles = tracks->size();
287 
288  double SumPx = 0;
289  double SumPy = 0;
290  double SumPz = 0;
291 
292  double SumE = 0;
293 
294  for(unsigned int i=0; i<numOfParticles; i++)
295  {
296  SumPx += tracks->at(i).px();
297  SumPy += tracks->at(i).py();
298  SumPz += tracks->at(i).pz();
299 
300  SumE += sqrt(pow(tracks->at(i).p(),2)+pow(MuonMass,2));
301  }
302 
303  double invmass = sqrt(pow(SumE,2)-pow(SumPx,2)-pow(SumPy,2)-pow(SumPz,2));
304 
305  return invmass;
306 }
307 
308 double Tau3MuReco::getDeltaR(const reco::Track& track1, const reco::Track& track2)
309 {
310  double dEta = track1.eta() - track2.eta();
311  double dPhi = track1.phi() - track2.phi();
312 
313  while(dPhi >= TMath::Pi()) dPhi -= (2.0*TMath::Pi());
314  while(dPhi < (-1.0*TMath::Pi())) dPhi += (2.0*TMath::Pi());
315 
316  return sqrt(pow(dEta,2)+pow(dPhi,2));
317 }
#define LogDebug(id)
Tau3MuReco(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
Definition: Tau3MuReco.cc:9
const double Pi
double getDeltaR(const reco::Track &track1, const reco::Track &track2)
Definition: Tau3MuReco.cc:308
std::vector< std::vector< UInt_t > > GetCombinations()
double getInvariantMass(const reco::TrackCollection *tracks, const double MuonMass=0.106)
Definition: Tau3MuReco.cc:284
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
const double m_kTauMass
Definition: Tau3MuReco.h:33
bool find3rdTrack(const edm::Event &iEvent, const edm::EventSetup &iSetup, const reco::TrackCollection &Tracks)
Definition: Tau3MuReco.cc:108
const double m_kMatchingPt
Definition: Tau3MuReco.h:31
reco::TrackCollection * m_TrackCollection
Definition: Tau3MuReco.h:40
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:645
bool findCorrectPairing()
Definition: Tau3MuReco.cc:216
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
int iEvent
Definition: GenABIO.cc:230
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:651
const double m_kMatchingDeltaR
Definition: Tau3MuReco.h:30
T sqrt(T t)
Definition: SSEVec.h:18
double pt() const
track transverse momentum
Definition: TrackBase.h:621
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool doTau3MuReco(const edm::Event &iEvent, const edm::EventSetup &iSetup, reco::MuonCollection *muonCollection, reco::TrackCollection *trackCollection)
Definition: Tau3MuReco.cc:24
T const * product() const
Definition: Handle.h:81
const double m_kTauMassCut
Definition: Tau3MuReco.h:32
const double MuonMass
const double m_kMuonMass
Definition: Tau3MuReco.h:34
fixed size matrix
bool removeIncorrectMuon()
Definition: Tau3MuReco.cc:259
HLT enums.
const edm::EDGetTokenT< reco::MuonCollection > m_kMuonSourceToken
Definition: Tau3MuReco.h:36
const edm::EDGetTokenT< reco::TrackCollection > m_kTrackSourceToken
Definition: Tau3MuReco.h:37
reco::MuonCollection * m_MuonCollection
Definition: Tau3MuReco.h:39
int charge() const
track electric charge
Definition: TrackBase.h:567
bool check4MuonTrack(const reco::Track &track)
Definition: Tau3MuReco.cc:86
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40