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  QWdebug_ = iPSet.getUntrackedParameter<bool>("QWdebug",false);
19 
20  hepmcCollectionToken_=consumes<HepMCProduct>(hepmcCollection_);
21 }
22 
24 
26 
28  i.setCurrentFolder("Generator/HeavyIon");
29 
30  // Number of analyzed events
31  nEvt = i.book1D("nEvt", "n analyzed Events", 1, 0., 1.);
32 
34  Ncoll_hard = i.book1D("Ncoll_hard", "Ncoll_hard", 700, 0, 700);
35  Npart_proj = i.book1D("Npart_proj", "Npart_proj", 250, 0, 250);
36  Npart_targ = i.book1D("Npart_targ", "Npart_targ", 250, 0, 250);
37  Ncoll = i.book1D("Ncoll", "Ncoll", 700, 0, 700);
38  N_Nwounded_collisions = i.book1D("N_Nwounded_collisions", "N_Nwounded_collisions", 250, 0, 250);
39  Nwounded_N_collisions = i.book1D("Nwounded_N_collisions", "Nwounded_N_collisions", 250, 0, 250);
40  Nwounded_Nwounded_collisions = i.book1D("Nwounded_Nwounded_collisions", "Nwounded_Nwounded_collisions", 250, 0, 250);
41  spectator_neutrons = i.book1D("spectator_neutrons", "spectator_neutrons", 250, 0, 250);
42  spectator_protons = i.book1D("spectator_protons", "spectator_protons", 250, 0, 250);
43  impact_parameter = i.book1D("impact_parameter", "impact_parameter", 50, 0, 50);
44  event_plane_angle = i.book1D("event_plane_angle", "event_plane_angle", 200, -CLHEP::pi, CLHEP::pi);
45  eccentricity = i.book1D("eccentricity", "eccentricity", 200, 0, 1.0);
46  sigma_inel_NN = i.book1D("sigma_inel_NN", "sigma_inel_NN", 200, 0, 10.0);
47 
48  return;
49 }
50 
53 
56  iEvent.getByToken(hepmcCollectionToken_, evt);
57 
58  //Get EVENT
59  //HepMC::GenEvent *myGenEvent = new HepMC::GenEvent(*(evt->GetEvent()));
60 
61 
62 
63  const HepMC::HeavyIon* ion = evt->GetEvent()->heavy_ion();
64 
65  if (!ion) {
66  if ( QWdebug_ ) std::cout << "!!QW!! HeavyIon == null" << std::endl;
67  return;
68  }
69 
70  double weight = wmanager_.weight(iEvent);
71  nEvt->Fill(0.5,weight);
72 
73  Ncoll_hard->Fill(ion->Ncoll_hard(), weight);
74  Npart_proj->Fill(ion->Npart_proj(), weight);
75  Npart_targ->Fill(ion->Npart_targ(), weight);
76  Ncoll->Fill(ion->Ncoll(), weight);
77  N_Nwounded_collisions->Fill(ion->N_Nwounded_collisions(), weight);
78  Nwounded_N_collisions->Fill(ion->Nwounded_N_collisions(), weight);
79  Nwounded_Nwounded_collisions->Fill(ion->Nwounded_Nwounded_collisions(), weight);
80  spectator_neutrons->Fill(ion->spectator_neutrons(), weight);
81  spectator_protons->Fill(ion->spectator_protons(), weight);
82  impact_parameter->Fill(ion->impact_parameter(), weight);
83  event_plane_angle->Fill(ion->event_plane_angle(), weight);
84  eccentricity->Fill(ion->eccentricity(), weight);
85  sigma_inel_NN->Fill(ion->sigma_inel_NN(), weight);
86 
87 
88  //delete myGenEvent;
89 }//analyze
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
virtual void analyze(edm::Event const &, edm::EventSetup const &) override
void Fill(long long x)
const Double_t pi
int iEvent
Definition: GenABIO.cc:230
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:113
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
virtual void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
BasicHepMCHeavyIonValidation(const edm::ParameterSet &)
tuple cout
Definition: gather_cfg.py:121
int weight
Definition: histoStyle.py:50
double weight(const edm::Event &)
Definition: Run.h:41