CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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
 
 EDAnalyzer (const EDAnalyzer &)=delete
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
const EDAnalyzeroperator= (const EDAnalyzer &)=delete
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
bool wantsInputProcessBlocks () const final
 
bool wantsProcessBlocks () 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
 
std::vector< ESProxyIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex >
const & 
esGetTokenRecordIndicesVector (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::array< std::vector< ModuleDescription const * > *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, 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
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
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

const std::string detectorBeam_
 
const 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::HepMCProduct
tok_hepMC_
 
edm::EDGetTokenT
< edm::PCaloHitContainer
tok_hitsBeam_
 
edm::EDGetTokenT
< edm::PCaloHitContainer
tok_hitsEE_
 
edm::EDGetTokenT
< edm::SimTrackContainer
tok_simTk_
 
edm::EDGetTokenT
< edm::SimVertexContainer
tok_simVtx_
 
const edm::ESGetToken
< HGCalDDDConstants,
IdealGeometryRecord
tokDDD_
 
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)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< B > consumes (edm::InputTag tag) noexcept
 
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<Transition Tr = Transition::Event>
constexpr auto esConsumes () noexcept
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag) noexcept
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
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)
 
void resetItemsToGetFrom (BranchType iType)
 

Detailed Description

Definition at line 44 of file HGCalTimingAnalyzer.cc.

Constructor & Destructor Documentation

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

Definition at line 78 of file HGCalTimingAnalyzer.cc.

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

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

Definition at line 124 of file HGCalTimingAnalyzer.cc.

124 {}

Member Function Documentation

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

Implements edm::one::EDAnalyzerBase.

Definition at line 167 of file HGCalTimingAnalyzer.cc.

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

167  {
168 #ifdef EDM_ML_DEBUG
169  // Generator input
171  iEvent.getByToken(tok_hepMC_, evtMC);
172  if (!evtMC.isValid()) {
173  edm::LogWarning("HGCal") << "no HepMCProduct found";
174  } else {
175  const HepMC::GenEvent* myGenEvent = evtMC->GetEvent();
176  unsigned int k(0);
177  for (HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end();
178  ++p, ++k) {
179  edm::LogVerbatim("HGCSim") << "Particle[" << k << "] with p " << (*p)->momentum().rho() << " theta "
180  << (*p)->momentum().theta() << " phi " << (*p)->momentum().phi();
181  }
182  }
183 #endif
184 
185  // Now the Simhits
187  iEvent.getByToken(tok_simTk_, SimTk);
189  iEvent.getByToken(tok_simVtx_, SimVtx);
190  analyzeSimTracks(SimTk, SimVtx);
191 
192  simHitCellIdEE_.clear();
193  simHitCellIdBeam_.clear();
194  simHitCellEnEE_.clear();
195  simHitCellEnBeam_.clear();
196  simHitCellTmEE_.clear();
197  simHitCellTmBeam_.clear();
198 
199  edm::Handle<edm::PCaloHitContainer> theCaloHitContainers;
200  std::vector<PCaloHit> caloHits;
201  iEvent.getByToken(tok_hitsEE_, theCaloHitContainers);
202  if (theCaloHitContainers.isValid()) {
203 #ifdef EDM_ML_DEBUG
204  edm::LogVerbatim("HGCSim") << "PcalohitContainer for " << detectorEE_ << " has " << theCaloHitContainers->size()
205  << " hits";
206 #endif
207  caloHits.clear();
208  caloHits.insert(caloHits.end(), theCaloHitContainers->begin(), theCaloHitContainers->end());
209  analyzeSimHits(0, caloHits);
210  } else {
211 #ifdef EDM_ML_DEBUG
212  edm::LogVerbatim("HGCSim") << "PCaloHitContainer does not exist for " << detectorEE_ << " !!!";
213 #endif
214  }
215 
216  iEvent.getByToken(tok_hitsBeam_, theCaloHitContainers);
217  if (theCaloHitContainers.isValid()) {
218 #ifdef EDM_ML_DEBUG
219  edm::LogVerbatim("HGCSim") << "PcalohitContainer for " << detectorBeam_ << " has " << theCaloHitContainers->size()
220  << " hits";
221 #endif
222  caloHits.clear();
223  caloHits.insert(caloHits.end(), theCaloHitContainers->begin(), theCaloHitContainers->end());
224  analyzeSimHits(1, caloHits);
225  } else {
226 #ifdef EDM_ML_DEBUG
227  edm::LogVerbatim("HGCSim") << "PCaloHitContainer does not exist for " << detectorBeam_ << " !!!";
228 #endif
229  }
230  if (doTree_)
231  tree_->Fill();
232 }
Log< level::Info, true > LogVerbatim
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_
const std::string detectorBeam_
bool isValid() const
Definition: HandleBase.h:70
std::vector< float > simHitCellEnBeam_
std::vector< uint32_t > simHitCellIdBeam_
std::vector< float > simHitCellTmEE_
std::vector< uint32_t > simHitCellIdEE_
edm::EDGetTokenT< edm::SimVertexContainer > tok_simVtx_
const std::string detectorEE_
Log< level::Warning, false > LogWarning
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 234 of file HGCalTimingAnalyzer.cc.

