CMS 3D CMS Logo

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