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
00018 _Monte=false;
00019 _Debug=false;
00020 }
00021
00022
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
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
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
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 if((simTracks.isValid())&&(simVertices.isValid())){
00100 for (unsigned int j=0; j<simTracks->size(); j++) {
00101 int pdgid = simTracks->at(j).type();
00102 if (abs(pdgid)!=13) continue;
00103 double pt = simTracks->at(j).momentum().pt();
00104 if (pt<5.) continue;
00105 double eta = simTracks->at(j).momentum().eta();
00106 if (abs(eta)>2.5) continue;
00107 if (simTracks->at(j).noVertex()) continue;
00108 int vertIndex = simTracks->at(j).vertIndex();
00109 double x = simVertices->at(vertIndex).position().x();
00110 double y = simVertices->at(vertIndex).position().y();
00111 double r = sqrt(x*x+y*y);
00112 if (r>150.) continue;
00113 double z = simVertices->at(vertIndex).position().z();
00114 if (abs(z)>300.) continue;
00115 mu3 += 1;
00116 break;
00117 }
00118
00119
00120 std::vector<PileupSummaryInfo>::const_iterator PVI;
00121 for(PVI = PupInfo->begin(); PVI != PupInfo->end(); ++PVI) {
00122
00123 int BX = PVI->getBunchCrossing();
00124 npvtrue = PVI->getTrueNumInteractions();
00125 npuvert = PVI->getPU_NumInteractions();
00126 if(BX == 0)
00127 {
00128 npubx0+=npvtrue;
00129 npuvertbx0+=npuvert;
00130 }
00131 }
00132
00133 }
00134
00135 if (mctruth.isValid()){
00136
00137 for (size_t i = 0; i < mctruth->size(); ++ i) {
00138 const reco::Candidate & p = (*mctruth)[i];
00139
00140 mcpid[nmc] = p.pdgId();
00141 mcstatus[nmc] = p.status();
00142 mcpt[nmc] = p.pt();
00143 mceta[nmc] = p.eta();
00144 mcphi[nmc] = p.phi();
00145 mcvx[nmc] = p.vx();
00146 mcvy[nmc] = p.vy();
00147 mcvz[nmc] = p.vz();
00148
00149 if ((mcpid[nmc]==24)||(mcpid[nmc]==-24)) {
00150 size_t idg = p.numberOfDaughters();
00151 for (size_t j=0; j != idg; ++j){
00152 const reco::Candidate & d = *p.daughter(j);
00153 if ((d.pdgId()==11)||(d.pdgId()==-11)){wel += 1;}
00154 if ((d.pdgId()==13)||(d.pdgId()==-13)){wmu += 1;}
00155
00156
00157 }
00158 }
00159 if (mcpid[nmc]==23) {
00160 size_t idg = p.numberOfDaughters();
00161 for (size_t j=0; j != idg; ++j){
00162 const reco::Candidate & d = *p.daughter(j);
00163 if (d.pdgId()==11){zee += 1;}
00164 if (d.pdgId()==-11){zee += 2;}
00165 if (d.pdgId()==13){zmumu += 1;}
00166 if (d.pdgId()==-13){zmumu += 2;}
00167 }
00168 }
00169
00170
00171
00172 if (((mcpid[nmc]==13)||(mcpid[nmc]==-13))&&(mcpt[nmc]>5.0)) {mu3 += 1;}
00173 if (((mcpid[nmc]==11)||(mcpid[nmc]==-11))&&(mcpt[nmc]>5.0)) {el3 += 1;}
00174
00175 if (mcpid[nmc]==-5) {mab += 1;}
00176 if (mcpid[nmc]==5) {mbb += 1;}
00177
00178 if ((mcpid[nmc]==13)||(mcpid[nmc]==-13))
00179 {if (p.pt()>ptMuMax) {ptMuMax=p.pt();} }
00180 if ((mcpid[nmc]==11)||(mcpid[nmc]==-11))
00181 {if (p.pt() > ptEleMax) ptEleMax=p.pt();}
00182
00183 nmc++;
00184 }
00185
00186 }
00187
00188
00189 nmcpart = nmc;
00190 nmu3 = mu3;
00191 nel3 = el3;
00192 nbb = mbb;
00193 nab = mab;
00194 nwenu = wel;
00195 nwmunu = wmu;
00196 if((zee%3)==0){nzee = zee/3;}
00197
00198 if ((zmumu%3)==0){nzmumu = zmumu/3;}
00199
00200
00201 }
00202
00203 }