CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
WValidation Class Reference

#include <WValidation.h>

Inheritance diagram for WValidation:
edm::EDAnalyzer

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
 
virtual void beginJob ()
 
virtual void beginRun (const edm::Run &, const edm::EventSetup &)
 
virtual void endJob ()
 
virtual void endRun (const edm::Run &, const edm::EventSetup &)
 
 WValidation (const edm::ParameterSet &)
 
virtual ~WValidation ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Private Attributes

int _flavor
 decay flavor More...
 
std::string _name
 decay flavor name More...
 
WeightManager _wmanager
 
MonitorElementcos_theta_gamma_lepton
 
DQMStoredbe
 ME's "container". More...
 
edm::ESHandle
< HepPDT::ParticleDataTable
fPDGTable
 PDT table. More...
 
MonitorElementgamma_energy
 
edm::InputTag hepmcCollection_
 
MonitorElementlepeta
 
MonitorElementlepmet_mT
 
MonitorElementlepmet_mTPeak
 
MonitorElementlepmet_pt
 
MonitorElementlepmet_ptLog
 
MonitorElementlepmet_rap
 
MonitorElementleppt
 
MonitorElementmet
 
MonitorElementnEvt
 
MonitorElementWdaughters
 
MonitorElementWmass
 
MonitorElementWmassPeak
 
MonitorElementWpt
 
MonitorElementWptLog
 
MonitorElementWrap
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

Definition at line 36 of file WValidation.h.

Constructor & Destructor Documentation

WValidation::WValidation ( const edm::ParameterSet iPSet)
explicit

Definition at line 20 of file WValidation.cc.

References dbe, and cppFunctionSkipper::operator.

20  :
21  _wmanager(iPSet),
22  hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")),
23  _flavor(iPSet.getParameter<int>("decaysTo")),
24  _name(iPSet.getParameter<std::string>("name"))
25 {
26  dbe = 0;
28 }
edm::InputTag hepmcCollection_
Definition: WValidation.h:51
T getParameter(std::string const &) const
std::string _name
decay flavor name
Definition: WValidation.h:68
DQMStore * dbe
ME&#39;s &quot;container&quot;.
Definition: WValidation.h:57
WeightManager _wmanager
Definition: WValidation.h:49
int _flavor
decay flavor
Definition: WValidation.h:66
WValidation::~WValidation ( )
virtual

Definition at line 30 of file WValidation.cc.

30 {}

Member Function Documentation

void WValidation::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
virtual

Gathering the HepMCProduct information

Implements edm::EDAnalyzer.

Definition at line 75 of file WValidation.cc.

References _flavor, _wmanager, abs, funct::cos(), cos_theta_gamma_lepton, MonitorElement::Fill(), HepMCValidationHelper::findFSRPhotons(), gamma_energy, HepMCValidationHelper::genMet(), edm::Event::getByLabel(), hepmcCollection_, lepeta, lepmet_mT, lepmet_mTPeak, lepmet_pt, lepmet_ptLog, leppt, met, nEvt, edm::es::products(), python.multivaluedict::sort(), HepMCValidationHelper::sortByPt(), lumiQTWidget::t, Wdaughters, WeightManager::weight(), CommonMethods::weight(), Wmass, WmassPeak, Wpt, WptLog, Wrap, vdt::x, detailsBasic3DVector::y, and detailsBasic3DVector::z.

