CMS 3D CMS Logo

/data/git/CMSSW_5_3_11_patch5/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 
00066 }
00067 
00068 /* **Analyze the event** */
00069 void HLTMCtruth::analyze(const edm::Handle<reco::CandidateView> & mctruth,
00070                          const double        & pthat,
00071                          const edm::Handle<std::vector<SimTrack> > & simTracks,
00072                          const edm::Handle<std::vector<SimVertex> > & simVertices,
00073                          TTree* HltTree) {
00074 
00075   //std::cout << " Beginning HLTMCtruth " << std::endl;
00076 
00077   if (_Monte) {
00078     int nmc = 0;
00079     int mu3 = 0;
00080     int el3 = 0;
00081     int mab = 0;
00082     int mbb = 0;
00083     int wel = 0;
00084     int wmu = 0;
00085     int zee = 0;
00086     int zmumu = 0;
00087 
00088     ptEleMax = -999.0;
00089     ptMuMax  = -999.0;    
00090     pthatf   = pthat;
00091 
00092     if((simTracks.isValid())&&(simVertices.isValid())){
00093       for (unsigned int j=0; j<simTracks->size(); j++) {
00094         int pdgid = simTracks->at(j).type();
00095         if (abs(pdgid)!=13) continue;
00096         double pt = simTracks->at(j).momentum().pt();
00097         if (pt<2.5) continue;
00098         double eta = simTracks->at(j).momentum().eta();
00099         if (abs(eta)>2.5) continue;
00100         if (simTracks->at(j).noVertex()) continue;
00101         int vertIndex = simTracks->at(j).vertIndex();
00102         double x = simVertices->at(vertIndex).position().x();
00103         double y = simVertices->at(vertIndex).position().y();
00104         double r = sqrt(x*x+y*y);
00105         if (r>150.) continue; // I think units are cm here
00106         double z = simVertices->at(vertIndex).position().z();
00107         if (abs(z)>300.) continue; // I think units are cm here
00108         mu3 += 1;
00109         break;
00110       }
00111     }
00112 
00113     if (mctruth.isValid()){
00114 
00115       for (size_t i = 0; i < mctruth->size(); ++ i) {
00116         const reco::Candidate & p = (*mctruth)[i];
00117 
00118         mcpid[nmc] = p.pdgId();
00119         mcstatus[nmc] = p.status();
00120         mcpt[nmc] = p.pt();
00121         mceta[nmc] = p.eta();
00122         mcphi[nmc] = p.phi();
00123         mcvx[nmc] = p.vx();
00124         mcvy[nmc] = p.vy();
00125         mcvz[nmc] = p.vz();
00126 
00127         if ((mcpid[nmc]==24)||(mcpid[nmc]==-24)) { // Checking W -> e/mu nu
00128           size_t idg = p.numberOfDaughters();
00129           for (size_t j=0; j != idg; ++j){
00130             const reco::Candidate & d = *p.daughter(j);
00131             if ((d.pdgId()==11)||(d.pdgId()==-11)){wel += 1;}
00132             if ((d.pdgId()==13)||(d.pdgId()==-13)){wmu += 1;}
00133 //          if ( (abs(d.pdgId())!=24) && ((mcpid[nmc])*(d.pdgId())>0) ) 
00134 //            {std::cout << "Wrong sign between mother-W and daughter !" << std::endl;}
00135           }
00136         }
00137         if (mcpid[nmc]==23) { // Checking Z -> 2 e/mu
00138           size_t idg = p.numberOfDaughters();
00139           for (size_t j=0; j != idg; ++j){
00140             const reco::Candidate & d = *p.daughter(j);
00141             if (d.pdgId()==11){zee += 1;}
00142             if (d.pdgId()==-11){zee += 2;}
00143             if (d.pdgId()==13){zmumu += 1;}
00144             if (d.pdgId()==-13){zmumu += 2;}
00145           }
00146         }
00147 
00148         // Set-up flags, based on Pythia-generator information, for avoiding double-counting events when
00149         // using both pp->{e,mu}X AND QCD samples
00150 //      if (((mcpid[nmc]==13)||(mcpid[nmc]==-13))&&(mcpt[nmc]>2.5)) {mu3 += 1;} // Flag for muons with pT > 2.5 GeV/c
00151         if (((mcpid[nmc]==11)||(mcpid[nmc]==-11))&&(mcpt[nmc]>2.5)) {el3 += 1;} // Flag for electrons with pT > 2.5 GeV/c
00152 
00153         if (mcpid[nmc]==-5) {mab += 1;} // Flag for bbar
00154         if (mcpid[nmc]==5) {mbb += 1;} // Flag for b
00155 
00156         if ((mcpid[nmc]==13)||(mcpid[nmc]==-13))
00157           {if (p.pt()>ptMuMax) {ptMuMax=p.pt();} } // Save max pt of generated Muons
00158         if ((mcpid[nmc]==11)||(mcpid[nmc]==-11))
00159           {if (p.pt() > ptEleMax) ptEleMax=p.pt();} // Save max pt of generated Electrons
00160 
00161         nmc++;
00162       }
00163 
00164     }
00165     //    else {std::cout << "%HLTMCtruth -- No MC truth information" << std::endl;}
00166 
00167     nmcpart = nmc;
00168     nmu3 = mu3;
00169     nel3 = el3;
00170     nbb = mbb;
00171     nab = mab;
00172     nwenu = wel;
00173     nwmunu = wmu;
00174     if((zee%3)==0){nzee = zee/3;}
00175 //     else {std::cout << "Z does not decay in e+ e- !" << std::endl;}
00176     if ((zmumu%3)==0){nzmumu = zmumu/3;}
00177 //     else {std::cout << "Z does not decay in mu+ mu- !" << std::endl;}
00178 
00179   }
00180 
00181 }