References TauDecayModes::dec, relval_parameters_module::energy, first, hgcons_, mps_fire::i, hit::id, submitPVResolutionJobs::key, phase1PixelTopology::layer, DetId::rawId(), edm::second(), simHitCellEnBeam_, simHitCellEnEE_, simHitCellIdBeam_, simHitCellIdEE_, simHitCellTmBeam_, simHitCellTmEE_, HGCalDDDConstants::simToReco(), timeUnit_, HGCalTestNumbering::unpackHexagonIndex(), HcalTestBeamNumbering::unpackIndex(), x, y, and ecaldqm::zside().

Referenced by analyze().

234  {
235 #ifdef EDM_ML_DEBUG
236  unsigned int i(0);
237 #endif
238  std::map<std::pair<uint32_t, uint64_t>, std::pair<double, double>> map_hits;
239  for (const auto& hit : hits) {
240  double energy = hit.energy();
241  double time = hit.time();
242  uint32_t id = hit.id();
243  if (type == 0) {
244  int subdet, zside, layer, sector, subsector, cell;
245  HGCalTestNumbering::unpackHexagonIndex(id, subdet, zside, layer, sector, subsector, cell);
246  std::pair<int, int> recoLayerCell = hgcons_->simToReco(cell, layer, sector, true);
247  id = HGCalDetId((ForwardSubdetector)(subdet), zside, recoLayerCell.second, subsector, sector, recoLayerCell.first)
248  .rawId();
249 #ifdef EDM_ML_DEBUG
250  edm::LogVerbatim("HGCSim") << "SimHit:Hit[" << i << "] Id " << subdet << ":" << zside << ":" << layer << ":"
251  << sector << ":" << subsector << ":" << recoLayerCell.first << ":"
252  << recoLayerCell.second << " Energy " << energy << " Time " << time;
253 #endif
254  } else {
255 #ifdef EDM_ML_DEBUG
256  int subdet, layer, x, y;
257  HcalTestBeamNumbering::unpackIndex(id, subdet, layer, x, y);
258  edm::LogVerbatim("HGCSim") << "SimHit:Hit[" << i << "] Beam Subdet " << subdet << " Layer " << layer << " x|y "
259  << x << ":" << y << " Energy " << energy << " Time " << time;
260 #endif
261  }
262  uint64_t tid = (uint64_t)((time + 50.0) / timeUnit_);
263  std::pair<uint32_t, uint64_t> key(id, tid);
264  auto itr = map_hits.find(key);
265  if (itr == map_hits.end()) {
266  map_hits[key] = std::pair<double, double>(time, 0.0);
267  itr = map_hits.find(key);
268  }
269  energy += (itr->second).second;
270  map_hits[key] = std::pair<double, double>((itr->second).first, energy);
271 #ifdef EDM_ML_DEBUG
272  ++i;
273 #endif
274  }
275 
276 #ifdef EDM_ML_DEBUG
277  edm::LogVerbatim("HGCSim") << "analyzeSimHits: Finds " << map_hits.size() << " hits "
278  << " from the Hit Vector of size " << hits.size() << " for type " << type;
279 #endif
280  for (const auto& itr : map_hits) {
281  uint32_t id = (itr.first).first;
282  double time = (itr.second).first;
283  double energy = (itr.second).second;
284  if (type == 0) {
285  simHitCellIdEE_.push_back(id);
286  simHitCellEnEE_.push_back(energy);
287  simHitCellTmEE_.push_back(time);
288  } else {
289  simHitCellIdBeam_.push_back(id);
290  simHitCellEnBeam_.push_back(energy);
291  simHitCellTmBeam_.push_back(time);
292  }
293 #ifdef EDM_ML_DEBUG
294  edm::LogVerbatim("HGCSim") << "SimHit::ID: " << std::hex << id << std::dec << " T: " << time << " E: " << energy;
295 #endif
296  }
297 }
Log< level::Info, true > LogVerbatim
const HGCalDDDConstants * hgcons_
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
int zside(DetId const &)
ForwardSubdetector
constexpr std::array< uint8_t, layerIndexSize > layer
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_
tuple key
prepare the HTCondor submission files and eventually submit them
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:13
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 299 of file HGCalTimingAnalyzer.cc.