76 {
77 
78  // we *DO NOT* rely on a Z entry in the particle listings!
79 
82  iEvent.getByLabel(hepmcCollection_, evt);
83 
84  //Get EVENT
85  const HepMC::GenEvent *myGenEvent = evt->GetEvent();
86 
87  double weight = _wmanager.weight(iEvent);
88 
89  nEvt->Fill(0.5,weight);
90 
91  std::vector<const HepMC::GenParticle*> allleptons;
92  std::vector<const HepMC::GenParticle*> allneutrinos;
93 
94  //requires status 1 for leptons and neutrinos (except tau)
95  int requiredstatus = (abs(_flavor) == 11 || abs(_flavor) ==13 ) ? 1 : 3;
96 
97  bool vetotau = true; //(abs(_flavor) == 11 || abs(_flavor) == 12 || abs(_flavor) ==13 || abs(_flavor) ==14 || abs(_flavor) ==16) ? true : false;
98 
99  for(HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin(); iter != myGenEvent->particles_end(); ++iter) {
100  if (vetotau) {
101  if ((*iter)->status()==3 && abs((*iter)->pdg_id() == 15) ) return;
102  }
103  if((*iter)->status()==requiredstatus) {
104  //@todo: improve this selection
105  if((*iter)->pdg_id()==_flavor)
106  allleptons.push_back(*iter);
107  else if (abs((*iter)->pdg_id()) == abs(_flavor)+1)
108  allneutrinos.push_back(*iter);
109  }
110  }
111 
112  //nothing to do if we don't have 2 particles
113  if (allleptons.empty() || allneutrinos.empty()) return;
114 
115  //sort them in pt
116  std::sort(allleptons.begin(), allleptons.end(), HepMCValidationHelper::sortByPt);
117  std::sort(allneutrinos.begin(), allneutrinos.end(), HepMCValidationHelper::sortByPt);
118 
119  //get the first lepton and the first neutrino, and check that one is particle one is antiparticle (product of pdgids < 0)
120  std::vector<const HepMC::GenParticle*> products;
121  if (allleptons.front()->pdg_id() * allneutrinos.front()->pdg_id() > 0) return;
122 
123  //require at least 20 GeV on the lepton
124  if (allleptons.front()->momentum().perp() < 20. || allneutrinos.front()->momentum().perp() < 20. ) return;
125 
126  //find possible qed fsr photons
127  std::vector<const HepMC::GenParticle*> selectedLepton;
128  selectedLepton.push_back(allleptons.front());
129  std::vector<const HepMC::GenParticle*> fsrphotons;
130  HepMCValidationHelper::findFSRPhotons(selectedLepton, myGenEvent, 0.1, fsrphotons);
131 
132  Wdaughters->Fill(allleptons.front()->pdg_id(),weight);
133  Wdaughters->Fill(allneutrinos.front()->pdg_id(),weight);
134 
135  //assemble FourMomenta
136  TLorentzVector lep1(allleptons[0]->momentum().x(), allleptons[0]->momentum().y(), allleptons[0]->momentum().z(), allleptons[0]->momentum().t());
137  TLorentzVector lep2(allneutrinos[0]->momentum().x(), allneutrinos[0]->momentum().y(), allneutrinos[0]->momentum().z(), allneutrinos[0]->momentum().t());
138  TLorentzVector dilepton_mom = lep1 + lep2;
139  TLorentzVector dilepton_andphoton_mom = dilepton_mom;
140  std::vector<TLorentzVector> gammasMomenta;
141  for (unsigned int ipho = 0; ipho < fsrphotons.size(); ++ipho){
142  TLorentzVector phomom(fsrphotons[ipho]->momentum().x(), fsrphotons[ipho]->momentum().y(), fsrphotons[ipho]->momentum().z(), fsrphotons[ipho]->momentum().t());
143  dilepton_andphoton_mom += phomom;
144  Wdaughters->Fill(fsrphotons[ipho]->pdg_id(),weight);
145  gammasMomenta.push_back(phomom);
146  }
147  //Fill "true" W histograms
148  Wmass->Fill(dilepton_andphoton_mom.M(),weight);
149  WmassPeak->Fill(dilepton_andphoton_mom.M(),weight);
150  Wpt->Fill(dilepton_andphoton_mom.Pt(),weight);
151  WptLog->Fill(log10(dilepton_andphoton_mom.Pt()),weight);
152  Wrap->Fill(dilepton_andphoton_mom.Rapidity(),weight);
153 
154  TLorentzVector met_mom = HepMCValidationHelper::genMet(myGenEvent, -3., 3.);
155  TLorentzVector lep1T(lep1.Px(), lep1.Py(), 0., lep1.Et());
156  TLorentzVector lepmet_mom = lep1T + met_mom;
157  //Fill lepmet histograms
158  lepmet_mT->Fill(lepmet_mom.M(),weight);
159  lepmet_mTPeak->Fill(lepmet_mom.M(),weight);
160  lepmet_pt->Fill(lepmet_mom.Pt(),weight);
161  lepmet_ptLog->Fill(log10(lepmet_mom.Pt()),weight);
162 
163  //Fill lepton histograms
164  leppt->Fill(lep1.Pt(),weight);
165  lepeta->Fill(lep1.Eta(),weight);
166  met->Fill(met_mom.Pt(),weight);
167 
168  //boost everything in the W frame
169  TVector3 boost = dilepton_andphoton_mom.BoostVector();
170  boost*=-1.;
171  lep1.Boost(boost);
172  lep2.Boost(boost);
173  for (unsigned int ipho = 0; ipho < gammasMomenta.size(); ++ipho){
174  gammasMomenta[ipho].Boost(boost);
175  }
176  std::sort(gammasMomenta.begin(), gammasMomenta.end(), HepMCValidationHelper::GreaterByE<TLorentzVector>);
177 
178  //fill gamma histograms
179  if (gammasMomenta.size() != 0 && dilepton_andphoton_mom.M() > 50.) {
180  gamma_energy->Fill(gammasMomenta.front().E(),weight);
181  double dphi = lep1.DeltaR(gammasMomenta.front());
182  cos_theta_gamma_lepton->Fill(cos(dphi),weight);
183  }
184 
185 
186 }//analyze
edm::InputTag hepmcCollection_
Definition: WValidation.h:51
MonitorElement * leppt
Definition: WValidation.h:62
void findFSRPhotons(const std::vector< const HepMC::GenParticle * > &leptons, const std::vector< const HepMC::GenParticle * > &all, double deltaR, std::vector< const HepMC::GenParticle * > &photons)
MonitorElement * Wrap
Definition: WValidation.h:60
MonitorElement * nEvt
Definition: WValidation.h:59
bool sortByPt(const HepMC::GenParticle *a, const HepMC::GenParticle *b)
#define abs(x)
Definition: mlp_lapack.h:159
MonitorElement * lepeta
Definition: WValidation.h:62
double double double z
MonitorElement * lepmet_mT
Definition: WValidation.h:61
void Fill(long long x)
MonitorElement * lepmet_pt
Definition: WValidation.h:61
ESProducts< T, S > products(const T &i1, const S &i2)
Definition: ESProducts.h:189
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
MonitorElement * gamma_energy
Definition: WValidation.h:63
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
MonitorElement * WptLog
Definition: WValidation.h:60
TLorentzVector genMet(const HepMC::GenEvent *all, double etamin=-9999., double etamax=9999.)
MonitorElement * lepmet_mTPeak
Definition: WValidation.h:61
MonitorElement * WmassPeak
Definition: WValidation.h:60
WeightManager _wmanager
Definition: WValidation.h:49
MonitorElement * met
Definition: WValidation.h:62
MonitorElement * Wpt
Definition: WValidation.h:60
MonitorElement * lepmet_ptLog
Definition: WValidation.h:61
x
Definition: VDTMath.h:216
MonitorElement * Wdaughters
Definition: WValidation.h:60
double weight(const edm::Event &)
MonitorElement * Wmass
Definition: WValidation.h:60
MonitorElement * cos_theta_gamma_lepton
Definition: WValidation.h:63
int _flavor
decay flavor
Definition: WValidation.h:66
void WValidation::beginJob ( void  )
virtual

