CMS 3D CMS Logo

List of all members | Public Member Functions | Private Types | Private Member Functions | Private Attributes | Static Private Attributes
CaloSteppingAction Class Reference
Inheritance diagram for CaloSteppingAction:
SimProducer Observer< const BeginOfRun *> Observer< const BeginOfEvent *> Observer< const EndOfEvent *> Observer< const G4Step *> SimWatcher

Public Member Functions

void beginRun (edm::EventSetup const &) override
 
 CaloSteppingAction (const edm::ParameterSet &p)
 
void produce (edm::Event &, const edm::EventSetup &) override
 
void registerConsumes (edm::ConsumesCollector) override
 
 ~CaloSteppingAction () override
 
- Public Member Functions inherited from SimProducer
const SimProduceroperator= (const SimProducer &)=delete
 
void registerProducts (edm::ProducesCollector producesCollector)
 
 SimProducer ()
 
 SimProducer (const SimProducer &)=delete
 
- Public Member Functions inherited from SimWatcher
bool isMT () const
 
const SimWatcheroperator= (const SimWatcher &)=delete
 
 SimWatcher ()
 
 SimWatcher (const SimWatcher &)=delete
 
virtual ~SimWatcher ()
 
- Public Member Functions inherited from Observer< const BeginOfRun *>
 Observer ()
 
void slotForUpdate (const BeginOfRun * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const BeginOfEvent *>
 Observer ()
 
void slotForUpdate (const BeginOfEvent * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const EndOfEvent *>
 Observer ()
 
void slotForUpdate (const EndOfEvent * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const G4Step *>
 Observer ()
 
void slotForUpdate (const G4Step * iT)
 
virtual ~Observer ()
 

Private Types

typedef std::tuple< const G4LogicalVolume *, uint32_t, int, int, double, double, double, double, double, double, double > PassiveData
 

Private Member Functions

double curve_LY (double crystalLength, double crystalDepth) const
 
void fillHit (uint32_t id, double dE, double time, int primID, uint16_t depth, double em, int flag)
 
void fillHits (edm::PCaloHitContainer &cc, int type)
 
void fillPassiveHits (edm::PassiveHitContainer &cc)
 
double getBirkHC (double dE, double step, double chg, double dens) const
 
double getBirkL3 (double dE, double step, double chg, double dens) const
 
uint16_t getDepth (bool flag, double crystalDepth, double radl) const
 
uint32_t getDetIDHC (int det, int lay, int depth, const math::XYZVectorD &pos) const
 
void NaNTrap (const G4Step *) const
 
void saveHits (int flag)
 
void update (const BeginOfRun *run) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const BeginOfEvent *evt) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const G4Step *step) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const EndOfEvent *evt) override
 This routine will be called when the appropriate signal arrives. More...
 

Private Attributes

const int allSteps_
 
const double birkC1EC_
 
const double birkC1HC_
 
const double birkC2HC_
 
const double birkC3HC_
 
const double birkCutEC_
 
const double birkSlopeEC_
 
int count_
 
std::unique_ptr< EcalBarrelNumberingSchemeebNumberingScheme_
 
std::unique_ptr< EcalEndcapNumberingSchemeeeNumberingScheme_
 
int eventID_
 
std::unique_ptr< HcalNumberingFromPShcNumberingPS_
 
std::unique_ptr< HcalNumberingSchemehcNumberingScheme_
 
std::map< std::pair< int, CaloHitID >, CaloGVHithitMap_ [nSD_]
 
const edm::ParameterSet iC_
 
std::map< const G4LogicalVolume *, std::string > mapLV_
 
const std::vector< std::string > nameEBSD_
 
const std::vector< std::string > nameEESD_
 
const std::vector< std::string > nameHCSD_
 
const std::vector< std::string > nameHitC_
 
std::unique_ptr< CaloSlaveSDslave_ [nSD_]
 
const double slopeLY_
 
std::vector< PassiveDatastore_
 
const double timeSliceUnit_
 
std::vector< const G4LogicalVolume * > volEBSD_
 
std::vector< const G4LogicalVolume * > volEESD_
 
std::vector< const G4LogicalVolume * > volHCSD_
 
std::map< const G4LogicalVolume *, double > xtalMap_
 

Static Private Attributes

static const int nSD_ = 3
 

Additional Inherited Members

- Protected Member Functions inherited from SimProducer
template<class T >
void produces ()
 
template<class T >
void produces (const std::string &instanceName)
 
- Protected Member Functions inherited from SimWatcher
void setMT (bool val)
 

Detailed Description

Definition at line 67 of file CaloSteppingAction.cc.

Member Typedef Documentation

◆ PassiveData

typedef std::tuple<const G4LogicalVolume*, uint32_t, int, int, double, double, double, double, double, double, double> CaloSteppingAction::PassiveData
private

Definition at line 122 of file CaloSteppingAction.cc.

Constructor & Destructor Documentation

◆ CaloSteppingAction()

CaloSteppingAction::CaloSteppingAction ( const edm::ParameterSet p)

Definition at line 126 of file CaloSteppingAction.cc.

References allSteps_, birkC1EC_, birkC1HC_, birkC2HC_, birkC3HC_, birkCutEC_, birkSlopeEC_, ebNumberingScheme_, eeNumberingScheme_, hcNumberingPS_, hcNumberingScheme_, dqmdumpme::k, nameEBSD_, nameEESD_, nameHCSD_, nameHitC_, nSD_, AlCaHLTBitMon_ParallelJobs::p, slave_, slopeLY_, and timeSliceUnit_.

