CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BasicHepMCHeavyIonValidation.cc
Go to the documentation of this file.
1 /*class BasicHepMCHeavyIonValidation
2  *
3  * Class to fill dqm monitor elements from existing EDM file
4  * Quan Wang - 04/2013
5  */
6 
8 
9 #include "CLHEP/Units/defs.h"
10 #include "CLHEP/Units/PhysicalConstants.h"
11 
12 using namespace edm;
13 
15  wmanager_(iPSet,consumesCollector()),
16  hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection"))
17 {
18  dbe = 0;
20  QWdebug_ = iPSet.getUntrackedParameter<bool>("QWdebug",false);
21 
22  hepmcCollectionToken_=consumes<HepMCProduct>(hepmcCollection_);
23 }
24 
26 
28 {
29  if(dbe){
31  dbe->setCurrentFolder("Generator/HeavyIon");
32 
33  // Number of analyzed events
34  nEvt = dbe->book1D("nEvt", "n analyzed Events", 1, 0., 1.);
35 
37  Ncoll_hard = dbe->book1D("Ncoll_hard", "Ncoll_hard", 700, 0, 700);
38  Npart_proj = dbe->book1D("Npart_proj", "Npart_proj", 250, 0, 250);
39  Npart_targ = dbe->book1D("Npart_targ", "Npart_targ", 250, 0, 250);
40  Ncoll = dbe->book1D("Ncoll", "Ncoll", 700, 0, 700);
41  N_Nwounded_collisions = dbe->book1D("N_Nwounded_collisions", "N_Nwounded_collisions", 250, 0, 250);
42  Nwounded_N_collisions = dbe->book1D("Nwounded_N_collisions", "Nwounded_N_collisions", 250, 0, 250);
43  Nwounded_Nwounded_collisions = dbe->book1D("Nwounded_Nwounded_collisions", "Nwounded_Nwounded_collisions", 250, 0, 250);
44  spectator_neutrons = dbe->book1D("spectator_neutrons", "spectator_neutrons", 250, 0, 250);
45  spectator_protons = dbe->book1D("spectator_protons", "spectator_protons", 250, 0, 250);
46  impact_parameter = dbe->book1D("impact_parameter", "impact_parameter", 50, 0, 50);
47  event_plane_angle = dbe->book1D("event_plane_angle", "event_plane_angle", 200, -CLHEP::pi, CLHEP::pi);
48  eccentricity = dbe->book1D("eccentricity", "eccentricity", 200, 0, 1.0);
49  sigma_inel_NN = dbe->book1D("sigma_inel_NN", "sigma_inel_NN", 200, 0, 10.0);
50 
51  }
52  return;
53 }
54 
57 {
59  //iSetup.getData( fPDGTable );
60  return;
61 }
62 void BasicHepMCHeavyIonValidation::endRun(const edm::Run& iRun,const edm::EventSetup& iSetup){return;}
64 {
66 
69  iEvent.getByToken(hepmcCollectionToken_, evt);
70 
71  //Get EVENT
72  //HepMC::GenEvent *myGenEvent = new HepMC::GenEvent(*(evt->GetEvent()));
73 
74 
75 
76  const HepMC::HeavyIon* ion = evt->GetEvent()->heavy_ion();
77 
78  if (!ion) {
79  if ( QWdebug_ ) std::cout << "!!QW!! HeavyIon == null" << std::endl;
80  return;
81  }
82 
83  double weight = wmanager_.weight(iEvent);
84  nEvt->Fill(0.5,weight);
85 
86  Ncoll_hard->Fill(ion->Ncoll_hard(), weight);
87  Npart_proj->Fill(ion->Npart_proj(), weight);
88  Npart_targ->Fill(ion->Npart_targ(), weight);
89  Ncoll->Fill(ion->Ncoll(), weight);
90  N_Nwounded_collisions->Fill(ion->N_Nwounded_collisions(), weight);
91  Nwounded_N_collisions->Fill(ion->Nwounded_N_collisions(), weight);
92  Nwounded_Nwounded_collisions->Fill(ion->Nwounded_Nwounded_collisions(), weight);
93  spectator_neutrons->Fill(ion->spectator_neutrons(), weight);
94  spectator_protons->Fill(ion->spectator_protons(), weight);
95  impact_parameter->Fill(ion->impact_parameter(), weight);
96  event_plane_angle->Fill(ion->event_plane_angle(), weight);
97  eccentricity->Fill(ion->eccentricity(), weight);
98  sigma_inel_NN->Fill(ion->sigma_inel_NN(), weight);
99 
100 
101  //delete myGenEvent;
102 }//analyze
virtual void endRun(const edm::Run &, const edm::EventSetup &)
T getUntrackedParameter(std::string const &, T const &) const
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:872
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:243
virtual void analyze(const edm::Event &, const edm::EventSetup &)
virtual void beginRun(const edm::Run &, const edm::EventSetup &)
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
BasicHepMCHeavyIonValidation(const edm::ParameterSet &)
tuple cout
Definition: gather_cfg.py:121
double pi
int weight
Definition: histoStyle.py:50
double weight(const edm::Event &)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:584
Definition: Run.h:41