CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/HLTrigger/HLTanalyzers/src/HLTMCtruth.cc

Go to the documentation of this file.
00001 #include <iostream>
00002 #include <sstream>
00003 #include <istream>
00004 #include <fstream>
00005 #include <iomanip>
00006 #include <string>
00007 #include <cmath>
00008 #include <functional>
00009 #include <stdlib.h>
00010 #include <string.h>
00011 
00012 #include "HLTrigger/HLTanalyzers/interface/HLTMCtruth.h"
00013 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00014 
00015 HLTMCtruth::HLTMCtruth() {
00016 
00017   //set parameter defaults 
00018   _Monte=false;
00019   _Debug=false;
00020 }
00021 
00022 /*  Setup the analysis to put the branch-variables into the tree. */
00023 void HLTMCtruth::setup(const edm::ParameterSet& pSet, TTree* HltTree) {
00024 
00025   edm::ParameterSet myMCParams = pSet.getParameter<edm::ParameterSet>("RunParameters") ;
00026   std::vector<std::string> parameterNames = myMCParams.getParameterNames() ;
00027   
00028   for ( std::vector<std::string>::iterator iParam = parameterNames.begin();
00029         iParam != parameterNames.end(); iParam++ ){
00030     if  ( (*iParam) == "Monte" ) _Monte =  myMCParams.getParameter<bool>( *iParam );
00031     else if ( (*iParam) == "Debug" ) _Debug =  myMCParams.getParameter<bool>( *iParam );
00032   }
00033 
00034   const int kMaxMcTruth = 10000;
00035   mcpid = new int[kMaxMcTruth];
00036   mcstatus = new int[kMaxMcTruth];
00037   mcvx = new float[kMaxMcTruth];
00038   mcvy = new float[kMaxMcTruth];
00039   mcvz = new float[kMaxMcTruth];
00040   mcpt = new float[kMaxMcTruth];
00041   mceta = new float[kMaxMcTruth];
00042   mcphi = new float[kMaxMcTruth];
00043 
00044   // MCtruth-specific branches of the tree 
00045   HltTree->Branch("NMCpart",&nmcpart,"NMCpart/I");
00046   HltTree->Branch("MCpid",mcpid,"MCpid[NMCpart]/I");
00047   HltTree->Branch("MCstatus",mcstatus,"MCstatus[NMCpart]/I");
00048   HltTree->Branch("MCvtxX",mcvx,"MCvtxX[NMCpart]/F");
00049   HltTree->Branch("MCvtxY",mcvy,"MCvtxY[NMCpart]/F");
00050   HltTree->Branch("MCvtxZ",mcvz,"MCvtxZ[NMCpart]/F");
00051   HltTree->Branch("MCpt",mcpt,"MCpt[NMCpart]/F");
00052   HltTree->Branch("MCeta",mceta,"MCeta[NMCpart]/F");
00053   HltTree->Branch("MCphi",mcphi,"MCphi[NMCpart]/F");
00054   HltTree->Branch("MCPtHat",&pthatf,"MCPtHat/F");
00055   HltTree->Branch("MCmu3",&nmu3,"MCmu3/I");
00056   HltTree->Branch("MCel3",&nel3,"MCel3/I");
00057   HltTree->Branch("MCbb",&nbb,"MCbb/I");
00058   HltTree->Branch("MCab",&nab,"MCab/I");
00059   HltTree->Branch("MCWenu",&nwenu,"MCWenu/I");
00060   HltTree->Branch("MCWmunu",&nwmunu,"MCmunu/I");
00061   HltTree->Branch("MCZee",&nzee,"MCZee/I");
00062   HltTree->Branch("MCZmumu",&nzmumu,"MCZmumu/I");
00063   HltTree->Branch("MCptEleMax",&ptEleMax,"MCptEleMax/F");
00064   HltTree->Branch("MCptMuMax",&ptMuMax,"MCptMuMax/F");
00065   HltTree->Branch("NPUTrueBX0",&npubx0, "NPUTrueBX0/I");
00066   HltTree->Branch("NPUgenBX0",&npuvertbx0, "NPUgenBX0/I");
00067 }
00068 
00069 /* **Analyze the event** */
00070 void HLTMCtruth::analyze(const edm::Handle<reco::CandidateView> & mctruth,
00071                          const double        & pthat,
00072                          const edm::Handle<std::vector<SimTrack> > & simTracks,
00073                          const edm::Handle<std::vector<SimVertex> > & simVertices,
00074                          const edm::Handle<std::vector< PileupSummaryInfo > > & PupInfo,
00075                          TTree* HltTree) {
00076 
00077   //std::cout << " Beginning HLTMCtruth " << std::endl;
00078 
00079   if (_Monte) {
00080     int nmc = 0;
00081     int mu3 = 0;
00082     int el3 = 0;
00083     int mab = 0;
00084     int mbb = 0;
00085     int wel = 0;
00086     int wmu = 0;
00087     int zee = 0;
00088     int zmumu = 0;
00089 
00090     ptEleMax = -999.0;
00091     ptMuMax  = -999.0;    
00092     pthatf   = pthat;
00093     npubx0  = 0.0;
00094     npuvertbx0 = 0;
00095 
00096     int npvtrue = 0; 
00097     int npuvert = 0;
00098     
00099 
00100     if((simTracks.isValid())&&(simVertices.isValid())){
00101       for (unsigned int j=0; j<simTracks->size(); j++) {
00102         int pdgid = simTracks->at(j).type();
00103         if (abs(pdgid)!=13) continue;
00104         double pt = simTracks->at(j).momentum().pt();
00105         if (pt<5.0) continue;
00106         double eta = simTracks->at(j).momentum().eta();
00107         if (abs(eta)>2.5) continue;
00108         if (simTracks->at(j).noVertex()) continue;
00109         int vertIndex = simTracks->at(j).vertIndex();
00110         double x = simVertices->at(vertIndex).position().x();
00111         double y = simVertices->at(vertIndex).position().y();
00112         double r = sqrt(x*x+y*y);
00113         if (r>150.) continue; // I think units are cm here
00114         double z = simVertices->at(vertIndex).position().z();
00115         if (abs(z)>300.) continue; // I think units are cm here
00116         mu3 += 1;
00117         break;
00118       }
00119 
00120 
00121       std::vector<PileupSummaryInfo>::const_iterator PVI;  
00122       for(PVI = PupInfo->begin(); PVI != PupInfo->end(); ++PVI) {  
00123         
00124         int BX = PVI->getBunchCrossing();  
00125         npvtrue = PVI->getTrueNumInteractions();  
00126         npuvert = PVI->getPU_NumInteractions();
00127         
00128         if(BX == 0)  
00129           {  
00130             npubx0+=npvtrue; 
00131             npuvertbx0+=npuvert;
00132           } 
00133       }  
00134       
00135     }
00136 
00137     if (mctruth.isValid()){
00138 
00139       for (size_t i = 0; i < mctruth->size(); ++ i) {
00140         const reco::Candidate & p = (*mctruth)[i];
00141 
00142         mcpid[nmc] = p.pdgId();
00143         mcstatus[nmc] = p.status();
00144         mcpt[nmc] = p.pt();
00145         mceta[nmc] = p.eta();
00146         mcphi[nmc] = p.phi();
00147         mcvx[nmc] = p.vx();
00148         mcvy[nmc] = p.vy();
00149         mcvz[nmc] = p.vz();
00150 
00151         if ((mcpid[nmc]==24)||(mcpid[nmc]==-24)) { // Checking W -> e/mu nu
00152           size_t idg = p.numberOfDaughters();
00153           for (size_t j=0; j != idg; ++j){
00154             const reco::Candidate & d = *p.daughter(j);
00155             if ((d.pdgId()==11)||(d.pdgId()==-11)){wel += 1;}
00156             if ((d.pdgId()==13)||(d.pdgId()==-13)){wmu += 1;}
00157 //          if ( (abs(d.pdgId())!=24) && ((mcpid[nmc])*(d.pdgId())>0) ) 
00158 //            {std::cout << "Wrong sign between mother-W and daughter !" << std::endl;}
00159           }
00160         }
00161         if (mcpid[nmc]==23) { // Checking Z -> 2 e/mu
00162           size_t idg = p.numberOfDaughters();
00163           for (size_t j=0; j != idg; ++j){
00164             const reco::Candidate & d = *p.daughter(j);
00165             if (d.pdgId()==11){zee += 1;}
00166             if (d.pdgId()==-11){zee += 2;}
00167             if (d.pdgId()==13){zmumu += 1;}
00168             if (d.pdgId()==-13){zmumu += 2;}
00169           }
00170         }
00171 
00172         // Set-up flags, based on Pythia-generator information, for avoiding double-counting events when
00173         // using both pp->{e,mu}X AND QCD samples
00174         //if (((mcpid[nmc]==13)||(mcpid[nmc]==-13))&&(mcpt[nmc]>5.0)) {mu3 += 1;} // Flag for muons with pT > 2.5 GeV/c
00175         if (((mcpid[nmc]==11)||(mcpid[nmc]==-11))&&(mcpt[nmc]>5.0)) {el3 += 1;} // Flag for electrons with pT > 2.5 GeV/c
00176 
00177         if (mcpid[nmc]==-5) {mab += 1;} // Flag for bbar
00178         if (mcpid[nmc]==5) {mbb += 1;} // Flag for b
00179 
00180         if ((mcpid[nmc]==13)||(mcpid[nmc]==-13))
00181           {if (p.pt()>ptMuMax) {ptMuMax=p.pt();} } // Save max pt of generated Muons
00182         if ((mcpid[nmc]==11)||(mcpid[nmc]==-11))
00183           {if (p.pt() > ptEleMax) ptEleMax=p.pt();} // Save max pt of generated Electrons
00184 
00185         nmc++;
00186       }
00187 
00188     }
00189     //    else {std::cout << "%HLTMCtruth -- No MC truth information" << std::endl;}
00190 
00191     nmcpart = nmc;
00192     nmu3 = mu3;
00193     nel3 = el3;
00194     nbb = mbb;
00195     nab = mab;
00196     nwenu = wel;
00197     nwmunu = wmu;
00198     if((zee%3)==0){nzee = zee/3;}
00199 //     else {std::cout << "Z does not decay in e+ e- !" << std::endl;}
00200     if ((zmumu%3)==0){nzmumu = zmumu/3;}
00201 //     else {std::cout << "Z does not decay in mu+ mu- !" << std::endl;}
00202 
00203   }
00204 
00205 }