127  : iC_(p.getParameter<edm::ParameterSet>("CaloSteppingAction")),
128  nameEBSD_(iC_.getParameter<std::vector<std::string> >("EBSDNames")),
129  nameEESD_(iC_.getParameter<std::vector<std::string> >("EESDNames")),
130  nameHCSD_(iC_.getParameter<std::vector<std::string> >("HCSDNames")),
131  nameHitC_(iC_.getParameter<std::vector<std::string> >("HitCollNames")),
132  allSteps_(iC_.getParameter<int>("AllSteps")),
133  slopeLY_(iC_.getParameter<double>("SlopeLightYield")),
134  birkC1EC_(iC_.getParameter<double>("BirkC1EC")),
135  birkSlopeEC_(iC_.getParameter<double>("BirkSlopeEC")),
136  birkCutEC_(iC_.getParameter<double>("BirkCutEC")),
137  birkC1HC_(iC_.getParameter<double>("BirkC1HC")),
138  birkC2HC_(iC_.getParameter<double>("BirkC2HC")),
139  birkC3HC_(iC_.getParameter<double>("BirkC3HC")),
140  timeSliceUnit_(iC_.getUntrackedParameter<double>("TimeSliceUnit", 1.0)),
141  count_(0) {
142  edm::ParameterSet iC = p.getParameter<edm::ParameterSet>("CaloSteppingAction");
143 
144  edm::LogVerbatim("Step") << "CaloSteppingAction:: " << nameEBSD_.size() << " names for EB SD's";
145  for (unsigned int k = 0; k < nameEBSD_.size(); ++k)
146  edm::LogVerbatim("Step") << "[" << k << "] " << nameEBSD_[k];
147  edm::LogVerbatim("Step") << "CaloSteppingAction:: " << nameEESD_.size() << " names for EE SD's";
148  for (unsigned int k = 0; k < nameEESD_.size(); ++k)
149  edm::LogVerbatim("Step") << "[" << k << "] " << nameEESD_[k];
150  edm::LogVerbatim("Step") << "CaloSteppingAction:: " << nameHCSD_.size() << " names for HC SD's";
151  for (unsigned int k = 0; k < nameHCSD_.size(); ++k)
152  edm::LogVerbatim("Step") << "[" << k << "] " << nameHCSD_[k];
153  edm::LogVerbatim("Step") << "CaloSteppingAction::Constants for ECAL: slope " << slopeLY_ << " Birk constants "
154  << birkC1EC_ << ":" << birkSlopeEC_ << ":" << birkCutEC_;
155  edm::LogVerbatim("Step") << "CaloSteppingAction::Constants for HCAL: Birk "
156  << "constants " << birkC1HC_ << ":" << birkC2HC_ << ":" << birkC3HC_;
157  edm::LogVerbatim("Step") << "CaloSteppingAction::Constant for time slice " << timeSliceUnit_;
158  edm::LogVerbatim("Step") << "CaloSteppingAction:: " << nameHitC_.size() << " hit collections";
159  for (unsigned int k = 0; k < nameHitC_.size(); ++k)
160  edm::LogVerbatim("Step") << "[" << k << "] " << nameHitC_[k];
161 
162  ebNumberingScheme_ = std::make_unique<EcalBarrelNumberingScheme>();
163  eeNumberingScheme_ = std::make_unique<EcalEndcapNumberingScheme>();
164  hcNumberingPS_ = std::make_unique<HcalNumberingFromPS>(iC);
165  hcNumberingScheme_ = std::make_unique<HcalNumberingScheme>();
166 #ifdef HcalNumberingTest
167  hcNumbering_.reset(nullptr);
168 #endif
169  for (int k = 0; k < CaloSteppingAction::nSD_; ++k) {
170  slave_[k] = std::make_unique<CaloSlaveSD>(nameHitC_[k]);
171  produces<edm::PCaloHitContainer>(nameHitC_[k]);
172  }
173  if (allSteps_ != 0)
174  produces<edm::PassiveHitContainer>("AllPassiveHits");
175  edm::LogVerbatim("Step") << "CaloSteppingAction:: All Steps Flag " << allSteps_ << " for passive hits";
176 }
Log< level::Info, true > LogVerbatim
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
const std::vector< std::string > nameHitC_
const std::vector< std::string > nameEESD_
const std::vector< std::string > nameEBSD_
T getUntrackedParameter(std::string const &, T const &) const
const edm::ParameterSet iC_
std::unique_ptr< CaloSlaveSD > slave_[nSD_]
std::unique_ptr< HcalNumberingFromPS > hcNumberingPS_
std::unique_ptr< HcalNumberingScheme > hcNumberingScheme_
std::unique_ptr< EcalBarrelNumberingScheme > ebNumberingScheme_
std::unique_ptr< EcalEndcapNumberingScheme > eeNumberingScheme_
const std::vector< std::string > nameHCSD_
static const int nSD_

◆ ~CaloSteppingAction()

CaloSteppingAction::~CaloSteppingAction ( )
override

Definition at line 178 of file CaloSteppingAction.cc.

References count_.

178  {
179  edm::LogVerbatim("Step") << "CaloSteppingAction: --------> Total number of "
180  << "selected entries : " << count_;
181 }
Log< level::Info, true > LogVerbatim

Member Function Documentation

◆ beginRun()

void CaloSteppingAction::beginRun ( edm::EventSetup const &  es)
overridevirtual

Reimplemented from SimWatcher.

Definition at line 233 of file CaloSteppingAction.cc.

References edm::EventSetup::getData().

233  {
234  edm::LogVerbatim("Step") << "CaloSteppingAction:: Enter BeginOfRun";
235 
236 #ifdef HcalNumberingTest
237  // Numbering From DDD
238  const HcalDDDSimConstants* hcons_ = &es.getData(ddconsToken_);
239  edm::LogVerbatim("Step") << "CaloSteppingAction:: Initialise "
240  << "HcalNumberingFromDDD";
241  hcNumbering_ = std::make_unique<HcalNumberingFromDDD>(hcons_);
242 #endif
243 }
Log< level::Info, true > LogVerbatim

◆ curve_LY()

double CaloSteppingAction::curve_LY ( double  crystalLength,
double  crystalDepth 
) const
private

