CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Member Functions | Private Member Functions | Private Attributes
CaloTrkProcessing Class Reference

#include <CaloTrkProcessing.h>

Inheritance diagram for CaloTrkProcessing:
SensitiveCaloDetector Observer< const BeginOfEvent * > Observer< const G4Step * > SensitiveDetector

Classes

struct  Detector
 

Public Member Functions

 CaloTrkProcessing (G4String aSDname, const DDCompactView &cpv, SensitiveDetectorCatalog &clg, edm::ParameterSet const &p, const SimTrackManager *)
 
virtual void clearHits ()
 
virtual void EndOfEvent (G4HCofThisEvent *)
 
void fillHits (edm::PCaloHitContainer &, std::string)
 
virtual void Initialize (G4HCofThisEvent *)
 
virtual bool ProcessHits (G4Step *, G4TouchableHistory *)
 
virtual uint32_t setDetUnitId (G4Step *step)
 
virtual ~CaloTrkProcessing ()
 
- Public Member Functions inherited from SensitiveCaloDetector
 SensitiveCaloDetector (std::string &iname, const DDCompactView &cpv, SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
 
- Public Member Functions inherited from SensitiveDetector
virtual void AssignSD (std::string &vname)
 
Local3DPoint ConvertToLocal3DPoint (const G4ThreeVector &point)
 
Local3DPoint FinalStepPosition (G4Step *s, coordinates)
 
virtual std::vector< std::string > getNames ()
 
Local3DPoint InitialStepPosition (G4Step *s, coordinates)
 
std::string nameOfSD ()
 
void NaNTrap (G4Step *step)
 
void Register ()
 
 SensitiveDetector (std::string &iname, const DDCompactView &cpv, SensitiveDetectorCatalog &, edm::ParameterSet const &p)
 
virtual ~SensitiveDetector ()
 
- Public Member Functions inherited from Observer< const BeginOfEvent * >
 Observer ()
 
void slotForUpdate (const BeginOfEvent *iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const G4Step * >
 Observer ()
 
void slotForUpdate (const G4Step *iT)
 
virtual ~Observer ()
 

Private Member Functions

void detectorLevel (const G4VTouchable *, int &, int *, G4String *) const
 
int detLevels (const G4VTouchable *) const
 
G4LogicalVolume * detLV (const G4VTouchable *, int) const
 
std::vector< std::string > getNames (G4String, const DDsvalues_type &)
 
std::vector< double > getNumbers (G4String, const DDsvalues_type &)
 
int isItCalo (const G4VTouchable *)
 
int isItInside (const G4VTouchable *, int, int)
 
void update (const BeginOfEvent *evt)
 This routine will be called when the appropriate signal arrives. More...
 
void update (const G4Step *)
 This routine will be called when the appropriate signal arrives. More...
 

Private Attributes

std::vector< Detectordetectors
 
double eMin
 
int lastTrackID
 
const SimTrackManagerm_trackManager
 
bool putHistory
 
bool testBeam
 

Additional Inherited Members

- Public Types inherited from SensitiveDetector
enum  coordinates { WorldCoordinates, LocalCoordinates }
 

Detailed Description

Definition at line 22 of file CaloTrkProcessing.h.

Constructor & Destructor Documentation

CaloTrkProcessing::CaloTrkProcessing ( G4String  aSDname,
const DDCompactView cpv,
SensitiveDetectorCatalog clg,
edm::ParameterSet const &  p,
const SimTrackManager manager 
)

Definition at line 24 of file CaloTrkProcessing.cc.

References DDFilteredView::addFilter(), detectors, eMin, DDSpecificsFilter::equals, edm::hlt::Exception, alcazmumu_cfi::filter, DDFilteredView::firstChild(), CaloTrkProcessing::Detector::fromDetL, CaloTrkProcessing::Detector::fromDets, CaloTrkProcessing::Detector::fromLevels, SensitiveDetector::getNames(), getNumbers(), edm::ParameterSet::getParameter(), i, gen::k, CaloTrkProcessing::Detector::level, LogDebug, CaloTrkProcessing::Detector::lv, DDFilteredView::mergedSpecifics(), SensitiveDetector::name, CaloTrkProcessing::Detector::name, putHistory, DDSpecificsFilter::setCriteria(), testBeam, and relativeConstraints::value.

28  :
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();
50  DDsvalues_type sv(fv.mergedSpecifics());
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 }
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
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()
std::vector< std::string > fromDets
int k[5][pyjets_maxn]
SensitiveCaloDetector(std::string &iname, const DDCompactView &cpv, SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
void setCriteria(const DDValue &nameVal, comp_op, log_op l=AND, bool asString=true, bool merged=true)
Definition: DDFilter.cc:285
const SimTrackManager * m_trackManager
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
Definition: DDFilter.h:37
CaloTrkProcessing::~CaloTrkProcessing ( )
virtual

Definition at line 164 of file CaloTrkProcessing.cc.

164  {
165  edm::LogInfo("CaloSim") << "CaloTrkProcessing: Deleted";
166 }

Member Function Documentation

virtual void CaloTrkProcessing::clearHits ( )
inlinevirtual

Implements SensitiveDetector.

Definition at line 33 of file CaloTrkProcessing.h.

33 {}
void CaloTrkProcessing::detectorLevel ( const G4VTouchable *  touch,
int &  level,
int *  copyno,
G4String *  name 
) const
private

Definition at line 417 of file CaloTrkProcessing.cc.

References i, cuy::ii, testEve_cfg::level, AlCaHLTBitMon_QueryRunRegistry::string, and susybsm::HSCParticleType::unknown.

Referenced by isItCalo(), and isItInside().

418  {
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 }
int i
Definition: DBlmapReader.cc:9
int ii
Definition: cuy.py:588
tuple level
Definition: testEve_cfg.py:34
int CaloTrkProcessing::detLevels ( const G4VTouchable *  touch) const
private

Definition at line 393 of file CaloTrkProcessing.cc.

Referenced by isItCalo(), and isItInside().

393  {
394 
395  //Return number of levels
396  if (touch)
397  return ((touch->GetHistoryDepth())+1);
398  else
399  return 0;
400 }
G4LogicalVolume * CaloTrkProcessing::detLV ( const G4VTouchable *  touch,
int  currentlevel 
) const
private

Definition at line 402 of file CaloTrkProcessing.cc.

References cuy::ii, and testEve_cfg::level.

Referenced by isItCalo(), and isItInside().

403  {
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 }
int ii
Definition: cuy.py:588
tuple level
Definition: testEve_cfg.py:34
virtual void CaloTrkProcessing::EndOfEvent ( G4HCofThisEvent *  )
inlinevirtual

Reimplemented from SensitiveDetector.

Definition at line 38 of file CaloTrkProcessing.h.

38 {}
void CaloTrkProcessing::fillHits ( edm::PCaloHitContainer ,
std::string   
)
inlinevirtual

Implements SensitiveCaloDetector.

Definition at line 39 of file CaloTrkProcessing.h.

39 {}
std::vector< std::string > CaloTrkProcessing::getNames ( G4String  str,
const DDsvalues_type sv 
)
private

Definition at line 252 of file CaloTrkProcessing.cc.

References DDfetch(), edm::hlt::Exception, LogDebug, DDValue::strings(), and relativeConstraints::value.

253  {
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 }
#define LogDebug(id)
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:102
std::vector< double > CaloTrkProcessing::getNumbers ( G4String  str,
const DDsvalues_type sv 
)
private

Definition at line 281 of file CaloTrkProcessing.cc.

References DDfetch(), DDValue::doubles(), edm::hlt::Exception, LogDebug, and relativeConstraints::value.

Referenced by CaloTrkProcessing().

282  {
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 }
#define LogDebug(id)
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:102
virtual void CaloTrkProcessing::Initialize ( G4HCofThisEvent *  )
inlinevirtual

Reimplemented from SensitiveDetector.

Definition at line 32 of file CaloTrkProcessing.h.

32 {}
int CaloTrkProcessing::isItCalo ( const G4VTouchable *  touch)
private

Definition at line 309 of file CaloTrkProcessing.cc.

References detectorLevel(), detectors, detLevels(), detLV(), testEve_cfg::level, LogDebug, convertSQLiteXML::ok, and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by update().

309  {
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 }
#define LogDebug(id)
G4LogicalVolume * detLV(const G4VTouchable *, int) const
int detLevels(const G4VTouchable *) const
void detectorLevel(const G4VTouchable *, int &, int *, G4String *) const
std::vector< Detector > detectors
tuple level
Definition: testEve_cfg.py:34
int CaloTrkProcessing::isItInside ( const G4VTouchable *  touch,
int  idcal,
int  idin 
)
private

Definition at line 338 of file CaloTrkProcessing.cc.

References detectorLevel(), detectors, detLevels(), detLV(), LogDebug, convertSQLiteXML::ok, and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by update().

339  {
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 }
#define LogDebug(id)
G4LogicalVolume * detLV(const G4VTouchable *, int) const
int detLevels(const G4VTouchable *) const
void detectorLevel(const G4VTouchable *, int &, int *, G4String *) const
std::vector< Detector > detectors
virtual bool CaloTrkProcessing::ProcessHits ( G4Step *  ,
G4TouchableHistory *   
)
inlinevirtual

Implements SensitiveDetector.

Definition at line 34 of file CaloTrkProcessing.h.

34  {
35  return true;
36  }
virtual uint32_t CaloTrkProcessing::setDetUnitId ( G4Step *  step)
inlinevirtual

Implements SensitiveDetector.

Definition at line 37 of file CaloTrkProcessing.h.

37 {return 0;}
void CaloTrkProcessing::update ( const BeginOfEvent )
privatevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfEvent * >.

Definition at line 168 of file CaloTrkProcessing.cc.

References lastTrackID.

Referenced by progressbar.ProgressBar::__next__(), relval_steps.Matrix::__setitem__(), relval_steps.Steps::__setitem__(), Vispa.Gui.VispaWidget.VispaWidget::autosize(), Vispa.Views.LineDecayView.LineDecayContainer::createObject(), Vispa.Views.LineDecayView.LineDecayContainer::deselectAllObjects(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::deselectAllWidgets(), Vispa.Gui.VispaWidget.VispaWidget::enableAutosizing(), progressbar.ProgressBar::finish(), Vispa.Gui.MenuWidget.MenuWidget::leaveEvent(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::mouseMoveEvent(), Vispa.Gui.MenuWidget.MenuWidget::mouseMoveEvent(), Vispa.Views.LineDecayView.LineDecayContainer::mouseMoveEvent(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::mouseReleaseEvent(), Vispa.Views.LineDecayView.LineDecayContainer::objectMoved(), relval_steps.Steps::overwrite(), Vispa.Views.LineDecayView.LineDecayContainer::removeObject(), Vispa.Gui.ConnectableWidget.ConnectableWidget::removePorts(), Vispa.Gui.FindDialog.FindDialog::reset(), Vispa.Gui.PortConnection.PointToPointConnection::select(), Vispa.Gui.VispaWidget.VispaWidget::select(), Vispa.Views.LineDecayView.LineDecayContainer::select(), Vispa.Gui.VispaWidget.VispaWidget::setText(), Vispa.Gui.VispaWidget.VispaWidget::setTitle(), Vispa.Gui.ZoomableWidget.ZoomableWidget::setZoom(), Vispa.Views.LineDecayView.LineDecayContainer::setZoom(), and Vispa.Gui.PortConnection.PointToPointConnection::updateConnection().

168  {
169  lastTrackID = -1;
170 }
void CaloTrkProcessing::update ( const G4Step *  )
privatevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const G4Step * >.

Definition at line 172 of file CaloTrkProcessing.cc.

References TrackInformation::caloIDChecked(), eMin, edm::hlt::Exception, TrackInformation::getIDCaloVolume(), TrackInformation::getIDLastVolume(), TrackInformation::getIDonCaloSurface(), isItCalo(), isItInside(), lastTrackID, LogDebug, putHistory, TrackInformation::putInHistory(), TrackInformation::setCaloIDChecked(), TrackInformation::setIDonCaloSurface(), and testBeam.

Referenced by progressbar.ProgressBar::__next__(), relval_steps.Matrix::__setitem__(), relval_steps.Steps::__setitem__(), Vispa.Gui.VispaWidget.VispaWidget::autosize(), Vispa.Views.LineDecayView.LineDecayContainer::createObject(), Vispa.Views.LineDecayView.LineDecayContainer::deselectAllObjects(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::deselectAllWidgets(), Vispa.Gui.VispaWidget.VispaWidget::enableAutosizing(), progressbar.ProgressBar::finish(), Vispa.Gui.MenuWidget.MenuWidget::leaveEvent(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::mouseMoveEvent(), Vispa.Gui.MenuWidget.MenuWidget::mouseMoveEvent(), Vispa.Views.LineDecayView.LineDecayContainer::mouseMoveEvent(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::mouseReleaseEvent(), Vispa.Views.LineDecayView.LineDecayContainer::objectMoved(), relval_steps.Steps::overwrite(), Vispa.Views.LineDecayView.LineDecayContainer::removeObject(), Vispa.Gui.ConnectableWidget.ConnectableWidget::removePorts(), Vispa.Gui.FindDialog.FindDialog::reset(), Vispa.Gui.PortConnection.PointToPointConnection::select(), Vispa.Gui.VispaWidget.VispaWidget::select(), Vispa.Views.LineDecayView.LineDecayContainer::select(), Vispa.Gui.VispaWidget.VispaWidget::setText(), Vispa.Gui.VispaWidget.VispaWidget::setTitle(), Vispa.Gui.ZoomableWidget.ZoomableWidget::setZoom(), Vispa.Views.LineDecayView.LineDecayContainer::setZoom(), and Vispa.Gui.PortConnection.PointToPointConnection::updateConnection().

172  {
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 }
#define LogDebug(id)
int getIDonCaloSurface() const
int isItCalo(const G4VTouchable *)
int getIDCaloVolume() const
int isItInside(const G4VTouchable *, int, int)
int getIDLastVolume() const
bool caloIDChecked() const
void setCaloIDChecked(bool f)
void setIDonCaloSurface(int id, int ical, int last, int pdgID, double p)

Member Data Documentation

std::vector<Detector> CaloTrkProcessing::detectors
private

Definition at line 70 of file CaloTrkProcessing.h.

Referenced by CaloTrkProcessing(), isItCalo(), and isItInside().

double CaloTrkProcessing::eMin
private

Definition at line 68 of file CaloTrkProcessing.h.

Referenced by CaloTrkProcessing(), and update().

int CaloTrkProcessing::lastTrackID
private

Definition at line 69 of file CaloTrkProcessing.h.

Referenced by update().

const SimTrackManager* CaloTrkProcessing::m_trackManager
private

Definition at line 71 of file CaloTrkProcessing.h.

bool CaloTrkProcessing::putHistory
private

Definition at line 67 of file CaloTrkProcessing.h.

Referenced by CaloTrkProcessing(), and update().

bool CaloTrkProcessing::testBeam
private

Definition at line 67 of file CaloTrkProcessing.h.

Referenced by CaloTrkProcessing(), and update().