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 DDCompactView& cpv,
27  const SensitiveDetectorCatalog& clg,
28  edm::ParameterSet const& p,
29  const SimTrackManager* manager)
30  : SensitiveCaloDetector(name, cpv, clg, p), lastTrackID(-1), m_trackManager(manager) {
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 
37  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: Initailised with TestBeam = " << testBeam << " Emin = " << eMin
38  << " MeV and"
39  << " History flag " << putHistory;
40 
41  //Get the names
42  G4String attribute = "ReadOutName";
43  DDSpecificsMatchesValueFilter filter{DDValue(attribute, name, 0)};
44  DDFilteredView fv(cpv, filter);
45  fv.firstChild();
47 
48  G4String value = "Calorimeter";
49  std::vector<std::string> caloNames = getNames(value, sv);
50 #ifdef EDM_ML_DEBUG
51  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: Names for " << value << ":";
52  for (unsigned int i = 0; i < caloNames.size(); i++)
53  edm::LogVerbatim("CaloSim") << " (" << i << ") " << caloNames[i];
54 #endif
55 
56  value = "Levels";
57  std::vector<double> levels = getNumbers(value, sv);
58 #ifdef EDM_ML_DEBUG
59  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: Names for " << value << ":";
60  for (unsigned int i = 0; i < levels.size(); i++)
61  edm::LogVerbatim("CaloSim") << " (" << i << ") " << levels[i];
62 #endif
63 
64  value = "Neighbours";
65  std::vector<double> neighbours = getNumbers(value, sv);
66 #ifdef EDM_ML_DEBUG
67  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: Names for " << value << ":";
68  for (unsigned int i = 0; i < neighbours.size(); i++)
69  edm::LogVerbatim("CaloSim") << " (" << i << ") " << neighbours[i];
70 #endif
71 
72  value = "Inside";
73  std::vector<std::string> insideNames = getNames(value, sv);
74 #ifdef EDM_ML_DEBUG
75  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: Names for " << value << ":";
76  for (unsigned int i = 0; i < insideNames.size(); i++)
77  edm::LogVerbatim("CaloSim") << " (" << i << ") " << insideNames[i];
78 #endif
79 
80  value = "InsideLevel";
81  std::vector<double> insideLevel = getNumbers(value, sv);
82 #ifdef EDM_ML_DEBUG
83  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: Names for " << value << ":";
84  for (unsigned int i = 0; i < insideLevel.size(); i++)
85  edm::LogVerbatim("CaloSim") << " (" << i << ") " << insideLevel[i];
86 #endif
87 
88  if (caloNames.size() < neighbours.size()) {
89  edm::LogError("CaloSim") << "CaloTrkProcessing: # of Calorimeter bins " << caloNames.size()
90  << " does not match with " << neighbours.size() << " ==> illegal ";
91  throw cms::Exception("Unknown", "CaloTrkProcessing")
92  << "Calorimeter array size does not match with size of neighbours\n";
93  }
94 
95  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
96  std::vector<G4LogicalVolume*>::const_iterator lvcite;
97  int istart = 0;
98  for (unsigned int i = 0; i < caloNames.size(); i++) {
99  G4LogicalVolume* lv = nullptr;
100  G4String name = caloNames[i];
101  int number = static_cast<int>(neighbours[i]);
102  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
103  if ((*lvcite)->GetName() == name) {
104  lv = (*lvcite);
105  break;
106  }
107  }
108  if (lv != nullptr) {
110  detector.name = name;
111  detector.lv = lv;
112  detector.level = static_cast<int>(levels[i]);
113  if (istart + number > (int)(insideNames.size())) {
114  edm::LogError("CaloSim") << "CaloTrkProcessing: # of InsideNames bins " << insideNames.size()
115  << " too few compaerd to " << istart + number << " requested ==> illegal ";
116  throw cms::Exception("Unknown", "CaloTrkProcessing")
117  << "InsideNames array size does not match with list of neighbours\n";
118  }
119  std::vector<std::string> inside;
120  std::vector<G4LogicalVolume*> insideLV;
121  std::vector<int> insideLevels;
122  for (int k = 0; k < number; k++) {
123  lv = nullptr;
124  name = insideNames[istart + k];
125  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++)
126  if ((*lvcite)->GetName() == name) {
127  lv = (*lvcite);
128  break;
129  }
130  inside.push_back(name);
131  insideLV.push_back(lv);
132  insideLevels.push_back(static_cast<int>(insideLevel[istart + k]));
133  }
134  detector.fromDets = inside;
135  detector.fromDetL = insideLV;
136  detector.fromLevels = insideLevels;
137  detectors.push_back(detector);
138  }
139  istart += number;
140  }
141 
142  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: with " << detectors.size() << " calorimetric volumes";
143  for (unsigned int i = 0; i < detectors.size(); i++) {
144  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: Calorimeter volume " << i << " " << detectors[i].name << " LV "
145  << detectors[i].lv << " at level " << detectors[i].level << " with "
146  << detectors[i].fromDets.size() << " neighbours";
147  for (unsigned int k = 0; k < detectors[i].fromDets.size(); k++)
148  edm::LogVerbatim("CaloSim") << " Element " << k << " " << detectors[i].fromDets[k] << " LV "
149  << detectors[i].fromDetL[k] << " at level " << detectors[i].fromLevels[k];
150  }
151 }
152 
153 CaloTrkProcessing::~CaloTrkProcessing() { edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: Deleted"; }
154 
156 
157 void CaloTrkProcessing::update(const G4Step* aStep) {
158  // define if you are at the surface of CALO
159 
160  G4Track* theTrack = aStep->GetTrack();
161  int id = theTrack->GetTrackID();
162 
163  TrackInformation* trkInfo = dynamic_cast<TrackInformation*>(theTrack->GetUserInformation());
164 
165  if (trkInfo == nullptr) {
166  edm::LogError("CaloSim") << "CaloTrkProcessing: No trk info !!!! abort ";
167  throw cms::Exception("Unknown", "CaloTrkProcessing") << "cannot get trkInfo for Track " << id << "\n";
168  }
169 
170  if (testBeam) {
171  if (trkInfo->getIDonCaloSurface() == 0) {
172 #ifdef EDM_ML_DEBUG
173  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing set IDonCaloSurface to " << id << " at step Number "
174  << theTrack->GetCurrentStepNumber();
175 #endif
176  trkInfo->setIDonCaloSurface(id, 0, 0, theTrack->GetDefinition()->GetPDGEncoding(), theTrack->GetMomentum().mag());
177  lastTrackID = id;
178  if (theTrack->GetKineticEnergy() / MeV > eMin)
179  trkInfo->putInHistory();
180  }
181  } else {
182  if (putHistory) {
183  trkInfo->putInHistory();
184  // trkInfo->setAncestor();
185  }
186 #ifdef EDM_ML_DEBUG
187  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing Entered for " << id << " at stepNumber "
188  << theTrack->GetCurrentStepNumber() << " IDonCaloSur.. "
189  << trkInfo->getIDonCaloSurface() << " CaloCheck " << trkInfo->caloIDChecked();
190 #endif
191  if (trkInfo->getIDonCaloSurface() != 0) {
192  if (trkInfo->caloIDChecked() == false) {
193  G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
194  const G4VTouchable* post_touch = postStepPoint->GetTouchable();
195 
196  if (isItInside(post_touch, trkInfo->getIDCaloVolume(), trkInfo->getIDLastVolume()) > 0) {
197  trkInfo->setIDonCaloSurface(0, -1, -1, 0, 0);
198  } else {
199  trkInfo->setCaloIDChecked(true);
200  }
201  }
202  } else {
203  G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
204  const G4VTouchable* post_touch = postStepPoint->GetTouchable();
205  int ical = isItCalo(post_touch);
206  if (ical >= 0) {
207  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
208  const G4VTouchable* pre_touch = preStepPoint->GetTouchable();
209  int inside = isItInside(pre_touch, ical, -1);
210  if (inside >= 0 || (theTrack->GetCurrentStepNumber() == 1)) {
211  trkInfo->setIDonCaloSurface(
212  id, ical, inside, theTrack->GetDefinition()->GetPDGEncoding(), theTrack->GetMomentum().mag());
213  trkInfo->setCaloIDChecked(true);
214  lastTrackID = id;
215  if (theTrack->GetKineticEnergy() / MeV > eMin)
216  trkInfo->putInHistory();
217 #ifdef EDM_ML_DEBUG
218  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: set ID on Calo " << ical << " surface (Inside " << inside
219  << ") to " << id << " of a Track with Kinetic Energy "
220  << theTrack->GetKineticEnergy() / MeV << " MeV";
221 #endif
222  }
223  }
224  }
225  }
226 }
227 
228 std::vector<std::string> CaloTrkProcessing::getNames(const G4String str, const DDsvalues_type& sv) {
229 #ifdef EDM_ML_DEBUG
230  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing::getNames called for " << str;
231 #endif
232  DDValue value(str);
233  if (DDfetch(&sv, value)) {
234 #ifdef EDM_ML_DEBUG
235  edm::LogVerbatim("CaloSim") << value;
236 #endif
237  const std::vector<std::string>& fvec = value.strings();
238  int nval = fvec.size();
239  if (nval < 1) {
240  edm::LogError("CaloSim") << "CaloTrkProcessing: # of " << str << " bins " << nval << " < 1 ==> illegal ";
241  throw cms::Exception("Unknown", "CaloTrkProcessing") << "nval < 2 for array " << str << "\n";
242  }
243 
244  return fvec;
245  } else {
246  edm::LogError("CaloSim") << "CaloTrkProcessing: cannot get array " << str;
247  throw cms::Exception("Unknown", "CaloTrkProcessing") << "cannot get array " << str << "\n";
248  }
249 }
250 
251 std::vector<double> CaloTrkProcessing::getNumbers(const G4String str, const DDsvalues_type& sv) {
252 #ifdef EDM_ML_DEBUG
253  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing::getNumbers called for " << str;
254 #endif
255  DDValue value(str);
256  if (DDfetch(&sv, value)) {
257 #ifdef EDM_ML_DEBUG
258  edm::LogVerbatim("CaloSim") << value;
259 #endif
260  const std::vector<double>& fvec = value.doubles();
261  int nval = fvec.size();
262  if (nval < 1) {
263  edm::LogError("CaloSim") << "CaloTrkProcessing: # of " << str << " bins " << nval << " < 1 ==> illegal ";
264  throw cms::Exception("Unknown", "CaloTrkProcessing") << "nval < 2 for array " << str << "\n";
265  }
266 
267  return fvec;
268  } else {
269  edm::LogError("CaloSim") << "CaloTrkProcessing: cannot get array " << str;
270  throw cms::Exception("Unknown", "CaloTrkProcessing") << "cannot get array " << str << "\n";
271  }
272 }
273 
274 int CaloTrkProcessing::isItCalo(const G4VTouchable* touch) {
275  int lastLevel = -1;
276  G4LogicalVolume* lv = nullptr;
277  for (unsigned int it = 0; it < detectors.size(); it++) {
278  if (lastLevel != detectors[it].level) {
279  lastLevel = detectors[it].level;
280  lv = detLV(touch, lastLevel);
281 #ifdef EDM_ML_DEBUG
282  std::string name1 = "Unknown";
283  if (lv != 0)
284  name1 = lv->GetName();
285  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: volume " << name1 << " at Level " << lastLevel;
286  int levels = detLevels(touch);
287  if (levels > 0) {
288  G4String name2[20];
289  int copyno2[20];
290  detectorLevel(touch, levels, copyno2, name2);
291  for (int i2 = 0; i2 < levels; i2++)
292  edm::LogVerbatim("CaloSim") << " " << i2 << " " << name2[i2] << " " << copyno2[i2];
293  }
294 #endif
295  }
296  bool ok = (lv == detectors[it].lv);
297  if (ok)
298  return it;
299  }
300  return -1;
301 }
302 
303 int CaloTrkProcessing::isItInside(const G4VTouchable* touch, int idcal, int idin) {
304  int lastLevel = -1;
305  G4LogicalVolume* lv = nullptr;
306  int id1, id2;
307  if (idcal < 0) {
308  id1 = 0;
309  id2 = static_cast<int>(detectors.size());
310  } else {
311  id1 = idcal;
312  id2 = id1 + 1;
313  }
314  for (int it1 = id1; it1 < id2; it1++) {
315  if (idin < 0) {
316  for (unsigned int it2 = 0; it2 < detectors[it1].fromDets.size(); it2++) {
317  if (lastLevel != detectors[it1].fromLevels[it2]) {
318  lastLevel = detectors[it1].fromLevels[it2];
319  lv = detLV(touch, lastLevel);
320 #ifdef EDM_ML_DEBUG
321  std::string name1 = "Unknown";
322  if (lv != 0)
323  name1 = lv->GetName();
324  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: volume " << name1 << " at Level " << lastLevel;
325  int levels = detLevels(touch);
326  if (levels > 0) {
327  G4String name2[20];
328  int copyno2[20];
329  detectorLevel(touch, levels, copyno2, name2);
330  for (int i2 = 0; i2 < levels; i2++)
331  edm::LogVerbatim("CaloSim") << " " << i2 << " " << name2[i2] << " " << copyno2[i2];
332  }
333 #endif
334  }
335  bool ok = (lv == detectors[it1].fromDetL[it2]);
336  if (ok)
337  return it2;
338  }
339  } else {
340  lastLevel = detectors[it1].fromLevels[idin];
341  lv = detLV(touch, lastLevel);
342 #ifdef EDM_ML_DEBUG
343  std::string name1 = "Unknown";
344  if (lv != 0)
345  name1 = lv->GetName();
346  edm::LogVerbatim("CaloSim") << "CaloTrkProcessing: volume " << name1 << " at Level " << lastLevel;
347  int levels = detLevels(touch);
348  if (levels > 0) {
349  G4String name2[20];
350  int copyno2[20];
351  detectorLevel(touch, levels, copyno2, name2);
352  for (int i2 = 0; i2 < levels; i2++)
353  edm::LogVerbatim("CaloSim") << " " << i2 << " " << name2[i2] << " " << copyno2[i2];
354  }
355 #endif
356  bool ok = (lv == detectors[it1].fromDetL[idin]);
357  if (ok)
358  return idin;
359  }
360  }
361  return -1;
362 }
363 
364 int CaloTrkProcessing::detLevels(const G4VTouchable* touch) const {
365  //Return number of levels
366  if (touch)
367  return ((touch->GetHistoryDepth()) + 1);
368  else
369  return 0;
370 }
371 
372 G4LogicalVolume* CaloTrkProcessing::detLV(const G4VTouchable* touch, int currentlevel) const {
373  G4LogicalVolume* lv = nullptr;
374  if (touch) {
375  int level = ((touch->GetHistoryDepth()) + 1);
376  if (level > 0 && level >= currentlevel) {
377  int ii = level - currentlevel;
378  lv = touch->GetVolume(ii)->GetLogicalVolume();
379  return lv;
380  }
381  }
382  return lv;
383 }
384 
385 void CaloTrkProcessing::detectorLevel(const G4VTouchable* touch, int& level, int* copyno, G4String* name) const {
386  static const std::string unknown("Unknown");
387  //Get name and copy numbers
388  if (level > 0) {
389  for (int ii = 0; ii < level; ii++) {
390  int i = level - ii - 1;
391  G4VPhysicalVolume* pv = touch->GetVolume(i);
392  if (pv != nullptr)
393  name[ii] = pv->GetName();
394  else
395  name[ii] = unknown;
396  copyno[ii] = touch->GetReplicaNumber(i);
397  }
398  }
399 }
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
const std::vector< double > & doubles() const
a reference to the double-valued values stored in the given instance of DDValue
Definition: DDValue.cc:140
CaloTrkProcessing(const std::string &aSDname, const DDCompactView &cpv, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p, const SimTrackManager *)
int getIDonCaloSurface() const
int detLevels(const G4VTouchable *) const
int isItCalo(const G4VTouchable *)
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:80
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:81
void detectorLevel(const G4VTouchable *, int &, int *, G4String *) const
const double MeV
int getIDCaloVolume() const
int isItInside(const G4VTouchable *, int, int)
const std::vector< std::string > & getNames() const
std::vector< int > fromLevels
std::vector< Detector > detectors
def pv(vc)
Definition: MetAnalyzer.py:7
std::vector< double > getNumbers(G4String, const DDsvalues_type &)
std::vector< G4LogicalVolume * > fromDetL
Definition: value.py:1
int getIDLastVolume() const
std::vector< std::string > fromDets
ii
Definition: cuy.py:590
const std::vector< std::string > & strings() const
a reference to the std::string-valued values stored in the given instance of DDValue ...
Definition: DDValue.h:61
int k[5][pyjets_maxn]
~CaloTrkProcessing() override
bool caloIDChecked() const
DDsvalues_type mergedSpecifics() const
void setCaloIDChecked(bool f)
void setIDonCaloSurface(int id, int ical, int last, int pdgID, double p)
std::vector< std::pair< unsigned int, DDValue > > DDsvalues_type
Definition: DDsvalues.h:12
bool firstChild()
set the current node to the first child ...
#define str(s)