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
22 
28 
35 
39 
40 #include "G4Step.hh"
41 #include "G4Track.hh"
42 #include "G4NavigationHistory.hh"
43 #include "G4TouchableHistory.hh"
44 
45 #include "CLHEP/Units/GlobalSystemOfUnits.h"
46 #include "CLHEP/Units/GlobalPhysicalConstants.h"
47 
48 #include <cmath>
49 #include <iostream>
50 #include <iomanip>
51 #include <memory>
52 #include <vector>
53 #include <string>
54 
56  public Observer<const BeginOfJob *>,
57  public Observer<const BeginOfEvent *>,
58  public Observer<const G4Step *> {
59 
60 public:
62  virtual ~SimG4HGCalValidation();
63 
64  void produce(edm::Event&, const edm::EventSetup&);
65 
66 private:
67  SimG4HGCalValidation(const SimG4HGCalValidation&); // stop default
69 
70  void init();
71 
72  // observer classes
73  void update(const BeginOfJob * job);
74  void update(const BeginOfEvent * evt);
75  void update(const G4Step * step);
76 
77  // analysis related class
79  void clear();
80 
81 private:
82  //Keep reference to instantiate HcalNumberingFromDDD later
84 
85  //HGCal numbering scheme
86  std::vector<HGCNumberingScheme*> hgcNumbering_;
87 
88  // to read from parameter set
89  std::vector<std::string> names_;
90  std::vector<int> types_, subdet_;
92 
93  // some private members for ananlysis
94  unsigned int count_;
96  std::vector<double> hgcEEedep_, hgcHEFedep_, hgcHEBedep_;
97  std::vector<unsigned int> dets_, hgchitDets_, hgchitIndex_;
98  std::vector<double> hgchitX_, hgchitY_, hgchitZ_;
99 };
100 
102  numberingFromDDD_(0), count_(0) {
103 
104  edm::ParameterSet m_Anal = p.getParameter<edm::ParameterSet>("SimG4HGCalValidation");
105  names_ = m_Anal.getParameter<std::vector<std::string> >("Names");
106  types_ = m_Anal.getParameter<std::vector<int> >("Types");
107  labelLayer_ = m_Anal.getParameter<std::string>("LabelLayerInfo");
108 
109  produces<PHGCalValidInfo>(labelLayer_);
110 
111  edm::LogInfo("ValidHGCal") << "HGCalTestAnalysis:: Initialised as observer "
112  << "of begin events and of G4step with Parameter "
113  << "values: \n\tLabel : " << labelLayer_
114  << " and with " << names_.size() << " detectors"
115  << std::endl;
116  for (unsigned int k=0; k<names_.size(); ++k)
117  edm::LogInfo("ValidHGCal") << " [" << k << "] " << names_[k] << " Type "
118  << types_[k] << std::endl;
119 }
120 
123  for (unsigned int k=0; k<hgcNumbering_.size(); ++k)
124  if (hgcNumbering_[k]) delete hgcNumbering_[k];
125 }
126 
128 
129  std::unique_ptr<PHGCalValidInfo> productLayer(new PHGCalValidInfo);
130  layerAnalysis(*productLayer);
131  e.put(std::move(productLayer),labelLayer_);
132 }
133 
135 
136  const edm::EventSetup* es = (*job)();
137  for (unsigned int type=0; type<types_.size(); ++type) {
138  int layers(0);
139  G4String nameX = "HGCal";
140  if (types_[type] <= 1) {
141  dets_.push_back((unsigned int)(DetId::Forward));
142  if (type == 0) {
143  subdet_.push_back((int)(ForwardSubdetector::HGCEE));
144  nameX = "HGCalEESensitive";
145  } else if (type == 1) {
146  subdet_.push_back((int)(ForwardSubdetector::HGCHEF));
147  nameX = "HGCalHESiliconSensitive";
148  } else {
149  subdet_.push_back((int)(ForwardSubdetector::HGCHEB));
150  nameX = "HGCalHEScintillatorSensitive";
151  }
153  es->get<IdealGeometryRecord>().get(nameX,hdc);
154  if (hdc.isValid()) {
156  hgcNumbering_.push_back(new HGCNumberingScheme(*hdc,nameX));
157  if (m_mode == HGCalGeometryMode::Square) types_[type] = 0;
158  else types_[type] = 1;
159  layers = hdc->layers(false);
160  } else {
161  edm::LogError("ValidHGCal") << "Cannot find HGCalDDDConstants for "
162  << nameX;
163  throw cms::Exception("Unknown", "ValidHGCal") << "Cannot find HGCalDDDConstants for " << nameX << "\n";
164  }
165  } else {
166  nameX = "HcalEndcap";
167  dets_.push_back((unsigned int)(DetId::Hcal));
168  subdet_.push_back((int)(HcalSubdetector::HcalEndcap));
170  es->get<HcalSimNumberingRecord>().get(hdc);
171  if (hdc.isValid()) {
172  HcalDDDSimConstants* hcalConstants = (HcalDDDSimConstants*)(&(*hdc));
173  numberingFromDDD_ = new HcalNumberingFromDDD(hcalConstants);
174  layers = 18;
175  } else {
176  edm::LogError("ValidHGCal") << "Cannot find HcalDDDSimConstant";
177  throw cms::Exception("Unknown", "ValidHGCal") << "Cannot find HcalDDDSimConstant\n";
178  }
179  }
180  if (type == 0) {
181  for (int i=0; i<layers; ++i) hgcEEedep_.push_back(0);
182  } else if (type == 1) {
183  for (int i=0; i<layers; ++i) hgcHEFedep_.push_back(0);
184  } else {
185  for (int i=0; i<layers; ++i) hgcHEBedep_.push_back(0);
186  }
187  edm::LogInfo("ValidHGCal") << "[" << type << "]: " << nameX << " det "
188  << dets_[type] << " subdet " << subdet_[type]
189  << " with " << layers << " layers" << std::endl;
190  }
191 }
192 
193 //=================================================================== per EVENT
195 
196  int iev = (*evt)()->GetEventID();
197  edm::LogInfo("ValidHGCal") << "SimG4HGCalValidation: =====> Begin event = "
198  << iev << std::endl;
199 
200  ++count_;
201  edepEE_ = edepHEF_ = edepHEB_ = 0.;
202 
203  //HGCal variables
204  for (unsigned int i = 0; i<hgcEEedep_.size(); i++) hgcEEedep_[i] = 0.;
205  for (unsigned int i = 0; i<hgcHEFedep_.size(); i++) hgcHEFedep_[i] = 0.;
206  for (unsigned int i = 0; i<hgcHEBedep_.size(); i++) hgcHEBedep_[i] = 0.;
207 
208  //Cache reset
209  clear();
210 }
211 
212 //=================================================================== each STEP
213 void SimG4HGCalValidation::update(const G4Step * aStep) {
214 
215  if (aStep != NULL) {
216  G4VPhysicalVolume* curPV = aStep->GetPreStepPoint()->GetPhysicalVolume();
217  G4VSensitiveDetector* curSD = aStep->GetPreStepPoint()->GetSensitiveDetector();
218 
219  // Only for Sensitive detector
220  if (curSD != 0) {
221  G4String name = curPV->GetName();
222  int type(-1);
223  for (unsigned int k=0; k<names_.size(); ++k) {
224  if (name.find(names_[k].c_str()) != std::string::npos) {
225  type = k; break;
226  }
227  }
228  edm::LogInfo("ValidHGCal") << "Step in " << name << " type " << type;
229  // Right type of SD
230  if (type >= 0) {
231  //Get the 32-bit index of the hit cell
232  G4TouchableHistory* touchable = (G4TouchableHistory*)aStep->GetPreStepPoint()->GetTouchable();
233  unsigned int index(0);
234  int layer(0);
235  G4ThreeVector hitPoint = aStep->GetPreStepPoint()->GetPosition();
236  if (types_[type] <= 1) {
237  // HGCal
238  G4ThreeVector localpos = touchable->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
239  float globalZ=touchable->GetTranslation(0).z();
240  int iz( globalZ>0 ? 1 : -1);
241  int module(0), cell(0);
242  if (types_[type] == 0) {
243  layer = touchable->GetReplicaNumber(0);
244  module = touchable->GetReplicaNumber(1);
245  } else {
246  layer = touchable->GetReplicaNumber(2);
247  module = touchable->GetReplicaNumber(1);
248  cell = touchable->GetReplicaNumber(0);
249  }
250  index = hgcNumbering_[type]->getUnitID((ForwardSubdetector)(subdet_[type]), layer, module, cell, iz, localpos);
251  } else {
252  // Hcal
253  int depth = (touchable->GetReplicaNumber(0))%10 + 1;
254  int lay = (touchable->GetReplicaNumber(0)/10)%100 + 1;
255  int det = (touchable->GetReplicaNumber(1))/1000;
256  HcalNumberingFromDDD::HcalID tmp = numberingFromDDD_->unitID(det, hitPoint, depth, lay);
257  index = HcalTestNumbering::packHcalIndex(tmp.subdet,tmp.zside,tmp.depth,tmp.etaR,tmp.phis,tmp.lay);
258  layer = tmp.lay;
259  edm::LogInfo("ValidHGCal") << "HCAL: " << det << ":" << depth << ":"
260  << lay << " o/p " << tmp.subdet << ":"
261  << tmp.zside << ":" << tmp.depth << ":"
262  << tmp.etaR << ":" << tmp.phis << ":"
263  << tmp.lay << " point " << hitPoint << " "
264  << hitPoint.rho() << ":" << hitPoint.eta()
265  << ":" << hitPoint.phi();
266  }
267 
268  double edeposit = aStep->GetTotalEnergyDeposit();
269  edm::LogInfo("ValidHGCal") << "Layer " << layer << " Index "
270  << std::hex << index << std::dec
271  << " Edep " << edeposit << " hit at "
272  << hitPoint;
273  if (type == 0) {
274  edepEE_ += edeposit;
275  if (layer < (int)(hgcEEedep_.size())) hgcEEedep_[layer] += edeposit;
276  } else if (type == 1) {
277  edepHEF_ += edeposit;
278  if (layer < (int)(hgcHEFedep_.size())) hgcHEFedep_[layer] += edeposit;
279  } else {
280  edepHEB_ += edeposit;
281  if (layer < (int)(hgcHEBedep_.size())) hgcHEBedep_[layer] += edeposit;
282  }
283  G4String nextVolume("XXX");
284  if (aStep->GetTrack()->GetNextVolume()!=0)
285  nextVolume = aStep->GetTrack()->GetNextVolume()->GetName();
286 
287  if (nextVolume.c_str()!=name.c_str()) { //save hit when it exits cell
288  if (std::find(hgchitIndex_.begin(),hgchitIndex_.end(),index) == hgchitIndex_.end()) {
289  hgchitDets_.push_back(dets_[type]);
290  hgchitIndex_.push_back(index);
291  hgchitX_.push_back(hitPoint.x());
292  hgchitY_.push_back(hitPoint.y());
293  hgchitZ_.push_back(hitPoint.z());
294  }
295  }
296  } // it is right type of SD
297  } // it is in a SD
298  }//if aStep!=NULL
299 }//end update aStep
300 
301 //================================================================ End of EVENT
302 
304 
305  edm::LogInfo("ValidHGCal") << "\n ===>>> SimG4HGCalValidation: Energy deposit"
306  << "\n at EE : " << std::setw(6) << edepEE_/MeV
307  << "\n at HEF: " << std::setw(6) << edepHEF_/MeV
308  << "\n at HEB: " << std::setw(6) << edepHEB_/MeV
309  << "\n";
310 
311 
312  //Fill HGC Variables
315 }
316 
317 //---------------------------------------------------
319 
320  hgchitDets_.erase(hgchitDets_.begin(),hgchitDets_.end());
321  hgchitIndex_.erase(hgchitIndex_.begin(),hgchitIndex_.end());
322  hgchitX_.erase(hgchitX_.begin(),hgchitX_.end());
323  hgchitY_.erase(hgchitY_.begin(),hgchitY_.end());
324  hgchitZ_.erase(hgchitZ_.begin(),hgchitZ_.end());
325 }
326 
std::vector< unsigned int > hgchitIndex_
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
std::vector< double > hgchitY_
std::vector< int > subdet_
#define DEFINE_SIMWATCHER(type)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
SimG4HGCalValidation(const edm::ParameterSet &p)
std::vector< double > hgchitZ_
#define NULL
Definition: scimark2.h:8
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
ForwardSubdetector
void produce(edm::Event &, const edm::EventSetup &)
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
void update(const BeginOfJob *job)
This routine will be called when the appropriate signal arrives.
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_
const T & get() const
Definition: EventSetup.h:55
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)
step
bool isValid() const
Definition: ESHandle.h:47
HcalID unitID(int det, const CLHEP::Hep3Vector &pos, int depth, int lay=-1) const
std::vector< double > hgchitX_
Definition: vlib.h:208
def move(src, dest)
Definition: eostools.py:510