CMS 3D CMS Logo

CaloTrkProcessing.cc
Go to the documentation of this file.
5 
7 
13 
14 #include "G4EventManager.hh"
15 
16 #include "G4LogicalVolumeStore.hh"
17 #include "G4LogicalVolume.hh"
18 #include "G4Step.hh"
19 #include "G4Track.hh"
20 
21 #include "G4SystemOfUnits.hh"
22 
23 //#define EDM_ML_DEBUG
24 
26  const edm::EventSetup& es,
27  const SensitiveDetectorCatalog& clg,
28  edm::ParameterSet const& p,
29  const SimTrackManager*)
30  : SensitiveCaloDetector(name, es, clg, p), lastTrackID_(-1) {
31  //Initialise the parameter set
32  edm::ParameterSet m_p = p.getParameter<edm::ParameterSet>("CaloTrkProcessing");
33  testBeam_ = m_p.getParameter<bool>("TestBeam");
34  eMin_ = m_p.getParameter<double>("EminTrack") * MeV;
35  putHistory_ = m_p.getParameter<bool>("PutHistory");
36  doFineCalo_ = m_p.getParameter<bool>("DoFineCalo");
37  eMinFine_ = m_p.getParameter<double>("EminFineTrack") * MeV;
38  eMinFinePhoton_ = m_p.getParameter<double>("EminFinePhoton") * MeV;
39 
40  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: Initailised with TestBeam = " << testBeam_ << " Emin = " << eMin_
41  << ":" << eMinFine_ << ":" << eMinFinePhoton_ << " MeV and Flags " << putHistory_
42  << " (History), " << doFineCalo_ << " (Special Calorimeter)";
43 
44  // Get pointers to HcalDDDConstant and HcalSimulationParameters
46  es.get<HcalParametersRcd>().get(csps);
47  if (csps.isValid()) {
48  const CaloSimulationParameters* csp = csps.product();
49 #ifdef EDM_ML_DEBUG
50  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: " << csp->caloNames_.size() << " entries for caloNames:";
51  for (unsigned int i = 0; i < csp->caloNames_.size(); i++)
52  edm::LogVerbatim("CaloSim") << " (" << i << ") " << csp->caloNames_[i];
53  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: " << csp->levels_.size() << " entries for levels:";
54  for (unsigned int i = 0; i < csp->levels_.size(); i++)
55  edm::LogVerbatim("CaloSim") << " (" << i << ") " << csp->levels_[i];
56  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: " << csp->neighbours_.size() << " entries for neighbours:";
57  for (unsigned int i = 0; i < csp->neighbours_.size(); i++)
58  edm::LogVerbatim("CaloSim") << " (" << i << ") " << csp->neighbours_[i];
59  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: " << csp->insideNames_.size() << " entries for insideNames:";
60  for (unsigned int i = 0; i < csp->insideNames_.size(); i++)
61  edm::LogVerbatim("CaloSim") << " (" << i << ") " << csp->insideNames_[i];
62  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: " << csp->insideLevel_.size() << " entries for insideLevel:";
63  for (unsigned int i = 0; i < csp->insideLevel_.size(); i++)
64  edm::LogVerbatim("CaloSim") << " (" << i << ") " << csp->insideLevel_[i];
65  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: " << csp->fCaloNames_.size()
66  << " entries for fineCalorimeterNames:";
67  for (unsigned int i = 0; i < csp->fCaloNames_.size(); i++)
68  edm::LogVerbatim("CaloSim") << " (" << i << ") " << csp->fCaloNames_[i];
69  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: " << csp->fLevels_.size()
70  << " entries for fineCalorimeterLevels:";
71  for (unsigned int i = 0; i < csp->fLevels_.size(); i++)
72  edm::LogVerbatim("CaloSim") << " (" << i << ") " << csp->fLevels_[i];
73 #endif
74 
75  if (csp->caloNames_.size() < csp->neighbours_.size()) {
76  edm::LogError("CaloSim") << "CaloTrkProcessing: # of Calorimeter bins " << csp->caloNames_.size()
77  << " does not match with " << csp->neighbours_.size() << " ==> illegal ";
78  throw cms::Exception("Unknown", "CaloTrkProcessing")
79  << "Calorimeter array size does not match with size of neighbours\n";
80  }
81 
82  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
83  std::vector<G4LogicalVolume*>::const_iterator lvcite;
84  int istart = 0;
85  for (unsigned int i = 0; i < csp->caloNames_.size(); i++) {
86  G4LogicalVolume* lv = nullptr;
87  G4String name = static_cast<G4String>(csp->caloNames_[i]);
88  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
89  if ((*lvcite)->GetName() == name) {
90  lv = (*lvcite);
91  break;
92  }
93  }
94  if (lv != nullptr) {
96  detector.name = name;
97  detector.lv = lv;
98  detector.level = csp->levels_[i];
99  if (istart + csp->neighbours_[i] > static_cast<int>(csp->insideNames_.size())) {
100  edm::LogError("CaloSim") << "CaloTrkProcessing: # of InsideNames bins " << csp->insideNames_.size()
101  << " too few compaerd to " << istart + csp->neighbours_[i]
102  << " requested ==> illegal ";
103  throw cms::Exception("Unknown", "CaloTrkProcessing")
104  << "InsideNames array size does not match with list of neighbours\n";
105  }
106  std::vector<std::string> inside;
107  std::vector<G4LogicalVolume*> insideLV;
108  std::vector<int> insideLevels;
109  for (int k = 0; k < csp->neighbours_[i]; k++) {
110  lv = nullptr;
111  name = static_cast<G4String>(csp->insideNames_[istart + k]);
112  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
113  if ((*lvcite)->GetName() == name) {
114  lv = (*lvcite);
115  break;
116  }
117  }
118  inside.push_back(name);
119  insideLV.push_back(lv);
120  insideLevels.push_back(csp->insideLevel_[istart + k]);
121  }
122  detector.fromDets = inside;
123  detector.fromDetL = insideLV;
124  detector.fromLevels = insideLevels;
125  detectors_.emplace_back(detector);
126  }
127  istart += csp->neighbours_[i];
128  }
129 
130  for (unsigned int i = 0; i < csp->fCaloNames_.size(); i++) {
131  G4LogicalVolume* lv = nullptr;
132  G4String name = static_cast<G4String>(csp->fCaloNames_[i]);
133  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
134  if ((*lvcite)->GetName() == name) {
135  lv = (*lvcite);
136  break;
137  }
138  }
139  if (lv != nullptr) {
141  detector.name = name;
142  detector.lv = lv;
143  detector.level = csp->fLevels_[i];
144  detector.fromDets.clear();
145  detector.fromDetL.clear();
146  detector.fromLevels.clear();
147  fineDetectors_.emplace_back(detector);
148  }
149  }
150  } else {
151  edm::LogError("HcalSim") << "CaloTrkProcessing: Cannot find CaloSimulationParameters";
152  throw cms::Exception("Unknown", "CaloSD") << "Cannot find CaloSimulationParameters\n";
153  }
154 
155  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: with " << detectors_.size() << " calorimetric volumes";
156  for (unsigned int i = 0; i < detectors_.size(); i++) {
157  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: Calorimeter volume " << i << " " << detectors_[i].name << " LV "
158  << detectors_[i].lv << " at level " << detectors_[i].level << " with "
159  << detectors_[i].fromDets.size() << " neighbours";
160  for (unsigned int k = 0; k < detectors_[i].fromDets.size(); k++)
161  edm::LogVerbatim("CaloSim") << " Element " << k << " " << detectors_[i].fromDets[k] << " LV "
162  << detectors_[i].fromDetL[k] << " at level " << detectors_[i].fromLevels[k];
163  }
164 
165  doFineCalo_ = !(fineDetectors_.empty());
166  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: with " << fineDetectors_.size() << " special calorimetric volumes";
167  for (unsigned int i = 0; i < detectors_.size(); i++)
168  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: Calorimeter volume " << i << " " << detectors_[i].name << " LV "
169  << detectors_[i].lv << " at level " << detectors_[i].level;
170 }
171 
173 
175 
176 void CaloTrkProcessing::update(const G4Step* aStep) {
177  // define if you are at the surface of CALO
178 
179  G4Track* theTrack = aStep->GetTrack();
180  int id = theTrack->GetTrackID();
181 
182  TrackInformation* trkInfo = dynamic_cast<TrackInformation*>(theTrack->GetUserInformation());
183 
184  if (trkInfo == nullptr) {
185  edm::LogError("CaloSim") << "CaloTrkProcessing: No trk info !!!! abort ";
186  throw cms::Exception("Unknown", "CaloTrkProcessing") << "cannot get trkInfo for Track " << id << "\n";
187  }
188 
189  if (testBeam_) {
190  if (trkInfo->getIDonCaloSurface() == 0) {
191 #ifdef EDM_ML_DEBUG
192  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing set IDonCaloSurface to " << id << " at step Number "
193  << theTrack->GetCurrentStepNumber();
194 #endif
195  trkInfo->setIDonCaloSurface(id, 0, 0, theTrack->GetDefinition()->GetPDGEncoding(), theTrack->GetMomentum().mag());
196  lastTrackID_ = id;
197  if (theTrack->GetKineticEnergy() / MeV > eMin_)
198  trkInfo->putInHistory();
199  }
200  } else {
201  if (putHistory_) {
202  trkInfo->putInHistory();
203  // trkInfo->setAncestor();
204  }
205 #ifdef EDM_ML_DEBUG
206  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing Entered for " << id << " at stepNumber "
207  << theTrack->GetCurrentStepNumber() << " IDonCaloSur.. "
208  << trkInfo->getIDonCaloSurface() << " CaloCheck " << trkInfo->caloIDChecked();
209 #endif
210  if (trkInfo->getIDonCaloSurface() != 0) {
211  if (trkInfo->caloIDChecked() == false) {
212  G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
213  const G4VTouchable* post_touch = postStepPoint->GetTouchable();
214 
215  if (isItInside(post_touch, trkInfo->getIDCaloVolume(), trkInfo->getIDLastVolume()) > 0) {
216  trkInfo->setIDonCaloSurface(0, -1, -1, 0, 0);
217  } else {
218  trkInfo->setCaloIDChecked(true);
219  }
220  }
221  } else {
222  G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
223  const G4VTouchable* post_touch = postStepPoint->GetTouchable();
224  int ical = isItCalo(post_touch, detectors_);
225  if (ical >= 0) {
226  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
227  const G4VTouchable* pre_touch = preStepPoint->GetTouchable();
228  int inside = isItInside(pre_touch, ical, -1);
229  if (inside >= 0 || (theTrack->GetCurrentStepNumber() == 1)) {
230  trkInfo->setIDonCaloSurface(
231  id, ical, inside, theTrack->GetDefinition()->GetPDGEncoding(), theTrack->GetMomentum().mag());
232  trkInfo->setCaloIDChecked(true);
233  lastTrackID_ = id;
234  if (theTrack->GetKineticEnergy() / MeV > eMin_)
235  trkInfo->putInHistory();
236 #ifdef EDM_ML_DEBUG
237  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: set ID on Calo " << ical << " surface (Inside " << inside
238  << ") to " << id << " of a Track with Kinetic Energy "
239  << theTrack->GetKineticEnergy() / MeV << " MeV";
240 #endif
241  }
242  }
243  }
244  }
245  if (doFineCalo_ && (!trkInfo->isInHistory())) {
246  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
247  if (isItCalo(touch, fineDetectors_) >= 0) {
248  int pdg = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
249  double cut = (pdg == 22) ? eMinFinePhoton_ : eMinFine_;
250  if (aStep->GetTrack()->GetKineticEnergy() / MeV > cut)
251  trkInfo->putInHistory();
252 #ifdef EDM_ML_DEBUG
253  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: the track with PDGID " << pdg << " and kinetic energy "
254  << aStep->GetTrack()->GetKineticEnergy() / MeV << " is tested against " << cut
255  << " to be put in history";
256 #endif
257  }
258  }
259 }
260 
261 int CaloTrkProcessing::isItCalo(const G4VTouchable* touch, const std::vector<Detector>& detectors) {
262  int lastLevel = -1;
263  G4LogicalVolume* lv = nullptr;
264  for (unsigned int it = 0; it < detectors.size(); it++) {
265  if (lastLevel != detectors[it].level) {
266  lastLevel = detectors[it].level;
267  lv = detLV(touch, lastLevel);
268 #ifdef EDM_ML_DEBUG
269  std::string name1 = "Unknown";
270  if (lv != 0)
271  name1 = lv->GetName();
272  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: volume " << name1 << " at Level " << lastLevel;
273  int levels = detLevels(touch);
274  if (levels > 0) {
275  G4String name2[20];
276  int copyno2[20];
277  detectorLevel(touch, levels, copyno2, name2);
278  for (int i2 = 0; i2 < levels; i2++)
279  edm::LogVerbatim("CaloSim") << " " << i2 << " " << name2[i2] << " " << copyno2[i2];
280  }
281 #endif
282  }
283  bool ok = (lv == detectors[it].lv);
284  if (ok)
285  return it;
286  }
287  return -1;
288 }
289 
290 int CaloTrkProcessing::isItInside(const G4VTouchable* touch, int idcal, int idin) {
291  int lastLevel = -1;
292  G4LogicalVolume* lv = nullptr;
293  int id1, id2;
294  if (idcal < 0) {
295  id1 = 0;
296  id2 = static_cast<int>(detectors_.size());
297  } else {
298  id1 = idcal;
299  id2 = id1 + 1;
300  }
301  for (int it1 = id1; it1 < id2; it1++) {
302  if (idin < 0) {
303  for (unsigned int it2 = 0; it2 < detectors_[it1].fromDets.size(); it2++) {
304  if (lastLevel != detectors_[it1].fromLevels[it2]) {
305  lastLevel = detectors_[it1].fromLevels[it2];
306  lv = detLV(touch, lastLevel);
307 #ifdef EDM_ML_DEBUG
308  std::string name1 = "Unknown";
309  if (lv != 0)
310  name1 = lv->GetName();
311  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: volume " << name1 << " at Level " << lastLevel;
312  int levels = detLevels(touch);
313  if (levels > 0) {
314  G4String name2[20];
315  int copyno2[20];
316  detectorLevel(touch, levels, copyno2, name2);
317  for (int i2 = 0; i2 < levels; i2++)
318  edm::LogVerbatim("CaloSim") << " " << i2 << " " << name2[i2] << " " << copyno2[i2];
319  }
320 #endif
321  }
322  bool ok = (lv == detectors_[it1].fromDetL[it2]);
323  if (ok)
324  return it2;
325  }
326  } else {
327  lastLevel = detectors_[it1].fromLevels[idin];
328  lv = detLV(touch, lastLevel);
329 #ifdef EDM_ML_DEBUG
330  std::string name1 = "Unknown";
331  if (lv != 0)
332  name1 = lv->GetName();
333  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: volume " << name1 << " at Level " << lastLevel;
334  int levels = detLevels(touch);
335  if (levels > 0) {
336  G4String name2[20];
337  int copyno2[20];
338  detectorLevel(touch, levels, copyno2, name2);
339  for (int i2 = 0; i2 < levels; i2++)
340  edm::LogVerbatim("CaloSim") << " " << i2 << " " << name2[i2] << " " << copyno2[i2];
341  }
342 #endif
343  bool ok = (lv == detectors_[it1].fromDetL[idin]);
344  if (ok)
345  return idin;
346  }
347  }
348  return -1;
349 }
350 
351 int CaloTrkProcessing::detLevels(const G4VTouchable* touch) const {
352  //Return number of levels
353  if (touch)
354  return ((touch->GetHistoryDepth()) + 1);
355  else
356  return 0;
357 }
358 
359 G4LogicalVolume* CaloTrkProcessing::detLV(const G4VTouchable* touch, int currentlevel) const {
360  G4LogicalVolume* lv = nullptr;
361  if (touch) {
362  int level = ((touch->GetHistoryDepth()) + 1);
363  if (level > 0 && level >= currentlevel) {
364  int ii = level - currentlevel;
365  lv = touch->GetVolume(ii)->GetLogicalVolume();
366  return lv;
367  }
368  }
369  return lv;
370 }
371 
372 void CaloTrkProcessing::detectorLevel(const G4VTouchable* touch, int& level, int* copyno, G4String* name) const {
373  static const std::string unknown("Unknown");
374  //Get name and copy numbers
375  if (level > 0) {
376  for (int ii = 0; ii < level; ii++) {
377  int i = level - ii - 1;
378  G4VPhysicalVolume* pv = touch->GetVolume(i);
379  if (pv != nullptr)
380  name[ii] = pv->GetName();
381  else
382  name[ii] = unknown;
383  copyno[ii] = touch->GetReplicaNumber(i);
384  }
385  }
386 }
void update(const BeginOfEvent *evt) override
This routine will be called when the appropriate signal arrives.
T getParameter(std::string const &) const
G4LogicalVolume * detLV(const G4VTouchable *, int) const
std::vector< std::string > insideNames_
int getIDonCaloSurface() const
int detLevels(const G4VTouchable *) const
void detectorLevel(const G4VTouchable *, int &, int *, G4String *) const
const double MeV
int getIDCaloVolume() const
std::vector< std::string > fCaloNames_
int isItInside(const G4VTouchable *, int, int)
CaloTrkProcessing(const std::string &aSDname, const edm::EventSetup &es, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p, const SimTrackManager *)
def detectors(dt=True, csc=True, me42=False, chambers=True, superlayers=False, layers=False)
std::vector< int > fromLevels
def pv(vc)
Definition: MetAnalyzer.py:7
std::vector< Detector > fineDetectors_
std::vector< G4LogicalVolume * > fromDetL
int getIDLastVolume() const
std::vector< std::string > fromDets
bool isInHistory() const
ii
Definition: cuy.py:590
int isItCalo(const G4VTouchable *, const std::vector< Detector > &)
~CaloTrkProcessing() override
bool caloIDChecked() const
void setCaloIDChecked(bool f)
void setIDonCaloSurface(int id, int ical, int last, int pdgID, double p)
T get() const
Definition: EventSetup.h:73
bool isValid() const
Definition: ESHandle.h:44
T const * product() const
Definition: ESHandle.h:86
std::vector< std::string > caloNames_
std::vector< Detector > detectors_