CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
CaloSteppingAction.cc
Go to the documentation of this file.
1 //#define EDM_ML_DEBUG
2 //#define HcalNumberingTest
3 
4 // to make hits in EB/EE/HC
5 
9 
16 
20 #ifdef HcalNumberingTest
24 #endif
25 
30 
34 
42 
43 #include "G4LogicalVolume.hh"
44 #include "G4LogicalVolumeStore.hh"
45 #include "G4NavigationHistory.hh"
46 #include "G4ParticleTable.hh"
47 #include "G4PhysicalVolumeStore.hh"
48 #include "G4Region.hh"
49 #include "G4RegionStore.hh"
50 #include "G4Step.hh"
51 #include "G4SystemOfUnits.hh"
52 #include "G4Track.hh"
53 #include "G4Trap.hh"
54 #include "G4UnitsTable.hh"
55 #include "G4UserSteppingAction.hh"
56 #include "G4VPhysicalVolume.hh"
57 #include "G4VTouchable.hh"
58 
59 #include <algorithm>
60 #include <cmath>
61 #include <iostream>
62 #include <iomanip>
63 #include <string>
64 #include <utility>
65 #include <vector>
66 
68  public Observer<const BeginOfRun*>,
69  public Observer<const BeginOfEvent*>,
70  public Observer<const EndOfEvent*>,
71  public Observer<const G4Step*> {
72 public:
74  ~CaloSteppingAction() override;
75 
77  void produce(edm::Event&, const edm::EventSetup&) override;
78  void beginRun(edm::EventSetup const&) override;
79 
80 private:
81  void fillHits(edm::PCaloHitContainer& cc, int type);
83  // observer classes
84  void update(const BeginOfRun* run) override;
85  void update(const BeginOfEvent* evt) override;
86  void update(const G4Step* step) override;
87  void update(const EndOfEvent* evt) override;
88 
89  void NaNTrap(const G4Step*) const;
90  uint32_t getDetIDHC(int det, int lay, int depth, const math::XYZVectorD& pos) const;
91  void fillHit(uint32_t id, double dE, double time, int primID, uint16_t depth, double em, int flag);
92  uint16_t getDepth(bool flag, double crystalDepth, double radl) const;
93  double curve_LY(double crystalLength, double crystalDepth) const;
94  double getBirkL3(double dE, double step, double chg, double dens) const;
95  double getBirkHC(double dE, double step, double chg, double dens) const;
96  void saveHits(int flag);
97 
98  static const int nSD_ = 3;
99  std::unique_ptr<EcalBarrelNumberingScheme> ebNumberingScheme_;
100  std::unique_ptr<EcalEndcapNumberingScheme> eeNumberingScheme_;
101  std::unique_ptr<HcalNumberingFromPS> hcNumberingPS_;
102 #ifdef HcalNumberingTest
104  std::unique_ptr<HcalNumberingFromDDD> hcNumbering_;
105 #endif
106  std::unique_ptr<HcalNumberingScheme> hcNumberingScheme_;
107  std::unique_ptr<CaloSlaveSD> slave_[nSD_];
108 
109  std::vector<std::string> nameEBSD_, nameEESD_, nameHCSD_;
110  std::vector<std::string> nameHitC_;
111  std::vector<const G4LogicalVolume*> volEBSD_, volEESD_, volHCSD_;
112  std::map<const G4LogicalVolume*, double> xtalMap_;
113  std::map<const G4LogicalVolume*, std::string> mapLV_;
118  std::map<std::pair<int, CaloHitID>, CaloGVHit> hitMap_[nSD_];
119  typedef std::tuple<const G4LogicalVolume*, uint32_t, int, int, double, double, double, double, double, double, double>
121  std::vector<PassiveData> store_;
122 };
123 
125  edm::ParameterSet iC = p.getParameter<edm::ParameterSet>("CaloSteppingAction");
126  nameEBSD_ = iC.getParameter<std::vector<std::string> >("EBSDNames");
127  nameEESD_ = iC.getParameter<std::vector<std::string> >("EESDNames");
128  nameHCSD_ = iC.getParameter<std::vector<std::string> >("HCSDNames");
129  nameHitC_ = iC.getParameter<std::vector<std::string> >("HitCollNames");
130  allSteps_ = iC.getParameter<int>("AllSteps");
131  slopeLY_ = iC.getParameter<double>("SlopeLightYield");
132  birkC1EC_ = iC.getParameter<double>("BirkC1EC");
133  birkSlopeEC_ = iC.getParameter<double>("BirkSlopeEC");
134  birkCutEC_ = iC.getParameter<double>("BirkCutEC");
135  birkC1HC_ = iC.getParameter<double>("BirkC1HC");
136  birkC2HC_ = iC.getParameter<double>("BirkC2HC");
137  birkC3HC_ = iC.getParameter<double>("BirkC3HC");
138  timeSliceUnit_ = iC.getUntrackedParameter<double>("TimeSliceUnit", 1.0);
139 
140  edm::LogVerbatim("Step") << "CaloSteppingAction:: " << nameEBSD_.size() << " names for EB SD's";
141  for (unsigned int k = 0; k < nameEBSD_.size(); ++k)
142  edm::LogVerbatim("Step") << "[" << k << "] " << nameEBSD_[k];
143  edm::LogVerbatim("Step") << "CaloSteppingAction:: " << nameEESD_.size() << " names for EE SD's";
144  for (unsigned int k = 0; k < nameEESD_.size(); ++k)
145  edm::LogVerbatim("Step") << "[" << k << "] " << nameEESD_[k];
146  edm::LogVerbatim("Step") << "CaloSteppingAction:: " << nameHCSD_.size() << " names for HC SD's";
147  for (unsigned int k = 0; k < nameHCSD_.size(); ++k)
148  edm::LogVerbatim("Step") << "[" << k << "] " << nameHCSD_[k];
149  edm::LogVerbatim("Step") << "CaloSteppingAction::Constants for ECAL: slope " << slopeLY_ << " Birk constants "
150  << birkC1EC_ << ":" << birkSlopeEC_ << ":" << birkCutEC_;
151  edm::LogVerbatim("Step") << "CaloSteppingAction::Constants for HCAL: Birk "
152  << "constants " << birkC1HC_ << ":" << birkC2HC_ << ":" << birkC3HC_;
153  edm::LogVerbatim("Step") << "CaloSteppingAction::Constant for time slice " << timeSliceUnit_;
154  edm::LogVerbatim("Step") << "CaloSteppingAction:: " << nameHitC_.size() << " hit collections";
155  for (unsigned int k = 0; k < nameHitC_.size(); ++k)
156  edm::LogVerbatim("Step") << "[" << k << "] " << nameHitC_[k];
157 
158  ebNumberingScheme_ = std::make_unique<EcalBarrelNumberingScheme>();
159  eeNumberingScheme_ = std::make_unique<EcalEndcapNumberingScheme>();
160  hcNumberingPS_ = std::make_unique<HcalNumberingFromPS>(iC);
161  hcNumberingScheme_ = std::make_unique<HcalNumberingScheme>();
162 #ifdef HcalNumberingTest
163  hcNumbering_.reset(nullptr);
164 #endif
165  for (int k = 0; k < CaloSteppingAction::nSD_; ++k) {
166  slave_[k] = std::make_unique<CaloSlaveSD>(nameHitC_[k]);
167  produces<edm::PCaloHitContainer>(nameHitC_[k]);
168  }
169  if (allSteps_ != 0)
170  produces<edm::PassiveHitContainer>("AllPassiveHits");
171  edm::LogVerbatim("Step") << "CaloSteppingAction:: All Steps Flag " << allSteps_ << " for passive hits";
172 }
173 
175  edm::LogVerbatim("Step") << "CaloSteppingAction: --------> Total number of "
176  << "selected entries : " << count_;
177 }
178 
180 #ifdef HcalNumberingTest
182  edm::LogVerbatim("Step") << "CaloSteppingAction::Initialize ESGetToken for HcalDDDSimConstants";
183 #endif
184 }
185 
187  for (int k = 0; k < CaloSteppingAction::nSD_; ++k) {
188  saveHits(k);
189  auto product = std::make_unique<edm::PCaloHitContainer>();
190  fillHits(*product, k);
191  e.put(std::move(product), nameHitC_[k]);
192  }
193  if (allSteps_ != 0) {
194  std::unique_ptr<edm::PassiveHitContainer> hgcPH(new edm::PassiveHitContainer);
195  fillPassiveHits(*hgcPH);
196  e.put(std::move(hgcPH), "AllPassiveHits");
197  }
198 }
199 
201  edm::LogVerbatim("Step") << "CaloSteppingAction::fillHits for type " << type << " with "
202  << slave_[type].get()->hits().size() << " hits";
203  cc = slave_[type].get()->hits();
204  slave_[type].get()->Clean();
205 }
206 
208  edm::LogVerbatim("Step") << "CaloSteppingAction::fillPassiveHits with " << store_.size() << " hits";
209  for (const auto& element : store_) {
210  auto lv = std::get<0>(element);
211  auto it = mapLV_.find(lv);
212  if (it != mapLV_.end()) {
213  PassiveHit hit(it->second,
214  std::get<1>(element),
215  std::get<5>(element),
216  std::get<6>(element),
217  std::get<4>(element),
218  std::get<2>(element),
219  std::get<3>(element),
220  std::get<7>(element),
221  std::get<8>(element),
222  std::get<9>(element),
223  std::get<10>(element));
224  cc.emplace_back(hit);
225  }
226  }
227 }
228 
230  edm::LogVerbatim("Step") << "CaloSteppingAction:: Enter BeginOfRun";
231 
232 #ifdef HcalNumberingTest
233  // Numbering From DDD
234  const HcalDDDSimConstants* hcons_ = &es.getData(ddconsToken_);
235  edm::LogVerbatim("Step") << "CaloSteppingAction:: Initialise "
236  << "HcalNumberingFromDDD";
237  hcNumbering_ = std::make_unique<HcalNumberingFromDDD>(hcons_);
238 #endif
239 }
240 
241 //==================================================================== per RUN
243  int irun = (*run)()->GetRunID();
244  edm::LogVerbatim("Step") << "CaloSteppingAction:: Begin of Run = " << irun;
245 
246  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
247  if (lvs) {
248  std::map<const std::string, const G4LogicalVolume*> nameMap;
249  std::map<const std::string, const G4LogicalVolume*>::const_iterator itr;
250  for (auto lvi = lvs->begin(), lve = lvs->end(); lvi != lve; ++lvi) {
251  nameMap.emplace((*lvi)->GetName(), *lvi);
252  if (allSteps_ < 0)
253  mapLV_[*lvi] = (*lvi)->GetName();
254  }
255 
256  for (auto const& name : nameEBSD_) {
257  for (itr = nameMap.begin(); itr != nameMap.end(); ++itr) {
258  const std::string& lvname = itr->first;
259  if (lvname.find(name) != std::string::npos) {
260  volEBSD_.emplace_back(itr->second);
261  int type = (lvname.find("refl") == std::string::npos) ? -1 : 1;
262  G4Trap* solid = static_cast<G4Trap*>(itr->second->GetSolid());
263  double dz = 2 * solid->GetZHalfLength() / CLHEP::mm;
264  xtalMap_.insert(std::pair<const G4LogicalVolume*, double>(itr->second, dz * type));
265  if ((allSteps_ > 0) && ((allSteps_ % 10) > 0))
266  mapLV_[itr->second] = itr->first;
267  }
268  }
269  }
270  for (auto const& name : nameEESD_) {
271  for (itr = nameMap.begin(); itr != nameMap.end(); ++itr) {
272  const std::string& lvname = itr->first;
273  if (lvname.find(name) != std::string::npos) {
274  volEESD_.emplace_back(itr->second);
275  int type = (lvname.find("refl") == std::string::npos) ? 1 : -1;
276  G4Trap* solid = static_cast<G4Trap*>(itr->second->GetSolid());
277  double dz = 2 * solid->GetZHalfLength() / CLHEP::mm;
278  xtalMap_.insert(std::pair<const G4LogicalVolume*, double>(itr->second, dz * type));
279  if ((allSteps_ > 0) && (((allSteps_ / 10) % 10) > 0))
280  mapLV_[itr->second] = itr->first;
281  }
282  }
283  }
284 
285  for (auto const& name : nameHCSD_) {
286  for (itr = nameMap.begin(); itr != nameMap.end(); ++itr) {
287  const std::string& lvname = itr->first;
288  if (lvname.find(name) != std::string::npos) {
289  volHCSD_.emplace_back(itr->second);
290  if ((allSteps_ > 0) && (((allSteps_ / 100) % 10) > 0))
291  mapLV_[itr->second] = itr->first;
292  }
293  }
294  }
295  }
296 #ifdef EDM_ML_DEBUG
297  edm::LogVerbatim("Step") << volEBSD_.size() << " logical volumes for EB SD";
298  for (unsigned int k = 0; k < volEBSD_.size(); ++k)
299  edm::LogVerbatim("Step") << "[" << k << "] " << volEBSD_[k];
300  edm::LogVerbatim("Step") << volEESD_.size() << " logical volumes for EE SD";
301  for (unsigned int k = 0; k < volEESD_.size(); ++k)
302  edm::LogVerbatim("Step") << "[" << k << "] " << volEESD_[k];
303  edm::LogVerbatim("Step") << volHCSD_.size() << " logical volumes for HC SD";
304  for (unsigned int k = 0; k < volHCSD_.size(); ++k)
305  edm::LogVerbatim("Step") << "[" << k << "] " << volHCSD_[k];
306  edm::LogVerbatim("Step") << mapLV_.size() << " logical volumes for Passive hits";
307  unsigned int k(0);
308  for (auto itr = mapLV_.begin(); itr != mapLV_.end(); ++itr) {
309  edm::LogVerbatim("Step") << "[" << k << "] " << itr->second << ":" << itr->first;
310  ++k;
311  }
312 #endif
313 }
314 
315 //=================================================================== per EVENT
317  eventID_ = (*evt)()->GetEventID();
318  edm::LogVerbatim("Step") << "CaloSteppingAction: Begin of event = " << eventID_;
319  for (int k = 0; k < CaloSteppingAction::nSD_; ++k) {
320  hitMap_[k].erase(hitMap_[k].begin(), hitMap_[k].end());
321  slave_[k].get()->Initialize();
322  }
323  if (allSteps_ != 0)
324  store_.clear();
325 }
326 
327 //=================================================================== each STEP
328 void CaloSteppingAction::update(const G4Step* aStep) {
329  // edm::LogVerbatim("Step") <<"CaloSteppingAction: At each Step";
330  NaNTrap(aStep);
331  auto lv = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
332  bool hc = (std::find(volHCSD_.begin(), volHCSD_.end(), lv) != volHCSD_.end());
333  bool eb = (std::find(volEBSD_.begin(), volEBSD_.end(), lv) != volEBSD_.end());
334  bool ee = (std::find(volEESD_.begin(), volEESD_.end(), lv) != volEESD_.end());
335  uint32_t unitID(0);
336  if (hc || eb || ee) {
337  double dEStep = aStep->GetTotalEnergyDeposit() / CLHEP::MeV;
338  auto const theTrack = aStep->GetTrack();
339  double time = theTrack->GetGlobalTime() / CLHEP::nanosecond;
340  int primID = theTrack->GetTrackID();
342  auto const touch = aStep->GetPreStepPoint()->GetTouchable();
343  auto const& hitPoint = aStep->GetPreStepPoint()->GetPosition();
344  if (hc) {
345  int depth = (touch->GetReplicaNumber(0)) % 10 + 1;
346  int lay = (touch->GetReplicaNumber(0) / 10) % 100 + 1;
347  int det = (touch->GetReplicaNumber(1)) / 1000;
348  unitID = getDetIDHC(det, lay, depth, math::XYZVectorD(hitPoint.x(), hitPoint.y(), hitPoint.z()));
349  if (unitID > 0 && dEStep > 0.0) {
350  dEStep *= getBirkHC(dEStep,
351  (aStep->GetStepLength() / CLHEP::cm),
352  aStep->GetPreStepPoint()->GetCharge(),
353  (aStep->GetPreStepPoint()->GetMaterial()->GetDensity() / (CLHEP::g / CLHEP::cm3)));
354  fillHit(unitID, dEStep, time, primID, 0, em, 2);
355  }
356  } else {
357  EcalBaseNumber theBaseNumber;
358  int size = touch->GetHistoryDepth() + 1;
359  if (theBaseNumber.getCapacity() < size)
360  theBaseNumber.setSize(size);
361  //Get name and copy numbers
362  if (size > 1) {
363  for (int ii = 0; ii < size; ii++) {
364  theBaseNumber.addLevel(touch->GetVolume(ii)->GetName(), touch->GetReplicaNumber(ii));
365  }
366  }
367  unitID = (eb ? (ebNumberingScheme_->getUnitID(theBaseNumber)) : (eeNumberingScheme_->getUnitID(theBaseNumber)));
368  if (unitID > 0 && dEStep > 0.0) {
369  auto local = touch->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
370  auto ite = xtalMap_.find(lv);
371  double crystalLength = ((ite == xtalMap_.end()) ? 230.0 : std::abs(ite->second));
372  double crystalDepth =
373  ((ite == xtalMap_.end()) ? 0.0 : (std::abs(0.5 * (ite->second) + (local.z() / CLHEP::mm))));
374  double radl = aStep->GetPreStepPoint()->GetMaterial()->GetRadlen() / CLHEP::mm;
375  bool flag = ((ite == xtalMap_.end()) ? true : (((ite->second) >= 0) ? true : false));
376  auto depth = getDepth(flag, crystalDepth, radl);
377  dEStep *= (getBirkL3(dEStep,
378  (aStep->GetStepLength() / CLHEP::cm),
379  aStep->GetPreStepPoint()->GetCharge(),
380  (aStep->GetPreStepPoint()->GetMaterial()->GetDensity() / (CLHEP::g / CLHEP::cm3))) *
381  curve_LY(crystalLength, crystalDepth));
382  fillHit(unitID, dEStep, time, primID, depth, em, (eb ? 0 : 1));
383  }
384  }
385  }
386 
387  if (allSteps_ != 0) {
388  auto it = mapLV_.find(lv);
389  if (it != mapLV_.end()) {
390  double energy = aStep->GetTotalEnergyDeposit() / CLHEP::MeV;
391  auto const touch = aStep->GetPreStepPoint()->GetTouchable();
392  double time = aStep->GetTrack()->GetGlobalTime() / CLHEP::nanosecond;
393  int trackId = aStep->GetTrack()->GetTrackID();
394  int pdg = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
395  double stepl = (aStep->GetStepLength() / CLHEP::cm);
396  double xp = aStep->GetPreStepPoint()->GetPosition().x() / CLHEP::cm;
397  double yp = aStep->GetPreStepPoint()->GetPosition().y() / CLHEP::cm;
398  double zp = aStep->GetPreStepPoint()->GetPosition().z() / CLHEP::cm;
399 #ifdef EDM_ML_DEBUG
400  edm::LogVerbatim("Step") << "CaloSteppingAction: Volume " << lv->GetName() << " History "
401  << touch->GetHistoryDepth() << " Pointers " << aStep->GetPostStepPoint() << ":"
402  << aStep->GetTrack()->GetNextVolume() << ":" << aStep->IsLastStepInVolume() << " E "
403  << energy << " T " << time << " PDG " << pdg << " step " << stepl << " Position (" << xp
404  << ", " << yp << ", " << zp << ")";
405 #endif
406  uint32_t copy = (allSteps_ < 0) ? 0 : unitID;
407  if (((aStep->GetPostStepPoint() == nullptr) || (aStep->GetTrack()->GetNextVolume() == nullptr)) &&
408  (aStep->IsLastStepInVolume())) {
409  energy += (aStep->GetPreStepPoint()->GetKineticEnergy() / CLHEP::MeV);
410  } else {
411  time = aStep->GetPostStepPoint()->GetGlobalTime() / CLHEP::nanosecond;
412  if (copy == 0)
413  copy = (touch->GetHistoryDepth() < 1)
414  ? static_cast<uint32_t>(touch->GetReplicaNumber(0))
415  : static_cast<uint32_t>(touch->GetReplicaNumber(0) + 1000 * touch->GetReplicaNumber(1));
416  }
417  PassiveData key(std::make_tuple(lv, copy, trackId, pdg, time, energy, energy, stepl, xp, yp, zp));
418  store_.push_back(key);
419  }
420  }
421 }
422 
423 //================================================================ End of EVENT
425  ++count_;
426  // Fill event input
427  edm::LogVerbatim("Step") << "CaloSteppingAction: EndOfEvent " << (*evt)()->GetEventID();
428 }
429 
430 void CaloSteppingAction::NaNTrap(const G4Step* aStep) const {
431  auto currentPos = aStep->GetTrack()->GetPosition();
432  double xyz = currentPos.x() + currentPos.y() + currentPos.z();
433  auto currentMom = aStep->GetTrack()->GetMomentum();
434  xyz += currentMom.x() + currentMom.y() + currentMom.z();
435 
436  if (edm::isNotFinite(xyz)) {
437  auto pCurrentVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
438  auto& nameOfVol = pCurrentVol->GetName();
439  throw cms::Exception("Unknown", "CaloSteppingAction")
440  << " Corrupted Event - NaN detected in volume " << nameOfVol << "\n";
441  }
442 }
443 
444 uint32_t CaloSteppingAction::getDetIDHC(int det, int lay, int depth, const math::XYZVectorD& pos) const {
445  HcalNumberingFromDDD::HcalID tmp = hcNumberingPS_.get()->unitID(det, lay, depth, pos);
446 #ifdef HcalNumberingTest
447  auto id0 = hcNumberingScheme_.get()->getUnitID(tmp);
448  HcalNumberingFromDDD::HcalID tmpO = hcNumbering_.get()->unitID(det, pos, depth, lay);
449  auto idO = hcNumberingScheme_.get()->getUnitID(tmpO);
450  std::string error = (id0 == idO) ? " ** OK **" : " ** ERROR **";
451  edm::LogVerbatim("Step") << "Det ID " << HcalDetId(id0) << " Original " << HcalDetId(idO) << error;
452 #endif
453  return (hcNumberingScheme_.get()->getUnitID(tmp));
454 }
455 
456 void CaloSteppingAction::fillHit(uint32_t id, double dE, double time, int primID, uint16_t depth, double em, int flag) {
457  CaloHitID currentID(id, time, primID, depth, timeSliceUnit_);
458  double edepEM = (em ? dE : 0);
459  double edepHAD = (em ? 0 : dE);
460  std::pair<int, CaloHitID> evID = std::make_pair(eventID_, currentID);
461  auto it = hitMap_[flag].find(evID);
462  if (it != hitMap_[flag].end()) {
463  (it->second).addEnergyDeposit(edepEM, edepHAD);
464  } else {
465  CaloGVHit aHit;
466  aHit.setEventID(eventID_);
467  aHit.setID(currentID);
468  aHit.addEnergyDeposit(edepEM, edepHAD);
469  hitMap_[flag][evID] = aHit;
470  }
471 }
472 
473 uint16_t CaloSteppingAction::getDepth(bool flag, double crystalDepth, double radl) const {
474  uint16_t depth1 = (flag ? 0 : PCaloHit::kEcalDepthRefz);
475  uint16_t depth2 = (uint16_t)floor(crystalDepth / radl);
476  uint16_t depth = (((depth2 & PCaloHit::kEcalDepthMask) << PCaloHit::kEcalDepthOffset) | depth1);
477 #ifdef EDM_ML_DEBUG
478  edm::LogVerbatim("Step") << "CaloSteppingAction::getDepth radl " << radl << ":" << crystalDepth << " depth " << depth;
479 #endif
480  return depth;
481 }
482 
483 double CaloSteppingAction::curve_LY(double crystalLength, double crystalDepth) const {
484  double weight = 1.;
485  double dapd = crystalLength - crystalDepth;
486  if (dapd >= -0.1 || dapd <= crystalLength + 0.1) {
487  if (dapd <= 100.)
488  weight = 1.0 + slopeLY_ - dapd * 0.01 * slopeLY_;
489 #ifdef EDM_ML_DEBUG
490  edm::LogVerbatim("Step") << "CaloSteppingAction::curve_LY " << crystalDepth << ":" << crystalLength << ":" << dapd
491  << ":" << weight;
492 #endif
493  } else {
494  edm::LogWarning("Step") << "CaloSteppingAction: light coll curve : wrong "
495  << "distance to APD " << dapd << " crlength = " << crystalLength
496  << " crystal Depth = " << crystalDepth << " weight = " << weight;
497  }
498  return weight;
499 }
500 
501 double CaloSteppingAction::getBirkL3(double dEStep, double step, double charge, double density) const {
502  double weight = 1.;
503  if (charge != 0. && step > 0.) {
504  double dedx = dEStep / step;
505  double rkb = birkC1EC_ / density;
506  if (dedx > 0) {
507  weight = 1. - birkSlopeEC_ * log(rkb * dedx);
508  if (weight < birkCutEC_)
509  weight = birkCutEC_;
510  else if (weight > 1.)
511  weight = 1.;
512  }
513 #ifdef EDM_ML_DEBUG
514  edm::LogVerbatim("Step") << "CaloSteppingAction::getBirkL3 Charge " << charge << " dE/dx " << dedx << " Birk Const "
515  << rkb << " Weight = " << weight << " dE " << dEStep << " step " << step;
516 #endif
517  }
518  return weight;
519 }
520 
521 double CaloSteppingAction::getBirkHC(double dEStep, double step, double charge, double density) const {
522  double weight = 1.;
523  if (charge != 0. && step > 0.) {
524  double dedx = dEStep / step;
525  double rkb = birkC1HC_ / density;
526  double c = birkC2HC_ * rkb * rkb;
527  if (std::abs(charge) >= 2.)
528  rkb /= birkC3HC_;
529  weight = 1. / (1. + rkb * dedx + c * dedx * dedx);
530 #ifdef EDM_ML_DEBUG
531  edm::LogVerbatim("Step") << "CaloSteppingAction::getBirkHC Charge " << charge << " dE/dx " << dedx << " Birk Const "
532  << rkb << ", " << c << " Weight = " << weight << " dE " << dEStep;
533 #endif
534  }
535  return weight;
536 }
537 
539  edm::LogVerbatim("Step") << "CaloSteppingAction:: saveHits for type " << type << " with " << hitMap_[type].size()
540  << " hits";
541  slave_[type].get()->ReserveMemory(hitMap_[type].size());
542  for (auto const& hit : hitMap_[type]) {
543  slave_[type].get()->processHits(hit.second.getUnitID(),
544  0.001 * hit.second.getEM(),
545  0.001 * hit.second.getHadr(),
546  hit.second.getTimeSlice(),
547  hit.second.getTrackID(),
548  hit.second.getDepth());
549  }
550 }
551 
554 
std::map< std::pair< int, CaloHitID >, CaloGVHit > hitMap_[nSD_]
Log< level::Info, true > LogVerbatim
T getUntrackedParameter(std::string const &, T const &) const
#define DEFINE_SIMWATCHER(type)
std::vector< const G4LogicalVolume * > volEBSD_
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:133
static std::vector< std::string > checklist log
const edm::EventSetup & c
double getBirkL3(double dE, double step, double chg, double dens) const
static const int kEcalDepthRefz
Definition: PCaloHit.h:63
const float chg[109]
Definition: CoreSimTrack.cc:5
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
void produce(edm::Event &, const edm::EventSetup &) override
void registerConsumes(edm::ConsumesCollector) override
uint32_t getDetIDHC(int det, int lay, int depth, const math::XYZVectorD &pos) const
std::map< const G4LogicalVolume *, double > xtalMap_
std::vector< const G4LogicalVolume * > volHCSD_
double getBirkHC(double dE, double step, double chg, double dens) const
double curve_LY(double crystalLength, double crystalDepth) const
int ii
Definition: cuy.py:589
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
bool getData(T &iHolder) const
Definition: EventSetup.h:128
const double MeV
void fillHits(edm::PCaloHitContainer &cc, int type)
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
std::vector< std::string > nameHitC_
CaloSteppingAction(const edm::ParameterSet &p)
void update(const BeginOfRun *run) override
This routine will be called when the appropriate signal arrives.
def move
Definition: eostools.py:511
static const int kEcalDepthMask
Definition: PCaloHit.h:61
std::vector< std::string > nameHCSD_
tuple key
prepare the HTCondor submission files and eventually submit them
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static const int kEcalDepthOffset
Definition: PCaloHit.h:62
std::map< const G4LogicalVolume *, std::string > mapLV_
void fillPassiveHits(edm::PassiveHitContainer &cc)
std::vector< PassiveHit > PassiveHitContainer
Definition: PassiveHit.h:100
void addEnergyDeposit(double em, double hd)
Definition: CaloGVHit.cc:29
std::unique_ptr< CaloSlaveSD > slave_[nSD_]
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_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void beginRun(edm::EventSetup const &) override
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)
string end
Definition: dataset.py:937
std::unique_ptr< EcalEndcapNumberingScheme > eeNumberingScheme_
std::vector< std::string > nameEESD_
std::vector< PassiveData > store_
step
Definition: StallMonitor.cc:94
Log< level::Warning, false > LogWarning
int weight
Definition: histoStyle.py:51
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
tuple size
Write out results.