CMS 3D CMS Logo

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:
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 
110  const std::vector<std::string> nameEBSD_, nameEESD_, nameHCSD_;
111  const std::vector<std::string> nameHitC_;
112  std::vector<const G4LogicalVolume*> volEBSD_, volEESD_, volHCSD_;
113  std::map<const G4LogicalVolume*, double> xtalMap_;
114  std::map<const G4LogicalVolume*, std::string> mapLV_;
115  const int allSteps_;
118  const double birkC3HC_, timeSliceUnit_;
120  std::map<std::pair<int, CaloHitID>, CaloGVHit> hitMap_[nSD_];
121  typedef std::tuple<const G4LogicalVolume*, uint32_t, int, int, double, double, double, double, double, double, double>
123  std::vector<PassiveData> store_;
124 };
125 
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 }
177 
179  edm::LogVerbatim("Step") << "CaloSteppingAction: --------> Total number of "
180  << "selected entries : " << count_;
181 }
182 
184 #ifdef HcalNumberingTest
186  edm::LogVerbatim("Step") << "CaloSteppingAction::Initialize ESGetToken for HcalDDDSimConstants";
187 #endif
188 }
189 
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 }
203 
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 }
210 
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 }
232 
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 }
244 
245 //==================================================================== per RUN
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 }
318 
319 //=================================================================== per EVENT
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 }
330 
331 //=================================================================== each STEP
332 void CaloSteppingAction::update(const G4Step* aStep) {
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 }
426 
427 //================================================================ End of EVENT
429  ++count_;
430  // Fill event input
431  edm::LogVerbatim("Step") << "CaloSteppingAction: EndOfEvent " << (*evt)()->GetEventID();
432 }
433 
434 void CaloSteppingAction::NaNTrap(const G4Step* aStep) const {
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 }
447 
448 uint32_t CaloSteppingAction::getDetIDHC(int det, int lay, int depth, const math::XYZVectorD& pos) const {
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 }
459 
460 void CaloSteppingAction::fillHit(uint32_t id, double dE, double time, int primID, uint16_t depth, double em, int flag) {
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 }
476 
477 uint16_t CaloSteppingAction::getDepth(bool flag, double crystalDepth, double radl) const {
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 }
486 
487 double CaloSteppingAction::curve_LY(double crystalLength, double crystalDepth) const {
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 }
504 
505 double CaloSteppingAction::getBirkL3(double dEStep, double step, double charge, double density) const {
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 }
524 
525 double CaloSteppingAction::getBirkHC(double dEStep, double step, double charge, double density) const {
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 }
541 
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 }
555 
558 
size
Write out results.
std::map< std::pair< int, CaloHitID >, CaloGVHit > hitMap_[nSD_]
Log< level::Info, true > LogVerbatim
#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
const std::vector< std::string > nameHitC_
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const std::vector< std::string > nameEESD_
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
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
void registerConsumes(edm::ConsumesCollector) override
std::map< const G4LogicalVolume *, double > xtalMap_
std::vector< const G4LogicalVolume * > volHCSD_
const std::vector< std::string > nameEBSD_
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
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 fillHits(edm::PCaloHitContainer &cc, int type)
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
CaloSteppingAction(const edm::ParameterSet &p)
void update(const BeginOfRun *run) override
This routine will be called when the appropriate signal arrives.
static const int kEcalDepthMask
Definition: PCaloHit.h:61
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_
const edm::ParameterSet iC_
uint32_t getDetIDHC(int det, int lay, int depth, const math::XYZVectorD &pos) const
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_]
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::unique_ptr< HcalNumberingFromPS > hcNumberingPS_
std::unique_ptr< HcalNumberingScheme > hcNumberingScheme_
void beginRun(edm::EventSetup const &) override
std::vector< const G4LogicalVolume * > volEESD_
std::unique_ptr< EcalBarrelNumberingScheme > ebNumberingScheme_
HLT enums.
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_
step
Definition: StallMonitor.cc:98
const std::vector< std::string > nameHCSD_
Log< level::Warning, false > LogWarning
static const int nSD_
void setSize(const int &size)
tmp
align.sh
Definition: createJobs.py:716
double curve_LY(double crystalLength, double crystalDepth) const
def move(src, dest)
Definition: eostools.py:511