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  secThrBool = cfg.getParameter<bool> ( "applySecThreshold" );
34  thesecThr = cfg.getParameter<double>( "secondThreshold" );
35  LogDebug("Alignment") << "> Massrange min,max : " << theMinMass << "," << theMaxMass
36  << "\n> Mass of daughter Particle : " << theDaughterMass;
37 
38  }else{
39  theMinMass = 0;
40  theMaxMass = 0;
41  theDaughterMass = 0;
42  }
43  theChargeSwitch = cfg.getParameter<bool>( "applyChargeFilter" );
44  if(theChargeSwitch){
45  theCharge = cfg.getParameter<int>( "charge" );
46  theUnsignedSwitch = cfg.getParameter<bool>( "useUnsignedCharge" );
47  if(theUnsignedSwitch)
49  LogDebug("Alignment") << "> Desired Charge, unsigned: "<<theCharge<<" , "<<theUnsignedSwitch;
50  }else{
51  theCharge =0;
52  theUnsignedSwitch = true;
53  }
54  theMissingETSwitch = cfg.getParameter<bool>( "applyMissingETFilter" );
56  theMissingETSource = cfg.getParameter<InputTag>( "missingETSource" );
57  LogDebug("Alignment") << "> missing Et Source: "<< theMissingETSource;
58  }
59  theAcoplanarityFilterSwitch = cfg.getParameter<bool>( "applyAcoplanarityFilter" );
61  theAcoplanarDistance = cfg.getParameter<double>( "acoplanarDistance" );
62  LogDebug("Alignment") << "> Acoplanar Distance: "<<theAcoplanarDistance;
63  }
64 
65 }
66 
67 // destructor -----------------------------------------------------------------
68 
70 {}
71 
72 
75 {
77 }
78 
79 // do selection ---------------------------------------------------------------
80 
83 {
85 
86  if (theMassrangeSwitch) {
88  result = checkMETMass(result,iEvent);
89  else
90  result = checkMass(result);
91  }
92 
93  LogDebug("Alignment") << "> TwoBodyDecay tracks all,kept: " << tracks.size() << "," << result.size();
94  return result;
95 }
96 
97 template<class T>
98 struct higherTwoBodyDecayPt : public std::binary_function<T,T,bool>
99 {
100  bool operator()( const T& a, const T& b )
101  {
102  return a.first > b.first ;
103  }
104 };
105 
109 {
110  Tracks result;
111 
112  LogDebug("Alignment") <<"> cands size : "<< cands.size();
113 
114  if (cands.size()<2) return result;
115 
116  TLorentzVector track0;
117  TLorentzVector track1;
118  TLorentzVector mother;
119  typedef pair<const reco::Track*,const reco::Track*> constTrackPair;
120  typedef pair<double,constTrackPair> candCollectionItem;
121  vector<candCollectionItem> candCollection;
122 
123  for (unsigned int iCand = 0; iCand < cands.size(); iCand++) {
124 
125  track0.SetXYZT(cands.at(iCand)->px(),
126  cands.at(iCand)->py(),
127  cands.at(iCand)->pz(),
128  sqrt( cands.at(iCand)->p()*cands.at(iCand)->p() + theDaughterMass*theDaughterMass ));
129 
130  for (unsigned int jCand = iCand+1; jCand < cands.size(); jCand++) {
131 
132  track1.SetXYZT(cands.at(jCand)->px(),
133  cands.at(jCand)->py(),
134  cands.at(jCand)->pz(),
135  sqrt( cands.at(jCand)->p()*cands.at(jCand)->p() + theDaughterMass*theDaughterMass ));
136  if (secThrBool==true && track1.Pt() < thesecThr && track0.Pt()< thesecThr) continue;
137  mother = track0 + track1;
138 
139  const reco::Track *trk1 = cands.at(iCand);
140  const reco::Track *trk2 = cands.at(jCand);
141 
142  bool correctCharge = true;
143  if (theChargeSwitch) correctCharge = this->checkCharge(trk1, trk2);
144 
145  bool acoplanarTracks = true;
146  if (theAcoplanarityFilterSwitch) acoplanarTracks = this->checkAcoplanarity(trk1, trk2);
147 
148  if (mother.M() > theMinMass &&
149  mother.M() < theMaxMass &&
150  correctCharge &&
151  acoplanarTracks) {
152  candCollection.push_back(candCollectionItem(mother.Pt(),
153  constTrackPair(trk1, trk2)));
154  }
155  }
156  }
157 
158  if (candCollection.size()==0) return result;
159 
160  sort(candCollection.begin(), candCollection.end(),
162 
163  std::map<const reco::Track*,unsigned int> uniqueTrackIndex;
164  std::map<const reco::Track*,unsigned int>::iterator it;
165  for (unsigned int i=0;
166  i<candCollection.size() && i<theCandNumber;
167  i++) {
168  constTrackPair & trackPair = candCollection[i].second;
169 
170  it = uniqueTrackIndex.find(trackPair.first);
171  if (it==uniqueTrackIndex.end()) {
172  result.push_back(trackPair.first);
173  uniqueTrackIndex[trackPair.first] = i;
174  }
175 
176  it = uniqueTrackIndex.find(trackPair.second);
177  if (it==uniqueTrackIndex.end()) {
178  result.push_back(trackPair.second);
179  uniqueTrackIndex[trackPair.second] = i;
180  }
181  }
182 
183  return result;
184 }
185 
189 {
190  Tracks result;
191 
192  LogDebug("Alignment") <<"> cands size : "<< cands.size();
193 
194  if (cands.size()==0) return result;
195 
196  TLorentzVector track;
197  TLorentzVector met4;
198  TLorentzVector mother;
199 
201  iEvent.getByLabel(theMissingETSource ,missingET);
202  if (!missingET.isValid()) {
203  LogError("Alignment")<< "@SUB=AlignmentTwoBodyDecayTrackSelector::checkMETMass"
204  << "> could not optain missingET Collection!";
205  return result;
206  }
207 
208  typedef pair<double,const reco::Track*> candCollectionItem;
209  vector<candCollectionItem> candCollection;
210 
211  for (reco::CaloMETCollection::const_iterator itMET = missingET->begin();
212  itMET != missingET->end();
213  ++itMET) {
214 
215  met4.SetXYZT((*itMET).px(),
216  (*itMET).py(),
217  (*itMET).pz(),
218  (*itMET).p());
219 
220  for (unsigned int iCand = 0; iCand < cands.size(); iCand++) {
221 
222  track.SetXYZT(cands.at(iCand)->px(),
223  cands.at(iCand)->py(),
224  cands.at(iCand)->pz(),
225  sqrt( cands.at(iCand)->p()*cands.at(iCand)->p() + theDaughterMass*theDaughterMass ));
226 
227  mother = track + met4;
228 
229  const reco::Track *trk = cands.at(iCand);
230  const reco::CaloMET *met = &(*itMET);
231 
232  bool correctCharge = true;
233  if (theChargeSwitch) correctCharge = this->checkCharge(trk);
234 
235  bool acoplanarTracks = true;
236  if (theAcoplanarityFilterSwitch) acoplanarTracks = this->checkMETAcoplanarity(trk, met);
237 
238  if (mother.M() > theMinMass &&
239  mother.M() < theMaxMass &&
240  correctCharge &&
241  acoplanarTracks) {
242  candCollection.push_back(candCollectionItem(mother.Pt(), trk));
243  }
244  }
245  }
246 
247  if (candCollection.size()==0) return result;
248 
249  sort(candCollection.begin(), candCollection.end(),
251 
252  std::map<const reco::Track*,unsigned int> uniqueTrackIndex;
253  std::map<const reco::Track*,unsigned int>::iterator it;
254  for (unsigned int i=0;
255  i<candCollection.size() && i<theCandNumber;
256  i++) {
257  it = uniqueTrackIndex.find(candCollection[i].second);
258  if (it==uniqueTrackIndex.end()) {
259  result.push_back(candCollection[i].second);
260  uniqueTrackIndex[candCollection[i].second] = i;
261  }
262  }
263 
264  return result;
265 }
266 
268 bool
270 {
271  int sumCharge = trk1->charge();
272  if (trk2) sumCharge += trk2->charge();
273  if (theUnsignedSwitch) sumCharge = std::abs(sumCharge);
274  if (sumCharge == theCharge) return true;
275  return false;
276 }
277 
279 bool
281 {
282  if (fabs(deltaPhi(trk1->phi(),trk2->phi()-M_PI)) < theAcoplanarDistance) return true;
283  return false;
284 }
285 
287 bool
289 {
290  if (fabs(deltaPhi(trk1->phi(),met->phi()-M_PI)) < theAcoplanarDistance) return true;
291  return false;
292 }
293 
294 //===================HELPERS===================
295 
298 {
299  int count = 0;
300  LogDebug("Alignment") << ">......................................";
301  for(Tracks::const_iterator it = col.begin();it < col.end();++it,++count){
302  LogDebug("Alignment")
303  <<"> Track No. "<< count <<": p = ("<<(*it)->px()<<","<<(*it)->py()<<","<<(*it)->pz()<<")\n"
304  <<"> pT = "<<(*it)->pt()<<" eta = "<<(*it)->eta()<<" charge = "<<(*it)->charge();
305  }
306  LogDebug("Alignment") << ">......................................";
307 }
308 
309 
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
tuple met
____________________________________________________________________________||
Definition: CaloMET_cfi.py:7
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:137
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)
virtual float phi() const GCC11_FINAL
momentum azimuthal angle
int iEvent
Definition: GenABIO.cc:243
T sqrt(T t)
Definition: SSEVec.h:48
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
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Tracks select(const Tracks &tracks, const edm::Event &iEvent, const edm::EventSetup &iSetup)
select tracks
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
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:111
int col
Definition: cuy.py:1008
long double T
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