CMS 3D CMS Logo

SimG4HGCalValidation.cc
Go to the documentation of this file.
1 // File: SimG4HGCalValidation.cc
3 // Description: Main analysis class for HGCal Validation of G4 Hits
5 
9 
10 // to retreive hits
15 
21 
24 
28 
29 #include "G4Step.hh"
30 #include "G4Track.hh"
31 #include "G4NavigationHistory.hh"
32 #include "G4TouchableHistory.hh"
33 
34 #include "CLHEP/Units/GlobalSystemOfUnits.h"
35 #include "CLHEP/Units/GlobalPhysicalConstants.h"
36 
37 #include <cmath>
38 #include <iostream>
39 #include <iomanip>
40 #include <memory>
41 #include <vector>
42 #include <string>
43 
44 class SimG4HGCalValidation : public SimProducer, public Observer<const BeginOfEvent*>, public Observer<const G4Step*> {
45 public:
47  ~SimG4HGCalValidation() override = default;
48 
50  void produce(edm::Event&, const edm::EventSetup&) override;
51  void beginRun(edm::EventSetup const&) override;
52 
53 private:
54  SimG4HGCalValidation(const SimG4HGCalValidation&); // stop default
56 
57  void init();
58 
59  // observer classes
60  void update(const BeginOfEvent* evt) override;
61  void update(const G4Step* step) override;
62 
63  // analysis related class
65  void clear();
66 
67 private:
68  //HGCal numbering scheme
69  std::vector<edm::ESGetToken<HGCalDDDConstants, IdealGeometryRecord>> ddconsToken_;
70  std::vector<std::unique_ptr<HGCalNumberingScheme>> hgcalNumbering_;
71 
72  // to read from parameter set
74  const std::vector<std::string> names_;
75  const std::vector<int> types_, detTypes_;
77  const int verbosity_;
78  std::vector<int> subdet_;
79  std::vector<std::string> nameXs_;
80 
81  // parameters from geometry
83 
84  // some private members for ananlysis
85  unsigned int count_;
87  std::vector<double> hgcEEedep_, hgcHEFedep_, hgcHEBedep_;
88  std::vector<unsigned int> dets_, hgchitDets_, hgchitIndex_;
89  std::vector<double> hgchitX_, hgchitY_, hgchitZ_;
90 };
91 
93  : m_Anal(p.getParameter<edm::ParameterSet>("SimG4HGCalValidation")),
94  names_(m_Anal.getParameter<std::vector<std::string>>("Names")),
95  types_(m_Anal.getParameter<std::vector<int>>("Types")),
96  detTypes_(m_Anal.getParameter<std::vector<int>>("DetTypes")),
97  labelLayer_(m_Anal.getParameter<std::string>("LabelLayerInfo")),
98  verbosity_(m_Anal.getUntrackedParameter<int>("Verbosity", 0)),
99  levelT1_(999),
100  levelT2_(999),
101  count_(0) {
102  produces<PHGCalValidInfo>(labelLayer_);
103 
104  if (verbosity_ > 0) {
105  edm::LogVerbatim("ValidHGCal") << "HGCalTestAnalysis:: Initialised as "
106  << "observer of begin events and of G4step "
107  << "with Parameter values: \n\tLabel : " << labelLayer_ << " and with "
108  << names_.size() << " detectors";
109  for (unsigned int k = 0; k < names_.size(); ++k)
110  edm::LogVerbatim("ValidHGCal") << " [" << k << "] " << names_[k] << " Type " << types_[k] << " DetType "
111  << detTypes_[k];
112  }
113 }
114 
116  for (unsigned int type = 0; type < types_.size(); ++type) {
117  int detType = detTypes_[type];
118  G4String nameX = "HGCal";
119  if (types_[type] <= 1) {
120  dets_.emplace_back(static_cast<unsigned int>(DetId::Forward));
121  if (detType == 0)
122  subdet_.emplace_back(static_cast<int>(ForwardSubdetector::HGCEE));
123  else if (detType == 1)
124  subdet_.emplace_back(static_cast<int>(ForwardSubdetector::HGCHEF));
125  else
126  subdet_.emplace_back(static_cast<int>(ForwardSubdetector::HGCHEB));
127  if (detType == 0)
128  nameX = "HGCalEESensitive";
129  else if (detType == 1)
130  nameX = "HGCalHESiliconSensitive";
131  else
132  nameX = "HGCalHEScintillatorSensitive";
133  nameXs_.emplace_back(nameX);
134  ddconsToken_.emplace_back(
136  edm::LogVerbatim("ValidHGCal") << "Defines ESGetToken for type " << type << ":" << dets_.back() << ":"
137  << subdet_.back() << ":" << nameX;
138  } else {
139  edm::LogError("ValidHGCal") << "Wrong Type " << types_[type];
140  throw cms::Exception("Unknown", "ValidHGCal") << "Wrong Type " << types_[type] << "\n";
141  }
142  }
143 }
144 
146  std::unique_ptr<PHGCalValidInfo> productLayer(new PHGCalValidInfo);
147  layerAnalysis(*productLayer);
148  e.put(std::move(productLayer), labelLayer_);
149 }
150 
152  for (unsigned int type = 0; type < types_.size(); ++type) {
153  int layers(0);
154  int detType = detTypes_[type];
156  if (hdc.isValid()) {
157  levelT1_ = hdc->levelTop(0);
158  levelT2_ = hdc->levelTop(1);
159  hgcalNumbering_.emplace_back(
160  std::make_unique<HGCalNumberingScheme>(*hdc, static_cast<DetId::Detector>(dets_[type]), nameXs_[type], ""));
161  layers = hdc->layers(false);
162  } else {
163  edm::LogError("ValidHGCal") << "Cannot find HGCalDDDConstants for " << nameXs_[type];
164  throw cms::Exception("Unknown", "ValidHGCal") << "Cannot find HGCalDDDConstants for " << nameXs_[type] << "\n";
165  }
166  if (detType == 0) {
167  for (int i = 0; i < layers; ++i)
168  hgcEEedep_.emplace_back(0);
169  } else if (detType == 1) {
170  for (int i = 0; i < layers; ++i)
171  hgcHEFedep_.emplace_back(0);
172  } else {
173  for (int i = 0; i < layers; ++i)
174  hgcHEBedep_.emplace_back(0);
175  }
176  if (verbosity_ > 0)
177  edm::LogVerbatim("ValidHGCal") << "[" << type << "]: " << nameXs_[type] << " det " << dets_[type] << " subdet "
178  << subdet_[type] << " with " << layers << " layers";
179  }
180 }
181 
182 //=================================================================== per EVENT
184  int iev = (*evt)()->GetEventID();
185  if (verbosity_ > 0)
186  edm::LogVerbatim("ValidHGCal") << "SimG4HGCalValidation: =====> Begin "
187  << "event = " << iev;
188 
189  ++count_;
190  edepEE_ = edepHEF_ = edepHEB_ = 0.;
191 
192  //HGCal variables
193  for (unsigned int i = 0; i < hgcEEedep_.size(); i++)
194  hgcEEedep_[i] = 0.;
195  for (unsigned int i = 0; i < hgcHEFedep_.size(); i++)
196  hgcHEFedep_[i] = 0.;
197  for (unsigned int i = 0; i < hgcHEBedep_.size(); i++)
198  hgcHEBedep_[i] = 0.;
199 
200  //Cache reset
201  clear();
202 }
203 
204 //=================================================================== each STEP
205 void SimG4HGCalValidation::update(const G4Step* aStep) {
206  if (aStep != nullptr) {
207  G4VPhysicalVolume* curPV = aStep->GetPreStepPoint()->GetPhysicalVolume();
208  const G4String& name = curPV->GetName();
209  G4VSensitiveDetector* curSD = aStep->GetPreStepPoint()->GetSensitiveDetector();
210  if (verbosity_ > 1)
211  edm::LogVerbatim("ValidHGCal") << "ValidHGCal::Step in " << name << " at "
212  << aStep->GetPreStepPoint()->GetPosition();
213 
214  // Only for Sensitive detector
215  if (curSD != nullptr) {
216  int type(-1);
217  for (unsigned int k = 0; k < names_.size(); ++k) {
218  if (name.find(names_[k].c_str()) != std::string::npos) {
219  type = k;
220  break;
221  }
222  }
223  int detType = (type >= 0) ? detTypes_[type] : -1;
224  if (verbosity_ > 0)
225  edm::LogVerbatim("ValidHGCal") << "ValidHGCal::In SD " << name << " type " << type << ":" << detType;
226 
227  // Right type of SD
228  if (type >= 0) {
229  //Get the 32-bit index of the hit cell
230  const G4TouchableHistory* touchable =
231  static_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
232  unsigned int index(0);
233  int layer(0);
234  G4ThreeVector hitPoint = aStep->GetPreStepPoint()->GetPosition();
235  // HGCal
236  G4ThreeVector localpos = touchable->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
237  float globalZ = touchable->GetTranslation(0).z();
238  int iz(globalZ > 0 ? 1 : -1);
239  int module(-1), cell(-1);
240  if ((touchable->GetHistoryDepth() == levelT1_) || (touchable->GetHistoryDepth() == levelT2_)) {
241  layer = touchable->GetReplicaNumber(0);
242  } else {
243  layer = touchable->GetReplicaNumber(3);
244  module = touchable->GetReplicaNumber(2);
245  cell = touchable->GetReplicaNumber(1);
246  }
247  double weight(0);
248  index = hgcalNumbering_[type]->getUnitID(layer, module, cell, iz, hitPoint, weight);
249  if (verbosity_ > 1)
250  edm::LogVerbatim("ValidHGCal") << "HGCal: " << name << " Layer " << layer << " Module " << module << " Cell "
251  << cell;
252 
253  double edeposit = aStep->GetTotalEnergyDeposit();
254  if (verbosity_ > 0)
255  edm::LogVerbatim("ValidHGCal") << "Layer " << layer << " Index " << std::hex << index << std::dec << " Edep "
256  << edeposit << " hit at " << hitPoint;
257  if (detType == 0) {
258  edepEE_ += edeposit;
259  if (layer < static_cast<int>(hgcEEedep_.size()))
260  hgcEEedep_[layer] += edeposit;
261  } else if (detType == 1) {
262  edepHEF_ += edeposit;
263  if (layer < static_cast<int>(hgcHEFedep_.size()))
264  hgcHEFedep_[layer] += edeposit;
265  } else {
266  edepHEB_ += edeposit;
267  if (layer < static_cast<int>(hgcHEBedep_.size()))
268  hgcHEBedep_[layer] += edeposit;
269  }
270  G4String nextVolume("XXX");
271  if (aStep->GetTrack()->GetNextVolume() != nullptr)
272  nextVolume = aStep->GetTrack()->GetNextVolume()->GetName();
273 
274  if (nextVolume.c_str() != name.c_str()) { //save hit when it exits cell
275  if (std::find(hgchitIndex_.begin(), hgchitIndex_.end(), index) == hgchitIndex_.end()) {
276  hgchitDets_.emplace_back(dets_[type]);
277  hgchitIndex_.emplace_back(index);
278  hgchitX_.emplace_back(hitPoint.x());
279  hgchitY_.emplace_back(hitPoint.y());
280  hgchitZ_.emplace_back(hitPoint.z());
281  }
282  }
283  } // it is right type of SD
284  } // it is in a SD
285  } //if aStep!=NULL
286 } //end update aStep
287 
288 //================================================================ End of EVENT
289 
291  if (verbosity_ > 0)
292  edm::LogVerbatim("ValidHGCal") << "\n ===>>> SimG4HGCalValidation: Energy "
293  << "deposit\n at EE : " << std::setw(6) << edepEE_ / CLHEP::MeV
294  << "\n at HEF: " << std::setw(6) << edepHEF_ / CLHEP::MeV
295  << "\n at HEB: " << std::setw(6) << edepHEB_ / CLHEP::MeV;
296 
297  //Fill HGC Variables
300 }
301 
302 //---------------------------------------------------
304  hgchitDets_.erase(hgchitDets_.begin(), hgchitDets_.end());
305  hgchitIndex_.erase(hgchitIndex_.begin(), hgchitIndex_.end());
306  hgchitX_.erase(hgchitX_.begin(), hgchitX_.end());
307  hgchitY_.erase(hgchitY_.begin(), hgchitY_.end());
308  hgchitZ_.erase(hgchitZ_.begin(), hgchitZ_.end());
309 }
310 
std::vector< unsigned int > hgchitIndex_
Log< level::Info, true > LogVerbatim
std::vector< double > hgchitY_
std::vector< int > subdet_
#define DEFINE_SIMWATCHER(type)
~SimG4HGCalValidation() override=default
void produce(edm::Event &, const edm::EventSetup &) override
SimG4HGCalValidation(const edm::ParameterSet &p)
std::vector< std::string > nameXs_
const std::vector< int > detTypes_
std::vector< double > hgchitZ_
const std::vector< int > types_
Definition: weight.py:1
Log< level::Error, false > LogError
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
constexpr std::array< uint8_t, layerIndexSize > layer
void beginRun(edm::EventSetup const &) override
const std::vector< std::string > names_
void update(const BeginOfEvent *evt) override
This routine will be called when the appropriate signal arrives.
void fillhgcLayers(const double edepEE, const double edepHEF, const double edepHEB, const std::vector< double > &eedep, const std::vector< double > &hefdep, const std::vector< double > &hebdep)
void registerConsumes(edm::ConsumesCollector) override
const SimG4HGCalValidation & operator=(const SimG4HGCalValidation &)
std::vector< unsigned int > dets_
void layerAnalysis(PHGCalValidInfo &)
unsigned int layers(bool reco) const
std::vector< unsigned int > hgchitDets_
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:151
std::vector< double > hgcHEFedep_
bool isValid() const
Definition: ESHandle.h:44
std::vector< edm::ESGetToken< HGCalDDDConstants, IdealGeometryRecord > > ddconsToken_
std::vector< double > hgcHEBedep_
HitContainer const *__restrict__ TkSoA const *__restrict__ Quality const *__restrict__ CAHitNtupletGeneratorKernelsGPU::HitToTuple const *__restrict__ int32_t int32_t int iev
std::vector< double > hgcEEedep_
void fillhgcHits(const std::vector< unsigned int > &hitdets, const std::vector< unsigned int > &hitindex, const std::vector< double > &hitvtxX, const std::vector< double > &hitvtxY, const std::vector< double > &hitvtxZ)
HLT enums.
step
Definition: StallMonitor.cc:98
const std::string labelLayer_
const edm::ParameterSet m_Anal
int levelTop(int ind=0) const
std::vector< double > hgchitX_
def move(src, dest)
Definition: eostools.py:511
std::vector< std::unique_ptr< HGCalNumberingScheme > > hgcalNumbering_