CMS 3D CMS Logo

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"
12 
13 using namespace edm;
14 
16  : wmanager_(iPSet, consumesCollector()), hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")) {
17  QWdebug_ = iPSet.getUntrackedParameter<bool>("QWdebug", false);
18 
19  hepmcCollectionToken_ = consumes<HepMCProduct>(hepmcCollection_);
20 }
21 
23 
26  DQMHelper dqm(&i);
27  i.setCurrentFolder("Generator/HeavyIon");
28 
29  // Number of analyzed events
30  nEvt = dqm.book1dHisto("nEvt", "n analyzed Events", 1, 0., 1., "", "Number of Events");
31 
33  Ncoll_hard =
34  dqm.book1dHisto("Ncoll_hard", "Ncoll_hard", 700, 0, 700, "Number of hard scatterings", "Number of Events");
35  Npart_proj =
36  dqm.book1dHisto("Npart_proj", "Npart_proj", 250, 0, 250, "Number of projectile participants", "Number of Events");
37  Npart_targ =
38  dqm.book1dHisto("Npart_targ", "Npart_targ", 250, 0, 250, "Number of target participants", "Number of Events");
39  Ncoll = dqm.book1dHisto("Ncoll", "Ncoll", 700, 0, 700, "Number of N-N collisions", "Number of Events");
40  N_Nwounded_collisions = dqm.book1dHisto("N_Nwounded_collisions",
41  "N_Nwounded_collisions",
42  250,
43  0,
44  250,
45  "Number of N-N wounded collisions",
46  "Number of Events");
47  Nwounded_N_collisions = dqm.book1dHisto("Nwounded_N_collisions",
48  "Nwounded_N_collisions",
49  250,
50  0,
51  250,
52  "Number of N wounded-N collisions",
53  "Number of Events");
54  Nwounded_Nwounded_collisions = dqm.book1dHisto("Nwounded_Nwounded_collisions",
55  "Nwounded_Nwounded_collisions",
56  250,
57  0,
58  250,
59  "Number of N wounded-N wounded collisions",
60  "Number of Events");
61  spectator_neutrons = dqm.book1dHisto(
62  "spectator_neutrons", "spectator_neutrons", 250, 0, 250, "Number of spectator neutrons", "Number of Events");
63  spectator_protons = dqm.book1dHisto(
64  "spectator_protons", "spectator_protons", 250, 0, 250, "Number of spectator protons", "Number of Events");
65  impact_parameter = dqm.book1dHisto(
66  "impact_parameter", "impact_parameter", 50, 0, 50, "Empact parameter of collision (fm)", "Number of Events");
67  event_plane_angle = dqm.book1dHisto(
68  "event_plane_angle", "event_plane_angle", 200, 0, 2 * CLHEP::pi, "#phi_{event plane} (rad)", "Number of Events");
69  eccentricity = dqm.book1dHisto("eccentricity", "eccentricity", 200, 0, 1.0, "Eccentricity", "Number of Events");
70  sigma_inel_NN = dqm.book1dHisto("sigma_inel_NN",
71  "sigma_inel_NN",
72  200,
73  0,
74  10.0,
75  "#sigma{nucleon-nucleon inelastic cross-section}",
76  "Number of Events");
77 
78  return;
79 }
80 
83 
86  iEvent.getByToken(hepmcCollectionToken_, evt);
87 
88  //Get EVENT
89  //HepMC::GenEvent *myGenEvent = new HepMC::GenEvent(*(evt->GetEvent()));
90 
91  const HepMC::HeavyIon* ion = evt->GetEvent()->heavy_ion();
92 
93  if (!ion) {
94  if (QWdebug_)
95  std::cout << "!!QW!! HeavyIon == null" << std::endl;
96  return;
97  }
98 
99  double weight = wmanager_.weight(iEvent);
100  nEvt->Fill(0.5, weight);
101 
102  Ncoll_hard->Fill(ion->Ncoll_hard(), weight);
103  Npart_proj->Fill(ion->Npart_proj(), weight);
104  Npart_targ->Fill(ion->Npart_targ(), weight);
105  Ncoll->Fill(ion->Ncoll(), weight);
106  N_Nwounded_collisions->Fill(ion->N_Nwounded_collisions(), weight);
107  Nwounded_N_collisions->Fill(ion->Nwounded_N_collisions(), weight);
108  Nwounded_Nwounded_collisions->Fill(ion->Nwounded_Nwounded_collisions(), weight);
109  spectator_neutrons->Fill(ion->spectator_neutrons(), weight);
110  spectator_protons->Fill(ion->spectator_protons(), weight);
111  impact_parameter->Fill(ion->impact_parameter(), weight);
112  event_plane_angle->Fill(ion->event_plane_angle(), weight);
113  eccentricity->Fill(ion->eccentricity(), weight);
114  sigma_inel_NN->Fill(ion->sigma_inel_NN(), weight);
115 
116  //delete myGenEvent;
117 } //analyze
mps_fire.i
i
Definition: mps_fire.py:428
edm::Run
Definition: Run.h:45
BasicHepMCHeavyIonValidation::analyze
void analyze(edm::Event const &, edm::EventSetup const &) override
Definition: BasicHepMCHeavyIonValidation.cc:81
edm
HLT enums.
Definition: AlignableModifier.h:19
mps_merge.weight
weight
Definition: mps_merge.py:88
gather_cfg.cout
cout
Definition: gather_cfg.py:144
BasicHepMCHeavyIonValidation::Ncoll
MonitorElement * Ncoll
Definition: BasicHepMCHeavyIonValidation.h:54
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
BasicHepMCHeavyIonValidation::impact_parameter
MonitorElement * impact_parameter
Definition: BasicHepMCHeavyIonValidation.h:60
edm::Handle
Definition: AssociativeIterator.h:50
BasicHepMCHeavyIonValidation::N_Nwounded_collisions
MonitorElement * N_Nwounded_collisions
Definition: BasicHepMCHeavyIonValidation.h:55
BasicHepMCHeavyIonValidation::Nwounded_N_collisions
MonitorElement * Nwounded_N_collisions
Definition: BasicHepMCHeavyIonValidation.h:56
BasicHepMCHeavyIonValidation::spectator_protons
MonitorElement * spectator_protons
Definition: BasicHepMCHeavyIonValidation.h:59
DQMHelper.h
BasicHepMCHeavyIonValidation::QWdebug_
bool QWdebug_
Definition: BasicHepMCHeavyIonValidation.h:43
BasicHepMCHeavyIonValidation::hepmcCollectionToken_
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
Definition: BasicHepMCHeavyIonValidation.h:68
BasicHepMCHeavyIonValidation::nEvt
MonitorElement * nEvt
PDT table.
Definition: BasicHepMCHeavyIonValidation.h:48
BasicHepMCHeavyIonValidation::spectator_neutrons
MonitorElement * spectator_neutrons
Definition: BasicHepMCHeavyIonValidation.h:58
dqm::impl::MonitorElement::Fill
void Fill(long long x)
Definition: MonitorElement.h:290
BasicHepMCHeavyIonValidation::Ncoll_hard
MonitorElement * Ncoll_hard
Definition: BasicHepMCHeavyIonValidation.h:51
BasicHepMCHeavyIonValidation::hepmcCollection_
edm::InputTag hepmcCollection_
Definition: BasicHepMCHeavyIonValidation.h:42
BasicHepMCHeavyIonValidation.h
edm::ParameterSet
Definition: ParameterSet.h:47
iEvent
int iEvent
Definition: GenABIO.cc:224
BasicHepMCHeavyIonValidation::wmanager_
WeightManager wmanager_
Definition: BasicHepMCHeavyIonValidation.h:41
BasicHepMCHeavyIonValidation::BasicHepMCHeavyIonValidation
BasicHepMCHeavyIonValidation(const edm::ParameterSet &)
Definition: BasicHepMCHeavyIonValidation.cc:15
edm::EventSetup
Definition: EventSetup.h:58
edm::HepMCProduct::GetEvent
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
BasicHepMCHeavyIonValidation::Npart_proj
MonitorElement * Npart_proj
Definition: BasicHepMCHeavyIonValidation.h:52
WeightManager::weight
double weight(const edm::Event &)
Definition: WeightManager.cc:24
BasicHepMCHeavyIonValidation::bookHistograms
void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
Definition: BasicHepMCHeavyIonValidation.cc:24
DQMHelper
Definition: DQMHelper.h:15
BasicHepMCHeavyIonValidation::Nwounded_Nwounded_collisions
MonitorElement * Nwounded_Nwounded_collisions
Definition: BasicHepMCHeavyIonValidation.h:57
BasicHepMCHeavyIonValidation::sigma_inel_NN
MonitorElement * sigma_inel_NN
Definition: BasicHepMCHeavyIonValidation.h:65
BasicHepMCHeavyIonValidation::eccentricity
MonitorElement * eccentricity
Definition: BasicHepMCHeavyIonValidation.h:62
dqm::implementation::IBooker
Definition: DQMStore.h:43
BasicHepMCHeavyIonValidation::~BasicHepMCHeavyIonValidation
~BasicHepMCHeavyIonValidation() override
Definition: BasicHepMCHeavyIonValidation.cc:22
pi
const Double_t pi
Definition: trackSplitPlot.h:36
dqm
Definition: DQMStore.h:18
BasicHepMCHeavyIonValidation::Npart_targ
MonitorElement * Npart_targ
Definition: BasicHepMCHeavyIonValidation.h:53
BasicHepMCHeavyIonValidation::event_plane_angle
MonitorElement * event_plane_angle
Definition: BasicHepMCHeavyIonValidation.h:61
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
weight
Definition: weight.py:1