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 
12 
13 // to retreive hits
23 
29 
36 
40 
41 #include "G4Step.hh"
42 #include "G4Track.hh"
43 #include "G4NavigationHistory.hh"
44 #include "G4TouchableHistory.hh"
45 
46 #include "CLHEP/Units/GlobalSystemOfUnits.h"
47 #include "CLHEP/Units/GlobalPhysicalConstants.h"
48 
49 #include <cmath>
50 #include <iostream>
51 #include <iomanip>
52 #include <memory>
53 #include <vector>
54 #include <string>
55 
57  public Observer<const BeginOfJob *>,
58  public Observer<const BeginOfEvent *>,
59  public Observer<const G4Step *> {
60 
61 public:
63  ~SimG4HGCalValidation() override;
64 
65  void produce(edm::Event&, const edm::EventSetup&) override;
66 
67 private:
68  SimG4HGCalValidation(const SimG4HGCalValidation&); // stop default
70 
71  void init();
72 
73  // observer classes
74  void update(const BeginOfJob * job) override;
75  void update(const BeginOfEvent * evt) override;
76  void update(const G4Step * step) override;
77 
78  // analysis related class
80  void clear();
81 
82 private:
83  //Keep reference to instantiate HcalNumberingFromDDD later
85 
86  //HGCal numbering scheme
87  std::vector<HGCNumberingScheme*> hgcNumbering_;
88  std::vector<HGCalNumberingScheme*> hgcalNumbering_;
89 
90  // to read from parameter set
91  std::vector<std::string> names_;
92  std::vector<int> types_, detTypes_, subdet_;
94 
95  // parameters from geometry
97 
98  // some private members for ananlysis
99  unsigned int count_;
102  std::vector<double> hgcEEedep_, hgcHEFedep_, hgcHEBedep_;
103  std::vector<unsigned int> dets_, hgchitDets_, hgchitIndex_;
104  std::vector<double> hgchitX_, hgchitY_, hgchitZ_;
105 };
106 
109 
110  edm::ParameterSet m_Anal = p.getParameter<edm::ParameterSet>("SimG4HGCalValidation");
111  names_ = m_Anal.getParameter<std::vector<std::string> >("Names");
112  types_ = m_Anal.getParameter<std::vector<int> >("Types");
113  detTypes_ = m_Anal.getParameter<std::vector<int> >("DetTypes");
114  labelLayer_ = m_Anal.getParameter<std::string>("LabelLayerInfo");
115  verbosity_ = m_Anal.getUntrackedParameter<int>("Verbosity",0);
116 
117  produces<PHGCalValidInfo>(labelLayer_);
118 
119  edm::LogVerbatim("ValidHGCal") << "HGCalTestAnalysis:: Initialised as "
120  << "observer of begin events and of G4step "
121  << "with Parameter values: \n\tLabel : "
122  << labelLayer_ << " and with "
123  << names_.size() << " detectors";
124  for (unsigned int k=0; k<names_.size(); ++k)
125  edm::LogVerbatim("ValidHGCal") << " [" << k << "] " << names_[k]
126  << " Type " << types_[k] << " DetType "
127  << detTypes_[k];
128 }
129 
131  delete numberingFromDDD_;
132  for (auto number : hgcNumbering_) delete number;
133  for (auto number : hgcalNumbering_) delete number;
134 }
135 
137 
138  std::unique_ptr<PHGCalValidInfo> productLayer(new PHGCalValidInfo);
139  layerAnalysis(*productLayer);
140  e.put(std::move(productLayer),labelLayer_);
141 }
142 
144 
145  const edm::EventSetup* es = (*job)();
146  for (unsigned int type=0; type<types_.size(); ++type) {
147  int layers(0);
148  int detType = detTypes_[type];
149  G4String nameX = "HGCal";
150  if (types_[type] <= 1) {
151  if (types_[type] == 0) {
153  if (detType == 0) dets_.emplace_back((int)(DetId::HGCalEE));
154  else if (detType == 1) dets_.emplace_back((int)(DetId::HGCalHSi));
155  else dets_.emplace_back((int)(DetId::HGCalHSc));
156  } else {
157  dets_.push_back((unsigned int)(DetId::Forward));
158  if (detType == 0) subdet_.emplace_back((int)(ForwardSubdetector::HGCEE));
159  else if (detType == 1) subdet_.emplace_back((int)(ForwardSubdetector::HGCHEF));
160  else subdet_.emplace_back((int)(ForwardSubdetector::HGCHEB));
161  }
162  if (detType == 0) nameX = "HGCalEESensitive";
163  else if (detType == 1) nameX = "HGCalHESiliconSensitive";
164  else nameX = "HGCalHEScintillatorSensitive";
166  es->get<IdealGeometryRecord>().get(nameX,hdc);
167  if (hdc.isValid()) {
169  levelT1_ = hdc->levelTop(0);
170  levelT2_ = hdc->levelTop(1);
171  if (m_mode == HGCalGeometryMode::Trapezoid) {
172  types_[type] =-1;
173  hgcalNumbering_.emplace_back(new HGCalNumberingScheme(*hdc,(DetId::Detector)(dets_[type]),nameX));
174  } else if ((m_mode == HGCalGeometryMode::Hexagon) ||
175  (m_mode == HGCalGeometryMode::HexagonFull)) {
176  types_[type] = 1;
177  hgcNumbering_.push_back(new HGCNumberingScheme(*hdc,nameX));
178  } else {
179  types_[type] = 0;
180  hgcalNumbering_.emplace_back(new HGCalNumberingScheme(*hdc,(DetId::Detector)(dets_[type]),nameX));
181  }
182  layers = hdc->layers(false);
183  } else {
184  edm::LogError("ValidHGCal") << "Cannot find HGCalDDDConstants for "
185  << nameX;
186  throw cms::Exception("Unknown", "ValidHGCal")
187  << "Cannot find HGCalDDDConstants for " << nameX << "\n";
188  }
189  } else {
190  nameX = "HcalEndcap";
191  dets_.push_back((unsigned int)(DetId::Hcal));
192  subdet_.push_back((int)(HcalSubdetector::HcalEndcap));
194  es->get<HcalSimNumberingRecord>().get(hdc);
195  if (hdc.isValid()) {
197  layers = 18;
198  } else {
199  edm::LogError("ValidHGCal") << "Cannot find HcalDDDSimConstant";
200  throw cms::Exception("Unknown", "ValidHGCal")
201  << "Cannot find HcalDDDSimConstant\n";
202  }
203  }
204  if (detType == 0) {
205  for (int i=0; i<layers; ++i) hgcEEedep_.push_back(0);
206  } else if (detType == 1) {
207  for (int i=0; i<layers; ++i) hgcHEFedep_.push_back(0);
208  } else {
209  for (int i=0; i<layers; ++i) hgcHEBedep_.push_back(0);
210  }
211  edm::LogVerbatim("ValidHGCal") << "[" << type << "]: " << nameX << " det "
212  << dets_[type] << " subdet "
213  << subdet_[type] << " with " << layers
214  << " layers";
215  }
216 }
217 
218 //=================================================================== per EVENT
220 
221  int iev = (*evt)()->GetEventID();
222  if (verbosity_ > 0)
223  edm::LogVerbatim("ValidHGCal") << "SimG4HGCalValidation: =====> Begin "
224  << "event = " << iev;
225 
226  ++count_;
227  edepEE_ = edepHEF_ = edepHEB_ = 0.;
228 
229  //HGCal variables
230  for (unsigned int i = 0; i<hgcEEedep_.size(); i++) hgcEEedep_[i] = 0.;
231  for (unsigned int i = 0; i<hgcHEFedep_.size(); i++) hgcHEFedep_[i] = 0.;
232  for (unsigned int i = 0; i<hgcHEBedep_.size(); i++) hgcHEBedep_[i] = 0.;
233 
234  //Cache reset
235  clear();
236 }
237 
238 //=================================================================== each STEP
239 void SimG4HGCalValidation::update(const G4Step * aStep) {
240 
241  if (aStep != nullptr) {
242  G4VPhysicalVolume* curPV = aStep->GetPreStepPoint()->GetPhysicalVolume();
243  const G4String& name = curPV->GetName();
244  G4VSensitiveDetector* curSD = aStep->GetPreStepPoint()->GetSensitiveDetector();
245  if (verbosity_ > 1)
246  edm::LogVerbatim("ValidHGCal") << "ValidHGCal::Step in " << name <<" at "
247  << aStep->GetPreStepPoint()->GetPosition();
248 
249  // Only for Sensitive detector
250  if (curSD != nullptr) {
251  int type(-1);
252  for (unsigned int k=0; k<names_.size(); ++k) {
253  if (name.find(names_[k].c_str()) != std::string::npos) {
254  type = k; break;
255  }
256  }
257  int detType = (type >= 0) ? detTypes_[type] : -1;
258  if (verbosity_ > 0)
259  edm::LogVerbatim("ValidHGCal") << "ValidHGCal::In SD " << name
260  << " type " << type << ":" << detType;
261 
262  // Right type of SD
263  if (type >= 0) {
264  //Get the 32-bit index of the hit cell
265  const G4TouchableHistory* touchable = static_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
266  unsigned int index(0);
267  int layer(0);
268  G4ThreeVector hitPoint = aStep->GetPreStepPoint()->GetPosition();
269  if (types_[type] <= 1) {
270  // HGCal
271  G4ThreeVector localpos = touchable->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
272  float globalZ = touchable->GetTranslation(0).z();
273  int iz(globalZ>0 ? 1 : -1);
274  int module(-1), cell(-1);
275  if (types_[type] == 1) {
276  if (touchable->GetHistoryDepth() == levelT1_) {
277  layer = touchable->GetReplicaNumber(0);
278  } else {
279  layer = touchable->GetReplicaNumber(2);
280  module = touchable->GetReplicaNumber(1);
281  cell = touchable->GetReplicaNumber(0);
282  }
283  index = hgcNumbering_[type]->getUnitID((ForwardSubdetector)(subdet_[type]), layer, module, cell, iz, localpos);
284  } else {
285  if ((touchable->GetHistoryDepth() == levelT1_) ||
286  (touchable->GetHistoryDepth() == levelT2_)) {
287  layer = touchable->GetReplicaNumber(0);
288  } else {
289  layer = touchable->GetReplicaNumber(3);
290  module = touchable->GetReplicaNumber(2);
291  cell = touchable->GetReplicaNumber(1);
292  }
293  double weight(0);
294  index = hgcalNumbering_[type]->getUnitID(layer, module, cell, iz,
295  hitPoint, weight);
296  }
297  if (verbosity_ > 1)
298  edm::LogVerbatim("ValidHGCal") << "HGCal: " << name << " Layer "
299  << layer << " Module " << module
300  << " Cell " << cell;
301  } else {
302  // Hcal
303  int depth = (touchable->GetReplicaNumber(0))%10 + 1;
304  int lay = (touchable->GetReplicaNumber(0)/10)%100 + 1;
305  int det = (touchable->GetReplicaNumber(1))/1000;
306  HcalNumberingFromDDD::HcalID tmp = numberingFromDDD_->unitID(det, hitPoint, depth, lay);
307  index = HcalTestNumbering::packHcalIndex(tmp.subdet,tmp.zside,tmp.depth,tmp.etaR,tmp.phis,tmp.lay);
308  layer = tmp.lay;
309  if (verbosity_ > 1)
310  edm::LogVerbatim("ValidHGCal") << "HCAL: " << det << ":" << depth
311  << ":" << lay << " o/p "
312  << tmp.subdet << ":" << tmp.zside
313  << ":" << tmp.depth << ":"
314  << tmp.etaR << ":" << tmp.phis
315  << ":" << tmp.lay << " point "
316  << hitPoint << " " << hitPoint.rho()
317  << ":" << hitPoint.eta() << ":"
318  << hitPoint.phi();
319  }
320 
321  double edeposit = aStep->GetTotalEnergyDeposit();
322  if (verbosity_ > 0)
323  edm::LogVerbatim("ValidHGCal") << "Layer " << layer << " Index "
324  << std::hex << index << std::dec
325  << " Edep " << edeposit << " hit at "
326  << hitPoint;
327  if (detType == 0) {
328  edepEE_ += edeposit;
329  if (layer < (int)(hgcEEedep_.size()))
330  hgcEEedep_[layer] += edeposit;
331  } else if (detType == 1) {
332  edepHEF_ += edeposit;
333  if (layer < (int)(hgcHEFedep_.size()))
334  hgcHEFedep_[layer] += edeposit;
335  } else {
336  edepHEB_ += edeposit;
337  if (layer < (int)(hgcHEBedep_.size()))
338  hgcHEBedep_[layer] += edeposit;
339  }
340  G4String nextVolume("XXX");
341  if (aStep->GetTrack()->GetNextVolume()!=nullptr)
342  nextVolume = aStep->GetTrack()->GetNextVolume()->GetName();
343 
344  if (nextVolume.c_str()!=name.c_str()) { //save hit when it exits cell
345  if (std::find(hgchitIndex_.begin(),hgchitIndex_.end(),index) == hgchitIndex_.end()) {
346  hgchitDets_.push_back(dets_[type]);
347  hgchitIndex_.push_back(index);
348  hgchitX_.push_back(hitPoint.x());
349  hgchitY_.push_back(hitPoint.y());
350  hgchitZ_.push_back(hitPoint.z());
351  }
352  }
353  } // it is right type of SD
354  } // it is in a SD
355  }//if aStep!=NULL
356 }//end update aStep
357 
358 //================================================================ End of EVENT
359 
361 
362  if (verbosity_ > 0)
363  edm::LogVerbatim("ValidHGCal") << "\n ===>>> SimG4HGCalValidation: Energy "
364  << "deposit\n at EE : " << std::setw(6)
365  << edepEE_/CLHEP::MeV << "\n at HEF: "
366  << std::setw(6) << edepHEF_/CLHEP::MeV
367  << "\n at HEB: " << std::setw(6)
368  << edepHEB_/CLHEP::MeV;
369 
370  //Fill HGC Variables
373 }
374 
375 //---------------------------------------------------
377 
378  hgchitDets_.erase(hgchitDets_.begin(),hgchitDets_.end());
379  hgchitIndex_.erase(hgchitIndex_.begin(),hgchitIndex_.end());
380  hgchitX_.erase(hgchitX_.begin(),hgchitX_.end());
381  hgchitY_.erase(hgchitY_.begin(),hgchitY_.end());
382  hgchitZ_.erase(hgchitZ_.begin(),hgchitZ_.end());
383 }
384 
std::vector< unsigned int > hgchitIndex_
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
std::vector< double > hgchitY_
T getUntrackedParameter(std::string const &, T const &) const
std::vector< int > subdet_
#define DEFINE_SIMWATCHER(type)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
std::vector< HGCalNumberingScheme * > hgcalNumbering_
void produce(edm::Event &, const edm::EventSetup &) override
SimG4HGCalValidation(const edm::ParameterSet &p)
std::vector< double > hgchitZ_
#define nullptr
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
ForwardSubdetector
std::vector< int > detTypes_
const double MeV
static uint32_t packHcalIndex(int det, int z, int depth, int eta, int phi, int lay)
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)
unsigned int layers(bool reco) const
const SimG4HGCalValidation & operator=(const SimG4HGCalValidation &)
std::vector< unsigned int > dets_
HcalNumberingFromDDD * numberingFromDDD_
void layerAnalysis(PHGCalValidInfo &)
HGCalGeometryMode::GeometryMode geomMode() const
std::vector< unsigned int > hgchitDets_
std::vector< double > hgcHEFedep_
int k[5][pyjets_maxn]
std::vector< HGCNumberingScheme * > hgcNumbering_
std::vector< std::string > names_
std::vector< double > hgcHEBedep_
Detector
Definition: DetId.h:26
std::vector< int > types_
std::vector< double > hgcEEedep_
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
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)
T get() const
Definition: EventSetup.h:68
step
int levelTop(int ind=0) const
void update(const BeginOfJob *job) override
This routine will be called when the appropriate signal arrives.
bool isValid() const
Definition: ESHandle.h:45
HcalID unitID(int det, const CLHEP::Hep3Vector &pos, int depth, int lay=-1) const
std::vector< double > hgchitX_
T const * product() const
Definition: ESHandle.h:84
Definition: vlib.h:208
def move(src, dest)
Definition: eostools.py:511