CMS 3D CMS Logo

CaloTrkProcessing.cc
Go to the documentation of this file.
5 
7 
9 
10 #include "G4EventManager.hh"
11 #include "G4LogicalVolumeStore.hh"
12 #include "G4LogicalVolume.hh"
13 #include "G4Step.hh"
14 #include "G4Track.hh"
15 #include "G4SystemOfUnits.hh"
16 #include "DD4hep/Filter.h"
17 
18 #include <sstream>
19 //#define EDM_ML_DEBUG
20 
22  const CaloSimulationParameters& csps,
23  const SensitiveDetectorCatalog& clg,
24  bool testBeam,
25  double eMin,
26  bool putHistory,
27  bool doFineCalo,
28  double eMinFine,
29  int addlevel,
30  const std::vector<std::string>& fineNames,
31  const std::vector<int>& fineLevels,
32  const std::vector<int>& useFines,
33  const SimTrackManager*)
35  testBeam_(testBeam),
36  eMin_(eMin),
37  putHistory_(putHistory),
38  doFineCalo_(doFineCalo),
39  eMinFine_(eMinFine),
40  addlevel_(addlevel),
41  lastTrackID_(-1) {
42  //Initialise the parameter set
43 
44  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: Initialised with TestBeam = " << testBeam_ << " Emin = " << eMin_
45  << " Flags " << putHistory_ << " (History), " << doFineCalo_ << " (Special Calorimeter)";
46  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: Have a possibility of " << fineNames.size()
47  << " fine calorimeters of which " << useFines.size() << " are selected";
48  for (unsigned int k = 0; k < fineNames.size(); ++k)
49  edm::LogVerbatim("CaloSim") << "[" << k << "] " << fineNames[k] << " at " << fineLevels[k];
50  std::ostringstream st1;
51  for (unsigned int k = 0; k < useFines.size(); ++k)
52  st1 << " [" << k << "] " << useFines[k] << ":" << fineNames[useFines[k]];
53  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing used calorimeters" << st1.str();
54 
55  // Debug prints
56  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: " << csps.caloNames_.size() << " entries for caloNames:";
57  for (unsigned int i = 0; i < csps.caloNames_.size(); i++)
58  edm::LogVerbatim("CaloSim") << " (" << i << ") " << csps.caloNames_[i];
59  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: " << csps.levels_.size() << " entries for levels:";
60  for (unsigned int i = 0; i < csps.levels_.size(); i++)
61  edm::LogVerbatim("CaloSim") << " (" << i << ") " << (csps.levels_[i] + addlevel_);
62  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: " << csps.neighbours_.size() << " entries for neighbours:";
63  for (unsigned int i = 0; i < csps.neighbours_.size(); i++)
64  edm::LogVerbatim("CaloSim") << " (" << i << ") " << csps.neighbours_[i];
65  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: " << csps.insideNames_.size() << " entries for insideNames:";
66  for (unsigned int i = 0; i < csps.insideNames_.size(); i++)
67  edm::LogVerbatim("CaloSim") << " (" << i << ") " << csps.insideNames_[i];
68  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: " << csps.insideLevel_.size() << " entries for insideLevel:";
69  for (unsigned int i = 0; i < csps.insideLevel_.size(); i++)
70  edm::LogVerbatim("CaloSim") << " (" << i << ") " << (csps.insideLevel_[i] + addlevel_);
71 
72  if (csps.caloNames_.size() < csps.neighbours_.size()) {
73  edm::LogError("CaloSim") << "CaloTrkProcessing: # of Calorimeter bins " << csps.caloNames_.size()
74  << " does not match with " << csps.neighbours_.size() << " ==> illegal ";
75  throw cms::Exception("Unknown", "CaloTrkProcessing")
76  << "Calorimeter array size does not match with size of neighbours\n";
77  }
78 
79  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
80  std::vector<G4LogicalVolume*>::const_iterator lvcite;
81  int istart = 0;
82  for (unsigned int i = 0; i < csps.caloNames_.size(); i++) {
83  G4LogicalVolume* lv = nullptr;
84  G4String name(csps.caloNames_[i]);
85  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
86  G4String namx(static_cast<std::string>(dd4hep::dd::noNamespace((*lvcite)->GetName())));
87  if (namx == name) {
88  lv = (*lvcite);
89  break;
90  }
91  }
92  if (lv != nullptr) {
94  detector.name = name;
95  detector.lv = lv;
96  detector.level = (csps.levels_[i] + addlevel_);
97  if (istart + csps.neighbours_[i] > static_cast<int>(csps.insideNames_.size())) {
98  edm::LogError("CaloSim") << "CaloTrkProcessing: # of InsideNames bins " << csps.insideNames_.size()
99  << " too few compaerd to " << istart + csps.neighbours_[i]
100  << " requested ==> illegal ";
101  throw cms::Exception("Unknown", "CaloTrkProcessing")
102  << "InsideNames array size does not match with list of neighbours\n";
103  }
104  std::vector<std::string> inside;
105  std::vector<G4LogicalVolume*> insideLV;
106  std::vector<int> insideLevels;
107  for (int k = 0; k < csps.neighbours_[i]; k++) {
108  lv = nullptr;
109  name = static_cast<G4String>(csps.insideNames_[istart + k]);
110  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
111  G4String namx(static_cast<std::string>(dd4hep::dd::noNamespace((*lvcite)->GetName())));
112  if (namx == name) {
113  lv = (*lvcite);
114  break;
115  }
116  }
117  inside.push_back(name);
118  insideLV.push_back(lv);
119  insideLevels.push_back(csps.insideLevel_[istart + k] + addlevel_);
120  }
121  detector.fromDets = inside;
122  detector.fromDetL = insideLV;
123  detector.fromLevels = insideLevels;
124  detectors_.emplace_back(detector);
125  }
126  istart += csps.neighbours_[i];
127  }
128 
129  for (unsigned int i = 0; i < useFines.size(); i++) {
130  G4LogicalVolume* lv = nullptr;
131  G4String name = static_cast<G4String>(fineNames[useFines[i]]);
132  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
133  G4String namx(static_cast<std::string>(dd4hep::dd::noNamespace((*lvcite)->GetName())));
134  if (namx == name) {
135  lv = (*lvcite);
136  break;
137  }
138  }
139  if (lv != nullptr) {
141  detector.name = name;
142  detector.lv = lv;
143  detector.level = fineLevels[useFines[i]];
144  detector.fromDets.clear();
145  detector.fromDetL.clear();
146  detector.fromLevels.clear();
147  fineDetectors_.emplace_back(detector);
148  }
149  }
150 
151  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: with " << detectors_.size() << " calorimetric volumes";
152  for (unsigned int i = 0; i < detectors_.size(); i++) {
153  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: Calorimeter volume " << i << " " << detectors_[i].name << " LV "
154  << detectors_[i].lv << " at level " << detectors_[i].level << " with "
155  << detectors_[i].fromDets.size() << " neighbours";
156  for (unsigned int k = 0; k < detectors_[i].fromDets.size(); k++)
157  edm::LogVerbatim("CaloSim") << " Element " << k << " " << detectors_[i].fromDets[k] << " LV "
158  << detectors_[i].fromDetL[k] << " at level " << detectors_[i].fromLevels[k];
159  }
160 
161  doFineCalo_ = doFineCalo_ && !(fineDetectors_.empty());
162  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: with " << fineDetectors_.size() << " special calorimetric volumes";
163  for (unsigned int i = 0; i < detectors_.size(); i++)
164  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: Calorimeter volume " << i << " " << detectors_[i].name << " LV "
165  << detectors_[i].lv << " at level " << detectors_[i].level;
166 }
167 
169 
171 
172 void CaloTrkProcessing::update(const G4Step* aStep) {
173  // define if you are at the surface of CALO
174 
175  G4Track* theTrack = aStep->GetTrack();
176  int id = theTrack->GetTrackID();
177 
178  TrackInformation* trkInfo = dynamic_cast<TrackInformation*>(theTrack->GetUserInformation());
179 
180  if (trkInfo == nullptr) {
181  edm::LogError("CaloSim") << "CaloTrkProcessing: No trk info !!!! abort ";
182  throw cms::Exception("Unknown", "CaloTrkProcessing") << "cannot get trkInfo for Track " << id << "\n";
183  }
184 
185  if (doFineCalo_) {
186  // Boundary-crossing logic
187  int prestepLV = isItCalo(aStep->GetPreStepPoint()->GetTouchable(), fineDetectors_);
188  int poststepLV = isItCalo(aStep->GetPostStepPoint()->GetTouchable(), fineDetectors_);
189  if (prestepLV < 0 && poststepLV >= 0
190  // Allow back-scattering and filter it out later; ensure consistency during the SIM step
191  // && std::abs(theTrack->GetStep()->GetPreStepPoint()->GetPosition().z()) < std::abs(theTrack->GetPosition().z())
192  ) {
193 #ifdef EDM_ML_DEBUG
194  edm::LogVerbatim("DoFineCalo") << "Entered fine volume " << poststepLV << ":"
195  << " Track " << id << " pdgid=" << theTrack->GetDefinition()->GetPDGEncoding()
196  << " prestepLV=" << prestepLV << " poststepLV=" << poststepLV
197  << " GetKineticEnergy[GeV]=" << theTrack->GetKineticEnergy() / CLHEP::GeV
198  << " GetVertexKineticEnergy[GeV]="
199  << theTrack->GetVertexKineticEnergy() / CLHEP::GeV << " prestepPosition[cm]=("
200  << theTrack->GetStep()->GetPreStepPoint()->GetPosition().x() / CLHEP::cm << ","
201  << theTrack->GetStep()->GetPreStepPoint()->GetPosition().y() / CLHEP::cm << ","
202  << theTrack->GetStep()->GetPreStepPoint()->GetPosition().z() / CLHEP::cm << ")"
203  << " poststepPosition[cm]=("
204  << theTrack->GetStep()->GetPostStepPoint()->GetPosition().x() / CLHEP::cm << ","
205  << theTrack->GetStep()->GetPostStepPoint()->GetPosition().y() / CLHEP::cm << ","
206  << theTrack->GetStep()->GetPostStepPoint()->GetPosition().z() / CLHEP::cm << ")"
207  << " position[cm]=(" << theTrack->GetPosition().x() / CLHEP::cm << ","
208  << theTrack->GetPosition().y() / CLHEP::cm << ","
209  << theTrack->GetPosition().z() / CLHEP::cm << ")"
210  << " vertex_position[cm]=(" << theTrack->GetVertexPosition().x() / CLHEP::cm << ","
211  << theTrack->GetVertexPosition().y() / CLHEP::cm << ","
212  << theTrack->GetVertexPosition().z() / CLHEP::cm << ")";
213 #endif
214  trkInfo->setCrossedBoundary(theTrack);
215  }
216 #ifdef EDM_ML_DEBUG
217  else if (prestepLV >= 0 && poststepLV < 0) {
218  edm::LogVerbatim("DoFineCalo") << "Exited fine volume " << prestepLV << ":"
219  << " Track " << id
220  << " GetKineticEnergy[GeV]=" << theTrack->GetKineticEnergy() / CLHEP::GeV
221  << " GetVertexKineticEnergy[GeV]="
222  << theTrack->GetVertexKineticEnergy() / CLHEP::GeV << " prestepPosition[cm]=("
223  << theTrack->GetStep()->GetPreStepPoint()->GetPosition().x() / CLHEP::cm << ","
224  << theTrack->GetStep()->GetPreStepPoint()->GetPosition().y() / CLHEP::cm << ","
225  << theTrack->GetStep()->GetPreStepPoint()->GetPosition().z() / CLHEP::cm << ")"
226  << " poststepPosition[cm]=("
227  << theTrack->GetStep()->GetPostStepPoint()->GetPosition().x() / CLHEP::cm << ","
228  << theTrack->GetStep()->GetPostStepPoint()->GetPosition().y() / CLHEP::cm << ","
229  << theTrack->GetStep()->GetPostStepPoint()->GetPosition().z() / CLHEP::cm << ")";
230  }
231 #endif
232  }
233 
234  if (testBeam_) {
235  if (trkInfo->getIDonCaloSurface() == 0) {
236 #ifdef EDM_ML_DEBUG
237  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing set IDonCaloSurface to " << id << " at step Number "
238  << theTrack->GetCurrentStepNumber();
239 #endif
240  trkInfo->setIDonCaloSurface(id, 0, 0, theTrack->GetDefinition()->GetPDGEncoding(), theTrack->GetMomentum().mag());
241  lastTrackID_ = id;
242  if (theTrack->GetKineticEnergy() / CLHEP::MeV > eMin_)
243  trkInfo->putInHistory();
244  }
245  } else {
246  if (putHistory_) {
247  trkInfo->putInHistory();
248  // trkInfo->setAncestor();
249  }
250 #ifdef EDM_ML_DEBUG
251  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing Entered for " << id << " at stepNumber "
252  << theTrack->GetCurrentStepNumber() << " IDonCaloSur.. "
253  << trkInfo->getIDonCaloSurface() << " CaloCheck " << trkInfo->caloIDChecked();
254 #endif
255  if (trkInfo->getIDonCaloSurface() != 0) {
256  if (trkInfo->caloIDChecked() == false) {
257  G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
258  const G4VTouchable* post_touch = postStepPoint->GetTouchable();
259 
260  if (isItInside(post_touch, trkInfo->getIDCaloVolume(), trkInfo->getIDLastVolume()) > 0) {
261  trkInfo->setIDonCaloSurface(0, -1, -1, 0, 0);
262  } else {
263  trkInfo->setCaloIDChecked(true);
264  }
265  }
266  } else {
267  G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
268  const G4VTouchable* post_touch = postStepPoint->GetTouchable();
269  int ical = isItCalo(post_touch, detectors_);
270  if (ical >= 0) {
271  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
272  const G4VTouchable* pre_touch = preStepPoint->GetTouchable();
273  int inside = isItInside(pre_touch, ical, -1);
274  if (inside >= 0 || (theTrack->GetCurrentStepNumber() == 1)) {
275  trkInfo->setIDonCaloSurface(
276  id, ical, inside, theTrack->GetDefinition()->GetPDGEncoding(), theTrack->GetMomentum().mag());
277  trkInfo->setCaloIDChecked(true);
278  trkInfo->setCrossedBoundary(theTrack);
279  lastTrackID_ = id;
280  if (theTrack->GetKineticEnergy() / CLHEP::MeV > eMin_)
281  trkInfo->putInHistory();
282 #ifdef EDM_ML_DEBUG
283  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: set ID on Calo " << ical << " surface (Inside " << inside
284  << ") to " << id << " of a Track with Kinetic Energy "
285  << theTrack->GetKineticEnergy() / CLHEP::MeV << " MeV";
286 #endif
287  }
288  }
289  }
290  }
291 }
292 
293 int CaloTrkProcessing::isItCalo(const G4VTouchable* touch, const std::vector<Detector>& detectors) {
294  int lastLevel = -1;
295  G4LogicalVolume* lv = nullptr;
296  for (unsigned int it = 0; it < detectors.size(); it++) {
297  if (lastLevel != detectors[it].level) {
298  lastLevel = detectors[it].level;
299  lv = detLV(touch, lastLevel);
300 #ifdef EDM_ML_DEBUG
301  std::string name1 = "Unknown";
302  if (lv != 0)
303  name1 = lv->GetName();
304  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: volume " << name1 << " at Level " << lastLevel;
305  int levels = detLevels(touch);
306  if (levels > 0) {
307  G4String name2[20];
308  int copyno2[20];
309  detectorLevel(touch, levels, copyno2, name2);
310  for (int i2 = 0; i2 < levels; i2++)
311  edm::LogVerbatim("CaloSim") << " " << i2 << " " << name2[i2] << " " << copyno2[i2];
312  }
313 #endif
314  }
315  bool ok = (lv == detectors[it].lv);
316  if (ok)
317  return it;
318  }
319  return -1;
320 }
321 
322 int CaloTrkProcessing::isItInside(const G4VTouchable* touch, int idcal, int idin) {
323  int lastLevel = -1;
324  G4LogicalVolume* lv = nullptr;
325  int id1, id2;
326  if (idcal < 0) {
327  id1 = 0;
328  id2 = static_cast<int>(detectors_.size());
329  } else {
330  id1 = idcal;
331  id2 = id1 + 1;
332  }
333  for (int it1 = id1; it1 < id2; it1++) {
334  if (idin < 0) {
335  for (unsigned int it2 = 0; it2 < detectors_[it1].fromDets.size(); it2++) {
336  if (lastLevel != detectors_[it1].fromLevels[it2]) {
337  lastLevel = detectors_[it1].fromLevels[it2];
338  lv = detLV(touch, lastLevel);
339 #ifdef EDM_ML_DEBUG
340  std::string name1 = "Unknown";
341  if (lv != 0)
342  name1 = lv->GetName();
343  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: volume " << name1 << " at Level " << lastLevel;
344  int levels = detLevels(touch);
345  if (levels > 0) {
346  G4String name2[20];
347  int copyno2[20];
348  detectorLevel(touch, levels, copyno2, name2);
349  for (int i2 = 0; i2 < levels; i2++)
350  edm::LogVerbatim("CaloSim") << " " << i2 << " " << name2[i2] << " " << copyno2[i2];
351  }
352 #endif
353  }
354  bool ok = (lv == detectors_[it1].fromDetL[it2]);
355  if (ok)
356  return it2;
357  }
358  } else {
359  lastLevel = detectors_[it1].fromLevels[idin];
360  lv = detLV(touch, lastLevel);
361 #ifdef EDM_ML_DEBUG
362  std::string name1 = "Unknown";
363  if (lv != 0)
364  name1 = lv->GetName();
365  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: volume " << name1 << " at Level " << lastLevel;
366  int levels = detLevels(touch);
367  if (levels > 0) {
368  G4String name2[20];
369  int copyno2[20];
370  detectorLevel(touch, levels, copyno2, name2);
371  for (int i2 = 0; i2 < levels; i2++)
372  edm::LogVerbatim("CaloSim") << " " << i2 << " " << name2[i2] << " " << copyno2[i2];
373  }
374 #endif
375  bool ok = (lv == detectors_[it1].fromDetL[idin]);
376  if (ok)
377  return idin;
378  }
379  }
380  return -1;
381 }
382 
383 int CaloTrkProcessing::detLevels(const G4VTouchable* touch) const {
384  //Return number of levels
385  if (touch)
386  return ((touch->GetHistoryDepth()) + 1);
387  else
388  return 0;
389 }
390 
391 G4LogicalVolume* CaloTrkProcessing::detLV(const G4VTouchable* touch, int currentlevel) const {
392  G4LogicalVolume* lv = nullptr;
393  if (touch) {
394  int level = ((touch->GetHistoryDepth()) + 1);
395  if (level > 0 && level >= currentlevel) {
396  int ii = level - currentlevel;
397  lv = touch->GetVolume(ii)->GetLogicalVolume();
398  return lv;
399  }
400  }
401  return lv;
402 }
403 
404 void CaloTrkProcessing::detectorLevel(const G4VTouchable* touch, int& level, int* copyno, G4String* name) const {
405  static const std::string unknown("Unknown");
406  //Get name and copy numbers
407  if (level > 0) {
408  for (int ii = 0; ii < level; ii++) {
409  int i = level - ii - 1;
410  G4VPhysicalVolume* pv = touch->GetVolume(i);
411  if (pv != nullptr)
412  name[ii] = pv->GetName();
413  else
414  name[ii] = unknown;
415  copyno[ii] = touch->GetReplicaNumber(i);
416  }
417  }
418 #ifdef EDM_ML_DEBUG
419  edm::LogVerbatim("CaloSimX") << "CaloTrkProcessing::detectorLevel "
420  << " with " << level << ":" << detLevels(touch) << " levels";
421  for (int ii = 0; ii < level; ii++)
422  edm::LogVerbatim("CaloSimX") << "[" << ii << "] " << name[ii] << ":" << copyno[ii];
423 #endif
424 }
personalPlayback.level
level
Definition: personalPlayback.py:22
MuonGeometrySanityCheck_cfi.detectors
def detectors(dt=True, csc=True, me42=False, chambers=True, superlayers=False, layers=False)
Definition: MuonGeometrySanityCheck_cfi.py:13
testProducerWithPsetDescEmpty_cfi.i2
i2
Definition: testProducerWithPsetDescEmpty_cfi.py:46
SimTrackManager
Definition: SimTrackManager.h:35
mps_fire.i
i
Definition: mps_fire.py:428
CaloSimulationParameters
Definition: CaloSimulationParameters.h:6
TrackInformation::getIDonCaloSurface
int getIDonCaloSurface() const
Definition: TrackInformation.h:43
TrackInformation::getIDLastVolume
int getIDLastVolume() const
Definition: TrackInformation.h:52
CaloSimulationParameters::neighbours_
std::vector< int > neighbours_
Definition: CaloSimulationParameters.h:13
CaloSimulationParameters::insideLevel_
std::vector< int > insideLevel_
Definition: CaloSimulationParameters.h:15
CaloTrkProcessing::detectorLevel
void detectorLevel(const G4VTouchable *, int &, int *, G4String *) const
Definition: CaloTrkProcessing.cc:404
globals_cff.id1
id1
Definition: globals_cff.py:33
TrackInformation::setIDonCaloSurface
void setIDonCaloSurface(int id, int ical, int last, int pdgID, double p)
Definition: TrackInformation.h:44
CaloSimulationParameters::levels_
std::vector< int > levels_
Definition: CaloSimulationParameters.h:12
MeV
const double MeV
TrackInformation::putInHistory
void putInHistory()
Definition: TrackInformation.h:37
TrackInformation::setCaloIDChecked
void setCaloIDChecked(bool f)
Definition: TrackInformation.h:54
convertSQLiteXML.ok
bool ok
Definition: convertSQLiteXML.py:98
CaloTrkProcessing.h
CaloTrkProcessing::lastTrackID_
int lastTrackID_
Definition: CaloTrkProcessing.h:79
CaloTrkProcessing::doFineCalo_
bool doFineCalo_
Definition: CaloTrkProcessing.h:76
CaloTrkProcessing::isItInside
int isItInside(const G4VTouchable *, int, int)
Definition: CaloTrkProcessing.cc:322
SensitiveCaloDetector
Definition: SensitiveCaloDetector.h:10
CaloSimulationParameters::caloNames_
std::vector< std::string > caloNames_
Definition: CaloSimulationParameters.h:11
CaloTrkProcessing::detLevels
int detLevels(const G4VTouchable *) const
Definition: CaloTrkProcessing.cc:383
CaloTrkProcessing::detLV
G4LogicalVolume * detLV(const G4VTouchable *, int) const
Definition: CaloTrkProcessing.cc:391
photonAnalyzer_cfi.eMin
eMin
Definition: photonAnalyzer_cfi.py:50
CaloTrkProcessing::update
void update(const BeginOfEvent *evt) override
This routine will be called when the appropriate signal arrives.
Definition: CaloTrkProcessing.cc:170
SensitiveDetectorCatalog
Definition: SensitiveDetectorCatalog.h:10
dqmdumpme.k
k
Definition: dqmdumpme.py:60
CaloTrkProcessing::fineDetectors_
std::vector< Detector > fineDetectors_
Definition: CaloTrkProcessing.h:80
TrackInformation::getIDCaloVolume
int getIDCaloVolume() const
Definition: TrackInformation.h:51
GeV
const double GeV
Definition: MathUtil.h:16
CaloTrkProcessing::addlevel_
const int addlevel_
Definition: CaloTrkProcessing.h:78
TrackInformation::setCrossedBoundary
void setCrossedBoundary(const G4Track *track)
Definition: TrackInformation.h:61
BeginOfEvent.h
TrackInformation
Definition: TrackInformation.h:12
MetAnalyzer.pv
def pv(vc)
Definition: MetAnalyzer.py:7
BeginOfEvent
Definition: BeginOfEvent.h:6
CaloSimulationParameters::insideNames_
std::vector< std::string > insideNames_
Definition: CaloSimulationParameters.h:14
TrackInformation.h
CaloTrkProcessing::testBeam_
const bool testBeam_
Definition: CaloTrkProcessing.h:73
edm::LogError
Log< level::Error, false > LogError
Definition: MessageLogger.h:123
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
CaloTrkProcessing::Detector
Definition: CaloTrkProcessing.h:53
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
triggerObjects_cff.id
id
Definition: triggerObjects_cff.py:29
SimTrackManager.h
Exception
Definition: hltDiff.cc:245
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
Exception.h
hgcalTestNeighbor_cfi.detector
detector
Definition: hgcalTestNeighbor_cfi.py:6
CaloTrkProcessing::isItCalo
int isItCalo(const G4VTouchable *, const std::vector< Detector > &)
Definition: CaloTrkProcessing.cc:293
TrackInformation::caloIDChecked
bool caloIDChecked() const
Definition: TrackInformation.h:53
CaloTrkProcessing::detectors_
std::vector< Detector > detectors_
Definition: CaloTrkProcessing.h:80
CaloTrkProcessing::putHistory_
const bool putHistory_
Definition: CaloTrkProcessing.h:75
globals_cff.id2
id2
Definition: globals_cff.py:34
hfnoseParametersInitialization_cfi.name2
name2
Definition: hfnoseParametersInitialization_cfi.py:8
CaloTrkProcessing::~CaloTrkProcessing
~CaloTrkProcessing() override
Definition: CaloTrkProcessing.cc:168
edm::Log
Definition: MessageLogger.h:70
cuy.ii
ii
Definition: cuy.py:589
CaloTrkProcessing::CaloTrkProcessing
CaloTrkProcessing(const std::string &aSDname, const CaloSimulationParameters &csps, const SensitiveDetectorCatalog &clg, bool testBeam, double eMin, bool putHistory, bool doFineCalo, double eMinFine, int addlevel, const std::vector< std::string > &fineNames, const std::vector< int > &fineLevels, const std::vector< int > &useFines, const SimTrackManager *)
Definition: CaloTrkProcessing.cc:21
jets_cff.levels
levels
Definition: jets_cff.py:21
CaloTrkProcessing::eMin_
const double eMin_
Definition: CaloTrkProcessing.h:74
TrackWithHistory.h
susybsm::HSCParticleType::unknown
Definition: HSCParticle.h:20