CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
HGCalTimingAnalyzer Class Reference
Inheritance diagram for HGCalTimingAnalyzer:
edm::one::EDAnalyzer< edm::one::WatchRuns, edm::one::SharedResources > edm::one::EDAnalyzerBase edm::EDConsumerBase

Public Member Functions

 HGCalTimingAnalyzer (edm::ParameterSet const &)
 
 ~HGCalTimingAnalyzer () override
 
- Public Member Functions inherited from edm::one::EDAnalyzer< edm::one::WatchRuns, edm::one::SharedResources >
 EDAnalyzer ()=default
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
- Public Member Functions inherited from edm::one::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDAnalyzerBase () 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
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
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)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
- Static Public Member Functions inherited from edm::one::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 

Private Member Functions

void analyze (edm::Event const &, edm::EventSetup const &) override
 
void analyzeSimHits (int type, std::vector< PCaloHit > const &hits)
 
void analyzeSimTracks (edm::Handle< edm::SimTrackContainer > const &SimTk, edm::Handle< edm::SimVertexContainer > const &SimVtx)
 
void beginJob () override
 
void beginRun (edm::Run const &, edm::EventSetup const &) override
 
void endRun (edm::Run const &, edm::EventSetup const &) override
 

Private Attributes

std::string detectorBeam_
 
std::string detectorEE_
 
bool doTree_
 
edm::Service< TFileServicefs_
 
bool groupHits_
 
const HGCalDDDConstantshgcons_
 
std::vector< int > idBeams_
 
double pBeam_
 
std::vector< float > simHitCellEnBeam_
 
std::vector< float > simHitCellEnEE_
 
std::vector< uint32_t > simHitCellIdBeam_
 
std::vector< uint32_t > simHitCellIdEE_
 
std::vector< float > simHitCellTmBeam_
 
std::vector< float > simHitCellTmEE_
 
double timeUnit_
 
edm::EDGetTokenT< edm::HepMCProducttok_hepMC_
 
edm::EDGetTokenT< edm::PCaloHitContainertok_hitsBeam_
 
edm::EDGetTokenT< edm::PCaloHitContainertok_hitsEE_
 
edm::EDGetTokenT< edm::SimTrackContainertok_simTk_
 
edm::EDGetTokenT< edm::SimVertexContainertok_simVtx_
 
TTree * tree_
 
double xBeam_
 
double yBeam_
 
double zBeam_
 

Additional Inherited Members

- Public Types inherited from edm::one::EDAnalyzerBase
typedef EDAnalyzerBase ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- 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 ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
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)
 

Detailed Description

Definition at line 46 of file HGCalTimingAnalyzer.cc.

Constructor & Destructor Documentation

HGCalTimingAnalyzer::HGCalTimingAnalyzer ( edm::ParameterSet const &  )
explicit

Definition at line 81 of file HGCalTimingAnalyzer.cc.

References gather_cfg::cout, detectorBeam_, detectorEE_, doTree_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), groupHits_, triggerObjects_cff::id, idBeams_, AlCaHLTBitMon_QueryRunRegistry::string, timeUnit_, tok_hepMC_, tok_hitsBeam_, tok_hitsEE_, tok_simTk_, and tok_simVtx_.

81  {
82  usesResource("TFileService");
83 
84  // now do whatever initialization is needed
85  detectorEE_ = iConfig.getParameter<std::string>("DetectorEE");
86  detectorBeam_ = iConfig.getParameter<std::string>("DetectorBeam");
87  // Group hits (if groupHits_ = true) if hits come within timeUnit_
88  groupHits_ = iConfig.getParameter<bool>("GroupHits");
89  timeUnit_ = iConfig.getParameter<double>("TimeUnit");
90  // Only look into the beam counters with ID's as in idBeams_
91  idBeams_ = iConfig.getParameter<std::vector<int>>("IDBeams");
92  doTree_ = iConfig.getUntrackedParameter<bool>("DoTree", false);
93  if (!groupHits_) timeUnit_ = 0.000001;
94 #ifdef EDM_ML_DEBUG
95  std::cout << "HGCalTimingAnalyzer:: Group Hits " << groupHits_ << " in "
96  << timeUnit_ << " IdBeam " << idBeams_.size() << ":";
97  for (const auto& id : idBeams_) std::cout << " " << id;
98  std::cout << std::endl;
99 #endif
100  if (idBeams_.empty()) idBeams_.push_back(1001);
101 
102  edm::InputTag tmp0 = iConfig.getParameter<edm::InputTag>("GeneratorSrc");
103  tok_hepMC_ = consumes<edm::HepMCProduct>(tmp0);
104 #ifdef EDM_ML_DEBUG
105  std::cout << "HGCalTimingAnalyzer:: GeneratorSource = " << tmp0 << std::endl;
106 #endif
107  tok_simTk_ = consumes<edm::SimTrackContainer>(edm::InputTag("g4SimHits"));
108  tok_simVtx_ = consumes<edm::SimVertexContainer>(edm::InputTag("g4SimHits"));
109  std::string tmp1 = iConfig.getParameter<std::string>("CaloHitSrcEE");
110  tok_hitsEE_ =
111  consumes<edm::PCaloHitContainer>(edm::InputTag("g4SimHits", tmp1));
112 #ifdef EDM_ML_DEBUG
113  std::cout << "HGCalTimingAnalyzer:: Detector " << detectorEE_ << " with tags "
114  << tmp1 << std::endl;
115 #endif
116  tmp1 = iConfig.getParameter<std::string>("CaloHitSrcBeam");
117  tok_hitsBeam_ =
118  consumes<edm::PCaloHitContainer>(edm::InputTag("g4SimHits", tmp1));
119 #ifdef EDM_ML_DEBUG
120  std::cout << "HGCalTimingAnalyzer:: Detector " << detectorBeam_
121  << " with tags " << tmp1 << std::endl;
122 #endif
123 }
edm::EDGetTokenT< edm::SimTrackContainer > tok_simTk_
edm::EDGetTokenT< edm::PCaloHitContainer > tok_hitsEE_
edm::EDGetTokenT< edm::HepMCProduct > tok_hepMC_
std::vector< int > idBeams_
edm::EDGetTokenT< edm::SimVertexContainer > tok_simVtx_
edm::EDGetTokenT< edm::PCaloHitContainer > tok_hitsBeam_
HGCalTimingAnalyzer::~HGCalTimingAnalyzer ( )
override

Definition at line 125 of file HGCalTimingAnalyzer.cc.

125 {}

Member Function Documentation

void HGCalTimingAnalyzer::analyze ( edm::Event const &  ,
edm::EventSetup const &   
)
overrideprivatevirtual

Implements edm::one::EDAnalyzerBase.

Definition at line 172 of file HGCalTimingAnalyzer.cc.

References analyzeSimHits(), analyzeSimTracks(), edm::HandleBase::clear(), gather_cfg::cout, detectorBeam_, detectorEE_, doTree_, edm::Event::getByToken(), edm::HepMCProduct::GetEvent(), edm::HandleBase::isValid(), gen::k, AlCaHLTBitMon_ParallelJobs::p, simHitCellEnBeam_, simHitCellEnEE_, simHitCellIdBeam_, simHitCellIdEE_, simHitCellTmBeam_, simHitCellTmEE_, tok_hepMC_, tok_hitsBeam_, tok_hitsEE_, tok_simTk_, tok_simVtx_, and tree_.

Referenced by endRun().

173  {
174 #ifdef EDM_ML_DEBUG
175  // Generator input
177  iEvent.getByToken(tok_hepMC_, evtMC);
178  if (!evtMC.isValid()) {
179  edm::LogWarning("HGCal") << "no HepMCProduct found";
180  } else {
181  const HepMC::GenEvent* myGenEvent = evtMC->GetEvent();
182  unsigned int k(0);
183  for (HepMC::GenEvent::particle_const_iterator p =
184  myGenEvent->particles_begin();
185  p != myGenEvent->particles_end(); ++p, ++k) {
186  std::cout << "Particle[" << k << "] with p " << (*p)->momentum().rho()
187  << " theta " << (*p)->momentum().theta() << " phi "
188  << (*p)->momentum().phi() << std::endl;
189  }
190  }
191 #endif
192 
193  // Now the Simhits
195  iEvent.getByToken(tok_simTk_, SimTk);
197  iEvent.getByToken(tok_simVtx_, SimVtx);
198  analyzeSimTracks(SimTk, SimVtx);
199 
200  simHitCellIdEE_.clear();
201  simHitCellIdBeam_.clear();
202  simHitCellEnEE_.clear();
203  simHitCellEnBeam_.clear();
204  simHitCellTmEE_.clear();
205  simHitCellTmBeam_.clear();
206 
207  edm::Handle<edm::PCaloHitContainer> theCaloHitContainers;
208  std::vector<PCaloHit> caloHits;
209  iEvent.getByToken(tok_hitsEE_, theCaloHitContainers);
210  if (theCaloHitContainers.isValid()) {
211 #ifdef EDM_ML_DEBUG
212  std::cout << "PcalohitContainer for " << detectorEE_ << " has "
213  << theCaloHitContainers->size() << " hits" << std::endl;
214 #endif
215  caloHits.clear();
216  caloHits.insert(caloHits.end(), theCaloHitContainers->begin(),
217  theCaloHitContainers->end());
218  analyzeSimHits(0, caloHits);
219  } else {
220 #ifdef EDM_ML_DEBUG
221  std::cout << "PCaloHitContainer does not exist for " << detectorEE_
222  << " !!!" << std::endl;
223 #endif
224  }
225 
226  iEvent.getByToken(tok_hitsBeam_, theCaloHitContainers);
227  if (theCaloHitContainers.isValid()) {
228 #ifdef EDM_ML_DEBUG
229  std::cout << "PcalohitContainer for " << detectorBeam_ << " has "
230  << theCaloHitContainers->size() << " hits" << std::endl;
231 #endif
232  caloHits.clear();
233  caloHits.insert(caloHits.end(), theCaloHitContainers->begin(),
234  theCaloHitContainers->end());
235  analyzeSimHits(1, caloHits);
236  } else {
237 #ifdef EDM_ML_DEBUG
238  std::cout << "PCaloHitContainer does not exist for " << detectorBeam_
239  << " !!!" << std::endl;
240 #endif
241  }
242  if (doTree_) tree_->Fill();
243 }
edm::EDGetTokenT< edm::SimTrackContainer > tok_simTk_
edm::EDGetTokenT< edm::PCaloHitContainer > tok_hitsEE_
edm::EDGetTokenT< edm::HepMCProduct > tok_hepMC_
int iEvent
Definition: GenABIO.cc:224
std::vector< float > simHitCellEnEE_
std::vector< float > simHitCellTmBeam_
bool isValid() const
Definition: HandleBase.h:74
std::vector< float > simHitCellEnBeam_
std::vector< uint32_t > simHitCellIdBeam_
std::vector< float > simHitCellTmEE_
std::vector< uint32_t > simHitCellIdEE_
int k[5][pyjets_maxn]
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
edm::EDGetTokenT< edm::SimVertexContainer > tok_simVtx_
void analyzeSimHits(int type, std::vector< PCaloHit > const &hits)
void analyzeSimTracks(edm::Handle< edm::SimTrackContainer > const &SimTk, edm::Handle< edm::SimVertexContainer > const &SimVtx)
edm::EDGetTokenT< edm::PCaloHitContainer > tok_hitsBeam_
void HGCalTimingAnalyzer::analyzeSimHits ( int  type,
std::vector< PCaloHit > const &  hits 
)
private

Definition at line 245 of file HGCalTimingAnalyzer.cc.

References gather_cfg::cout, TauDecayModes::dec, plotBeamSpotDB::first, hgcons_, mps_fire::i, hit::id, crabWrapper::key, DetId::rawId(), edm::second(), simHitCellEnBeam_, simHitCellEnEE_, simHitCellIdBeam_, simHitCellIdEE_, simHitCellTmBeam_, simHitCellTmEE_, HGCalDDDConstants::simToReco(), ntuplemaker::time, timeUnit_, HGCalTestNumbering::unpackHexagonIndex(), HcalTestBeamNumbering::unpackIndex(), x, y, and ecaldqm::zside().

Referenced by analyze(), and endRun().