References gpuVertexFinder::iv, pBeam_, xBeam_, yBeam_, and zBeam_.

Referenced by analyze().

300  {
301  xBeam_ = yBeam_ = zBeam_ = pBeam_ = -1000000;
302  int vertIndex(-1);
303  for (edm::SimTrackContainer::const_iterator simTrkItr = SimTk->begin(); simTrkItr != SimTk->end(); simTrkItr++) {
304 #ifdef EDM_ML_DEBUG
305  edm::LogVerbatim("HGCSim") << "Track " << simTrkItr->trackId() << " Vertex " << simTrkItr->vertIndex() << " Type "
306  << simTrkItr->type() << " Charge " << simTrkItr->charge() << " momentum "
307  << simTrkItr->momentum() << " " << simTrkItr->momentum().P();
308 #endif
309  if (vertIndex == -1) {
310  vertIndex = simTrkItr->vertIndex();
311  pBeam_ = simTrkItr->momentum().P();
312  }
313  }
314  if (vertIndex != -1 && vertIndex < (int)SimVtx->size()) {
315  edm::SimVertexContainer::const_iterator simVtxItr = SimVtx->begin();
316  for (int iv = 0; iv < vertIndex; iv++)
317  simVtxItr++;
318 #ifdef EDM_ML_DEBUG
319  edm::LogVerbatim("HGCSim") << "Vertex " << vertIndex << " position " << simVtxItr->position();
320 #endif
321  xBeam_ = simVtxItr->position().X();
322  yBeam_ = simVtxItr->position().Y();
323  zBeam_ = simVtxItr->position().Z();
324  }
325 }
Log< level::Info, true > LogVerbatim
int32_t *__restrict__ iv
void HGCalTimingAnalyzer::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 141 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_.

141  {
143  if (doTree_) {
144  tree_ = fs_->make<TTree>("HGCTB", "SimHitEnergy");
145  tree_->Branch("xBeam", &xBeam_, "xBeam/D");
146  tree_->Branch("yBeam", &yBeam_, "yBeam/D");
147  tree_->Branch("zBeam", &zBeam_, "zBeam/D");
148  tree_->Branch("pBeam", &pBeam_, "pBeam/D");
149  tree_->Branch("simHitCellIdEE_", &simHitCellIdEE_);
150  tree_->Branch("simHitCellEnEE_", &simHitCellEnEE_);
151  tree_->Branch("simHitCellTmEE_", &simHitCellTmEE_);
152  tree_->Branch("simHitCellIdBeam_", &simHitCellIdBeam_);
153  tree_->Branch("simHitCellEnBeam_", &simHitCellEnBeam_);
154  tree_->Branch("simHitCellTmBeam_", &simHitCellTmBeam_);
155  }
156 }
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_
const std::string detectorEE_
void HGCalTimingAnalyzer::beginRun ( edm::Run const &  ,
edm::EventSetup const &  iSetup 
)
overrideprivate

Definition at line 158 of file HGCalTimingAnalyzer.cc.

References detectorEE_, edm::EventSetup::getData(), hgcons_, HGCalDDDConstants::layers(), and tokDDD_.

