CMS 3D CMS Logo

CaloSteppingAction.cc
Go to the documentation of this file.
1 // to make hits in EB/EE/HC
4 
7 #ifdef HcalNumberingTest
10 #endif
11 
17 
18 #include "G4LogicalVolumeStore.hh"
19 #include "G4NavigationHistory.hh"
20 #include "G4ParticleTable.hh"
21 #include "G4PhysicalVolumeStore.hh"
22 #include "G4RegionStore.hh"
23 #include "G4Trap.hh"
24 #include "G4UnitsTable.hh"
25 #include "G4SystemOfUnits.hh"
26 
27 #include <cmath>
28 #include <iostream>
29 #include <iomanip>
30 
31 //#define EDM_ML_DEBUG
32 
34  edm::ParameterSet iC = p.getParameter<edm::ParameterSet>("CaloSteppingAction");
35  nameEBSD_ = iC.getParameter<std::vector<std::string> >("EBSDNames");
36  nameEESD_ = iC.getParameter<std::vector<std::string> >("EESDNames");
37  nameHCSD_ = iC.getParameter<std::vector<std::string> >("HCSDNames");
38  nameHitC_ = iC.getParameter<std::vector<std::string> >("HitCollNames");
39  allSteps_ = iC.getParameter<int>("AllSteps");
40  slopeLY_ = iC.getParameter<double>("SlopeLightYield");
41  birkC1EC_ = iC.getParameter<double>("BirkC1EC");
42  birkSlopeEC_ = iC.getParameter<double>("BirkSlopeEC");
43  birkCutEC_ = iC.getParameter<double>("BirkCutEC");
44  birkC1HC_ = iC.getParameter<double>("BirkC1HC");
45  birkC2HC_ = iC.getParameter<double>("BirkC2HC");
46  birkC3HC_ = iC.getParameter<double>("BirkC3HC");
47  timeSliceUnit_ = iC.getUntrackedParameter<double>("TimeSliceUnit", 1.0);
48 
49  edm::LogVerbatim("Step") << "CaloSteppingAction:: " << nameEBSD_.size() << " names for EB SD's";
50  for (unsigned int k = 0; k < nameEBSD_.size(); ++k)
51  edm::LogVerbatim("Step") << "[" << k << "] " << nameEBSD_[k];
52  edm::LogVerbatim("Step") << "CaloSteppingAction:: " << nameEESD_.size() << " names for EE SD's";
53  for (unsigned int k = 0; k < nameEESD_.size(); ++k)
54  edm::LogVerbatim("Step") << "[" << k << "] " << nameEESD_[k];
55  edm::LogVerbatim("Step") << "CaloSteppingAction:: " << nameHCSD_.size() << " names for HC SD's";
56  for (unsigned int k = 0; k < nameHCSD_.size(); ++k)
57  edm::LogVerbatim("Step") << "[" << k << "] " << nameHCSD_[k];
58  edm::LogVerbatim("Step") << "CaloSteppingAction::Constants for ECAL: slope " << slopeLY_ << " Birk constants "
59  << birkC1EC_ << ":" << birkSlopeEC_ << ":" << birkCutEC_;
60  edm::LogVerbatim("Step") << "CaloSteppingAction::Constants for HCAL: Birk "
61  << "constants " << birkC1HC_ << ":" << birkC2HC_ << ":" << birkC3HC_;
62  edm::LogVerbatim("Step") << "CaloSteppingAction::Constant for time slice " << timeSliceUnit_;
63  edm::LogVerbatim("Step") << "CaloSteppingAction:: " << nameHitC_.size() << " hit collections";
64  for (unsigned int k = 0; k < nameHitC_.size(); ++k)
65  edm::LogVerbatim("Step") << "[" << k << "] " << nameHitC_[k];
66 
67  ebNumberingScheme_ = std::make_unique<EcalBarrelNumberingScheme>();
68  eeNumberingScheme_ = std::make_unique<EcalEndcapNumberingScheme>();
69  hcNumberingPS_ = std::make_unique<HcalNumberingFromPS>(iC);
70  hcNumberingScheme_ = std::make_unique<HcalNumberingScheme>();
71 #ifdef HcalNumberingTest
72  hcNumbering_.reset(nullptr);
73 #endif
74  for (int k = 0; k < CaloSteppingAction::nSD_; ++k) {
75  slave_[k] = std::make_unique<CaloSlaveSD>(nameHitC_[k]);
76  produces<edm::PCaloHitContainer>(nameHitC_[k]);
77  }
78  if (allSteps_ > 0)
79  produces<edm::PassiveHitContainer>("AllPassiveHits");
80  edm::LogVerbatim("Step") << "CaloSteppingAction:: All Steps Flag " << allSteps_ << " for passive hits";
81 }
82 
84  edm::LogVerbatim("Step") << "CaloSteppingAction: --------> Total number of "
85  << "selected entries : " << count_;
86 }
87 
89  for (int k = 0; k < CaloSteppingAction::nSD_; ++k) {
90  saveHits(k);
91  auto product = std::make_unique<edm::PCaloHitContainer>();
92  fillHits(*product, k);
93  e.put(std::move(product), nameHitC_[k]);
94  }
95  if (allSteps_ > 0) {
96  std::unique_ptr<edm::PassiveHitContainer> hgcPH(new edm::PassiveHitContainer);
97  fillPassiveHits(*hgcPH);
98  e.put(std::move(hgcPH), "AllPassiveHits");
99  }
100 }
101 
103  edm::LogVerbatim("Step") << "CaloSteppingAction::fillHits for type " << type << " with "
104  << slave_[type].get()->hits().size() << " hits";
105  cc = slave_[type].get()->hits();
106  slave_[type].get()->Clean();
107 }
108 
110  edm::LogVerbatim("Step") << "CaloSteppingAction::fillPassiveHits with " << store_.size() << " hits";
111  for (const auto& element : store_) {
112  auto lv = std::get<0>(element);
113  auto it = mapLV_.find(lv);
114  if (it != mapLV_.end()) {
115  PassiveHit hit(it->second,
116  std::get<1>(element),
117  std::get<5>(element),
118  std::get<6>(element),
119  std::get<4>(element),
120  std::get<2>(element),
121  std::get<3>(element),
122  std::get<7>(element));
123  cc.emplace_back(hit);
124  }
125  }
126 }
127 
129  edm::LogVerbatim("Step") << "CaloSteppingAction:: Enter BeginOfJob";
130 
131 #ifdef HcalNumberingTest
132  // Numbering From DDD
134  (*job)()->get<HcalSimNumberingRecord>().get(hdc);
135  const HcalDDDSimConstants* hcons_ = hdc.product();
136  edm::LogVerbatim("Step") << "CaloSteppingAction:: Initialise "
137  << "HcalNumberingFromDDD";
138  hcNumbering_ = std::make_unique<HcalNumberingFromDDD>(hcons_);
139 #endif
140 }
141 
142 //==================================================================== per RUN
144  int irun = (*run)()->GetRunID();
145  edm::LogVerbatim("Step") << "CaloSteppingAction:: Begin of Run = " << irun;
146 
147  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
148  if (lvs) {
149  std::map<const std::string, const G4LogicalVolume*> nameMap;
150  std::map<const std::string, const G4LogicalVolume*>::const_iterator itr;
151  for (auto lvi = lvs->begin(), lve = lvs->end(); lvi != lve; ++lvi) {
152  nameMap.emplace((*lvi)->GetName(), *lvi);
153  if (allSteps_ > 0)
154  mapLV_[*lvi] = (*lvi)->GetName();
155  }
156  for (auto const& name : nameEBSD_) {
157  for (itr = nameMap.begin(); itr != nameMap.end(); ++itr) {
158  const std::string& lvname = itr->first;
159  if (lvname.find(name) != std::string::npos) {
160  volEBSD_.emplace_back(itr->second);
161  int type = (lvname.find("refl") == std::string::npos) ? -1 : 1;
162  G4Trap* solid = static_cast<G4Trap*>(itr->second->GetSolid());
163  double dz = 2 * solid->GetZHalfLength() / CLHEP::mm;
164  xtalMap_.insert(std::pair<const G4LogicalVolume*, double>(itr->second, dz * type));
165  }
166  }
167  }
168  for (auto const& name : nameEESD_) {
169  for (itr = nameMap.begin(); itr != nameMap.end(); ++itr) {
170  const std::string& lvname = itr->first;
171  if (lvname.find(name) != std::string::npos) {
172  volEESD_.emplace_back(itr->second);
173  int type = (lvname.find("refl") == std::string::npos) ? 1 : -1;
174  G4Trap* solid = static_cast<G4Trap*>(itr->second->GetSolid());
175  double dz = 2 * solid->GetZHalfLength() / CLHEP::mm;
176  xtalMap_.insert(std::pair<const G4LogicalVolume*, double>(itr->second, dz * type));
177  }
178  }
179  }
180  for (auto const& name : nameHCSD_) {
181  for (itr = nameMap.begin(); itr != nameMap.end(); ++itr) {
182  const std::string& lvname = itr->first;
183  if (lvname.find(name) != std::string::npos)
184  volHCSD_.emplace_back(itr->second);
185  }
186  }
187  }
188 #ifdef EDM_ML_DEBUG
189  edm::LogVerbatim("Step") << volEBSD_.size() << " logical volumes for EB SD";
190  for (unsigned int k = 0; k < volEBSD_.size(); ++k)
191  edm::LogVerbatim("Step") << "[" << k << "] " << volEBSD_[k];
192  edm::LogVerbatim("Step") << volEESD_.size() << " logical volumes for EE SD";
193  for (unsigned int k = 0; k < volEESD_.size(); ++k)
194  edm::LogVerbatim("Step") << "[" << k << "] " << volEESD_[k];
195  edm::LogVerbatim("Step") << volHCSD_.size() << " logical volumes for HC SD";
196  for (unsigned int k = 0; k < volHCSD_.size(); ++k)
197  edm::LogVerbatim("Step") << "[" << k << "] " << volHCSD_[k];
198 #endif
199 }
200 
201 //=================================================================== per EVENT
203  eventID_ = (*evt)()->GetEventID();
204  edm::LogVerbatim("Step") << "CaloSteppingAction: Begin of event = " << eventID_;
205  for (int k = 0; k < CaloSteppingAction::nSD_; ++k) {
206  hitMap_[k].erase(hitMap_[k].begin(), hitMap_[k].end());
207  slave_[k].get()->Initialize();
208  }
209  if (allSteps_ > 0)
210  store_.clear();
211 }
212 
213 //=================================================================== each STEP
214 void CaloSteppingAction::update(const G4Step* aStep) {
215  // edm::LogVerbatim("Step") <<"CaloSteppingAction: At each Step";
216  NaNTrap(aStep);
217  auto lv = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
218  bool hc = (std::find(volHCSD_.begin(), volHCSD_.end(), lv) != volHCSD_.end());
219  bool eb = (std::find(volEBSD_.begin(), volEBSD_.end(), lv) != volEBSD_.end());
220  bool ee = (std::find(volEESD_.begin(), volEESD_.end(), lv) != volEESD_.end());
221  if (hc || eb || ee) {
222  double dEStep = aStep->GetTotalEnergyDeposit() / CLHEP::MeV;
223  auto const theTrack = aStep->GetTrack();
224  double time = theTrack->GetGlobalTime() / CLHEP::nanosecond;
225  int primID = theTrack->GetTrackID();
227  auto const touch = aStep->GetPreStepPoint()->GetTouchable();
228  auto const& hitPoint = aStep->GetPreStepPoint()->GetPosition();
229  if (hc) {
230  int depth = (touch->GetReplicaNumber(0)) % 10 + 1;
231  int lay = (touch->GetReplicaNumber(0) / 10) % 100 + 1;
232  int det = (touch->GetReplicaNumber(1)) / 1000;
233  auto unitID = getDetIDHC(det, lay, depth, math::XYZVectorD(hitPoint.x(), hitPoint.y(), hitPoint.z()));
234  if (unitID > 0 && dEStep > 0.0) {
235  dEStep *= getBirkHC(dEStep,
236  (aStep->GetStepLength() / CLHEP::cm),
237  aStep->GetPreStepPoint()->GetCharge(),
238  (aStep->GetPreStepPoint()->GetMaterial()->GetDensity() / (CLHEP::g / CLHEP::cm3)));
239  fillHit(unitID, dEStep, time, primID, 0, em, 2);
240  }
241  } else {
242  EcalBaseNumber theBaseNumber;
243  int size = touch->GetHistoryDepth() + 1;
244  if (theBaseNumber.getCapacity() < size)
245  theBaseNumber.setSize(size);
246  //Get name and copy numbers
247  if (size > 1) {
248  for (int ii = 0; ii < size; ii++) {
249  theBaseNumber.addLevel(touch->GetVolume(ii)->GetName(), touch->GetReplicaNumber(ii));
250  }
251  }
252  auto unitID =
253  (eb ? (ebNumberingScheme_->getUnitID(theBaseNumber)) : (eeNumberingScheme_->getUnitID(theBaseNumber)));
254  if (unitID > 0 && dEStep > 0.0) {
255  auto local = touch->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
256  auto ite = xtalMap_.find(lv);
257  double crystalLength = ((ite == xtalMap_.end()) ? 230.0 : std::abs(ite->second));
258  double crystalDepth =
259  ((ite == xtalMap_.end()) ? 0.0 : (std::abs(0.5 * (ite->second) + (local.z() / CLHEP::mm))));
260  double radl = aStep->GetPreStepPoint()->GetMaterial()->GetRadlen() / CLHEP::mm;
261  bool flag = ((ite == xtalMap_.end()) ? true : (((ite->second) >= 0) ? true : false));
262  auto depth = getDepth(flag, crystalDepth, radl);
263  dEStep *= (getBirkL3(dEStep,
264  (aStep->GetStepLength() / CLHEP::cm),
265  aStep->GetPreStepPoint()->GetCharge(),
266  (aStep->GetPreStepPoint()->GetMaterial()->GetDensity() / (CLHEP::g / CLHEP::cm3))) *
267  curve_LY(crystalLength, crystalDepth));
268  fillHit(unitID, dEStep, time, primID, depth, em, (eb ? 0 : 1));
269  }
270  }
271  }
272 
273  if (allSteps_ > 0) {
274  auto it = mapLV_.find(lv);
275  double energy = aStep->GetTotalEnergyDeposit() / CLHEP::MeV;
276  auto const touch = aStep->GetPreStepPoint()->GetTouchable();
277  double time = aStep->GetTrack()->GetGlobalTime() / CLHEP::nanosecond;
278  int trackId = aStep->GetTrack()->GetTrackID();
279  int pdg = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
280  double stepl = (aStep->GetStepLength() / CLHEP::cm);
281 #ifdef EDM_ML_DEBUG
282  edm::LogVerbatim("Step") << "CaloSteppingAction: Volume " << lv->GetName() << " History "
283  << touch->GetHistoryDepth() << " Pointers " << aStep->GetPostStepPoint() << ":"
284  << aStep->GetTrack()->GetNextVolume() << ":" << aStep->IsLastStepInVolume() << " E "
285  << energy << " T " << time << " PDG " << pdg << " step " << stepl;
286 #endif
287  uint32_t copy = 0;
288  if (((aStep->GetPostStepPoint() == nullptr) || (aStep->GetTrack()->GetNextVolume() == nullptr)) &&
289  (aStep->IsLastStepInVolume())) {
290  energy += (aStep->GetPreStepPoint()->GetKineticEnergy() / CLHEP::MeV);
291  } else {
292  time = aStep->GetPostStepPoint()->GetGlobalTime() / CLHEP::nanosecond;
293  copy = (touch->GetHistoryDepth() < 1)
294  ? static_cast<uint32_t>(touch->GetReplicaNumber(0))
295  : static_cast<uint32_t>(touch->GetReplicaNumber(0) + 1000 * touch->GetReplicaNumber(1));
296  }
297  if (it != mapLV_.end()) {
298  PassiveData key(std::make_tuple(lv, copy, trackId, pdg, time, energy, energy, stepl));
299  store_.push_back(key);
300  }
301  }
302 }
303 
304 //================================================================ End of EVENT
306  ++count_;
307  // Fill event input
308  edm::LogVerbatim("Step") << "CaloSteppingAction: EndOfEvent " << (*evt)()->GetEventID();
309 }
310 
311 void CaloSteppingAction::NaNTrap(const G4Step* aStep) const {
312  auto currentPos = aStep->GetTrack()->GetPosition();
313  double xyz = currentPos.x() + currentPos.y() + currentPos.z();
314  auto currentMom = aStep->GetTrack()->GetMomentum();
315  xyz += currentMom.x() + currentMom.y() + currentMom.z();
316 
317  if (edm::isNotFinite(xyz)) {
318  auto pCurrentVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
319  auto& nameOfVol = pCurrentVol->GetName();
320  throw cms::Exception("Unknown", "CaloSteppingAction")
321  << " Corrupted Event - NaN detected in volume " << nameOfVol << "\n";
322  }
323 }
324 
325 uint32_t CaloSteppingAction::getDetIDHC(int det, int lay, int depth, const math::XYZVectorD& pos) const {
326  HcalNumberingFromDDD::HcalID tmp = hcNumberingPS_.get()->unitID(det, lay, depth, pos);
327 #ifdef HcalNumberingTest
328  auto id0 = hcNumberingScheme_.get()->getUnitID(tmp);
329  HcalNumberingFromDDD::HcalID tmpO = hcNumbering_.get()->unitID(det, pos, depth, lay);
330  auto idO = hcNumberingScheme_.get()->getUnitID(tmpO);
331  std::string error = (id0 == idO) ? " ** OK **" : " ** ERROR **";
332  edm::LogVerbatim("Step") << "Det ID " << HcalDetId(id0) << " Original " << HcalDetId(idO) << error;
333 #endif
334  return (hcNumberingScheme_.get()->getUnitID(tmp));
335 }
336 
337 void CaloSteppingAction::fillHit(uint32_t id, double dE, double time, int primID, uint16_t depth, double em, int flag) {
338  CaloHitID currentID(id, time, primID, depth, timeSliceUnit_);
339  double edepEM = (em ? dE : 0);
340  double edepHAD = (em ? 0 : dE);
341  std::pair<int, CaloHitID> evID = std::make_pair(eventID_, currentID);
342  auto it = hitMap_[flag].find(evID);
343  if (it != hitMap_[flag].end()) {
344  (it->second).addEnergyDeposit(edepEM, edepHAD);
345  } else {
346  CaloGVHit aHit;
347  aHit.setEventID(eventID_);
348  aHit.setID(currentID);
349  aHit.addEnergyDeposit(edepEM, edepHAD);
350  hitMap_[flag][evID] = aHit;
351  }
352 }
353 
354 uint16_t CaloSteppingAction::getDepth(bool flag, double crystalDepth, double radl) const {
355  uint16_t depth1 = (flag ? 0 : PCaloHit::kEcalDepthRefz);
356  uint16_t depth2 = (uint16_t)floor(crystalDepth / radl);
357  uint16_t depth = (((depth2 & PCaloHit::kEcalDepthMask) << PCaloHit::kEcalDepthOffset) | depth1);
358 #ifdef EDM_ML_DEBUG
359  edm::LogVerbatim("Step") << "CaloSteppingAction::getDepth radl " << radl << ":" << crystalDepth << " depth " << depth;
360 #endif
361  return depth;
362 }
363 
364 double CaloSteppingAction::curve_LY(double crystalLength, double crystalDepth) const {
365  double weight = 1.;
366  double dapd = crystalLength - crystalDepth;
367  if (dapd >= -0.1 || dapd <= crystalLength + 0.1) {
368  if (dapd <= 100.)
369  weight = 1.0 + slopeLY_ - dapd * 0.01 * slopeLY_;
370 #ifdef EDM_ML_DEBUG
371  edm::LogVerbatim("Step") << "CaloSteppingAction::curve_LY " << crystalDepth << ":" << crystalLength << ":" << dapd
372  << ":" << weight;
373 #endif
374  } else {
375  edm::LogWarning("Step") << "CaloSteppingAction: light coll curve : wrong "
376  << "distance to APD " << dapd << " crlength = " << crystalLength
377  << " crystal Depth = " << crystalDepth << " weight = " << weight;
378  }
379  return weight;
380 }
381 
382 double CaloSteppingAction::getBirkL3(double dEStep, double step, double charge, double density) const {
383  double weight = 1.;
384  if (charge != 0. && step > 0.) {
385  double dedx = dEStep / step;
386  double rkb = birkC1EC_ / density;
387  if (dedx > 0) {
388  weight = 1. - birkSlopeEC_ * log(rkb * dedx);
389  if (weight < birkCutEC_)
390  weight = birkCutEC_;
391  else if (weight > 1.)
392  weight = 1.;
393  }
394 #ifdef EDM_ML_DEBUG
395  edm::LogVerbatim("Step") << "CaloSteppingAction::getBirkL3 Charge " << charge << " dE/dx " << dedx << " Birk Const "
396  << rkb << " Weight = " << weight << " dE " << dEStep << " step " << step;
397 #endif
398  }
399  return weight;
400 }
401 
402 double CaloSteppingAction::getBirkHC(double dEStep, double step, double charge, double density) const {
403  double weight = 1.;
404  if (charge != 0. && step > 0.) {
405  double dedx = dEStep / step;
406  double rkb = birkC1HC_ / density;
407  double c = birkC2HC_ * rkb * rkb;
408  if (std::abs(charge) >= 2.)
409  rkb /= birkC3HC_;
410  weight = 1. / (1. + rkb * dedx + c * dedx * dedx);
411 #ifdef EDM_ML_DEBUG
412  edm::LogVerbatim("Step") << "CaloSteppingAction::getBirkHC Charge " << charge << " dE/dx " << dedx << " Birk Const "
413  << rkb << ", " << c << " Weight = " << weight << " dE " << dEStep;
414 #endif
415  }
416  return weight;
417 }
418 
420  edm::LogVerbatim("Step") << "CaloSteppingAction:: saveHits for type " << type << " with " << hitMap_[type].size()
421  << " hits";
422  slave_[type].get()->ReserveMemory(hitMap_[type].size());
423  for (auto const& hit : hitMap_[type]) {
424  slave_[type].get()->processHits(hit.second.getUnitID(),
425  0.001 * hit.second.getEM(),
426  0.001 * hit.second.getHadr(),
427  hit.second.getTimeSlice(),
428  hit.second.getTrackID(),
429  hit.second.getDepth());
430  }
431 }
size
Write out results.
std::map< std::pair< int, CaloHitID >, CaloGVHit > hitMap_[nSD_]
std::vector< const G4LogicalVolume * > volEBSD_
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
void setID(uint32_t i, double d, int j, uint16_t k=0)
Definition: CaloGVHit.h:55
std::vector< PCaloHit > PCaloHitContainer
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
double getBirkL3(double dE, double step, double chg, double dens) const
static const int kEcalDepthRefz
Definition: PCaloHit.h:63
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
void produce(edm::Event &, const edm::EventSetup &) override
uint32_t getDetIDHC(int det, int lay, int depth, const math::XYZVectorD &pos) const
double getBirkHC(double dE, double step, double chg, double dens) const
double curve_LY(double crystalLength, double crystalDepth) const
Definition: weight.py:1
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
void setEventID(int id)
Definition: CaloGVHit.h:46
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
const double MeV
void fillHits(edm::PCaloHitContainer &cc, int type)
void addLevel(const std::string &name, const int &copyNumber)
std::vector< const G4LogicalVolume * > volEESD_
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
std::vector< std::string > nameHitC_
std::map< const G4LogicalVolume *, std::string > mapLV_
CaloSteppingAction(const edm::ParameterSet &p)
std::vector< const G4LogicalVolume * > volHCSD_
static const int kEcalDepthMask
Definition: PCaloHit.h:61
std::vector< std::string > nameHCSD_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::map< const G4LogicalVolume *, double > xtalMap_
static const int kEcalDepthOffset
Definition: PCaloHit.h:62
#define end
Definition: vmac.h:39
void fillPassiveHits(edm::PassiveHitContainer &cc)
std::vector< PassiveHit > PassiveHitContainer
Definition: PassiveHit.h:75
void addEnergyDeposit(double em, double hd)
Definition: CaloGVHit.cc:29
std::unique_ptr< CaloSlaveSD > slave_[nSD_]
ii
Definition: cuy.py:590
std::vector< std::string > nameEBSD_
void fillHit(uint32_t id, double dE, double time, int primID, uint16_t depth, double em, int flag)
void NaNTrap(const G4Step *) const
std::unique_ptr< HcalNumberingFromPS > hcNumberingPS_
std::unique_ptr< HcalNumberingScheme > hcNumberingScheme_
void update(const BeginOfJob *job) override
This routine will be called when the appropriate signal arrives.
std::unique_ptr< EcalBarrelNumberingScheme > ebNumberingScheme_
#define begin
Definition: vmac.h:32
static bool isGammaElectronPositron(int pdgCode)
susybsm::HSCParticleCollection hc
Definition: classes.h:25
std::unique_ptr< EcalEndcapNumberingScheme > eeNumberingScheme_
std::vector< std::string > nameEESD_
std::vector< PassiveData > store_
step
Definition: StallMonitor.cc:94
static const int nSD_
void setSize(const int &size)
tmp
align.sh
Definition: createJobs.py:716
uint16_t getDepth(bool flag, double crystalDepth, double radl) const
T const * product() const
Definition: ESHandle.h:86
def move(src, dest)
Definition: eostools.py:511
~CaloSteppingAction() override
std::tuple< const G4LogicalVolume *, uint32_t, int, int, double, double, double, double > PassiveData