CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTMCtruth.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <sstream>
3 #include <istream>
4 #include <fstream>
5 #include <iomanip>
6 #include <string>
7 #include <cmath>
8 #include <functional>
9 #include <stdlib.h>
10 #include <string.h>
11 
14 
16 
17  //set parameter defaults
18  _Monte=false;
19  _Debug=false;
20 }
21 
22 /* Setup the analysis to put the branch-variables into the tree. */
23 void HLTMCtruth::setup(const edm::ParameterSet& pSet, TTree* HltTree) {
24 
25  edm::ParameterSet myMCParams = pSet.getParameter<edm::ParameterSet>("RunParameters") ;
26  std::vector<std::string> parameterNames = myMCParams.getParameterNames() ;
27 
28  for ( std::vector<std::string>::iterator iParam = parameterNames.begin();
29  iParam != parameterNames.end(); iParam++ ){
30  if ( (*iParam) == "Monte" ) _Monte = myMCParams.getParameter<bool>( *iParam );
31  else if ( (*iParam) == "Debug" ) _Debug = myMCParams.getParameter<bool>( *iParam );
32  }
33 
34  const int kMaxMcTruth = 10000;
35  mcpid = new int[kMaxMcTruth];
36  mcstatus = new int[kMaxMcTruth];
37  mcvx = new float[kMaxMcTruth];
38  mcvy = new float[kMaxMcTruth];
39  mcvz = new float[kMaxMcTruth];
40  mcpt = new float[kMaxMcTruth];
41  mceta = new float[kMaxMcTruth];
42  mcphi = new float[kMaxMcTruth];
43 
44  // MCtruth-specific branches of the tree
45  HltTree->Branch("NMCpart",&nmcpart,"NMCpart/I");
46  HltTree->Branch("MCpid",mcpid,"MCpid[NMCpart]/I");
47  HltTree->Branch("MCstatus",mcstatus,"MCstatus[NMCpart]/I");
48  HltTree->Branch("MCvtxX",mcvx,"MCvtxX[NMCpart]/F");
49  HltTree->Branch("MCvtxY",mcvy,"MCvtxY[NMCpart]/F");
50  HltTree->Branch("MCvtxZ",mcvz,"MCvtxZ[NMCpart]/F");
51  HltTree->Branch("MCpt",mcpt,"MCpt[NMCpart]/F");
52  HltTree->Branch("MCeta",mceta,"MCeta[NMCpart]/F");
53  HltTree->Branch("MCphi",mcphi,"MCphi[NMCpart]/F");
54  HltTree->Branch("MCPtHat",&pthatf,"MCPtHat/F");
55  HltTree->Branch("MCmu3",&nmu3,"MCmu3/I");
56  HltTree->Branch("MCel3",&nel3,"MCel3/I");
57  HltTree->Branch("MCbb",&nbb,"MCbb/I");
58  HltTree->Branch("MCab",&nab,"MCab/I");
59  HltTree->Branch("MCWenu",&nwenu,"MCWenu/I");
60  HltTree->Branch("MCWmunu",&nwmunu,"MCmunu/I");
61  HltTree->Branch("MCZee",&nzee,"MCZee/I");
62  HltTree->Branch("MCZmumu",&nzmumu,"MCZmumu/I");
63  HltTree->Branch("MCptEleMax",&ptEleMax,"MCptEleMax/F");
64  HltTree->Branch("MCptMuMax",&ptMuMax,"MCptMuMax/F");
65  HltTree->Branch("NPUTrueBX0",&npubx0, "NPUTrueBX0/I");
66  HltTree->Branch("NPUgenBX0",&npuvertbx0, "NPUgenBX0/I");
67 }
68 
69 /* **Analyze the event** */
71  const double & pthat,
72  const edm::Handle<std::vector<SimTrack> > & simTracks,
73  const edm::Handle<std::vector<SimVertex> > & simVertices,
74  const edm::Handle<std::vector< PileupSummaryInfo > > & PupInfo,
75  TTree* HltTree) {
76 
77  //std::cout << " Beginning HLTMCtruth " << std::endl;
78 
79  if (_Monte) {
80  int nmc = 0;
81  int mu3 = 0;
82  int el3 = 0;
83  int mab = 0;
84  int mbb = 0;
85  int wel = 0;
86  int wmu = 0;
87  int zee = 0;
88  int zmumu = 0;
89 
90  ptEleMax = -999.0;
91  ptMuMax = -999.0;
92  pthatf = pthat;
93  npubx0 = 0.0;
94  npuvertbx0 = 0;
95 
96  int npvtrue = 0;
97  int npuvert = 0;
98 
99 
100  if((simTracks.isValid())&&(simVertices.isValid())){
101  for (unsigned int j=0; j<simTracks->size(); j++) {
102  int pdgid = simTracks->at(j).type();
103  if (abs(pdgid)!=13) continue;
104  double pt = simTracks->at(j).momentum().pt();
105  if (pt<5.0) continue;
106  double eta = simTracks->at(j).momentum().eta();
107  if (abs(eta)>2.5) continue;
108  if (simTracks->at(j).noVertex()) continue;
109  int vertIndex = simTracks->at(j).vertIndex();
110  double x = simVertices->at(vertIndex).position().x();
111  double y = simVertices->at(vertIndex).position().y();
112  double r = sqrt(x*x+y*y);
113  if (r>150.) continue; // I think units are cm here
114  double z = simVertices->at(vertIndex).position().z();
115  if (abs(z)>300.) continue; // I think units are cm here
116  mu3 += 1;
117  break;
118  }
119 
120 
121  std::vector<PileupSummaryInfo>::const_iterator PVI;
122  for(PVI = PupInfo->begin(); PVI != PupInfo->end(); ++PVI) {
123 
124  int BX = PVI->getBunchCrossing();
125  npvtrue = PVI->getTrueNumInteractions();
126  npuvert = PVI->getPU_NumInteractions();
127 
128  if(BX == 0)
129  {
130  npubx0+=npvtrue;
131  npuvertbx0+=npuvert;
132  }
133  }
134 
135  }
136 
137  if (mctruth.isValid()){
138 
139  for (size_t i = 0; i < mctruth->size(); ++ i) {
140  const reco::Candidate & p = (*mctruth)[i];
141 
142  mcpid[nmc] = p.pdgId();
143  mcstatus[nmc] = p.status();
144  mcpt[nmc] = p.pt();
145  mceta[nmc] = p.eta();
146  mcphi[nmc] = p.phi();
147  mcvx[nmc] = p.vx();
148  mcvy[nmc] = p.vy();
149  mcvz[nmc] = p.vz();
150 
151  if ((mcpid[nmc]==24)||(mcpid[nmc]==-24)) { // Checking W -> e/mu nu
152  size_t idg = p.numberOfDaughters();
153  for (size_t j=0; j != idg; ++j){
154  const reco::Candidate & d = *p.daughter(j);
155  if ((d.pdgId()==11)||(d.pdgId()==-11)){wel += 1;}
156  if ((d.pdgId()==13)||(d.pdgId()==-13)){wmu += 1;}
157 // if ( (abs(d.pdgId())!=24) && ((mcpid[nmc])*(d.pdgId())>0) )
158 // {std::cout << "Wrong sign between mother-W and daughter !" << std::endl;}
159  }
160  }
161  if (mcpid[nmc]==23) { // Checking Z -> 2 e/mu
162  size_t idg = p.numberOfDaughters();
163  for (size_t j=0; j != idg; ++j){
164  const reco::Candidate & d = *p.daughter(j);
165  if (d.pdgId()==11){zee += 1;}
166  if (d.pdgId()==-11){zee += 2;}
167  if (d.pdgId()==13){zmumu += 1;}
168  if (d.pdgId()==-13){zmumu += 2;}
169  }
170  }
171 
172  // Set-up flags, based on Pythia-generator information, for avoiding double-counting events when
173  // using both pp->{e,mu}X AND QCD samples
174  //if (((mcpid[nmc]==13)||(mcpid[nmc]==-13))&&(mcpt[nmc]>5.0)) {mu3 += 1;} // Flag for muons with pT > 2.5 GeV/c
175  if (((mcpid[nmc]==11)||(mcpid[nmc]==-11))&&(mcpt[nmc]>5.0)) {el3 += 1;} // Flag for electrons with pT > 2.5 GeV/c
176 
177  if (mcpid[nmc]==-5) {mab += 1;} // Flag for bbar
178  if (mcpid[nmc]==5) {mbb += 1;} // Flag for b
179 
180  if ((mcpid[nmc]==13)||(mcpid[nmc]==-13))
181  {if (p.pt()>ptMuMax) {ptMuMax=p.pt();} } // Save max pt of generated Muons
182  if ((mcpid[nmc]==11)||(mcpid[nmc]==-11))
183  {if (p.pt() > ptEleMax) ptEleMax=p.pt();} // Save max pt of generated Electrons
184 
185  nmc++;
186  }
187 
188  }
189  // else {std::cout << "%HLTMCtruth -- No MC truth information" << std::endl;}
190 
191  nmcpart = nmc;
192  nmu3 = mu3;
193  nel3 = el3;
194  nbb = mbb;
195  nab = mab;
196  nwenu = wel;
197  nwmunu = wmu;
198  if((zee%3)==0){nzee = zee/3;}
199 // else {std::cout << "Z does not decay in e+ e- !" << std::endl;}
200  if ((zmumu%3)==0){nzmumu = zmumu/3;}
201 // else {std::cout << "Z does not decay in mu+ mu- !" << std::endl;}
202 
203  }
204 
205 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
int * mcpid
Definition: HLTMCtruth.h:48
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
float * mcpt
Definition: HLTMCtruth.h:47
virtual double pt() const =0
transverse momentum
float pthatf
Definition: HLTMCtruth.h:51
virtual int status() const =0
status word
virtual double vx() const =0
x coordinate of vertex position
#define abs(x)
Definition: mlp_lapack.h:159
T eta() const
bool _Debug
Definition: HLTMCtruth.h:54
double double double z
virtual double vy() const =0
y coordinate of vertex position
float * mceta
Definition: HLTMCtruth.h:47
int nmcpart
Definition: HLTMCtruth.h:49
int nzmumu
Definition: HLTMCtruth.h:49
virtual size_type numberOfDaughters() const =0
number of daughters
float ptEleMax
Definition: HLTMCtruth.h:52
float ptMuMax
Definition: HLTMCtruth.h:52
T sqrt(T t)
Definition: SSEVec.h:46
float * mcvy
Definition: HLTMCtruth.h:47
int npubx0
Definition: HLTMCtruth.h:50
int j
Definition: DBlmapReader.cc:9
bool _Monte
Definition: HLTMCtruth.h:54
std::vector< std::string > getParameterNames() const
bool isValid() const
Definition: HandleBase.h:76
virtual int pdgId() const =0
PDG identifier.
void analyze(const edm::Handle< reco::CandidateView > &mctruth, const double &pthat, const edm::Handle< std::vector< SimTrack > > &simTracks, const edm::Handle< std::vector< SimVertex > > &simVertices, const edm::Handle< std::vector< PileupSummaryInfo > > &PupInfo, TTree *tree)
Definition: HLTMCtruth.cc:70
float * mcvx
Definition: HLTMCtruth.h:47
float * mcphi
Definition: HLTMCtruth.h:47
virtual double vz() const =0
z coordinate of vertex position
float * mcvz
Definition: HLTMCtruth.h:47
void setup(const edm::ParameterSet &pSet, TTree *tree)
Definition: HLTMCtruth.cc:23
int nwmunu
Definition: HLTMCtruth.h:49
x
Definition: VDTMath.h:216
int npuvertbx0
Definition: HLTMCtruth.h:50
int * mcstatus
Definition: HLTMCtruth.h:48
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity