CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
HLTMCtruth Class Reference

#include <HLTMCtruth.h>

Public Member Functions

void analyze (const edm::Handle< reco::CandidateView > &mctruth, const double &pthat, const double &weight, const edm::Handle< std::vector< SimTrack > > &simTracks, const edm::Handle< std::vector< SimVertex > > &simVertices, const edm::Handle< std::vector< PileupSummaryInfo > > &PupInfo, TTree *tree)
 
 HLTMCtruth ()
 
void setup (const edm::ParameterSet &pSet, TTree *tree)
 

Private Attributes

bool _Debug
 
bool _Gen
 
bool _Monte
 
float * mceta
 
float * mcphi
 
int * mcpid
 
float * mcpt
 
int * mcstatus
 
float * mcvx
 
float * mcvy
 
float * mcvz
 
int nab
 
int nbb
 
int nel3
 
int nmcpart
 
int nmu3
 
int npubx0
 
int npuvertbx0
 
int nwenu
 
int nwmunu
 
int nzee
 
int nzmumu
 
float ptEleMax
 
float pthatf
 
float ptMuMax
 
float weightf
 
float weightsignf
 

Detailed Description

$Date: November 2006 $Revision:

Author
P. Bargassa - Rice U.

Definition at line 30 of file HLTMCtruth.h.

Constructor & Destructor Documentation

HLTMCtruth::HLTMCtruth ( )

Definition at line 15 of file HLTMCtruth.cc.

References _Debug, _Gen, and _Monte.

15  {
16 
17  //set parameter defaults
18  _Monte=false;
19  _Gen=false;
20  _Debug=false;
21 }
bool _Gen
Definition: HLTMCtruth.h:55
bool _Debug
Definition: HLTMCtruth.h:55
bool _Monte
Definition: HLTMCtruth.h:55

Member Function Documentation

void HLTMCtruth::analyze ( const edm::Handle< reco::CandidateView > &  mctruth,
const double &  pthat,
const double &  weight,
const edm::Handle< std::vector< SimTrack > > &  simTracks,
const edm::Handle< std::vector< SimVertex > > &  simVertices,
const edm::Handle< std::vector< PileupSummaryInfo > > &  PupInfo,
TTree *  tree 
)

Analyze the Data

Definition at line 75 of file HLTMCtruth.cc.

References _Gen, _Monte, funct::abs(), rpcdqm::BX, edmIntegrityCheck::d, reco::Candidate::daughter(), stringResolutionProvider_cfi::eta, reco::Candidate::eta(), i, edm::HandleBase::isValid(), j, mceta, mcphi, mcpid, mcpt, mcstatus, mcvx, mcvy, mcvz, nab, nbb, nel3, nmcpart, nmu3, npubx0, npuvertbx0, reco::Candidate::numberOfDaughters(), nwenu, nwmunu, nzee, nzmumu, AlCaHLTBitMon_ParallelJobs::p, BPhysicsValidation_cfi::pdgid, reco::Candidate::pdgId(), reco::Candidate::phi(), EnergyCorrector::pt, reco::Candidate::pt(), ptEleMax, pthatf, ptMuMax, alignCSCRings::r, tkConvValidator_cfi::simTracks, edm::View< T >::size(), mathSSE::sqrt(), reco::Candidate::status(), reco::Candidate::vx(), reco::Candidate::vy(), reco::Candidate::vz(), mps_merge::weight, weightf, weightsignf, x, y, and z.

Referenced by HLTBitAnalyzer::analyze().

81  {
82 
83  //std::cout << " Beginning HLTMCtruth " << std::endl;
84 
85  if (_Monte) {
86  int nmc = 0;
87  int mu3 = 0;
88  int el3 = 0;
89  int mab = 0;
90  int mbb = 0;
91  int wel = 0;
92  int wmu = 0;
93  int zee = 0;
94  int zmumu = 0;
95 
96  ptEleMax = -999.0;
97  ptMuMax = -999.0;
98  pthatf = pthat;
99  weightf = weight;
100  weightsignf = (weight > 0) ? 1. : -1.;
101  npubx0 = 0.0;
102  npuvertbx0 = 0;
103 
104  int npvtrue = 0;
105  int npuvert = 0;
106 
107  if(PupInfo.isValid()) {
108  std::vector<PileupSummaryInfo>::const_iterator PVI;
109  for(PVI = PupInfo->begin(); PVI != PupInfo->end(); ++PVI) {
110 
111  int BX = PVI->getBunchCrossing();
112  npvtrue = PVI->getTrueNumInteractions();
113  npuvert = PVI->getPU_NumInteractions();
114 
115  if(BX == 0)
116  {
117  npubx0+=npvtrue;
118  npuvertbx0+=npuvert;
119  }
120  }
121  }
122 
123  if((simTracks.isValid())&&(simVertices.isValid())){
124  for (unsigned int j=0; j<simTracks->size(); j++) {
125  int pdgid = simTracks->at(j).type();
126  if (abs(pdgid)!=13) continue;
127  double pt = simTracks->at(j).momentum().pt();
128  if (pt<5.0) continue;
129  double eta = simTracks->at(j).momentum().eta();
130  if (abs(eta)>2.5) continue;
131  if (simTracks->at(j).noVertex()) continue;
132  int vertIndex = simTracks->at(j).vertIndex();
133  double x = simVertices->at(vertIndex).position().x();
134  double y = simVertices->at(vertIndex).position().y();
135  double r = sqrt(x*x+y*y);
136  if (r>200.) continue; // I think units are cm here
137  double z = simVertices->at(vertIndex).position().z();
138  if (abs(z)>400.) continue; // I think units are cm here
139  mu3 += 1;
140  break;
141  }
142 
143  }
144 
145  if (_Gen && mctruth.isValid()){
146 
147  for (size_t i = 0; i < mctruth->size(); ++ i) {
148  const reco::Candidate & p = (*mctruth)[i];
149 
150  mcpid[nmc] = p.pdgId();
151  mcstatus[nmc] = p.status();
152  mcpt[nmc] = p.pt();
153  mceta[nmc] = p.eta();
154  mcphi[nmc] = p.phi();
155  mcvx[nmc] = p.vx();
156  mcvy[nmc] = p.vy();
157  mcvz[nmc] = p.vz();
158 
159  if ((mcpid[nmc]==24)||(mcpid[nmc]==-24)) { // Checking W -> e/mu nu
160  size_t idg = p.numberOfDaughters();
161  for (size_t j=0; j != idg; ++j){
162  const reco::Candidate & d = *p.daughter(j);
163  if ((d.pdgId()==11)||(d.pdgId()==-11)){wel += 1;}
164  if ((d.pdgId()==13)||(d.pdgId()==-13)){wmu += 1;}
165 // if ( (abs(d.pdgId())!=24) && ((mcpid[nmc])*(d.pdgId())>0) )
166 // {std::cout << "Wrong sign between mother-W and daughter !" << std::endl;}
167  }
168  }
169  if (mcpid[nmc]==23) { // Checking Z -> 2 e/mu
170  size_t idg = p.numberOfDaughters();
171  for (size_t j=0; j != idg; ++j){
172  const reco::Candidate & d = *p.daughter(j);
173  if (d.pdgId()==11){zee += 1;}
174  if (d.pdgId()==-11){zee += 2;}
175  if (d.pdgId()==13){zmumu += 1;}
176  if (d.pdgId()==-13){zmumu += 2;}
177  }
178  }
179 
180  // Set-up flags, based on Pythia-generator information, for avoiding double-counting events when
181  // using both pp->{e,mu}X AND QCD samples
182 // if (((mcpid[nmc]==13)||(mcpid[nmc]==-13))&&(mcpt[nmc]>2.5)) {mu3 += 1;} // Flag for muons with pT > 2.5 GeV/c
183  if (((mcpid[nmc]==11)||(mcpid[nmc]==-11))&&(mcpt[nmc]>2.5)) {el3 += 1;} // Flag for electrons with pT > 2.5 GeV/c
184 
185  if (mcpid[nmc]==-5) {mab += 1;} // Flag for bbar
186  if (mcpid[nmc]==5) {mbb += 1;} // Flag for b
187 
188  if ((mcpid[nmc]==13)||(mcpid[nmc]==-13))
189  {if (p.pt()>ptMuMax) {ptMuMax=p.pt();} } // Save max pt of generated Muons
190  if ((mcpid[nmc]==11)||(mcpid[nmc]==-11))
191  {if (p.pt() > ptEleMax) ptEleMax=p.pt();} // Save max pt of generated Electrons
192 
193  nmc++;
194  }
195 
196  }
197  // else {std::cout << "%HLTMCtruth -- No MC truth information" << std::endl;}
198 
199  nmcpart = nmc;
200  nmu3 = mu3;
201  nel3 = el3;
202  nbb = mbb;
203  nab = mab;
204  nwenu = wel;
205  nwmunu = wmu;
206  if((zee%3)==0){nzee = zee/3;}
207 // else {std::cout << "Z does not decay in e+ e- !" << std::endl;}
208  if ((zmumu%3)==0){nzmumu = zmumu/3;}
209 // else {std::cout << "Z does not decay in mu+ mu- !" << std::endl;}
210 
211  }
212 
213 }
int i
Definition: DBlmapReader.cc:9
int * mcpid
Definition: HLTMCtruth.h:49
bool _Gen
Definition: HLTMCtruth.h:55
float * mcpt
Definition: HLTMCtruth.h:48
virtual double vx() const =0
x coordinate of vertex position
float pthatf
Definition: HLTMCtruth.h:52
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
size_type size() const
Definition: weight.py:1
float weightf
Definition: HLTMCtruth.h:52
virtual double vy() const =0
y coordinate of vertex position
float * mceta
Definition: HLTMCtruth.h:48
int nmcpart
Definition: HLTMCtruth.h:50
virtual int status() const =0
status word
int nzmumu
Definition: HLTMCtruth.h:50
float ptEleMax
Definition: HLTMCtruth.h:53
float ptMuMax
Definition: HLTMCtruth.h:53
virtual int pdgId() const =0
PDG identifier.
T sqrt(T t)
Definition: SSEVec.h:18
float * mcvy
Definition: HLTMCtruth.h:48
int npubx0
Definition: HLTMCtruth.h:51
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
bool _Monte
Definition: HLTMCtruth.h:55
float weightsignf
Definition: HLTMCtruth.h:52
bool isValid() const
Definition: HandleBase.h:75
virtual double eta() const =0
momentum pseudorapidity
virtual double pt() const =0
transverse momentum
float * mcvx
Definition: HLTMCtruth.h:48
float * mcphi
Definition: HLTMCtruth.h:48
float * mcvz
Definition: HLTMCtruth.h:48
int nwmunu
Definition: HLTMCtruth.h:50
virtual double vz() const =0
z coordinate of vertex position
virtual size_type numberOfDaughters() const =0
number of daughters
virtual double phi() const =0
momentum azimuthal angle
int npuvertbx0
Definition: HLTMCtruth.h:51
int * mcstatus
Definition: HLTMCtruth.h:49
void HLTMCtruth::setup ( const edm::ParameterSet pSet,
TTree *  tree 
)