246  {
247 #ifdef EDM_ML_DEBUG
248  unsigned int i(0);
249 #endif
250  std::map<std::pair<uint32_t, uint64_t>, std::pair<double, double>> map_hits;
251  for (const auto& hit : hits) {
252  double energy = hit.energy();
253  double time = hit.time();
254  uint32_t id = hit.id();
255  if (type == 0) {
256  int subdet, zside, layer, sector, subsector, cell;
257  HGCalTestNumbering::unpackHexagonIndex(id, subdet, zside, layer, sector,
258  subsector, cell);
259  std::pair<int, int> recoLayerCell =
260  hgcons_->simToReco(cell, layer, sector, true);
261  id = HGCalDetId((ForwardSubdetector)(subdet), zside, recoLayerCell.second,
262  subsector, sector, recoLayerCell.first)
263  .rawId();
264 #ifdef EDM_ML_DEBUG
265  std::cout << "SimHit:Hit[" << i << "] Id " << subdet << ":" << zside
266  << ":" << layer << ":" << sector << ":" << subsector << ":"
267  << recoLayerCell.first << ":" << recoLayerCell.second
268  << " Energy " << energy << " Time " << time << std::endl;
269 #endif
270  } else {
271 #ifdef EDM_ML_DEBUG
272  int subdet, layer, x, y;
273  HcalTestBeamNumbering::unpackIndex(id, subdet, layer, x, y);
274  std::cout << "SimHit:Hit[" << i << "] Beam Subdet " << subdet << " Layer "
275  << layer << " x|y " << x << ":" << y << " Energy " << energy
276  << " Time " << time << std::endl;
277 #endif
278  }
279  uint64_t tid = (uint64_t)((time + 50.0) / timeUnit_);
280  std::pair<uint32_t, uint64_t> key(id, tid);
281  auto itr = map_hits.find(key);
282  if (itr == map_hits.end()) {
283  map_hits[key] = std::pair<double, double>(time, 0.0);
284  itr = map_hits.find(key);
285  }
286  energy += (itr->second).second;
287  map_hits[key] = std::pair<double, double>((itr->second).first, energy);
288 #ifdef EDM_ML_DEBUG
289  ++i;
290 #endif
291  }
292 
293 #ifdef EDM_ML_DEBUG
294  std::cout << "analyzeSimHits: Finds " << map_hits.size() << " hits "
295  << " from the Hit Vector of size " << hits.size() << " for type "
296  << type << std::endl;
297 #endif
298  for (const auto& itr : map_hits) {
299  uint32_t id = (itr.first).first;
300  double time = (itr.second).first;
301  double energy = (itr.second).second;
302  if (type == 0) {
303  simHitCellIdEE_.push_back(id);
304  simHitCellEnEE_.push_back(energy);
305  simHitCellTmEE_.push_back(time);
306  } else {
307  simHitCellIdBeam_.push_back(id);
308  simHitCellEnBeam_.push_back(energy);
309  simHitCellTmBeam_.push_back(time);
310  }
311 #ifdef EDM_ML_DEBUG
312  std::cout << "SimHit::ID: " << std::hex << id << std::dec << " T: " << time
313  << " E: " << energy << std::endl;
314 #endif
315  }
316 }
type
Definition: HCALResponse.h:21
const HGCalDDDConstants * hgcons_
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
int zside(DetId const &)
ForwardSubdetector
U second(std::pair< T, U > const &p)
std::pair< int, int > simToReco(int cell, int layer, int mod, bool half) const
std::vector< float > simHitCellEnEE_
std::vector< float > simHitCellTmBeam_
std::vector< float > simHitCellEnBeam_
std::vector< uint32_t > simHitCellIdBeam_
std::vector< float > simHitCellTmEE_
std::vector< uint32_t > simHitCellIdEE_
unsigned int id
unsigned long long uint64_t
Definition: Time.h:15
static void unpackIndex(const uint32_t &idx, int &det, int &lay, int &x, int &y)
static void unpackHexagonIndex(const uint32_t &idx, int &subdet, int &z, int &lay, int &wafer, int &celltyp, int &cell)
void HGCalTimingAnalyzer::analyzeSimTracks ( edm::Handle< edm::SimTrackContainer > const &  SimTk,
edm::Handle< edm::SimVertexContainer > const &  SimVtx 
)
private

Definition at line 318 of file HGCalTimingAnalyzer.cc.

References gather_cfg::cout, DEFINE_FWK_MODULE, pBeam_, xBeam_, yBeam_, and zBeam_.

Referenced by analyze(), and endRun().

