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