Definition at line 487 of file CaloSteppingAction.cc.

References slopeLY_, and mps_merge::weight.

Referenced by update().

487  {
488  double weight = 1.;
489  double dapd = crystalLength - crystalDepth;
490  if (dapd >= -0.1 || dapd <= crystalLength + 0.1) {
491  if (dapd <= 100.)
492  weight = 1.0 + slopeLY_ - dapd * 0.01 * slopeLY_;
493 #ifdef EDM_ML_DEBUG
494  edm::LogVerbatim("Step") << "CaloSteppingAction::curve_LY " << crystalDepth << ":" << crystalLength << ":" << dapd
495  << ":" << weight;
496 #endif
497  } else {
498  edm::LogWarning("Step") << "CaloSteppingAction: light coll curve : wrong "
499  << "distance to APD " << dapd << " crlength = " << crystalLength
500  << " crystal Depth = " << crystalDepth << " weight = " << weight;
501  }
502  return weight;
503 }
Log< level::Info, true > LogVerbatim
Definition: weight.py:1
Log< level::Warning, false > LogWarning

◆ fillHit()

void CaloSteppingAction::fillHit ( uint32_t  id,
double  dE,
double  time,
int  primID,
uint16_t  depth,
double  em,
int  flag 
)
private

Definition at line 460 of file CaloSteppingAction.cc.

References CaloGVHit::addEnergyDeposit(), LEDCalibrationChannels::depth, mps_fire::end, eventID_, RemoveAddSevLevel::flag, hitMap_, CaloGVHit::setEventID(), CaloGVHit::setID(), protons_cff::time, and timeSliceUnit_.

Referenced by update().

460  {
461  CaloHitID currentID(id, time, primID, depth, timeSliceUnit_);
462  double edepEM = (em ? dE : 0);
463  double edepHAD = (em ? 0 : dE);
464  std::pair<int, CaloHitID> evID = std::make_pair(eventID_, currentID);
465  auto it = hitMap_[flag].find(evID);
466  if (it != hitMap_[flag].end()) {
467  (it->second).addEnergyDeposit(edepEM, edepHAD);
468  } else {
469  CaloGVHit aHit;
470  aHit.setEventID(eventID_);
471  aHit.setID(currentID);
472  aHit.addEnergyDeposit(edepEM, edepHAD);
473  hitMap_[flag][evID] = aHit;
474  }
475 }
std::map< std::pair< int, CaloHitID >, CaloGVHit > hitMap_[nSD_]
void setID(uint32_t i, double d, int j, uint16_t k=0)
Definition: CaloGVHit.h:55
void setEventID(int id)
Definition: CaloGVHit.h:46
void addEnergyDeposit(double em, double hd)
Definition: CaloGVHit.cc:29

◆ fillHits()

void CaloSteppingAction::fillHits ( edm::PCaloHitContainer cc,
int  type 
)
private

Definition at line 204 of file CaloSteppingAction.cc.

References slave_.

Referenced by produce().

204  {
205  edm::LogVerbatim("Step") << "CaloSteppingAction::fillHits for type " << type << " with "
206  << slave_[type].get()->hits().size() << " hits";
207  cc = slave_[type].get()->hits();
208  slave_[type].get()->Clean();
209 }
Log< level::Info, true > LogVerbatim
std::unique_ptr< CaloSlaveSD > slave_[nSD_]

◆ fillPassiveHits()

void CaloSteppingAction::fillPassiveHits ( edm::PassiveHitContainer cc)
private

Definition at line 211 of file CaloSteppingAction.cc.

References mapLV_, and store_.

Referenced by produce().

211  {
212  edm::LogVerbatim("Step") << "CaloSteppingAction::fillPassiveHits with " << store_.size() << " hits";
213  for (const auto& element : store_) {
214  auto lv = std::get<0>(element);
215  auto it = mapLV_.find(lv);
216  if (it != mapLV_.end()) {
217  PassiveHit hit(it->second,
218  std::get<1>(element),
219  std::get<5>(element),
220  std::get<6>(element),
221  std::get<4>(element),
222  std::get<2>(element),
223  std::get<3>(element),
224  std::get<7>(element),
225  std::get<8>(element),
226  std::get<9>(element),
227  std::get<10>(element));
228  cc.emplace_back(hit);
229  }
230  }
231 }
Log< level::Info, true > LogVerbatim
std::map< const G4LogicalVolume *, std::string > mapLV_
std::vector< PassiveData > store_

◆ getBirkHC()

double CaloSteppingAction::getBirkHC ( double  dE,
double  step,
double  chg,
double  dens 
) const
private

Definition at line 525 of file CaloSteppingAction.cc.

References funct::abs(), birkC1HC_, birkC2HC_, birkC3HC_, c, ALCARECOTkAlJpsiMuMu_cff::charge, fastSimProducer_cff::density, and mps_merge::weight.

Referenced by update().

525  {
526  double weight = 1.;
527  if (charge != 0. && step > 0.) {
528  double dedx = dEStep / step;
529  double rkb = birkC1HC_ / density;
530  double c = birkC2HC_ * rkb * rkb;
531  if (std::abs(charge) >= 2.)
532  rkb /= birkC3HC_;
533  weight = 1. / (1. + rkb * dedx + c * dedx * dedx);
534 #ifdef EDM_ML_DEBUG
535  edm::LogVerbatim("Step") << "CaloSteppingAction::getBirkHC Charge " << charge << " dE/dx " << dedx << " Birk Const "
536  << rkb << ", " << c << " Weight = " << weight << " dE " << dEStep;
537 #endif
538  }
539  return weight;
540 }
Log< level::Info, true > LogVerbatim
Definition: weight.py:1
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
step
Definition: StallMonitor.cc:98

◆ getBirkL3()

double CaloSteppingAction::getBirkL3 ( double  dE,
double  step,
double  chg,
double  dens 
) const
private

