CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
BasicHepMCHeavyIonValidation Class Reference

#include <BasicHepMCHeavyIonValidation.h>

Inheritance diagram for BasicHepMCHeavyIonValidation:
DQMEDAnalyzer edm::one::EDProducer< edm::Accumulator, edm::EndLuminosityBlockProducer, edm::EndRunProducer, edm::one::WatchLuminosityBlocks, edm::one::WatchRuns > edm::one::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

void analyze (edm::Event const &, edm::EventSetup const &) override
 
 BasicHepMCHeavyIonValidation (const edm::ParameterSet &)
 
void bookHistograms (DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
 
 ~BasicHepMCHeavyIonValidation () override
 
- Public Member Functions inherited from DQMEDAnalyzer
void accumulate (edm::Event const &ev, edm::EventSetup const &es) final
 
void beginLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) override
 
void beginRun (edm::Run const &run, edm::EventSetup const &setup) final
 
virtual void dqmBeginRun (edm::Run const &, edm::EventSetup const &)
 
 DQMEDAnalyzer ()
 
 DQMEDAnalyzer (DQMEDAnalyzer const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer &&)=delete
 
void endLuminosityBlock (edm::LuminosityBlock const &, edm::EventSetup const &) override
 
void endLuminosityBlockProduce (edm::LuminosityBlock &lumi, edm::EventSetup const &setup) final
 
void endRun (edm::Run const &run, edm::EventSetup const &setup) override
 
void endRunProduce (edm::Run &run, edm::EventSetup const &setup) override
 
 ~DQMEDAnalyzer () override=default
 
- Public Member Functions inherited from edm::one::EDProducer< edm::Accumulator, edm::EndLuminosityBlockProducer, edm::EndRunProducer, edm::one::WatchLuminosityBlocks, edm::one::WatchRuns >
 EDProducer ()=default
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
bool hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () const final
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
- Public Member Functions inherited from edm::one::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDProducerBase () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
 ~ProducerBase () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Attributes

MonitorElementeccentricity
 
MonitorElementevent_plane_angle
 
edm::InputTag hepmcCollection_
 
edm::EDGetTokenT< edm::HepMCProducthepmcCollectionToken_
 
MonitorElementimpact_parameter
 
MonitorElementN_Nwounded_collisions
 
MonitorElementNcoll
 
MonitorElementNcoll_hard
 
MonitorElementnEvt
 PDT table. More...
 
MonitorElementNpart_proj
 
MonitorElementNpart_targ
 
MonitorElementNwounded_N_collisions
 
MonitorElementNwounded_Nwounded_collisions
 
bool QWdebug_
 
MonitorElementsigma_inel_NN
 
MonitorElementspectator_neutrons
 
MonitorElementspectator_protons
 
WeightManager wmanager_
 

Additional Inherited Members

- Public Types inherited from edm::one::EDProducerBase
typedef EDProducerBase ModuleType
 
- Public Types inherited from edm::ProducerBase
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex >>
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDProducerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
- Protected Attributes inherited from DQMEDAnalyzer
edm::EDPutTokenT< DQMTokenlumiToken_
 
edm::EDPutTokenT< DQMTokenrunToken_
 

Detailed Description

Definition at line 34 of file BasicHepMCHeavyIonValidation.h.

Constructor & Destructor Documentation

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

Definition at line 15 of file BasicHepMCHeavyIonValidation.cc.

References edm::ParameterSet::getUntrackedParameter(), hepmcCollection_, hepmcCollectionToken_, and QWdebug_.