Setting the DQM top directories

Reimplemented from edm::EDAnalyzer.

Definition at line 32 of file WValidation.cc.

References _name, DQMStore::book1D(), cos_theta_gamma_lepton, dbe, gamma_energy, lepeta, lepmet_mT, lepmet_mTPeak, lepmet_pt, lepmet_ptLog, leppt, met, nEvt, DQMStore::setCurrentFolder(), Wdaughters, Wmass, WmassPeak, Wpt, WptLog, and Wrap.

33 {
34  if(dbe){
36  std::string folderName = "Generator/W";
37  folderName+=_name;
38  dbe->setCurrentFolder(folderName.c_str());
39 
40  // Number of analyzed events
41  nEvt = dbe->book1D("nEvt", "n analyzed Events", 1, 0., 1.);
42 
43  //Kinematics
44  Wmass = dbe->book1D("Wmass","inv. Mass W", 70 ,0,140);
45  WmassPeak = dbe->book1D("WmassPeak","inv. Mass W", 80 ,80 ,100);
46  Wpt = dbe->book1D("Wpt","W pt",100,0,200);
47  WptLog = dbe->book1D("WptLog","log(W pt)",100,0.,5.);
48  Wrap = dbe->book1D("Wrap", "W y", 100, -5, 5);
49  Wdaughters = dbe->book1D("Wdaughters", "W daughters", 60, -30, 30);
50 
51  lepmet_mT = dbe->book1D("lepmet_mT","lepton-met transverse mass", 70 ,0,140);
52  lepmet_mTPeak = dbe->book1D("lepmet_mTPeak","lepton-met transverse mass", 80 ,80 ,100);
53  lepmet_pt = dbe->book1D("lepmet_pt","lepton-met",100,0,200);
54  lepmet_ptLog = dbe->book1D("lepmet_ptLog","log(lepton-met pt)",100,0.,5.);
55 
56  gamma_energy = dbe->book1D("gamma_energy", "photon energy in W rest frame", 200, 0., 100.);
57  cos_theta_gamma_lepton = dbe->book1D("cos_theta_gamma_lepton", "cos_theta_gamma_lepton in W rest frame", 200, -1, 1);
58 
59  leppt = dbe->book1D("leadpt","lepton pt", 200, 0., 200.);
60  met = dbe->book1D("met","met", 200, 0., 200.);
61  lepeta = dbe->book1D("leadeta","leading lepton eta", 100, -5., 5.);
62 
63  }
64  return;
65 }
MonitorElement * leppt
Definition: WValidation.h:62
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:717
MonitorElement * Wrap
Definition: WValidation.h:60
MonitorElement * nEvt
Definition: WValidation.h:59
std::string _name
decay flavor name
Definition: WValidation.h:68
DQMStore * dbe
ME&#39;s &quot;container&quot;.
Definition: WValidation.h:57
MonitorElement * lepeta
Definition: WValidation.h:62
MonitorElement * lepmet_mT
Definition: WValidation.h:61
MonitorElement * lepmet_pt
Definition: WValidation.h:61
MonitorElement * gamma_energy
Definition: WValidation.h:63
MonitorElement * WptLog
Definition: WValidation.h:60
MonitorElement * lepmet_mTPeak
Definition: WValidation.h:61
MonitorElement * WmassPeak
Definition: WValidation.h:60
MonitorElement * met
Definition: WValidation.h:62
MonitorElement * Wpt
Definition: WValidation.h:60
MonitorElement * lepmet_ptLog
Definition: WValidation.h:61
MonitorElement * Wdaughters
Definition: WValidation.h:60
MonitorElement * Wmass
Definition: WValidation.h:60
MonitorElement * cos_theta_gamma_lepton
Definition: WValidation.h:63
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:429
void WValidation::beginRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
)
virtual