Definition at line 505 of file CaloSteppingAction.cc.

References birkC1EC_, birkCutEC_, birkSlopeEC_, ALCARECOTkAlJpsiMuMu_cff::charge, fastSimProducer_cff::density, dqm-mbProfile::log, and mps_merge::weight.

Referenced by update().

505  {
506  double weight = 1.;
507  if (charge != 0. && step > 0.) {
508  double dedx = dEStep / step;
509  double rkb = birkC1EC_ / density;
510  if (dedx > 0) {
511  weight = 1. - birkSlopeEC_ * log(rkb * dedx);
512  if (weight < birkCutEC_)
513  weight = birkCutEC_;
514  else if (weight > 1.)
515  weight = 1.;
516  }
517 #ifdef EDM_ML_DEBUG
518  edm::LogVerbatim("Step") << "CaloSteppingAction::getBirkL3 Charge " << charge << " dE/dx " << dedx << " Birk Const "
519  << rkb << " Weight = " << weight << " dE " << dEStep << " step " << step;
520 #endif
521  }
522  return weight;
523 }
Log< level::Info, true > LogVerbatim
Definition: weight.py:1
step
Definition: StallMonitor.cc:98

◆ getDepth()

uint16_t CaloSteppingAction::getDepth ( bool  flag,
double  crystalDepth,
double  radl 
) const
private

Definition at line 477 of file CaloSteppingAction.cc.

References LEDCalibrationChannels::depth, RemoveAddSevLevel::flag, PCaloHit::kEcalDepthMask, PCaloHit::kEcalDepthOffset, and PCaloHit::kEcalDepthRefz.

Referenced by update().

477  {
478  uint16_t depth1 = (flag ? 0 : PCaloHit::kEcalDepthRefz);
479  uint16_t depth2 = (uint16_t)floor(crystalDepth / radl);
480  uint16_t depth = (((depth2 & PCaloHit::kEcalDepthMask) << PCaloHit::kEcalDepthOffset) | depth1);
481 #ifdef EDM_ML_DEBUG
482  edm::LogVerbatim("Step") << "CaloSteppingAction::getDepth radl " << radl << ":" << crystalDepth << " depth " << depth;
483 #endif
484  return depth;
485 }
Log< level::Info, true > LogVerbatim
static const int kEcalDepthRefz
Definition: PCaloHit.h:63
static const int kEcalDepthMask
Definition: PCaloHit.h:61
static const int kEcalDepthOffset
Definition: PCaloHit.h:62

◆ getDetIDHC()

uint32_t CaloSteppingAction::getDetIDHC ( int  det,
int  lay,
int  depth,
const math::XYZVectorD pos 
) const
private

Definition at line 448 of file CaloSteppingAction.cc.

References LEDCalibrationChannels::depth, relativeConstraints::error, hcNumberingPS_, hcNumberingScheme_, AlCaHLTBitMon_QueryRunRegistry::string, and createJobs::tmp.

Referenced by update().

448  {
449  HcalNumberingFromDDD::HcalID tmp = hcNumberingPS_.get()->unitID(det, lay, depth, pos);
450 #ifdef HcalNumberingTest
451  auto id0 = hcNumberingScheme_.get()->getUnitID(tmp);
452  HcalNumberingFromDDD::HcalID tmpO = hcNumbering_.get()->unitID(det, pos, depth, lay);
453  auto idO = hcNumberingScheme_.get()->getUnitID(tmpO);
454  std::string error = (id0 == idO) ? " ** OK **" : " ** ERROR **";
455  edm::LogVerbatim("Step") << "Det ID " << HcalDetId(id0) << " Original " << HcalDetId(idO) << error;
456 #endif
457  return (hcNumberingScheme_.get()->getUnitID(tmp));
458 }
Log< level::Info, true > LogVerbatim
std::unique_ptr< HcalNumberingFromPS > hcNumberingPS_
std::unique_ptr< HcalNumberingScheme > hcNumberingScheme_
tmp
align.sh
Definition: createJobs.py:716

◆ NaNTrap()

void CaloSteppingAction::NaNTrap ( const G4Step *  aStep) const
private

Definition at line 434 of file CaloSteppingAction.cc.

References Exception, and edm::isNotFinite().

Referenced by update().

434  {
435  auto currentPos = aStep->GetTrack()->GetPosition();
436  double xyz = currentPos.x() + currentPos.y() + currentPos.z();
437  auto currentMom = aStep->GetTrack()->GetMomentum();
438  xyz += currentMom.x() + currentMom.y() + currentMom.z();
439 
440  if (edm::isNotFinite(xyz)) {
441  auto pCurrentVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
442  auto& nameOfVol = pCurrentVol->GetName();
443  throw cms::Exception("Unknown", "CaloSteppingAction")
444  << " Corrupted Event - NaN detected in volume " << nameOfVol << "\n";
445  }
446 }
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9

◆ produce()

void CaloSteppingAction::produce ( edm::Event e,
const edm::EventSetup  
)
overridevirtual

Implements SimProducer.

Definition at line 190 of file CaloSteppingAction.cc.

References allSteps_, MillePedeFileConverter_cfg::e, fillHits(), fillPassiveHits(), dqmdumpme::k, eostools::move(), nameHitC_, nSD_, and saveHits().

190  {
191  for (int k = 0; k < CaloSteppingAction::nSD_; ++k) {
192  saveHits(k);
193  auto product = std::make_unique<edm::PCaloHitContainer>();
194  fillHits(*product, k);
195  e.put(std::move(product), nameHitC_[k]);
196  }
197  if (allSteps_ != 0) {
198  std::unique_ptr<edm::PassiveHitContainer> hgcPH(new edm::PassiveHitContainer);
199  fillPassiveHits(*hgcPH);
200  e.put(std::move(hgcPH), "AllPassiveHits");
201  }
202 }
const std::vector< std::string > nameHitC_
void fillHits(edm::PCaloHitContainer &cc, int type)
void fillPassiveHits(edm::PassiveHitContainer &cc)
std::vector< PassiveHit > PassiveHitContainer
Definition: PassiveHit.h:100
static const int nSD_
def move(src, dest)
Definition: eostools.py:511

