CMS 3D CMS Logo

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