Definition at line 24 of file HLTMCtruth.cc.

References _Debug, _Gen, _Monte, edm::ParameterSet::getParameter(), edm::ParameterSet::getParameterNames(), mceta, mcphi, mcpid, mcpt, mcstatus, mcvx, mcvy, mcvz, nab, nbb, nel3, nmcpart, nmu3, npubx0, npuvertbx0, nwenu, nwmunu, nzee, nzmumu, ptEleMax, pthatf, ptMuMax, weightf, and weightsignf.

Referenced by HLTBitAnalyzer::HLTBitAnalyzer().

24  {
25 
26  edm::ParameterSet myMCParams = pSet.getParameter<edm::ParameterSet>("RunParameters") ;
27  std::vector<std::string> parameterNames = myMCParams.getParameterNames() ;
28 
29  for ( std::vector<std::string>::iterator iParam = parameterNames.begin();
30  iParam != parameterNames.end(); iParam++ ){
31  if ( (*iParam) == "Monte" ) _Monte = myMCParams.getParameter<bool>( *iParam );
32  else if ( (*iParam) == "GenTracks" ) _Gen = myMCParams.getParameter<bool>( *iParam );
33  else if ( (*iParam) == "Debug" ) _Debug = myMCParams.getParameter<bool>( *iParam );
34  }
35 
36  const int kMaxMcTruth = 10000;
37  mcpid = new int[kMaxMcTruth];
38  mcstatus = new int[kMaxMcTruth];
39  mcvx = new float[kMaxMcTruth];
40  mcvy = new float[kMaxMcTruth];
41  mcvz = new float[kMaxMcTruth];
42  mcpt = new float[kMaxMcTruth];
43  mceta = new float[kMaxMcTruth];
44  mcphi = new float[kMaxMcTruth];
45 
46  // MCtruth-specific branches of the tree
47  HltTree->Branch("NMCpart",&nmcpart,"NMCpart/I");
48  HltTree->Branch("MCpid",mcpid,"MCpid[NMCpart]/I");
49  HltTree->Branch("MCstatus",mcstatus,"MCstatus[NMCpart]/I");
50  HltTree->Branch("MCvtxX",mcvx,"MCvtxX[NMCpart]/F");
51  HltTree->Branch("MCvtxY",mcvy,"MCvtxY[NMCpart]/F");
52  HltTree->Branch("MCvtxZ",mcvz,"MCvtxZ[NMCpart]/F");
53  HltTree->Branch("MCpt",mcpt,"MCpt[NMCpart]/F");
54  HltTree->Branch("MCeta",mceta,"MCeta[NMCpart]/F");
55  HltTree->Branch("MCphi",mcphi,"MCphi[NMCpart]/F");
56  HltTree->Branch("MCPtHat",&pthatf,"MCPtHat/F");
57  HltTree->Branch("MCWeight",&weightf,"MCWeight/F");
58  HltTree->Branch("MCWeightSign",&weightsignf,"MCWeightSign/F");
59  HltTree->Branch("MCmu3",&nmu3,"MCmu3/I");
60  HltTree->Branch("MCel3",&nel3,"MCel3/I");
61  HltTree->Branch("MCbb",&nbb,"MCbb/I");
62  HltTree->Branch("MCab",&nab,"MCab/I");
63  HltTree->Branch("MCWenu",&nwenu,"MCWenu/I");
64  HltTree->Branch("MCWmunu",&nwmunu,"MCmunu/I");
65  HltTree->Branch("MCZee",&nzee,"MCZee/I");
66  HltTree->Branch("MCZmumu",&nzmumu,"MCZmumu/I");
67  HltTree->Branch("MCptEleMax",&ptEleMax,"MCptEleMax/F");
68  HltTree->Branch("MCptMuMax",&ptMuMax,"MCptMuMax/F");
69  HltTree->Branch("NPUTrueBX0",&npubx0, "NPUTrueBX0/I");
70  HltTree->Branch("NPUgenBX0",&npuvertbx0, "NPUgenBX0/I");
71 
72 }
T getParameter(std::string const &) const
int * mcpid
Definition: HLTMCtruth.h:49
bool _Gen
Definition: HLTMCtruth.h:55
float * mcpt
Definition: HLTMCtruth.h:48
float pthatf
Definition: HLTMCtruth.h:52
float weightf
Definition: HLTMCtruth.h:52
bool _Debug
Definition: HLTMCtruth.h:55
float * mceta
Definition: HLTMCtruth.h:48
int nmcpart
Definition: HLTMCtruth.h:50
int nzmumu
Definition: HLTMCtruth.h:50
float ptEleMax
Definition: HLTMCtruth.h:53
float ptMuMax
Definition: HLTMCtruth.h:53
float * mcvy
Definition: HLTMCtruth.h:48
int npubx0
Definition: HLTMCtruth.h:51
bool _Monte
Definition: HLTMCtruth.h:55
float weightsignf
Definition: HLTMCtruth.h:52
std::vector< std::string > getParameterNames() const
float * mcvx
Definition: HLTMCtruth.h:48
float * mcphi
Definition: HLTMCtruth.h:48
float * mcvz
Definition: HLTMCtruth.h:48
int nwmunu
Definition: HLTMCtruth.h:50
int npuvertbx0
Definition: HLTMCtruth.h:51
int * mcstatus
Definition: HLTMCtruth.h:49