◆ registerConsumes()

void CaloSteppingAction::registerConsumes ( edm::ConsumesCollector  cc)
overridevirtual

Reimplemented from SimWatcher.

Definition at line 183 of file CaloSteppingAction.cc.

References edm::BeginRun, and edm::ConsumesCollector::esConsumes().

183  {
184 #ifdef HcalNumberingTest
186  edm::LogVerbatim("Step") << "CaloSteppingAction::Initialize ESGetToken for HcalDDDSimConstants";
187 #endif
188 }
Log< level::Info, true > LogVerbatim

◆ saveHits()

void CaloSteppingAction::saveHits ( int  flag)
private

Definition at line 542 of file CaloSteppingAction.cc.

References hitMap_, findQualityFiles::size, and slave_.

Referenced by produce().

542  {
543  edm::LogVerbatim("Step") << "CaloSteppingAction:: saveHits for type " << type << " with " << hitMap_[type].size()
544  << " hits";
545  slave_[type].get()->ReserveMemory(hitMap_[type].size());
546  for (auto const& hit : hitMap_[type]) {
547  slave_[type].get()->processHits(hit.second.getUnitID(),
548  0.001 * hit.second.getEM(),
549  0.001 * hit.second.getHadr(),
550  hit.second.getTimeSlice(),
551  hit.second.getTrackID(),
552  hit.second.getDepth());
553  }
554 }
size
Write out results.
std::map< std::pair< int, CaloHitID >, CaloGVHit > hitMap_[nSD_]
Log< level::Info, true > LogVerbatim
std::unique_ptr< CaloSlaveSD > slave_[nSD_]

◆ update() [1/4]

void CaloSteppingAction::update ( const BeginOfRun )
overrideprivatevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfRun *>.

Definition at line 246 of file CaloSteppingAction.cc.

References allSteps_, PVValHelper::dz, dqmdumpme::k, mapLV_, Skims_PA_cff::name, nameEBSD_, nameEESD_, nameHCSD_, AlCaHLTBitMon_QueryRunRegistry::string, volEBSD_, volEESD_, volHCSD_, and xtalMap_.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), progressbar.ProgressBar::finish(), and MatrixUtil.Steps::overwrite().

246  {
247  int irun = (*run)()->GetRunID();
248  edm::LogVerbatim("Step") << "CaloSteppingAction:: Begin of Run = " << irun;
249 
250  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
251  if (lvs) {
252  std::map<const std::string, const G4LogicalVolume*> nameMap;
253  std::map<const std::string, const G4LogicalVolume*>::const_iterator itr;
254  for (auto lvi = lvs->begin(), lve = lvs->end(); lvi != lve; ++lvi) {
255  nameMap.emplace((*lvi)->GetName(), *lvi);
256  if (allSteps_ < 0)
257  mapLV_[*lvi] = (*lvi)->GetName();
258  }
259 
260  for (auto const& name : nameEBSD_) {
261  for (itr = nameMap.begin(); itr != nameMap.end(); ++itr) {
262  const std::string& lvname = itr->first;
263  if (lvname.find(name) != std::string::npos) {
264  volEBSD_.emplace_back(itr->second);
265  int type = (lvname.find("refl") == std::string::npos) ? -1 : 1;
266  G4Trap* solid = static_cast<G4Trap*>(itr->second->GetSolid());
267  double dz = 2 * solid->GetZHalfLength() / CLHEP::mm;
268  xtalMap_.insert(std::pair<const G4LogicalVolume*, double>(itr->second, dz * type));
269  if ((allSteps_ > 0) && ((allSteps_ % 10) > 0))
270  mapLV_[itr->second] = itr->first;
271  }
272  }
273  }
274  for (auto const& name : nameEESD_) {
275  for (itr = nameMap.begin(); itr != nameMap.end(); ++itr) {
276  const std::string& lvname = itr->first;
277  if (lvname.find(name) != std::string::npos) {
278  volEESD_.emplace_back(itr->second);
279  int type = (lvname.find("refl") == std::string::npos) ? 1 : -1;
280  G4Trap* solid = static_cast<G4Trap*>(itr->second->GetSolid());
281  double dz = 2 * solid->GetZHalfLength() / CLHEP::mm;
282  xtalMap_.insert(std::pair<const G4LogicalVolume*, double>(itr->second, dz * type));
283  if ((allSteps_ > 0) && (((allSteps_ / 10) % 10) > 0))
284  mapLV_[itr->second] = itr->first;
285  }
286  }
287  }
288 
289  for (auto const& name : nameHCSD_) {
290  for (itr = nameMap.begin(); itr != nameMap.end(); ++itr) {
291  const std::string& lvname = itr->first;
292  if (lvname.find(name) != std::string::npos) {
293  volHCSD_.emplace_back(itr->second);
294  if ((allSteps_ > 0) && (((allSteps_ / 100) % 10) > 0))
295  mapLV_[itr->second] = itr->first;
296  }
297  }
298  }
299  }
300 #ifdef EDM_ML_DEBUG
301  edm::LogVerbatim("Step") << volEBSD_.size() << " logical volumes for EB SD";
302  for (unsigned int k = 0; k < volEBSD_.size(); ++k)
303  edm::LogVerbatim("Step") << "[" << k << "] " << volEBSD_[k];
304  edm::LogVerbatim("Step") << volEESD_.size() << " logical volumes for EE SD";
305  for (unsigned int k = 0; k < volEESD_.size(); ++k)
306  edm::LogVerbatim("Step") << "[" << k << "] " << volEESD_[k];
307  edm::LogVerbatim("Step") << volHCSD_.size() << " logical volumes for HC SD";
308  for (unsigned int k = 0; k < volHCSD_.size(); ++k)
309  edm::LogVerbatim("Step") << "[" << k << "] " << volHCSD_[k];
310  edm::LogVerbatim("Step") << mapLV_.size() << " logical volumes for Passive hits";
311  unsigned int k(0);
312  for (auto itr = mapLV_.begin(); itr != mapLV_.end(); ++itr) {
313  edm::LogVerbatim("Step") << "[" << k << "] " << itr->second << ":" << itr->first;
314  ++k;
315  }
316 #endif
317 }
Log< level::Info, true > LogVerbatim
std::vector< const G4LogicalVolume * > volEBSD_
const std::vector< std::string > nameEESD_
std::map< const G4LogicalVolume *, double > xtalMap_
std::vector< const G4LogicalVolume * > volHCSD_
const std::vector< std::string > nameEBSD_
std::map< const G4LogicalVolume *, std::string > mapLV_
std::vector< const G4LogicalVolume * > volEESD_
const std::vector< std::string > nameHCSD_