15  :
17  hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection"))
18 {
19  QWdebug_ = iPSet.getUntrackedParameter<bool>("QWdebug",false);
20 
21  hepmcCollectionToken_=consumes<HepMCProduct>(hepmcCollection_);
22 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
BasicHepMCHeavyIonValidation::~BasicHepMCHeavyIonValidation ( )
override

Definition at line 24 of file BasicHepMCHeavyIonValidation.cc.

24 {}

Member Function Documentation

void BasicHepMCHeavyIonValidation::analyze ( edm::Event const &  ,
edm::EventSetup const &   
)
overridevirtual

counters to zero for every event

Gathering the HepMCProduct information

Reimplemented from DQMEDAnalyzer.

Definition at line 52 of file BasicHepMCHeavyIonValidation.cc.

References gather_cfg::cout, eccentricity, event_plane_angle, MonitorElement::Fill(), edm::Event::getByToken(), edm::HepMCProduct::GetEvent(), hepmcCollectionToken_, impact_parameter, N_Nwounded_collisions, Ncoll, Ncoll_hard, nEvt, Npart_proj, Npart_targ, Nwounded_N_collisions, Nwounded_Nwounded_collisions, QWdebug_, sigma_inel_NN, spectator_neutrons, spectator_protons, WeightManager::weight(), mps_merge::weight, and wmanager_.

52  {
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
Definition: weight.py:1
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
double weight(const edm::Event &)
void BasicHepMCHeavyIonValidation::bookHistograms ( DQMStore::IBooker i,
edm::Run const &  ,
edm::EventSetup const &   
)
overridevirtual

Setting the DQM top directories

Booking the ME's

Implements DQMEDAnalyzer.

Definition at line 26 of file BasicHepMCHeavyIonValidation.cc.

References DQMHelper::book1dHisto(), eccentricity, event_plane_angle, impact_parameter, N_Nwounded_collisions, Ncoll, Ncoll_hard, nEvt, Npart_proj, Npart_targ, Nwounded_N_collisions, Nwounded_Nwounded_collisions, pi, DQMStore::IBooker::setCurrentFolder(), sigma_inel_NN, spectator_neutrons, and spectator_protons.

26  {
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 }
const Double_t pi
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:279

Member Data Documentation

MonitorElement* BasicHepMCHeavyIonValidation::eccentricity
private

Definition at line 64 of file BasicHepMCHeavyIonValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCHeavyIonValidation::event_plane_angle
private

Definition at line 63 of file BasicHepMCHeavyIonValidation.h.

Referenced by analyze(), and bookHistograms().

edm::InputTag BasicHepMCHeavyIonValidation::hepmcCollection_
private

Definition at line 44 of file BasicHepMCHeavyIonValidation.h.

Referenced by BasicHepMCHeavyIonValidation().

edm::EDGetTokenT<edm::HepMCProduct> BasicHepMCHeavyIonValidation::hepmcCollectionToken_
private

Definition at line 70 of file BasicHepMCHeavyIonValidation.h.

Referenced by analyze(), and BasicHepMCHeavyIonValidation().

MonitorElement* BasicHepMCHeavyIonValidation::impact_parameter
private

Definition at line 62 of file BasicHepMCHeavyIonValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCHeavyIonValidation::N_Nwounded_collisions
private

Definition at line 57 of file BasicHepMCHeavyIonValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCHeavyIonValidation::Ncoll
private

Definition at line 56 of file BasicHepMCHeavyIonValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCHeavyIonValidation::Ncoll_hard
private

Definition at line 53 of file BasicHepMCHeavyIonValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCHeavyIonValidation::nEvt
private

PDT table.

Definition at line 50 of file BasicHepMCHeavyIonValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCHeavyIonValidation::Npart_proj
private

Definition at line 54 of file BasicHepMCHeavyIonValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCHeavyIonValidation::Npart_targ
private

Definition at line 55 of file BasicHepMCHeavyIonValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCHeavyIonValidation::Nwounded_N_collisions
private

Definition at line 58 of file BasicHepMCHeavyIonValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCHeavyIonValidation::Nwounded_Nwounded_collisions
private

Definition at line 59 of file BasicHepMCHeavyIonValidation.h.

Referenced by analyze(), and bookHistograms().

bool BasicHepMCHeavyIonValidation::QWdebug_
private

Definition at line 45 of file BasicHepMCHeavyIonValidation.h.

Referenced by analyze(), and BasicHepMCHeavyIonValidation().

MonitorElement* BasicHepMCHeavyIonValidation::sigma_inel_NN
private

Definition at line 67 of file BasicHepMCHeavyIonValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCHeavyIonValidation::spectator_neutrons
private

Definition at line 60 of file BasicHepMCHeavyIonValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCHeavyIonValidation::spectator_protons
private

Definition at line 61 of file BasicHepMCHeavyIonValidation.h.

Referenced by analyze(), and bookHistograms().

WeightManager BasicHepMCHeavyIonValidation::wmanager_
private

Definition at line 43 of file BasicHepMCHeavyIonValidation.h.

Referenced by analyze().