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()),
17  hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection"))
18 {
19  QWdebug_ = iPSet.getUntrackedParameter<bool>("QWdebug",false);
20 
21  hepmcCollectionToken_=consumes<HepMCProduct>(hepmcCollection_);
22 }
23 
25 
27 
29  DQMHelper dqm(&i); i.setCurrentFolder("Generator/HeavyIon");
30 
31  // Number of analyzed events
32  nEvt = dqm.book1dHisto("nEvt", "n analyzed Events", 1, 0., 1.,"","Number of Events");
33 
35  Ncoll_hard = dqm.book1dHisto("Ncoll_hard", "Ncoll_hard", 700, 0, 700,"Number of hard scatterings","Number of Events");
36  Npart_proj = dqm.book1dHisto("Npart_proj", "Npart_proj", 250, 0, 250,"Number of projectile participants","Number of Events");
37  Npart_targ = dqm.book1dHisto("Npart_targ", "Npart_targ", 250, 0, 250,"Number of target participants","Number of Events");
38  Ncoll = dqm.book1dHisto("Ncoll", "Ncoll", 700, 0, 700,"Number of N-N collisions","Number of Events");
39  N_Nwounded_collisions = dqm.book1dHisto("N_Nwounded_collisions", "N_Nwounded_collisions", 250, 0, 250,"Number of N-N wounded collisions","Number of Events");
40  Nwounded_N_collisions = dqm.book1dHisto("Nwounded_N_collisions", "Nwounded_N_collisions", 250, 0, 250,"Number of N wounded-N collisions","Number of Events");
41  Nwounded_Nwounded_collisions = dqm.book1dHisto("Nwounded_Nwounded_collisions", "Nwounded_Nwounded_collisions", 250, 0, 250,"Number of N wounded-N wounded collisions","Number of Events");
42  spectator_neutrons = dqm.book1dHisto("spectator_neutrons", "spectator_neutrons", 250, 0, 250,"Number of spectator neutrons","Number of Events");
43  spectator_protons = dqm.book1dHisto("spectator_protons", "spectator_protons", 250, 0, 250,"Number of spectator protons","Number of Events");
44  impact_parameter = dqm.book1dHisto("impact_parameter", "impact_parameter", 50, 0, 50,"Empact parameter of collision (fm)","Number of Events");
45  event_plane_angle = dqm.book1dHisto("event_plane_angle", "event_plane_angle", 200, 0, 2*CLHEP::pi,"#phi_{event plane} (rad)","Number of Events");
46  eccentricity = dqm.book1dHisto("eccentricity", "eccentricity", 200, 0, 1.0,"Eccentricity","Number of Events");
47  sigma_inel_NN = dqm.book1dHisto("sigma_inel_NN", "sigma_inel_NN", 200, 0, 10.0,"#sigma{nucleon-nucleon inelastic cross-section}","Number of Events");
48 
49  return;
50 }
51 
54 
57  iEvent.getByToken(hepmcCollectionToken_, evt);
58 
59  //Get EVENT
60  //HepMC::GenEvent *myGenEvent = new HepMC::GenEvent(*(evt->GetEvent()));
61 
62 
63 
64  const HepMC::HeavyIon* ion = evt->GetEvent()->heavy_ion();
65 
66  if (!ion) {
67  if ( QWdebug_ ) std::cout << "!!QW!! HeavyIon == null" << std::endl;
68  return;
69  }
70 
71  double weight = wmanager_.weight(iEvent);
72  nEvt->Fill(0.5,weight);
73 
74  Ncoll_hard->Fill(ion->Ncoll_hard(), weight);
75  Npart_proj->Fill(ion->Npart_proj(), weight);
76  Npart_targ->Fill(ion->Npart_targ(), weight);
77  Ncoll->Fill(ion->Ncoll(), weight);
78  N_Nwounded_collisions->Fill(ion->N_Nwounded_collisions(), weight);
79  Nwounded_N_collisions->Fill(ion->Nwounded_N_collisions(), weight);
80  Nwounded_Nwounded_collisions->Fill(ion->Nwounded_Nwounded_collisions(), weight);
81  spectator_neutrons->Fill(ion->spectator_neutrons(), weight);
82  spectator_protons->Fill(ion->spectator_protons(), weight);
83  impact_parameter->Fill(ion->impact_parameter(), weight);
84  event_plane_angle->Fill(ion->event_plane_angle(), weight);
85  eccentricity->Fill(ion->eccentricity(), weight);
86  sigma_inel_NN->Fill(ion->sigma_inel_NN(), weight);
87 
88 
89  //delete myGenEvent;
90 }//analyze
T getUntrackedParameter(std::string const &, T const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
void analyze(edm::Event const &, edm::EventSetup const &) override
Definition: weight.py:1
MonitorElement * book1dHisto(std::string name, std::string title, int n, double xmin, double xmax, std::string xaxis, std::string yaxis)
Definition: DQMHelper.cc:12
void Fill(long long x)
const Double_t pi
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
BasicHepMCHeavyIonValidation(const edm::ParameterSet &)
HLT enums.
double weight(const edm::Event &)
Definition: Run.h:44