◆ update() [2/4]

void CaloSteppingAction::update ( const BeginOfEvent )
overrideprivatevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfEvent *>.

Definition at line 320 of file CaloSteppingAction.cc.

References allSteps_, mps_fire::end, eventID_, hitMap_, dqmdumpme::k, nSD_, slave_, and store_.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), progressbar.ProgressBar::finish(), and MatrixUtil.Steps::overwrite().

320  {
321  eventID_ = (*evt)()->GetEventID();
322  edm::LogVerbatim("Step") << "CaloSteppingAction: Begin of event = " << eventID_;
323  for (int k = 0; k < CaloSteppingAction::nSD_; ++k) {
324  hitMap_[k].erase(hitMap_[k].begin(), hitMap_[k].end());
325  slave_[k].get()->Initialize();
326  }
327  if (allSteps_ != 0)
328  store_.clear();
329 }
std::map< std::pair< int, CaloHitID >, CaloGVHit > hitMap_[nSD_]
Log< level::Info, true > LogVerbatim
std::unique_ptr< CaloSlaveSD > slave_[nSD_]
std::vector< PassiveData > store_
static const int nSD_

◆ update() [3/4]

void CaloSteppingAction::update ( const G4Step *  )
overrideprivatevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const G4Step *>.

Definition at line 332 of file CaloSteppingAction.cc.

References funct::abs(), EcalBaseNumber::addLevel(), allSteps_, filterCSVwithJSON::copy, curve_LY(), LEDCalibrationChannels::depth, ebNumberingScheme_, eeNumberingScheme_, HCALHighEnergyHPDFilter_cfi::energy, fillHit(), spr::find(), RemoveAddSevLevel::flag, g, getBirkHC(), getBirkL3(), EcalBaseNumber::getCapacity(), getDepth(), getDetIDHC(), cuy::ii, G4TrackToParticleID::isGammaElectronPositron(), crabWrapper::key, DTRecHitClients_cfi::local, mapLV_, NaNTrap(), EcalBaseNumber::setSize(), findQualityFiles::size, store_, protons_cff::time, funct::true, volEBSD_, volEESD_, volHCSD_, and xtalMap_.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), progressbar.ProgressBar::finish(), and MatrixUtil.Steps::overwrite().

