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  slopeLY_ = iC.getParameter<double>("SlopeLightYield");
40  birkC1EC_ = iC.getParameter<double>("BirkC1EC") * (g / (MeV * cm2));
41  birkSlopeEC_ = iC.getParameter<double>("BirkSlopeEC");
42  birkCutEC_ = iC.getParameter<double>("BirkCutEC");
43  birkC1HC_ = iC.getParameter<double>("BirkC1HC") * (g / (MeV * cm2));
44  birkC2HC_ = iC.getParameter<double>("BirkC2HC");
45  birkC3HC_ = iC.getParameter<double>("BirkC3HC");
46 
47  edm::LogVerbatim("Step") << "CaloSteppingAction:: " << nameEBSD_.size() << " names for EB SD's";
48  for (unsigned int k = 0; k < nameEBSD_.size(); ++k)
49  edm::LogVerbatim("Step") << "[" << k << "] " << nameEBSD_[k];
50  edm::LogVerbatim("Step") << "CaloSteppingAction:: " << nameEESD_.size() << " names for EE SD's";
51  for (unsigned int k = 0; k < nameEESD_.size(); ++k)
52  edm::LogVerbatim("Step") << "[" << k << "] " << nameEESD_[k];
53  edm::LogVerbatim("Step") << "CaloSteppingAction:: " << nameHCSD_.size() << " names for HC SD's";
54  for (unsigned int k = 0; k < nameHCSD_.size(); ++k)
55  edm::LogVerbatim("Step") << "[" << k << "] " << nameHCSD_[k];
56  edm::LogVerbatim("Step") << "CaloSteppingAction::Constants for ECAL: slope " << slopeLY_ << " Birk constants "
57  << birkC1EC_ << ":" << birkSlopeEC_ << ":" << birkCutEC_;
58  edm::LogVerbatim("Step") << "CaloSteppingAction::Constants for HCAL: Birk "
59  << "constants " << birkC1HC_ << ":" << birkC2HC_ << ":" << birkC3HC_;
60  edm::LogVerbatim("Step") << "CaloSteppingAction:: " << nameHitC_.size() << " hit collections";
61  for (unsigned int k = 0; k < nameHitC_.size(); ++k)
62  edm::LogVerbatim("Step") << "[" << k << "] " << nameHitC_[k];
63 
64  ebNumberingScheme_ = std::make_unique<EcalBarrelNumberingScheme>();
65  eeNumberingScheme_ = std::make_unique<EcalEndcapNumberingScheme>();
66  hcNumberingPS_ = std::make_unique<HcalNumberingFromPS>(iC);
67  hcNumberingScheme_ = std::make_unique<HcalNumberingScheme>();
68 #ifdef HcalNumberingTest
69  hcNumbering_.reset(nullptr);
70 #endif
71  for (int k = 0; k < CaloSteppingAction::nSD_; ++k) {
72  slave_[k] = std::make_unique<CaloSlaveSD>(nameHitC_[k]);
73  produces<edm::PCaloHitContainer>(nameHitC_[k]);
74  }
75 }
76 
78  edm::LogVerbatim("Step") << "CaloSteppingAction: --------> Total number of "
79  << "selected entries : " << count_;
80 }
81 
83  for (int k = 0; k < CaloSteppingAction::nSD_; ++k) {
84  saveHits(k);
85  auto product = std::make_unique<edm::PCaloHitContainer>();
86  fillHits(*product, k);
87  e.put(std::move(product), nameHitC_[k]);
88  }
89 }
90 
92  edm::LogVerbatim("Step") << "CaloSteppingAction::fillHits for type " << type << " with "
93  << slave_[type].get()->hits().size() << " hits";
94  cc = slave_[type].get()->hits();
95  slave_[type].get()->Clean();
96 }
97 
99  edm::LogVerbatim("Step") << "CaloSteppingAction:: Enter BeginOfJob";
100 
101 #ifdef HcalNumberingTest
102  // Numbering From DDD
104  (*job)()->get<HcalSimNumberingRecord>().get(hdc);
105  const HcalDDDSimConstants* hcons_ = hdc.product();
106  edm::LogVerbatim("Step") << "CaloSteppingAction:: Initialise "
107  << "HcalNumberingFromDDD";
108  hcNumbering_ = std::make_unique<HcalNumberingFromDDD>(hcons_);
109 #endif
110 }
111 
112 //==================================================================== per RUN
114  int irun = (*run)()->GetRunID();
115  edm::LogVerbatim("Step") << "CaloSteppingAction:: Begin of Run = " << irun;
116 
117  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
118  if (lvs) {
119  std::map<const std::string, const G4LogicalVolume*> nameMap;
120  std::map<const std::string, const G4LogicalVolume*>::const_iterator itr;
121  for (auto lvi = lvs->begin(), lve = lvs->end(); lvi != lve; ++lvi)
122  nameMap.emplace((*lvi)->GetName(), *lvi);
123  for (auto const& name : nameEBSD_) {
124  for (itr = nameMap.begin(); itr != nameMap.end(); ++itr) {
125  const std::string& lvname = itr->first;
126  if (lvname.find(name) != std::string::npos) {
127  volEBSD_.emplace_back(itr->second);
128  int type = (lvname.find("refl") == std::string::npos) ? -1 : 1;
129  G4Trap* solid = static_cast<G4Trap*>(itr->second->GetSolid());
130  double dz = 2 * solid->GetZHalfLength();
131  xtalMap_.insert(std::pair<const G4LogicalVolume*, double>(itr->second, dz * type));
132  }
133  }
134  }
135  for (auto const& name : nameEESD_) {
136  for (itr = nameMap.begin(); itr != nameMap.end(); ++itr) {
137  const std::string& lvname = itr->first;
138  if (lvname.find(name) != std::string::npos) {
139  volEESD_.emplace_back(itr->second);
140  int type = (lvname.find("refl") == std::string::npos) ? 1 : -1;
141  G4Trap* solid = static_cast<G4Trap*>(itr->second->GetSolid());
142  double dz = 2 * solid->GetZHalfLength();
143  xtalMap_.insert(std::pair<const G4LogicalVolume*, double>(itr->second, dz * type));
144  }
145  }
146  }
147  for (auto const& name : nameHCSD_) {
148  for (itr = nameMap.begin(); itr != nameMap.end(); ++itr) {
149  const std::string& lvname = itr->first;
150  if (lvname.find(name) != std::string::npos)
151  volHCSD_.emplace_back(itr->second);
152  }
153  }
154  }
155 #ifdef EDM_ML_DEBUG
156  edm::LogVerbatim("Step") << volEBSD_.size() << " logical volumes for EB SD";
157  for (unsigned int k = 0; k < volEBSD_.size(); ++k)
158  edm::LogVerbatim("Step") << "[" << k << "] " << volEBSD_[k];
159  edm::LogVerbatim("Step") << volEESD_.size() << " logical volumes for EE SD";
160  for (unsigned int k = 0; k < volEESD_.size(); ++k)
161  edm::LogVerbatim("Step") << "[" << k << "] " << volEESD_[k];
162  edm::LogVerbatim("Step") << volHCSD_.size() << " logical volumes for HC SD";
163  for (unsigned int k = 0; k < volHCSD_.size(); ++k)
164  edm::LogVerbatim("Step") << "[" << k << "] " << volHCSD_[k];
165 #endif
166 }
167 
168 //=================================================================== per EVENT
170  eventID_ = (*evt)()->GetEventID();
171  edm::LogVerbatim("Step") << "CaloSteppingAction: Begin of event = " << eventID_;
172  for (int k = 0; k < CaloSteppingAction::nSD_; ++k) {
173  hitMap_[k].erase(hitMap_[k].begin(), hitMap_[k].end());
174  slave_[k].get()->Initialize();
175  }
176 }
177 
178 //=================================================================== each STEP
179 void CaloSteppingAction::update(const G4Step* aStep) {
180  // edm::LogVerbatim("Step") <<"CaloSteppingAction: At each Step";
181  NaNTrap(aStep);
182  auto lv = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
183  bool hc = (std::find(volHCSD_.begin(), volHCSD_.end(), lv) != volHCSD_.end());
184  bool eb = (std::find(volEBSD_.begin(), volEBSD_.end(), lv) != volEBSD_.end());
185  bool ee = (std::find(volEESD_.begin(), volEESD_.end(), lv) != volEESD_.end());
186  if (hc || eb || ee) {
187  double dEStep = aStep->GetTotalEnergyDeposit();
188  auto const theTrack = aStep->GetTrack();
189  double time = theTrack->GetGlobalTime() / nanosecond;
190  int primID = theTrack->GetTrackID();
192  auto const touch = aStep->GetPreStepPoint()->GetTouchable();
193  auto const& hitPoint = aStep->GetPreStepPoint()->GetPosition();
194  if (hc) {
195  int depth = (touch->GetReplicaNumber(0)) % 10 + 1;
196  int lay = (touch->GetReplicaNumber(0) / 10) % 100 + 1;
197  int det = (touch->GetReplicaNumber(1)) / 1000;
198  auto unitID = getDetIDHC(det, lay, depth, math::XYZVectorD(hitPoint.x(), hitPoint.y(), hitPoint.z()));
199  if (unitID > 0 && dEStep > 0.0) {
200  dEStep *= getBirkHC(dEStep,
201  aStep->GetStepLength(),
202  aStep->GetPreStepPoint()->GetCharge(),
203  aStep->GetPreStepPoint()->GetMaterial()->GetDensity());
204  fillHit(unitID, dEStep, time, primID, 0, em, 2);
205  }
206  } else {
207  EcalBaseNumber theBaseNumber;
208  int size = touch->GetHistoryDepth() + 1;
209  if (theBaseNumber.getCapacity() < size)
210  theBaseNumber.setSize(size);
211  //Get name and copy numbers
212  if (size > 1) {
213  for (int ii = 0; ii < size; ii++) {
214  theBaseNumber.addLevel(touch->GetVolume(ii)->GetName(), touch->GetReplicaNumber(ii));
215  }
216  }
217  auto unitID =
218  (eb ? (ebNumberingScheme_->getUnitID(theBaseNumber)) : (eeNumberingScheme_->getUnitID(theBaseNumber)));
219  if (unitID > 0 && dEStep > 0.0) {
220  auto local = touch->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
221  auto ite = xtalMap_.find(lv);
222  double crystalLength = ((ite == xtalMap_.end()) ? 230.0 : std::abs(ite->second));
223  double crystalDepth = ((ite == xtalMap_.end()) ? 0.0 : (std::abs(0.5 * (ite->second) + local.z())));
224  double radl = aStep->GetPreStepPoint()->GetMaterial()->GetRadlen();
225  bool flag = ((ite == xtalMap_.end()) ? true : (((ite->second) >= 0) ? true : false));
226  auto depth = getDepth(flag, crystalDepth, radl);
227  dEStep *= (getBirkL3(dEStep,
228  aStep->GetStepLength(),
229  aStep->GetPreStepPoint()->GetCharge(),
230  aStep->GetPreStepPoint()->GetMaterial()->GetDensity()) *
231  curve_LY(crystalLength, crystalDepth));
232  fillHit(unitID, dEStep, time, primID, depth, em, (eb ? 0 : 1));
233  }
234  }
235  }
236 }
237 
238 //================================================================ End of EVENT
240  ++count_;
241  // Fill event input
242  edm::LogVerbatim("Step") << "CaloSteppingAction: EndOfEvent " << (*evt)()->GetEventID();
243 }
244 
245 void CaloSteppingAction::NaNTrap(const G4Step* aStep) const {
246  auto currentPos = aStep->GetTrack()->GetPosition();
247  double xyz = currentPos.x() + currentPos.y() + currentPos.z();
248  auto currentMom = aStep->GetTrack()->GetMomentum();
249  xyz += currentMom.x() + currentMom.y() + currentMom.z();
250 
251  if (edm::isNotFinite(xyz)) {
252  auto pCurrentVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
253  auto& nameOfVol = pCurrentVol->GetName();
254  throw cms::Exception("Unknown", "CaloSteppingAction")
255  << " Corrupted Event - NaN detected in volume " << nameOfVol << "\n";
256  }
257 }
258 
259 uint32_t CaloSteppingAction::getDetIDHC(int det, int lay, int depth, const math::XYZVectorD& pos) const {
260  HcalNumberingFromDDD::HcalID tmp = hcNumberingPS_.get()->unitID(det, lay, depth, pos);
261 #ifdef HcalNumberingTest
262  auto id0 = hcNumberingScheme_.get()->getUnitID(tmp);
263  HcalNumberingFromDDD::HcalID tmpO = hcNumbering_.get()->unitID(det, pos, depth, lay);
264  auto idO = hcNumberingScheme_.get()->getUnitID(tmpO);
265  std::string error = (id0 == idO) ? " ** OK **" : " ** ERROR **";
266  edm::LogVerbatim("Step") << "Det ID " << HcalDetId(id0) << " Original " << HcalDetId(idO) << error;
267 #endif
268  return (hcNumberingScheme_.get()->getUnitID(tmp));
269 }
270 
271 void CaloSteppingAction::fillHit(uint32_t id, double dE, double time, int primID, uint16_t depth, double em, int flag) {
272  CaloHitID currentID(id, time, primID, depth);
273  double edepEM = (em ? dE : 0);
274  double edepHAD = (em ? 0 : dE);
275  std::pair<int, CaloHitID> evID = std::make_pair(eventID_, currentID);
276  auto it = hitMap_[flag].find(evID);
277  if (it != hitMap_[flag].end()) {
278  (it->second).addEnergyDeposit(edepEM, edepHAD);
279  } else {
280  CaloGVHit aHit;
281  aHit.setEventID(eventID_);
282  aHit.setID(currentID);
283  aHit.addEnergyDeposit(edepEM, edepHAD);
284  hitMap_[flag][evID] = aHit;
285  }
286 }
287 
288 uint16_t CaloSteppingAction::getDepth(bool flag, double crystalDepth, double radl) const {
289  uint16_t depth1 = (flag ? 0 : PCaloHit::kEcalDepthRefz);
290  uint16_t depth2 = (uint16_t)floor(crystalDepth / radl);
291  uint16_t depth = (((depth2 & PCaloHit::kEcalDepthMask) << PCaloHit::kEcalDepthOffset) | depth1);
292 #ifdef EDM_ML_DEBUG
293  edm::LogVerbatim("Step") << "CaloSteppingAction::getDepth radl " << radl << ":" << crystalDepth << " depth " << depth;
294 #endif
295  return depth;
296 }
297 
298 double CaloSteppingAction::curve_LY(double crystalLength, double crystalDepth) const {
299  double weight = 1.;
300  double dapd = crystalLength - crystalDepth;
301  if (dapd >= -0.1 || dapd <= crystalLength + 0.1) {
302  if (dapd <= 100.)
303  weight = 1.0 + slopeLY_ - dapd * 0.01 * slopeLY_;
304 #ifdef EDM_ML_DEBUG
305  edm::LogVerbatim("Step") << "CaloSteppingAction::curve_LY " << crystalDepth << ":" << crystalLength << ":" << dapd
306  << ":" << weight;
307 #endif
308  } else {
309  edm::LogWarning("Step") << "CaloSteppingAction: light coll curve : wrong "
310  << "distance to APD " << dapd << " crlength = " << crystalLength
311  << " crystal Depth = " << crystalDepth << " weight = " << weight;
312  }
313  return weight;
314 }
315 
316 double CaloSteppingAction::getBirkL3(double dEStep, double step, double charge, double density) const {
317  double weight = 1.;
318  if (charge != 0. && step > 0.) {
319  double dedx = dEStep / step;
320  double rkb = birkC1EC_ / density;
321  if (dedx > 0) {
322  weight = 1. - birkSlopeEC_ * log(rkb * dedx);
323  if (weight < birkCutEC_)
324  weight = birkCutEC_;
325  else if (weight > 1.)
326  weight = 1.;
327  }
328 #ifdef EDM_ML_DEBUG
329  edm::LogVerbatim("Step") << "CaloSteppingAction::getBirkL3 Charge " << charge << " dE/dx " << dedx << " Birk Const "
330  << rkb << " Weight = " << weight << " dE " << dEStep << " step " << step;
331 #endif
332  }
333  return weight;
334 }
335 
336 double CaloSteppingAction::getBirkHC(double dEStep, double step, double charge, double density) const {
337  double weight = 1.;
338  if (charge != 0. && step > 0.) {
339  double dedx = dEStep / step;
340  double rkb = birkC1HC_ / density;
341  double c = birkC2HC_ * rkb * rkb;
342  if (std::abs(charge) >= 2.)
343  rkb /= birkC3HC_;
344  weight = 1. / (1. + rkb * dedx + c * dedx * dedx);
345 #ifdef EDM_ML_DEBUG
346  edm::LogVerbatim("Step") << "CaloSteppingAction::getBirkHC Charge " << charge << " dE/dx " << dedx << " Birk Const "
347  << rkb << ", " << c << " Weight = " << weight << " dE " << dEStep;
348 #endif
349  }
350  return weight;
351 }
352 
354  edm::LogVerbatim("Step") << "CaloSteppingAction:: saveHits for type " << type << " with " << hitMap_[type].size()
355  << " hits";
356  slave_[type].get()->ReserveMemory(hitMap_[type].size());
357  for (auto const& hit : hitMap_[type]) {
358  slave_[type].get()->processHits(hit.second.getUnitID(),
359  hit.second.getEM() / GeV,
360  hit.second.getHadr() / GeV,
361  hit.second.getTimeSlice(),
362  hit.second.getTrackID(),
363  hit.second.getDepth());
364  }
365 }
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
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:125
const double GeV
Definition: MathUtil.h:16
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:20
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_
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 addEnergyDeposit(double em, double hd)
Definition: CaloGVHit.cc:29
std::unique_ptr< CaloSlaveSD > slave_[nSD_]
ii
Definition: cuy.py:590
int k[5][pyjets_maxn]
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::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
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_
step
Definition: StallMonitor.cc:94
static const int nSD_
void setSize(const int &size)
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