Get PDT Table

Reimplemented from edm::EDAnalyzer.

Definition at line 68 of file WValidation.cc.

References fPDGTable, and edm::EventSetup::getData().

69 {
71  iSetup.getData( fPDGTable );
72  return;
73 }
void getData(T &iHolder) const
Definition: EventSetup.h:67
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
PDT table.
Definition: WValidation.h:54
void WValidation::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 67 of file WValidation.cc.

67 {return;}
void WValidation::endRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
)
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 74 of file WValidation.cc.

74 {return;}

Member Data Documentation

int WValidation::_flavor
private

decay flavor

Definition at line 66 of file WValidation.h.

Referenced by analyze().

std::string WValidation::_name
private

decay flavor name

Definition at line 68 of file WValidation.h.

Referenced by beginJob().

WeightManager WValidation::_wmanager
private

Definition at line 49 of file WValidation.h.

Referenced by analyze().

MonitorElement * WValidation::cos_theta_gamma_lepton
private

Definition at line 63 of file WValidation.h.

Referenced by analyze(), and beginJob().

DQMStore* WValidation::dbe
private

ME's "container".

Definition at line 57 of file WValidation.h.

Referenced by beginJob(), and WValidation().

edm::ESHandle<HepPDT::ParticleDataTable> WValidation::fPDGTable
private

PDT table.

Definition at line 54 of file WValidation.h.

Referenced by beginRun().

MonitorElement* WValidation::gamma_energy
private

Definition at line 63 of file WValidation.h.

Referenced by analyze(), and beginJob().

edm::InputTag WValidation::hepmcCollection_
private

Definition at line 51 of file WValidation.h.

Referenced by analyze().

MonitorElement * WValidation::lepeta
private

Definition at line 62 of file WValidation.h.

Referenced by analyze(), and beginJob().

MonitorElement* WValidation::lepmet_mT
private

Definition at line 61 of file WValidation.h.

Referenced by analyze(), and beginJob().

MonitorElement * WValidation::lepmet_mTPeak
private

Definition at line 61 of file WValidation.h.

Referenced by analyze(), and beginJob().

MonitorElement * WValidation::lepmet_pt
private

Definition at line 61 of file WValidation.h.

Referenced by analyze(), and beginJob().

MonitorElement * WValidation::lepmet_ptLog
private

Definition at line 61 of file WValidation.h.

Referenced by analyze(), and beginJob().

MonitorElement * WValidation::lepmet_rap
private

Definition at line 61 of file WValidation.h.

MonitorElement* WValidation::leppt
private

Definition at line 62 of file WValidation.h.

Referenced by analyze(), and beginJob().

MonitorElement * WValidation::met
private

Definition at line 62 of file WValidation.h.

Referenced by analyze(), and beginJob().

MonitorElement* WValidation::nEvt
private

Definition at line 59 of file WValidation.h.

Referenced by analyze(), and beginJob().

MonitorElement * WValidation::Wdaughters
private

Definition at line 60 of file WValidation.h.

Referenced by analyze(), and beginJob().

MonitorElement* WValidation::Wmass
private

Definition at line 60 of file WValidation.h.

Referenced by analyze(), and beginJob().

MonitorElement * WValidation::WmassPeak
private

Definition at line 60 of file WValidation.h.

Referenced by analyze(), and beginJob().

MonitorElement * WValidation::Wpt
private

Definition at line 60 of file WValidation.h.

Referenced by analyze(), and beginJob().

MonitorElement * WValidation::WptLog
private

Definition at line 60 of file WValidation.h.

Referenced by analyze(), and beginJob().

MonitorElement * WValidation::Wrap
private

Definition at line 60 of file WValidation.h.

Referenced by analyze(), and beginJob().