CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/Alignment/CommonAlignmentProducer/src/AlignmentTwoBoyDecayTrackSelector.cc

Go to the documentation of this file.
00001 //Framework
00002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00003 
00004 #include "FWCore/Framework/interface/Event.h"
00005 #include "FWCore/Utilities/interface/EDMException.h"
00006 
00007 //DataFormats
00008 #include <DataFormats/TrackReco/interface/Track.h>
00009 #include <DataFormats/METReco/interface/CaloMET.h>
00010 #include <DataFormats/Math/interface/deltaPhi.h>
00011 
00012 //STL
00013 #include <math.h>
00014 //ROOT
00015 #include "TLorentzVector.h"
00016 
00017 #include "Alignment/CommonAlignmentProducer/interface/AlignmentTwoBodyDecayTrackSelector.h"
00018 //TODO put those namespaces into functions?
00019 using namespace std;
00020 using namespace edm; 
00021 // constructor ----------------------------------------------------------------
00022 
00023 AlignmentTwoBodyDecayTrackSelector::AlignmentTwoBodyDecayTrackSelector(const edm::ParameterSet & cfg) :
00024   theMissingETSource("met")
00025 {
00026  LogDebug("Alignment")   << "> applying two body decay Trackfilter ...";
00027   theMassrangeSwitch = cfg.getParameter<bool>( "applyMassrangeFilter" );
00028   if (theMassrangeSwitch){
00029     theMinMass = cfg.getParameter<double>( "minXMass" );
00030     theMaxMass = cfg.getParameter<double>( "maxXMass" );
00031     theDaughterMass = cfg.getParameter<double>( "daughterMass" );
00032     theCandNumber = cfg.getParameter<unsigned int>( "numberOfCandidates" );//Number of candidates to keep
00033     LogDebug("Alignment") << ">  Massrange min,max         :   " << theMinMass   << "," << theMaxMass 
00034                          << "\n>  Mass of daughter Particle :   " << theDaughterMass;
00035 
00036   }else{
00037     theMinMass = 0;
00038     theMaxMass = 0;
00039     theDaughterMass = 0;
00040   }
00041   theChargeSwitch = cfg.getParameter<bool>( "applyChargeFilter" );
00042   if(theChargeSwitch){
00043     theCharge = cfg.getParameter<int>( "charge" );
00044     theUnsignedSwitch = cfg.getParameter<bool>( "useUnsignedCharge" );
00045     if(theUnsignedSwitch) 
00046       theCharge=std::abs(theCharge);
00047     LogDebug("Alignment") << ">  Desired Charge, unsigned: "<<theCharge<<" , "<<theUnsignedSwitch;
00048   }else{
00049     theCharge =0;
00050     theUnsignedSwitch = true;
00051   }
00052   theMissingETSwitch = cfg.getParameter<bool>( "applyMissingETFilter" );
00053   if(theMissingETSwitch){
00054     theMissingETSource = cfg.getParameter<InputTag>( "missingETSource" );
00055     LogDebug("Alignment") << ">  missing Et Source: "<< theMissingETSource;
00056   }
00057   theAcoplanarityFilterSwitch = cfg.getParameter<bool>( "applyAcoplanarityFilter" );
00058   if(theAcoplanarityFilterSwitch){
00059     theAcoplanarDistance = cfg.getParameter<double>( "acoplanarDistance" );
00060     LogDebug("Alignment") << ">  Acoplanar Distance: "<<theAcoplanarDistance;
00061   }
00062   
00063 }
00064 
00065 // destructor -----------------------------------------------------------------
00066 
00067 AlignmentTwoBodyDecayTrackSelector::~AlignmentTwoBodyDecayTrackSelector()
00068 {}
00069 
00070 
00072 bool AlignmentTwoBodyDecayTrackSelector::useThisFilter()
00073 {
00074   return theMassrangeSwitch || theChargeSwitch || theAcoplanarityFilterSwitch;
00075 }
00076 
00077 // do selection ---------------------------------------------------------------
00078 
00079 AlignmentTwoBodyDecayTrackSelector::Tracks 
00080 AlignmentTwoBodyDecayTrackSelector::select(const Tracks& tracks, const edm::Event& iEvent) 
00081 {
00082   Tracks result = tracks;
00083 
00084   if (theMassrangeSwitch) {  
00085     if (theMissingETSwitch)
00086       result = checkMETMass(result,iEvent); 
00087     else
00088       result = checkMass(result); 
00089   }
00090 
00091   LogDebug("Alignment") << ">  TwoBodyDecay tracks all,kept: " << tracks.size() << "," << result.size();
00092   return result;
00093 }
00094 
00095 template<class T>
00096 struct higherTwoBodyDecayPt : public std::binary_function<T,T,bool>
00097 {
00098   bool operator()( const T& a, const T& b ) 
00099   { 
00100     return a.first > b.first ; 
00101   }
00102 };
00103 
00105 AlignmentTwoBodyDecayTrackSelector::Tracks 
00106 AlignmentTwoBodyDecayTrackSelector::checkMass(const Tracks& cands) const
00107 {
00108   Tracks result;
00109   
00110   LogDebug("Alignment") <<">  cands size : "<< cands.size();
00111   
00112   if (cands.size()<2) return result;
00113 
00114   TLorentzVector track0;
00115   TLorentzVector track1;
00116   TLorentzVector mother;
00117   typedef pair<const reco::Track*,const reco::Track*> constTrackPair;
00118   typedef pair<double,constTrackPair> candCollectionItem;
00119   vector<candCollectionItem> candCollection;
00120   
00121   for (unsigned int iCand = 0; iCand < cands.size(); iCand++) {
00122     
00123     track0.SetXYZT(cands.at(iCand)->px(),
00124                    cands.at(iCand)->py(),
00125                    cands.at(iCand)->pz(),
00126                    sqrt( cands.at(iCand)->p()*cands.at(iCand)->p() + theDaughterMass*theDaughterMass ));
00127     
00128     for (unsigned int jCand = iCand+1; jCand < cands.size(); jCand++) {
00129       
00130       track1.SetXYZT(cands.at(jCand)->px(),
00131                      cands.at(jCand)->py(),
00132                      cands.at(jCand)->pz(),
00133                      sqrt( cands.at(jCand)->p()*cands.at(jCand)->p() + theDaughterMass*theDaughterMass ));
00134       
00135       mother = track0 + track1;
00136       
00137       const reco::Track *trk1 = cands.at(iCand);
00138       const reco::Track *trk2 = cands.at(jCand);
00139 
00140       bool correctCharge = true;
00141       if (theChargeSwitch) correctCharge = this->checkCharge(trk1, trk2);
00142 
00143       bool acoplanarTracks = true;
00144       if (theAcoplanarityFilterSwitch) acoplanarTracks = this->checkAcoplanarity(trk1, trk2);
00145 
00146       if (mother.M() > theMinMass &&
00147           mother.M() < theMaxMass &&
00148           correctCharge &&
00149           acoplanarTracks) {
00150         candCollection.push_back(candCollectionItem(mother.Pt(),
00151                                                     constTrackPair(trk1, trk2)));
00152       }
00153     }
00154   }
00155 
00156   if (candCollection.size()==0) return result;
00157 
00158   sort(candCollection.begin(), candCollection.end(), 
00159        higherTwoBodyDecayPt<candCollectionItem>());
00160 
00161   std::map<const reco::Track*,unsigned int> uniqueTrackIndex;
00162   std::map<const reco::Track*,unsigned int>::iterator it;
00163   for (unsigned int i=0;
00164        i<candCollection.size() && i<theCandNumber;
00165        i++) {
00166     constTrackPair & trackPair = candCollection[i].second;
00167     
00168     it = uniqueTrackIndex.find(trackPair.first);
00169     if (it==uniqueTrackIndex.end()) {
00170       result.push_back(trackPair.first);
00171       uniqueTrackIndex[trackPair.first] = i;
00172     }
00173     
00174     it = uniqueTrackIndex.find(trackPair.second);
00175     if (it==uniqueTrackIndex.end()) {
00176       result.push_back(trackPair.second);
00177       uniqueTrackIndex[trackPair.second] = i;
00178     }
00179   }
00180 
00181   return result;
00182 }
00183 
00185 AlignmentTwoBodyDecayTrackSelector::Tracks 
00186 AlignmentTwoBodyDecayTrackSelector::checkMETMass(const Tracks& cands,const edm::Event& iEvent) const
00187 {
00188   Tracks result;
00189   
00190   LogDebug("Alignment") <<">  cands size : "<< cands.size();
00191   
00192   if (cands.size()==0) return result;
00193 
00194   TLorentzVector track;
00195   TLorentzVector met4;
00196   TLorentzVector mother;
00197 
00198   Handle<reco::CaloMETCollection> missingET;
00199   iEvent.getByLabel(theMissingETSource ,missingET);
00200   if (!missingET.isValid()) {
00201     LogError("Alignment")<< "@SUB=AlignmentTwoBodyDecayTrackSelector::checkMETMass"
00202                          << ">  could not optain missingET Collection!";
00203     return result;
00204   }
00205 
00206   typedef pair<double,const reco::Track*> candCollectionItem;
00207   vector<candCollectionItem> candCollection;
00208 
00209   for (reco::CaloMETCollection::const_iterator itMET = missingET->begin();
00210        itMET != missingET->end();
00211        ++itMET) {
00212     
00213     met4.SetXYZT((*itMET).px(),
00214                  (*itMET).py(),
00215                  (*itMET).pz(),
00216                  (*itMET).p());
00217   
00218     for (unsigned int iCand = 0; iCand < cands.size(); iCand++) {
00219     
00220       track.SetXYZT(cands.at(iCand)->px(),
00221                     cands.at(iCand)->py(),
00222                     cands.at(iCand)->pz(),
00223                     sqrt( cands.at(iCand)->p()*cands.at(iCand)->p() + theDaughterMass*theDaughterMass ));
00224       
00225       mother = track + met4;
00226       
00227       const reco::Track *trk = cands.at(iCand);
00228       const reco::CaloMET *met = &(*itMET);
00229 
00230       bool correctCharge = true;
00231       if (theChargeSwitch) correctCharge = this->checkCharge(trk);
00232 
00233       bool acoplanarTracks = true;
00234       if (theAcoplanarityFilterSwitch) acoplanarTracks = this->checkMETAcoplanarity(trk, met);
00235 
00236       if (mother.M() > theMinMass &&
00237           mother.M() < theMaxMass &&
00238           correctCharge &&
00239           acoplanarTracks) {
00240         candCollection.push_back(candCollectionItem(mother.Pt(), trk));
00241       }
00242     }
00243   }
00244 
00245   if (candCollection.size()==0) return result;
00246 
00247   sort(candCollection.begin(), candCollection.end(), 
00248        higherTwoBodyDecayPt<candCollectionItem>());
00249   
00250   std::map<const reco::Track*,unsigned int> uniqueTrackIndex;
00251   std::map<const reco::Track*,unsigned int>::iterator it;
00252   for (unsigned int i=0;
00253        i<candCollection.size() && i<theCandNumber;
00254        i++) {
00255     it = uniqueTrackIndex.find(candCollection[i].second);
00256     if (it==uniqueTrackIndex.end()) {
00257       result.push_back(candCollection[i].second);
00258       uniqueTrackIndex[candCollection[i].second] = i;
00259     }
00260   }
00261 
00262   return result;
00263 }
00264 
00266 bool
00267 AlignmentTwoBodyDecayTrackSelector::checkCharge(const reco::Track* trk1, const reco::Track* trk2)const
00268 {
00269   int sumCharge = trk1->charge();
00270   if (trk2) sumCharge += trk2->charge();
00271   if (theUnsignedSwitch) sumCharge = std::abs(sumCharge);
00272   if (sumCharge == theCharge) return true;
00273   return false;
00274 }
00275 
00277 bool
00278 AlignmentTwoBodyDecayTrackSelector::checkAcoplanarity(const reco::Track* trk1, const reco::Track* trk2)const
00279 {
00280   if (fabs(deltaPhi(trk1->phi(),trk2->phi()-M_PI)) < theAcoplanarDistance) return true;
00281   return false;
00282 }
00283 
00285 bool
00286 AlignmentTwoBodyDecayTrackSelector::checkMETAcoplanarity(const reco::Track* trk1, const reco::CaloMET* met)const
00287 {
00288   if (fabs(deltaPhi(trk1->phi(),met->phi()-M_PI)) < theAcoplanarDistance) return true;
00289   return false;
00290 }
00291 
00292 //===================HELPERS===================
00293 
00295 void AlignmentTwoBodyDecayTrackSelector::printTracks(const Tracks& col) const
00296 {
00297   int count = 0;
00298   LogDebug("Alignment") << ">......................................";
00299   for(Tracks::const_iterator it = col.begin();it < col.end();++it,++count){
00300     LogDebug("Alignment") 
00301       <<">  Track No. "<< count <<": p = ("<<(*it)->px()<<","<<(*it)->py()<<","<<(*it)->pz()<<")\n"
00302       <<">                        pT = "<<(*it)->pt()<<" eta = "<<(*it)->eta()<<" charge = "<<(*it)->charge();    
00303   }
00304   LogDebug("Alignment") << ">......................................";
00305 }
00306 
00307