320  {
321  xBeam_ = yBeam_ = zBeam_ = pBeam_ = -1000000;
322  int vertIndex(-1);
323  for (edm::SimTrackContainer::const_iterator simTrkItr = SimTk->begin();
324  simTrkItr != SimTk->end(); simTrkItr++) {
325 #ifdef EDM_ML_DEBUG
326  std::cout << "Track " << simTrkItr->trackId() << " Vertex "
327  << simTrkItr->vertIndex() << " Type " << simTrkItr->type()
328  << " Charge " << simTrkItr->charge() << " momentum "
329  << simTrkItr->momentum() << " " << simTrkItr->momentum().P()
330  << std::endl;
331 #endif
332  if (vertIndex == -1) {
333  vertIndex = simTrkItr->vertIndex();
334  pBeam_ = simTrkItr->momentum().P();
335  }
336  }
337  if (vertIndex != -1 && vertIndex < (int)SimVtx->size()) {
338  edm::SimVertexContainer::const_iterator simVtxItr = SimVtx->begin();
339  for (int iv = 0; iv < vertIndex; iv++) simVtxItr++;
340 #ifdef EDM_ML_DEBUG
341  std::cout << "Vertex " << vertIndex << " position " << simVtxItr->position()
342  << std::endl;
343 #endif
344  xBeam_ = simVtxItr->position().X();
345  yBeam_ = simVtxItr->position().Y();
346  zBeam_ = simVtxItr->position().Z();
347  }
348 }
void HGCalTimingAnalyzer::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 143 of file HGCalTimingAnalyzer.cc.

References detectorEE_, doTree_, fs_, TFileService::make(), pBeam_, simHitCellEnBeam_, simHitCellEnEE_, simHitCellIdBeam_, simHitCellIdEE_, simHitCellTmBeam_, simHitCellTmEE_, AlCaHLTBitMon_QueryRunRegistry::string, tree_, xBeam_, yBeam_, and zBeam_.

143  {
145  if (doTree_) {
146  tree_ = fs_->make<TTree>("HGCTB", "SimHitEnergy");
147  tree_->Branch("xBeam", &xBeam_, "xBeam/D");
148  tree_->Branch("yBeam", &yBeam_, "yBeam/D");
149  tree_->Branch("zBeam", &zBeam_, "zBeam/D");
150  tree_->Branch("pBeam", &pBeam_, "pBeam/D");
151  tree_->Branch("simHitCellIdEE_", &simHitCellIdEE_);
152  tree_->Branch("simHitCellEnEE_", &simHitCellEnEE_);
153  tree_->Branch("simHitCellTmEE_", &simHitCellTmEE_);
154  tree_->Branch("simHitCellIdBeam_", &simHitCellIdBeam_);
155  tree_->Branch("simHitCellEnBeam_", &simHitCellEnBeam_);
156  tree_->Branch("simHitCellTmBeam_", &simHitCellTmBeam_);
157  }
158 }
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
std::vector< float > simHitCellEnEE_
std::vector< float > simHitCellTmBeam_
std::vector< float > simHitCellEnBeam_
std::vector< uint32_t > simHitCellIdBeam_
std::vector< float > simHitCellTmEE_
std::vector< uint32_t > simHitCellIdEE_
edm::Service< TFileService > fs_
void HGCalTimingAnalyzer::beginRun ( edm::Run const &  ,
edm::EventSetup const &   
)
overrideprivate

Definition at line 160 of file HGCalTimingAnalyzer.cc.

References gather_cfg::cout, detectorEE_, edm::EventSetup::get(), hgcons_, and HGCalDDDConstants::layers().

161  {
163  iSetup.get<IdealGeometryRecord>().get(detectorEE_, pHGDC);
164  hgcons_ = &(*pHGDC);
165 
166 #ifdef EDM_ML_DEBUG
167  std::cout << "HGCalTimingAnalyzer::" << detectorEE_ << " defined with "
168  << hgcons_->layers(false) << " layers" << std::endl;
169 #endif
170 }
const HGCalDDDConstants * hgcons_
unsigned int layers(bool reco) const
void HGCalTimingAnalyzer::endRun ( edm::Run const &  ,
edm::EventSetup const &   
)
inlineoverrideprivate
void HGCalTimingAnalyzer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 127 of file HGCalTimingAnalyzer.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), edm::ParameterSetDescription::addUntracked(), and AlCaHLTBitMon_QueryRunRegistry::string.

128  {
130  desc.add<std::string>("DetectorEE", "HGCalEESensitive");
131  desc.add<std::string>("DetectorBeam", "HcalTB06BeamDetector");
132  desc.add<bool>("GroupHits", false);
133  desc.add<double>("TimeUnit", 0.001);
134  std::vector<int> ids = {1001, 1002, 1003, 1004, 1005};
135  desc.add<std::vector<int>>("IDBeams", ids);
136  desc.addUntracked<bool>("DoTree", true);
137  desc.add<edm::InputTag>("GeneratorSrc", edm::InputTag("generatorSmeared"));
138  desc.add<std::string>("CaloHitSrcEE", "HGCHitsEE");
139  desc.add<std::string>("CaloHitSrcBeam", "HcalTB06BeamHits");
140  descriptions.add("HGCalTimingAnalyzer", desc);
141 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)