332  {
333  // edm::LogVerbatim("Step") <<"CaloSteppingAction: At each Step";
334  NaNTrap(aStep);
335  auto lv = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
336  bool hc = (std::find(volHCSD_.begin(), volHCSD_.end(), lv) != volHCSD_.end());
337  bool eb = (std::find(volEBSD_.begin(), volEBSD_.end(), lv) != volEBSD_.end());
338  bool ee = (std::find(volEESD_.begin(), volEESD_.end(), lv) != volEESD_.end());
339  uint32_t unitID(0);
340  if (hc || eb || ee) {
341  double dEStep = aStep->GetTotalEnergyDeposit() / CLHEP::MeV;
342  auto const theTrack = aStep->GetTrack();
343  double time = theTrack->GetGlobalTime() / CLHEP::nanosecond;
344  int primID = theTrack->GetTrackID();
346  auto const touch = aStep->GetPreStepPoint()->GetTouchable();
347  auto const& hitPoint = aStep->GetPreStepPoint()->GetPosition();
348  if (hc) {
349  int depth = (touch->GetReplicaNumber(0)) % 10 + 1;
350  int lay = (touch->GetReplicaNumber(0) / 10) % 100 + 1;
351  int det = (touch->GetReplicaNumber(1)) / 1000;
352  unitID = getDetIDHC(det, lay, depth, math::XYZVectorD(hitPoint.x(), hitPoint.y(), hitPoint.z()));
353  if (unitID > 0 && dEStep > 0.0) {
354  dEStep *= getBirkHC(dEStep,
355  (aStep->GetStepLength() / CLHEP::cm),
356  aStep->GetPreStepPoint()->GetCharge(),
357  (aStep->GetPreStepPoint()->GetMaterial()->GetDensity() / (CLHEP::g / CLHEP::cm3)));
358  fillHit(unitID, dEStep, time, primID, 0, em, 2);
359  }
360  } else {
361  EcalBaseNumber theBaseNumber;
362  int size = touch->GetHistoryDepth() + 1;
363  if (theBaseNumber.getCapacity() < size)
364  theBaseNumber.setSize(size);
365  //Get name and copy numbers
366  if (size > 1) {
367  for (int ii = 0; ii < size; ii++) {
368  theBaseNumber.addLevel(touch->GetVolume(ii)->GetName(), touch->GetReplicaNumber(ii));
369  }
370  }
371  unitID = (eb ? (ebNumberingScheme_->getUnitID(theBaseNumber)) : (eeNumberingScheme_->getUnitID(theBaseNumber)));
372  if (unitID > 0 && dEStep > 0.0) {
373  auto local = touch->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
374  auto ite = xtalMap_.find(lv);
375  double crystalLength = ((ite == xtalMap_.end()) ? 230.0 : std::abs(ite->second));
376  double crystalDepth =
377  ((ite == xtalMap_.end()) ? 0.0 : (std::abs(0.5 * (ite->second) + (local.z() / CLHEP::mm))));
378  double radl = aStep->GetPreStepPoint()->GetMaterial()->GetRadlen() / CLHEP::mm;
379  bool flag = ((ite == xtalMap_.end()) ? true : (((ite->second) >= 0) ? true : false));
380  auto depth = getDepth(flag, crystalDepth, radl);
381  dEStep *= (getBirkL3(dEStep,
382  (aStep->GetStepLength() / CLHEP::cm),
383  aStep->GetPreStepPoint()->GetCharge(),
384  (aStep->GetPreStepPoint()->GetMaterial()->GetDensity() / (CLHEP::g / CLHEP::cm3))) *
385  curve_LY(crystalLength, crystalDepth));
386  fillHit(unitID, dEStep, time, primID, depth, em, (eb ? 0 : 1));
387  }
388  }
389  }
390 
391  if (allSteps_ != 0) {
392  auto it = mapLV_.find(lv);
393  if (it != mapLV_.end()) {
394  double energy = aStep->GetTotalEnergyDeposit() / CLHEP::MeV;
395  auto const touch = aStep->GetPreStepPoint()->GetTouchable();
396  double time = aStep->GetTrack()->GetGlobalTime() / CLHEP::nanosecond;
397  int trackId = aStep->GetTrack()->GetTrackID();
398  int pdg = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
399  double stepl = (aStep->GetStepLength() / CLHEP::cm);
400  double xp = aStep->GetPreStepPoint()->GetPosition().x() / CLHEP::cm;
401  double yp = aStep->GetPreStepPoint()->GetPosition().y() / CLHEP::cm;
402  double zp = aStep->GetPreStepPoint()->GetPosition().z() / CLHEP::cm;
403 #ifdef EDM_ML_DEBUG
404  edm::LogVerbatim("Step") << "CaloSteppingAction: Volume " << lv->GetName() << " History "
405  << touch->GetHistoryDepth() << " Pointers " << aStep->GetPostStepPoint() << ":"
406  << aStep->GetTrack()->GetNextVolume() << ":" << aStep->IsLastStepInVolume() << " E "
407  << energy << " T " << time << " PDG " << pdg << " step " << stepl << " Position (" << xp
408  << ", " << yp << ", " << zp << ")";
409 #endif
410  uint32_t copy = (allSteps_ < 0) ? 0 : unitID;
411  if (((aStep->GetPostStepPoint() == nullptr) || (aStep->GetTrack()->GetNextVolume() == nullptr)) &&
412  (aStep->IsLastStepInVolume())) {
413  energy += (aStep->GetPreStepPoint()->GetKineticEnergy() / CLHEP::MeV);
414  } else {
415  time = aStep->GetPostStepPoint()->GetGlobalTime() / CLHEP::nanosecond;
416  if (copy == 0)
417  copy = (touch->GetHistoryDepth() < 1)
418  ? static_cast<uint32_t>(touch->GetReplicaNumber(0))
419  : static_cast<uint32_t>(touch->GetReplicaNumber(0) + 1000 * touch->GetReplicaNumber(1));
420  }
421  PassiveData key(std::make_tuple(lv, copy, trackId, pdg, time, energy, energy, stepl, xp, yp, zp));
422  store_.push_back(key);
423  }
424  }
425 }
size
Write out results.
Log< level::Info, true > LogVerbatim
std::vector< const G4LogicalVolume * > volEBSD_
std::map< const G4LogicalVolume *, double > xtalMap_
std::vector< const G4LogicalVolume * > volHCSD_
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
uint16_t getDepth(bool flag, double crystalDepth, double radl) const
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
void NaNTrap(const G4Step *) const
void addLevel(const std::string &name, const int &copyNumber)
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::map< const G4LogicalVolume *, std::string > mapLV_
uint32_t getDetIDHC(int det, int lay, int depth, const math::XYZVectorD &pos) const
ii
Definition: cuy.py:589
void fillHit(uint32_t id, double dE, double time, int primID, uint16_t depth, double em, int flag)
double getBirkHC(double dE, double step, double chg, double dens) const
std::vector< const G4LogicalVolume * > volEESD_
std::unique_ptr< EcalBarrelNumberingScheme > ebNumberingScheme_
std::tuple< const G4LogicalVolume *, uint32_t, int, int, double, double, double, double, double, double, double > PassiveData
static bool isGammaElectronPositron(int pdgCode)
std::unique_ptr< EcalEndcapNumberingScheme > eeNumberingScheme_
double getBirkL3(double dE, double step, double chg, double dens) const
std::vector< PassiveData > store_
void setSize(const int &size)
double curve_LY(double crystalLength, double crystalDepth) const

◆ update() [4/4]

void CaloSteppingAction::update ( const EndOfEvent )
overrideprivatevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const EndOfEvent *>.

Definition at line 428 of file CaloSteppingAction.cc.

References count_.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), progressbar.ProgressBar::finish(), and MatrixUtil.Steps::overwrite().

428  {
429  ++count_;
430  // Fill event input
431  edm::LogVerbatim("Step") << "CaloSteppingAction: EndOfEvent " << (*evt)()->GetEventID();
432 }
Log< level::Info, true > LogVerbatim

Member Data Documentation

◆ allSteps_

const int CaloSteppingAction::allSteps_
private

Definition at line 115 of file CaloSteppingAction.cc.

Referenced by CaloSteppingAction(), produce(), and update().

◆ birkC1EC_

const double CaloSteppingAction::birkC1EC_
private