158  {
159  hgcons_ = &iSetup.getData(tokDDD_);
160 
161 #ifdef EDM_ML_DEBUG
162  edm::LogVerbatim("HGCSim") << "HGCalTimingAnalyzer::" << detectorEE_ << " defined with " << hgcons_->layers(false)
163  << " layers";
164 #endif
165 }
Log< level::Info, true > LogVerbatim
const HGCalDDDConstants * hgcons_
const edm::ESGetToken< HGCalDDDConstants, IdealGeometryRecord > tokDDD_
unsigned int layers(bool reco) const
const std::string detectorEE_
void HGCalTimingAnalyzer::endRun ( edm::Run const &  ,
edm::EventSetup const &   
)
inlineoverrideprivate

Definition at line 54 of file HGCalTimingAnalyzer.cc.

54 {}
void HGCalTimingAnalyzer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 126 of file HGCalTimingAnalyzer.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), edm::ParameterSetDescription::addUntracked(), submitPVResolutionJobs::desc, HLT_FULL_cff::InputTag, and AlCaHLTBitMon_QueryRunRegistry::string.

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

Member Data Documentation

const std::string HGCalTimingAnalyzer::detectorBeam_
private

Definition at line 61 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), and HGCalTimingAnalyzer().

const std::string HGCalTimingAnalyzer::detectorEE_
private

Definition at line 61 of file HGCalTimingAnalyzer.cc.

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

bool HGCalTimingAnalyzer::doTree_
private

Definition at line 64 of file HGCalTimingAnalyzer.cc.

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

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

Definition at line 60 of file HGCalTimingAnalyzer.cc.

Referenced by beginJob().

bool HGCalTimingAnalyzer::groupHits_
private

Definition at line 64 of file HGCalTimingAnalyzer.cc.

Referenced by HGCalTimingAnalyzer().

const HGCalDDDConstants* HGCalTimingAnalyzer::hgcons_
private

Definition at line 63 of file HGCalTimingAnalyzer.cc.

Referenced by analyzeSimHits(), and beginRun().

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

Definition at line 66 of file HGCalTimingAnalyzer.cc.

Referenced by HGCalTimingAnalyzer().

double HGCalTimingAnalyzer::pBeam_
private

Definition at line 75 of file HGCalTimingAnalyzer.cc.

Referenced by analyzeSimTracks(), and beginJob().

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

Definition at line 73 of file HGCalTimingAnalyzer.cc.

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

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

Definition at line 73 of file HGCalTimingAnalyzer.cc.

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

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

Definition at line 72 of file HGCalTimingAnalyzer.cc.

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

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

Definition at line 72 of file HGCalTimingAnalyzer.cc.

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

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

Definition at line 74 of file HGCalTimingAnalyzer.cc.

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

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

Definition at line 74 of file HGCalTimingAnalyzer.cc.

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

double HGCalTimingAnalyzer::timeUnit_
private

Definition at line 65 of file HGCalTimingAnalyzer.cc.

Referenced by analyzeSimHits(), and HGCalTimingAnalyzer().

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

Definition at line 70 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), and HGCalTimingAnalyzer().

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

Definition at line 67 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), and HGCalTimingAnalyzer().

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

Definition at line 67 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), and HGCalTimingAnalyzer().

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

Definition at line 68 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), and HGCalTimingAnalyzer().

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

Definition at line 69 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), and HGCalTimingAnalyzer().

const edm::ESGetToken<HGCalDDDConstants, IdealGeometryRecord> HGCalTimingAnalyzer::tokDDD_
private

Definition at line 62 of file HGCalTimingAnalyzer.cc.

Referenced by beginRun().

TTree* HGCalTimingAnalyzer::tree_
private

Definition at line 71 of file HGCalTimingAnalyzer.cc.

Referenced by analyze(), and beginJob().

double HGCalTimingAnalyzer::xBeam_
private

Definition at line 75 of file HGCalTimingAnalyzer.cc.

Referenced by analyzeSimTracks(), and beginJob().

double HGCalTimingAnalyzer::yBeam_
private

Definition at line 75 of file HGCalTimingAnalyzer.cc.

Referenced by analyzeSimTracks(), and beginJob().

double HGCalTimingAnalyzer::zBeam_
private

Definition at line 75 of file HGCalTimingAnalyzer.cc.

Referenced by analyzeSimTracks(), and beginJob().