Member Data Documentation

std::string HGCalTimingAnalyzer::detectorBeam_
private

Definition at line 67 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), and HGCalTimingAnalyzer().

std::string HGCalTimingAnalyzer::detectorEE_
private

Definition at line 67 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), beginJob(), beginRun(), and HGCalTimingAnalyzer().

bool HGCalTimingAnalyzer::doTree_
private

Definition at line 66 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), beginJob(), and HGCalTimingAnalyzer().

edm::Service<TFileService> HGCalTimingAnalyzer::fs_
private

Definition at line 64 of file HGCalTimingAnalyzer.cc.

Referenced by beginJob().

bool HGCalTimingAnalyzer::groupHits_
private

Definition at line 66 of file HGCalTimingAnalyzer.cc.

Referenced by HGCalTimingAnalyzer().

const HGCalDDDConstants* HGCalTimingAnalyzer::hgcons_
private

Definition at line 65 of file HGCalTimingAnalyzer.cc.

Referenced by analyzeSimHits(), and beginRun().

std::vector<int> HGCalTimingAnalyzer::idBeams_
private

Definition at line 69 of file HGCalTimingAnalyzer.cc.

Referenced by HGCalTimingAnalyzer().

double HGCalTimingAnalyzer::pBeam_
private

Definition at line 78 of file HGCalTimingAnalyzer.cc.

Referenced by analyzeSimTracks(), and beginJob().

std::vector<float> HGCalTimingAnalyzer::simHitCellEnBeam_
private

Definition at line 76 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), analyzeSimHits(), and beginJob().

std::vector<float> HGCalTimingAnalyzer::simHitCellEnEE_
private

Definition at line 76 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), analyzeSimHits(), and beginJob().

std::vector<uint32_t> HGCalTimingAnalyzer::simHitCellIdBeam_
private

Definition at line 75 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), analyzeSimHits(), and beginJob().

std::vector<uint32_t> HGCalTimingAnalyzer::simHitCellIdEE_
private

Definition at line 75 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), analyzeSimHits(), and beginJob().

std::vector<float> HGCalTimingAnalyzer::simHitCellTmBeam_
private

Definition at line 77 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), analyzeSimHits(), and beginJob().

std::vector<float> HGCalTimingAnalyzer::simHitCellTmEE_
private

Definition at line 77 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), analyzeSimHits(), and beginJob().

double HGCalTimingAnalyzer::timeUnit_
private

Definition at line 68 of file HGCalTimingAnalyzer.cc.

Referenced by analyzeSimHits(), and HGCalTimingAnalyzer().

edm::EDGetTokenT<edm::HepMCProduct> HGCalTimingAnalyzer::tok_hepMC_
private

Definition at line 73 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), and HGCalTimingAnalyzer().

edm::EDGetTokenT<edm::PCaloHitContainer> HGCalTimingAnalyzer::tok_hitsBeam_
private

Definition at line 70 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), and HGCalTimingAnalyzer().

edm::EDGetTokenT<edm::PCaloHitContainer> HGCalTimingAnalyzer::tok_hitsEE_
private

Definition at line 70 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), and HGCalTimingAnalyzer().

edm::EDGetTokenT<edm::SimTrackContainer> HGCalTimingAnalyzer::tok_simTk_
private

Definition at line 71 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), and HGCalTimingAnalyzer().

edm::EDGetTokenT<edm::SimVertexContainer> HGCalTimingAnalyzer::tok_simVtx_
private

Definition at line 72 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), and HGCalTimingAnalyzer().

TTree* HGCalTimingAnalyzer::tree_
private

Definition at line 74 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), and beginJob().

double HGCalTimingAnalyzer::xBeam_
private

Definition at line 78 of file HGCalTimingAnalyzer.cc.

Referenced by analyzeSimTracks(), and beginJob().

double HGCalTimingAnalyzer::yBeam_
private

Definition at line 78 of file HGCalTimingAnalyzer.cc.

Referenced by analyzeSimTracks(), and beginJob().

double HGCalTimingAnalyzer::zBeam_
private

Definition at line 78 of file HGCalTimingAnalyzer.cc.

Referenced by analyzeSimTracks(), and beginJob().