Definition at line 116 of file CaloSteppingAction.cc.

Referenced by CaloSteppingAction(), and getBirkL3().

◆ birkC1HC_

const double CaloSteppingAction::birkC1HC_
private

Definition at line 117 of file CaloSteppingAction.cc.

Referenced by CaloSteppingAction(), and getBirkHC().

◆ birkC2HC_

const double CaloSteppingAction::birkC2HC_
private

Definition at line 117 of file CaloSteppingAction.cc.

Referenced by CaloSteppingAction(), and getBirkHC().

◆ birkC3HC_

const double CaloSteppingAction::birkC3HC_
private

Definition at line 118 of file CaloSteppingAction.cc.

Referenced by CaloSteppingAction(), and getBirkHC().

◆ birkCutEC_

const double CaloSteppingAction::birkCutEC_
private

Definition at line 117 of file CaloSteppingAction.cc.

Referenced by CaloSteppingAction(), and getBirkL3().

◆ birkSlopeEC_

const double CaloSteppingAction::birkSlopeEC_
private

Definition at line 116 of file CaloSteppingAction.cc.

Referenced by CaloSteppingAction(), and getBirkL3().

◆ count_

int CaloSteppingAction::count_
private

Definition at line 119 of file CaloSteppingAction.cc.

Referenced by update(), and ~CaloSteppingAction().

◆ ebNumberingScheme_

std::unique_ptr<EcalBarrelNumberingScheme> CaloSteppingAction::ebNumberingScheme_
private

Definition at line 99 of file CaloSteppingAction.cc.

Referenced by CaloSteppingAction(), and update().

◆ eeNumberingScheme_

std::unique_ptr<EcalEndcapNumberingScheme> CaloSteppingAction::eeNumberingScheme_
private

Definition at line 100 of file CaloSteppingAction.cc.

Referenced by CaloSteppingAction(), and update().

◆ eventID_

int CaloSteppingAction::eventID_
private

Definition at line 119 of file CaloSteppingAction.cc.

Referenced by fillHit(), and update().

◆ hcNumberingPS_

std::unique_ptr<HcalNumberingFromPS> CaloSteppingAction::hcNumberingPS_
private

Definition at line 101 of file CaloSteppingAction.cc.

Referenced by CaloSteppingAction(), and getDetIDHC().

◆ hcNumberingScheme_

std::unique_ptr<HcalNumberingScheme> CaloSteppingAction::hcNumberingScheme_
private

Definition at line 106 of file CaloSteppingAction.cc.

Referenced by CaloSteppingAction(), and getDetIDHC().

◆ hitMap_

std::map<std::pair<int, CaloHitID>, CaloGVHit> CaloSteppingAction::hitMap_[nSD_]
private

Definition at line 120 of file CaloSteppingAction.cc.

Referenced by fillHit(), saveHits(), and update().

◆ iC_

const edm::ParameterSet CaloSteppingAction::iC_
private

Definition at line 109 of file CaloSteppingAction.cc.

◆ mapLV_

std::map<const G4LogicalVolume*, std::string> CaloSteppingAction::mapLV_
private

Definition at line 114 of file CaloSteppingAction.cc.

Referenced by fillPassiveHits(), and update().

◆ nameEBSD_

const std::vector<std::string> CaloSteppingAction::nameEBSD_
private

Definition at line 110 of file CaloSteppingAction.cc.

Referenced by CaloSteppingAction(), and update().

◆ nameEESD_

const std::vector<std::string> CaloSteppingAction::nameEESD_
private

Definition at line 110 of file CaloSteppingAction.cc.

Referenced by CaloSteppingAction(), and update().

◆ nameHCSD_

const std::vector<std::string> CaloSteppingAction::nameHCSD_
private

Definition at line 110 of file CaloSteppingAction.cc.

Referenced by CaloSteppingAction(), and update().

◆ nameHitC_

const std::vector<std::string> CaloSteppingAction::nameHitC_
private

Definition at line 111 of file CaloSteppingAction.cc.

Referenced by CaloSteppingAction(), and produce().

◆ nSD_

const int CaloSteppingAction::nSD_ = 3
staticprivate

Definition at line 98 of file CaloSteppingAction.cc.

Referenced by CaloSteppingAction(), produce(), and update().

◆ slave_

std::unique_ptr<CaloSlaveSD> CaloSteppingAction::slave_[nSD_]
private

Definition at line 107 of file CaloSteppingAction.cc.

Referenced by CaloSteppingAction(), fillHits(), saveHits(), and update().

◆ slopeLY_

const double CaloSteppingAction::slopeLY_
private

Definition at line 116 of file CaloSteppingAction.cc.

Referenced by CaloSteppingAction(), and curve_LY().

◆ store_

std::vector<PassiveData> CaloSteppingAction::store_
private

Definition at line 123 of file CaloSteppingAction.cc.

Referenced by fillPassiveHits(), and update().

◆ timeSliceUnit_

const double CaloSteppingAction::timeSliceUnit_
private

Definition at line 118 of file CaloSteppingAction.cc.

Referenced by CaloSteppingAction(), and fillHit().

◆ volEBSD_

std::vector<const G4LogicalVolume*> CaloSteppingAction::volEBSD_
private

Definition at line 112 of file CaloSteppingAction.cc.

Referenced by update().

◆ volEESD_

std::vector<const G4LogicalVolume*> CaloSteppingAction::volEESD_
private

Definition at line 112 of file CaloSteppingAction.cc.

Referenced by update().

◆ volHCSD_

std::vector<const G4LogicalVolume*> CaloSteppingAction::volHCSD_
private

Definition at line 112 of file CaloSteppingAction.cc.

Referenced by update().

◆ xtalMap_

std::map<const G4LogicalVolume*, double> CaloSteppingAction::xtalMap_
private

Definition at line 113 of file CaloSteppingAction.cc.

Referenced by update().