CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AlignmentTwoBoyDecayTrackSelector.cc
Go to the documentation of this file.
1 //Framework
3 
6 
7 //DataFormats
11 
12 //STL
13 #include <math.h>
14 //ROOT
15 #include "TLorentzVector.h"
16 
18 //TODO put those namespaces into functions?
19 using namespace std;
20 using namespace edm;
21 // constructor ----------------------------------------------------------------
22 
24  theMissingETSource("met")
25 {
26  LogDebug("Alignment") << "> applying two body decay Trackfilter ...";
27  theMassrangeSwitch = cfg.getParameter<bool>( "applyMassrangeFilter" );
28  if (theMassrangeSwitch){
29  theMinMass = cfg.getParameter<double>( "minXMass" );
30  theMaxMass = cfg.getParameter<double>( "maxXMass" );
31  theDaughterMass = cfg.getParameter<double>( "daughterMass" );
32  theCandNumber = cfg.getParameter<unsigned int>( "numberOfCandidates" );//Number of candidates to keep
33  LogDebug("Alignment") << "> Massrange min,max : " << theMinMass << "," << theMaxMass
34  << "\n> Mass of daughter Particle : " << theDaughterMass;
35 
36  }else{
37  theMinMass = 0;
38  theMaxMass = 0;
39  theDaughterMass = 0;
40  }
41  theChargeSwitch = cfg.getParameter<bool>( "applyChargeFilter" );
42  if(theChargeSwitch){
43  theCharge = cfg.getParameter<int>( "charge" );
44  theUnsignedSwitch = cfg.getParameter<bool>( "useUnsignedCharge" );
45  if(theUnsignedSwitch)
47  LogDebug("Alignment") << "> Desired Charge, unsigned: "<<theCharge<<" , "<<theUnsignedSwitch;
48  }else{
49  theCharge =0;
50  theUnsignedSwitch = true;
51  }
52  theMissingETSwitch = cfg.getParameter<bool>( "applyMissingETFilter" );
54  theMissingETSource = cfg.getParameter<InputTag>( "missingETSource" );
55  LogDebug("Alignment") << "> missing Et Source: "<< theMissingETSource;
56  }
57  theAcoplanarityFilterSwitch = cfg.getParameter<bool>( "applyAcoplanarityFilter" );
59  theAcoplanarDistance = cfg.getParameter<double>( "acoplanarDistance" );
60  LogDebug("Alignment") << "> Acoplanar Distance: "<<theAcoplanarDistance;
61  }
62 
63 }
64 
65 // destructor -----------------------------------------------------------------
66 
68 {}
69 
70 
73 {
75 }
76 
77 // do selection ---------------------------------------------------------------
78 
81 {
83 
84  if (theMassrangeSwitch) {
86  result = checkMETMass(result,iEvent);
87  else
88  result = checkMass(result);
89  }
90 
91  LogDebug("Alignment") << "> TwoBodyDecay tracks all,kept: " << tracks.size() << "," << result.size();
92  return result;
93 }
94 
95 template<class T>
96 struct higherTwoBodyDecayPt : public std::binary_function<T,T,bool>
97 {
98  bool operator()( const T& a, const T& b )
99  {
100  return a.first > b.first ;
101  }
102 };
103 
107 {
108  Tracks result;
109 
110  LogDebug("Alignment") <<"> cands size : "<< cands.size();
111 
112  if (cands.size()<2) return result;
113 
114  TLorentzVector track0;
115  TLorentzVector track1;
116  TLorentzVector mother;
117  typedef pair<const reco::Track*,const reco::Track*> constTrackPair;
118  typedef pair<double,constTrackPair> candCollectionItem;
119  vector<candCollectionItem> candCollection;
120 
121  for (unsigned int iCand = 0; iCand < cands.size(); iCand++) {
122 
123  track0.SetXYZT(cands.at(iCand)->px(),
124  cands.at(iCand)->py(),
125  cands.at(iCand)->pz(),
126  sqrt( cands.at(iCand)->p()*cands.at(iCand)->p() + theDaughterMass*theDaughterMass ));
127 
128  for (unsigned int jCand = iCand+1; jCand < cands.size(); jCand++) {
129 
130  track1.SetXYZT(cands.at(jCand)->px(),
131  cands.at(jCand)->py(),
132  cands.at(jCand)->pz(),
133  sqrt( cands.at(jCand)->p()*cands.at(jCand)->p() + theDaughterMass*theDaughterMass ));
134 
135  mother = track0 + track1;
136 
137  const reco::Track *trk1 = cands.at(iCand);
138  const reco::Track *trk2 = cands.at(jCand);
139 
140  bool correctCharge = true;
141  if (theChargeSwitch) correctCharge = this->checkCharge(trk1, trk2);
142 
143  bool acoplanarTracks = true;
144  if (theAcoplanarityFilterSwitch) acoplanarTracks = this->checkAcoplanarity(trk1, trk2);
145 
146  if (mother.M() > theMinMass &&
147  mother.M() < theMaxMass &&
148  correctCharge &&
149  acoplanarTracks) {
150  candCollection.push_back(candCollectionItem(mother.Pt(),
151  constTrackPair(trk1, trk2)));
152  }
153  }
154  }
155 
156  if (candCollection.size()==0) return result;
157 
158  sort(candCollection.begin(), candCollection.end(),
160 
161  std::map<const reco::Track*,unsigned int> uniqueTrackIndex;
162  std::map<const reco::Track*,unsigned int>::iterator it;
163  for (unsigned int i=0;
164  i<candCollection.size() && i<theCandNumber;
165  i++) {
166  constTrackPair & trackPair = candCollection[i].second;
167 
168  it = uniqueTrackIndex.find(trackPair.first);
169  if (it==uniqueTrackIndex.end()) {
170  result.push_back(trackPair.first);
171  uniqueTrackIndex[trackPair.first] = i;
172  }
173 
174  it = uniqueTrackIndex.find(trackPair.second);
175  if (it==uniqueTrackIndex.end()) {
176  result.push_back(trackPair.second);
177  uniqueTrackIndex[trackPair.second] = i;
178  }
179  }
180 
181  return result;
182 }
183 
187 {
188  Tracks result;
189 
190  LogDebug("Alignment") <<"> cands size : "<< cands.size();
191 
192  if (cands.size()==0) return result;
193 
194  TLorentzVector track;
195  TLorentzVector met4;
196  TLorentzVector mother;
197 
199  iEvent.getByLabel(theMissingETSource ,missingET);
200  if (!missingET.isValid()) {
201  LogError("Alignment")<< "@SUB=AlignmentTwoBodyDecayTrackSelector::checkMETMass"
202  << "> could not optain missingET Collection!";
203  return result;
204  }
205 
206  typedef pair<double,const reco::Track*> candCollectionItem;
207  vector<candCollectionItem> candCollection;
208 
209  for (reco::CaloMETCollection::const_iterator itMET = missingET->begin();
210  itMET != missingET->end();
211  ++itMET) {
212 
213  met4.SetXYZT((*itMET).px(),
214  (*itMET).py(),
215  (*itMET).pz(),
216  (*itMET).p());
217 
218  for (unsigned int iCand = 0; iCand < cands.size(); iCand++) {
219 
220  track.SetXYZT(cands.at(iCand)->px(),
221  cands.at(iCand)->py(),
222  cands.at(iCand)->pz(),
223  sqrt( cands.at(iCand)->p()*cands.at(iCand)->p() + theDaughterMass*theDaughterMass ));
224 
225  mother = track + met4;
226 
227  const reco::Track *trk = cands.at(iCand);
228  const reco::CaloMET *met = &(*itMET);
229 
230  bool correctCharge = true;
231  if (theChargeSwitch) correctCharge = this->checkCharge(trk);
232 
233  bool acoplanarTracks = true;
234  if (theAcoplanarityFilterSwitch) acoplanarTracks = this->checkMETAcoplanarity(trk, met);
235 
236  if (mother.M() > theMinMass &&
237  mother.M() < theMaxMass &&
238  correctCharge &&
239  acoplanarTracks) {
240  candCollection.push_back(candCollectionItem(mother.Pt(), trk));
241  }
242  }
243  }
244 
245  if (candCollection.size()==0) return result;
246 
247  sort(candCollection.begin(), candCollection.end(),
249 
250  std::map<const reco::Track*,unsigned int> uniqueTrackIndex;
251  std::map<const reco::Track*,unsigned int>::iterator it;
252  for (unsigned int i=0;
253  i<candCollection.size() && i<theCandNumber;
254  i++) {
255  it = uniqueTrackIndex.find(candCollection[i].second);
256  if (it==uniqueTrackIndex.end()) {
257  result.push_back(candCollection[i].second);
258  uniqueTrackIndex[candCollection[i].second] = i;
259  }
260  }
261 
262  return result;
263 }
264 
266 bool
268 {
269  int sumCharge = trk1->charge();
270  if (trk2) sumCharge += trk2->charge();
271  if (theUnsignedSwitch) sumCharge = std::abs(sumCharge);
272  if (sumCharge == theCharge) return true;
273  return false;
274 }
275 
277 bool
279 {
280  if (fabs(deltaPhi(trk1->phi(),trk2->phi()-M_PI)) < theAcoplanarDistance) return true;
281  return false;
282 }
283 
285 bool
287 {
288  if (fabs(deltaPhi(trk1->phi(),met->phi()-M_PI)) < theAcoplanarDistance) return true;
289  return false;
290 }
291 
292 //===================HELPERS===================
293 
296 {
297  int count = 0;
298  LogDebug("Alignment") << ">......................................";
299  for(Tracks::const_iterator it = col.begin();it < col.end();++it,++count){
300  LogDebug("Alignment")
301  <<"> Track No. "<< count <<": p = ("<<(*it)->px()<<","<<(*it)->py()<<","<<(*it)->pz()<<")\n"
302  <<"> pT = "<<(*it)->pt()<<" eta = "<<(*it)->eta()<<" charge = "<<(*it)->charge();
303  }
304  LogDebug("Alignment") << ">......................................";
305 }
306 
307 
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
Tracks select(const Tracks &tracks, const edm::Event &iEvent)
select tracks
#define abs(x)
Definition: mlp_lapack.h:159
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:139
bool checkMETAcoplanarity(const reco::Track *trk, const reco::CaloMET *met) const
checks if [cands] contains a acoplanar track w.r.t missing ET (returns empty set if not) ...
bool checkAcoplanarity(const reco::Track *trk1, const reco::Track *trk2) const
checks if the [cands] are acoplanar (returns empty set if not)
U second(std::pair< T, U > const &p)
int iEvent
Definition: GenABIO.cc:243
T sqrt(T t)
Definition: SSEVec.h:46
bool checkCharge(const reco::Track *trk1, const reco::Track *trk2=0) const
checks if the mother has charge = [theCharge]
bool useThisFilter()
returns if any of the Filters is used.
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
bool operator()(const T &a, const T &b)
#define M_PI
Definition: BFit3D.cc:3
tuple tracks
Definition: testEve_cfg.py:39
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
AlignmentTwoBodyDecayTrackSelector(const edm::ParameterSet &cfg)
constructor
int charge() const
track electric charge
Definition: TrackBase.h:113
long double T
virtual double phi() const
momentum azimuthal angle
Tracks checkMETMass(const Tracks &cands, const edm::Event &iEvent) const
checks if the mass of the mother is in the mass region adding missing E_T
Tracks checkMass(const Tracks &cands) const
checks if the mass of the mother is in the mass region