Member Data Documentation

bool HLTMCtruth::_Debug
private

Definition at line 55 of file HLTMCtruth.h.

Referenced by HLTMCtruth(), and setup().

bool HLTMCtruth::_Gen
private

Definition at line 55 of file HLTMCtruth.h.

Referenced by analyze(), HLTMCtruth(), and setup().

bool HLTMCtruth::_Monte
private

Definition at line 55 of file HLTMCtruth.h.

Referenced by analyze(), HLTMCtruth(), and setup().

float * HLTMCtruth::mceta
private

Definition at line 48 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

float * HLTMCtruth::mcphi
private

Definition at line 48 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

int* HLTMCtruth::mcpid
private

Definition at line 49 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

float * HLTMCtruth::mcpt
private

Definition at line 48 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

int * HLTMCtruth::mcstatus
private

Definition at line 49 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

float* HLTMCtruth::mcvx
private

Definition at line 48 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

float * HLTMCtruth::mcvy
private

Definition at line 48 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

float * HLTMCtruth::mcvz
private

Definition at line 48 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

int HLTMCtruth::nab
private

Definition at line 50 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

int HLTMCtruth::nbb
private

Definition at line 50 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

int HLTMCtruth::nel3
private

Definition at line 50 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

int HLTMCtruth::nmcpart
private

Definition at line 50 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

int HLTMCtruth::nmu3
private

Definition at line 50 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

int HLTMCtruth::npubx0
private

Definition at line 51 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

int HLTMCtruth::npuvertbx0
private

Definition at line 51 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

int HLTMCtruth::nwenu
private

Definition at line 50 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

int HLTMCtruth::nwmunu
private

Definition at line 50 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

int HLTMCtruth::nzee
private

Definition at line 50 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

int HLTMCtruth::nzmumu
private

Definition at line 50 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

float HLTMCtruth::ptEleMax
private

Definition at line 53 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

float HLTMCtruth::pthatf
private

Definition at line 52 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

float HLTMCtruth::ptMuMax
private

Definition at line 53 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

float HLTMCtruth::weightf
private

Definition at line 52 of file HLTMCtruth.h.

Referenced by analyze(), and setup().

float HLTMCtruth::weightsignf
private

Definition at line 52 of file HLTMCtruth.h.

Referenced by analyze(), and setup().