CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/GeneratorInterface/GenFilters/src/XtoFFbarFilter.cc

Go to the documentation of this file.
00001 #include "GeneratorInterface/GenFilters/interface/XtoFFbarFilter.h"
00002 
00003 using namespace std;
00004 using namespace reco;
00005 
00006 XtoFFbarFilter::XtoFFbarFilter(const edm::ParameterSet& iConfig) :
00007   src_(iConfig.getParameter<edm::InputTag>("src")),
00008   idMotherX_(iConfig.getParameter<vector<int> >("idMotherX")),
00009   idDaughterF_(iConfig.getParameter<vector<int> >("idDaughterF")),
00010   idMotherY_(iConfig.getParameter<vector<int> >("idMotherY")),
00011   idDaughterG_(iConfig.getParameter<vector<int> >("idDaughterG")),
00012   xTotal_(0), xSumPt_(0.), xSumR_(0.), totalEvents_(0),rejectedEvents_(0)
00013 { 
00014   // Note if if not searching for Y --> g-gbar.
00015   requireY_ = (idMotherY_.size() > 0 && idDaughterG_.size() > 0);
00016 }
00017 
00018 
00019 bool XtoFFbarFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup)  {
00020 
00021   iEvent.getByLabel(src_, genParticles_);
00022 
00023   totalEvents_++;
00024 
00025   unsigned int numX = 0;
00026   unsigned int numY = 0;
00027   unsigned int numXorY = 0;
00028 
00029   for (unsigned int j = 0; j < genParticles_->size(); j++) {
00030     GenParticleRef moth(genParticles_, j);
00031 
00032     // Is it X -> f fbar ?
00033     bool isXtoFFbar = this->foundXtoFFbar(moth, idMotherX_, idDaughterF_);
00034     if (isXtoFFbar) numX++; 
00035 
00036     if (!requireY_) {
00037 
00038       // Has X been found already ?
00039       if (numX >= 1) return true;
00040 
00041     } else {
00042 
00043       // Is it Y -> g gbar ?
00044       bool isYtoGGbar = this->foundXtoFFbar(moth, idMotherY_, idDaughterG_);
00045       if (isYtoGGbar) numY++;
00046       if (isXtoFFbar || isYtoGGbar) numXorY++;
00047 
00048       // Have X and Y been found already ?
00049       if (numX >= 1 && numY >= 1 && numXorY >= 2) return true;
00050     }
00051   }
00052 
00053   rejectedEvents_++;
00054   //  cout<<"REJECTED "<<totalEvents_<<endl;
00055   return false;
00056 }
00057 
00058 
00059 bool XtoFFbarFilter::foundXtoFFbar(const GenParticleRef& moth, 
00060                                    const vector<int>& idMother, 
00061                                    const vector<int>& idDaughter) {
00062   // Check if given particle "moth" is X-->f fbar
00063   bool isXtoFFbar = false;
00064   int pdgIdMoth = moth->pdgId();
00065   double rho = -9.9e9;
00066 
00067   if (this->found(idMother, pdgIdMoth)) {
00068     bool foundF    = false;
00069     bool foundFbar = false;
00070     unsigned int nDau = moth->numberOfDaughters();
00071 
00072     for (unsigned int i = 0; i < nDau; i++) {
00073       GenParticleRef dau = moth->daughterRef(i);
00074       int pdgIdDau = dau->pdgId();
00075       if (this->found(idDaughter, -pdgIdDau))  foundFbar = true; 
00076       if (this->found(idDaughter,  pdgIdDau)) {foundF    = true; 
00077 
00078         // Just for statistics, get transverse decay length.
00079         // This is the normal case
00080         rho = dau->vertex().Rho();
00081         // Unfortunately, duplicate particles can appear in the event record. Handle this as follows:
00082         for (unsigned int j = 0; j < dau->numberOfDaughters(); j++) {
00083           GenParticleRef granddau = dau->daughterRef(j);        
00084           if (granddau->pdgId() == pdgIdDau) rho = granddau->vertex().Rho();
00085         }
00086       }
00087     }
00088     if (foundF && foundFbar) isXtoFFbar = true;
00089   }
00090 
00091   if (isXtoFFbar) {
00092     // Get statistics 
00093     xTotal_++;
00094     xSumPt_ += moth->pt();
00095     xSumR_  += rho;
00096   }
00097 
00098   return isXtoFFbar;
00099 }
00100 
00101 void XtoFFbarFilter::endJob() {
00102   cout<<endl;
00103   cout<<"=== XtoFFbarFilter statistics of selected X->ffbar or Y->ggbar"<<endl;
00104   if (xTotal_ > 0) {
00105     cout<<"===   mean X & Y Pt = "<<xSumPt_/xTotal_<<" and transverse decay length = "<<xSumR_/xTotal_<<endl;
00106   } else {
00107     cout<<"===   WARNING: NONE FOUND !"<<endl;
00108   }
00109   cout<<"===   events rejected = "<<rejectedEvents_<<"/"<<totalEvents_<<endl;                     
00110 }