CMS 3D CMS Logo

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 (const std::string &aSDname, const DDCompactView &cpv, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p, const SimTrackManager *)
 
void clearHits () override
 
void EndOfEvent (G4HCofThisEvent *) override
 
void fillHits (edm::PCaloHitContainer &, const std::string &) override
 
void Initialize (G4HCofThisEvent *) override
 
bool ProcessHits (G4Step *, G4TouchableHistory *) override
 
uint32_t setDetUnitId (const G4Step *step) override
 
 ~CaloTrkProcessing () override
 
- Public Member Functions inherited from SensitiveCaloDetector
 SensitiveCaloDetector (const std::string &iname, const DDCompactView &cpv, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
 
- Public Member Functions inherited from SensitiveDetector
void EndOfEvent (G4HCofThisEvent *eventHC) override
 
const std::vector< std::string > & getNames () const
 
void Initialize (G4HCofThisEvent *eventHC) override
 
bool isCaloSD () const
 
 SensitiveDetector (const std::string &iname, const DDCompactView &cpv, const SensitiveDetectorCatalog &, edm::ParameterSet const &p, bool calo)
 
 ~SensitiveDetector () override
 
- 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) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const G4Step *) override
 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

- Protected Types inherited from SensitiveDetector
enum  coordinates { WorldCoordinates, LocalCoordinates }
 
- Protected Member Functions inherited from SensitiveDetector
TrackInformationcmsTrackInformation (const G4Track *aTrack)
 
Local3DPoint ConvertToLocal3DPoint (const G4ThreeVector &point) const
 
Local3DPoint FinalStepPosition (const G4Step *step, coordinates) const
 
Local3DPoint InitialStepPosition (const G4Step *step, coordinates) const
 
Local3DPoint LocalPostStepPosition (const G4Step *step) const
 
Local3DPoint LocalPreStepPosition (const G4Step *step) const
 
void NaNTrap (const G4Step *step) const
 
void setNames (const std::vector< std::string > &)
 

Detailed Description

Definition at line 23 of file CaloTrkProcessing.h.

Constructor & Destructor Documentation

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

Definition at line 25 of file CaloTrkProcessing.cc.

References gamEcalExtractorBlocks_cff::detector, detectors, eMin, Exception, ALCARECOTkAlBeamHalo_cff::filter, DDFilteredView::firstChild(), CaloTrkProcessing::Detector::fromDetL, CaloTrkProcessing::Detector::fromDets, CaloTrkProcessing::Detector::fromLevels, SensitiveDetector::getNames(), getNumbers(), edm::ParameterSet::getParameter(), mps_fire::i, gen::k, CaloTrkProcessing::Detector::level, jets_cff::levels, LogDebug, CaloTrkProcessing::Detector::lv, DDFilteredView::mergedSpecifics(), MeV, dataset::name, CaloTrkProcessing::Detector::name, putHistory, pfDeepBoostedJetPreprocessParams_cfi::sv, and testBeam.

29  :
30  SensitiveCaloDetector(name, cpv, clg, p), lastTrackID(-1),
31  m_trackManager(manager) {
32 
33  //Initialise the parameter set
34  edm::ParameterSet m_p = p.getParameter<edm::ParameterSet>("CaloTrkProcessing");
35  testBeam = m_p.getParameter<bool>("TestBeam");
36  eMin = m_p.getParameter<double>("EminTrack")*MeV;
37  putHistory = m_p.getParameter<bool>("PutHistory");
38 
39  edm::LogInfo("CaloSim") << "CaloTrkProcessing: Initailised with TestBeam = "
40  << testBeam << " Emin = " << eMin << " MeV and"
41  << " History flag " << putHistory;
42 
43  //Get the names
44  G4String attribute = "ReadOutName";
46  DDFilteredView fv(cpv,filter);
47  fv.firstChild();
48  DDsvalues_type sv(fv.mergedSpecifics());
49 
50  G4String value = "Calorimeter";
51  std::vector<std::string> caloNames = getNames (value, sv);
52 #ifdef DebugLog
53  LogDebug("CaloSim") << "CaloTrkProcessing: Names for " << value << ":";
54  for (unsigned int i=0; i<caloNames.size(); i++)
55  LogDebug("CaloSim") << " (" << i << ") " << caloNames[i];
56 #endif
57 
58  value = "Levels";
59  std::vector<double> levels = getNumbers (value, sv);
60 #ifdef DebugLog
61  LogDebug("CaloSim") << "CaloTrkProcessing: Names for " << value << ":";
62  for (unsigned int i=0; i<levels.size(); i++)
63  LogDebug("CaloSim") << " (" << i << ") " << levels[i];
64 #endif
65 
66  value = "Neighbours";
67  std::vector<double> neighbours = getNumbers (value, sv);
68 #ifdef DebugLog
69  LogDebug("CaloSim") << "CaloTrkProcessing: Names for " << value << ":";
70  for (unsigned int i=0; i<neighbours.size(); i++)
71  LogDebug("CaloSim") << " (" << i << ") " << neighbours[i];
72 #endif
73 
74  value = "Inside";
75  std::vector<std::string> insideNames = getNames (value, sv);
76 #ifdef DebugLog
77  LogDebug("CaloSim") << "CaloTrkProcessing: Names for " << value << ":";
78  for (unsigned int i=0; i<insideNames.size(); i++)
79  LogDebug("CaloSim") << " (" << i << ") " << insideNames[i];
80 #endif
81 
82  value = "InsideLevel";
83  std::vector<double> insideLevel = getNumbers (value, sv);
84 #ifdef DebugLog
85  LogDebug("CaloSim") << "CaloTrkProcessing: Names for " << value << ":";
86  for (unsigned int i=0; i<insideLevel.size(); i++)
87  LogDebug("CaloSim") << " (" << i << ") " << insideLevel[i];
88 #endif
89 
90  if (caloNames.size() < neighbours.size()) {
91  edm::LogError("CaloSim") << "CaloTrkProcessing: # of Calorimeter bins "
92  << caloNames.size() << " does not match with "
93  << neighbours.size() << " ==> illegal ";
94  throw cms::Exception("Unknown", "CaloTrkProcessing")
95  << "Calorimeter array size does not match with size of neighbours\n";
96  }
97 
98  const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
99  std::vector<G4LogicalVolume *>::const_iterator lvcite;
100  int istart = 0;
101  for (unsigned int i=0; i<caloNames.size(); i++) {
102  G4LogicalVolume* lv = nullptr;
103  G4String name = caloNames[i];
104  int number = static_cast<int>(neighbours[i]);
105  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
106  if ((*lvcite)->GetName() == name) {
107  lv = (*lvcite);
108  break;
109  }
110  }
111  if (lv != nullptr) {
113  detector.name = name;
114  detector.lv = lv;
115  detector.level = static_cast<int>(levels[i]);
116  if (istart+number > (int)(insideNames.size())) {
117  edm::LogError("CaloSim") << "CaloTrkProcessing: # of InsideNames bins "
118  << insideNames.size() <<" too few compaerd to "
119  << istart+number << " requested ==> illegal ";
120  throw cms::Exception("Unknown", "CaloTrkProcessing")
121  << "InsideNames array size does not match with list of neighbours\n";
122  }
123  std::vector<std::string> inside;
124  std::vector<G4LogicalVolume*> insideLV;
125  std::vector<int> insideLevels;
126  for (int k = 0; k < number; k++) {
127  lv = nullptr;
128  name = insideNames[istart+k];
129  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++)
130  if ((*lvcite)->GetName() == name) {
131  lv = (*lvcite);
132  break;
133  }
134  inside.push_back(name);
135  insideLV.push_back(lv);
136  insideLevels.push_back(static_cast<int>(insideLevel[istart+k]));
137  }
138  detector.fromDets = inside;
139  detector.fromDetL = insideLV;
140  detector.fromLevels = insideLevels;
141  detectors.push_back(detector);
142  }
143  istart += number;
144  }
145 
146  edm::LogInfo("CaloSim") << "CaloTrkProcessing: with " << detectors.size()
147  << " calorimetric volumes";
148  for (unsigned int i=0; i<detectors.size(); i++) {
149  edm::LogInfo("CaloSim") << "CaloTrkProcessing: Calorimeter volume " << i
150  << " " << detectors[i].name << " LV "
151  << detectors[i].lv << " at level "
152  << detectors[i].level << " with "
153  << detectors[i].fromDets.size() << " neighbours";
154  for (unsigned int k=0; k<detectors[i].fromDets.size(); k++)
155  edm::LogInfo("CaloSim") << " Element " << k << " "
156  << detectors[i].fromDets[k] << " LV "
157  << detectors[i].fromDetL[k] << " at level "
158  << detectors[i].fromLevels[k];
159  }
160 }
#define LogDebug(id)
T getParameter(std::string const &) const
SensitiveCaloDetector(const std::string &iname, const DDCompactView &cpv, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
const double MeV
const std::vector< std::string > & getNames() const
std::vector< int > fromLevels
std::vector< Detector > detectors
std::vector< double > getNumbers(G4String, const DDsvalues_type &)
std::vector< G4LogicalVolume * > fromDetL
Definition: value.py:1
std::vector< std::string > fromDets
int k[5][pyjets_maxn]
std::vector< std::pair< unsigned int, DDValue > > DDsvalues_type
Definition: DDsvalues.h:12
const SimTrackManager * m_trackManager
CaloTrkProcessing::~CaloTrkProcessing ( )
override

Definition at line 162 of file CaloTrkProcessing.cc.

162  {
163  edm::LogInfo("CaloSim") << "CaloTrkProcessing: Deleted";
164 }

Member Function Documentation

void CaloTrkProcessing::clearHits ( )
inlineoverridevirtual

Implements SensitiveDetector.

Definition at line 34 of file CaloTrkProcessing.h.

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

Definition at line 415 of file CaloTrkProcessing.cc.

References mps_fire::i, cuy::ii, hcalDigis_cfi::level, MetAnalyzer::pv(), AlCaHLTBitMon_QueryRunRegistry::string, and susybsm::HSCParticleType::unknown.

Referenced by isItCalo(), and isItInside().

416  {
417 
418  static const std::string unknown("Unknown");
419  //Get name and copy numbers
420  if (level > 0) {
421  for (int ii = 0; ii < level; ii++) {
422  int i = level - ii - 1;
423  G4VPhysicalVolume* pv = touch->GetVolume(i);
424  if (pv != nullptr)
425  name[ii] = pv->GetName();
426  else
427  name[ii] = unknown;
428  copyno[ii] = touch->GetReplicaNumber(i);
429  }
430  }
431 }
def pv(vc)
Definition: MetAnalyzer.py:7
ii
Definition: cuy.py:590
int CaloTrkProcessing::detLevels ( const G4VTouchable *  touch) const
private

Definition at line 391 of file CaloTrkProcessing.cc.

Referenced by isItCalo(), and isItInside().

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

Definition at line 400 of file CaloTrkProcessing.cc.

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

Referenced by isItCalo(), and isItInside().

401  {
402 
403  G4LogicalVolume* lv=nullptr;
404  if (touch) {
405  int level = ((touch->GetHistoryDepth())+1);
406  if (level > 0 && level >= currentlevel) {
407  int ii = level - currentlevel;
408  lv = touch->GetVolume(ii)->GetLogicalVolume();
409  return lv;
410  }
411  }
412  return lv;
413 }
ii
Definition: cuy.py:590
void CaloTrkProcessing::EndOfEvent ( G4HCofThisEvent *  )
inlineoverride

Definition at line 39 of file CaloTrkProcessing.h.

39 {}
void CaloTrkProcessing::fillHits ( edm::PCaloHitContainer ,
const std::string &   
)
inlineoverridevirtual

Implements SensitiveCaloDetector.

Definition at line 40 of file CaloTrkProcessing.h.

References SensitiveDetector::getNames(), getNumbers(), isItCalo(), isItInside(), and update().

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

Definition at line 250 of file CaloTrkProcessing.cc.

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

251  {
252 
253 #ifdef DebugLog
254  LogDebug("CaloSim") << "CaloTrkProcessing::getNames called for " << str;
255 #endif
256  DDValue value(str);
257  if (DDfetch(&sv,value)) {
258 #ifdef DebugLog
259  LogDebug("CaloSim") << value;
260 #endif
261  const std::vector<std::string> & fvec = value.strings();
262  int nval = fvec.size();
263  if (nval < 1) {
264 
265  edm::LogError("CaloSim") << "CaloTrkProcessing: # of " << str
266  << " bins " << nval << " < 1 ==> illegal ";
267  throw cms::Exception("Unknown", "CaloTrkProcessing")
268  << "nval < 2 for array " << str << "\n";
269  }
270 
271  return fvec;
272  } else {
273  edm::LogError("CaloSim") << "CaloTrkProcessing: cannot get array " << str ;
274  throw cms::Exception("Unknown", "CaloTrkProcessing")
275  << "cannot get array " << str << "\n";
276  }
277 }
#define LogDebug(id)
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:81
Definition: value.py:1
#define str(s)
std::vector< double > CaloTrkProcessing::getNumbers ( G4String  str,
const DDsvalues_type sv 
)
private

Definition at line 279 of file CaloTrkProcessing.cc.

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

Referenced by CaloTrkProcessing(), and fillHits().

280  {
281 
282 #ifdef DebugLog
283  LogDebug("CaloSim") << "CaloTrkProcessing::getNumbers called for " << str;
284 #endif
285  DDValue value(str);
286  if (DDfetch(&sv,value)) {
287 #ifdef DebugLog
288  LogDebug("CaloSim") << value;
289 #endif
290  const std::vector<double> & fvec = value.doubles();
291  int nval = fvec.size();
292  if (nval < 1) {
293  edm::LogError("CaloSim") << "CaloTrkProcessing: # of " << str
294  << " bins " << nval << " < 1 ==> illegal ";
295  throw cms::Exception("Unknown", "CaloTrkProcessing")
296  << "nval < 2 for array " << str << "\n";
297  }
298 
299  return fvec;
300  } else {
301  edm::LogError("CaloSim") << "CaloTrkProcessing: cannot get array " << str ;
302  throw cms::Exception("Unknown", "CaloTrkProcessing")
303  << "cannot get array " << str << "\n";
304  }
305 }
#define LogDebug(id)
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:81
Definition: value.py:1
#define str(s)
void CaloTrkProcessing::Initialize ( G4HCofThisEvent *  )
inlineoverride

Definition at line 33 of file CaloTrkProcessing.h.

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

Definition at line 307 of file CaloTrkProcessing.cc.

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

Referenced by fillHits(), and update().

307  {
308 
309  int lastLevel = -1;
310  G4LogicalVolume* lv=nullptr;
311  for (unsigned int it=0; it < detectors.size(); it++) {
312  if (lastLevel != detectors[it].level) {
313  lastLevel = detectors[it].level;
314  lv = detLV(touch, lastLevel);
315 #ifdef DebugLog
316  std::string name1 = "Unknown";
317  if (lv != 0) name1 = lv->GetName();
318  LogDebug("CaloSim") << "CaloTrkProcessing: volume " << name1
319  << " at Level " << lastLevel;
320  int levels = detLevels(touch);
321  if (levels > 0) {
322  G4String name2[20]; int copyno2[20];
323  detectorLevel(touch, levels, copyno2, name2);
324  for (int i2=0; i2<levels; i2++)
325  LogDebug("CaloSim") << " " << i2 << " " << name2[i2] << " "
326  << copyno2[i2];
327  }
328 #endif
329  }
330  bool ok = (lv == detectors[it].lv);
331  if (ok) return it;
332  }
333  return -1;
334 }
#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
int CaloTrkProcessing::isItInside ( const G4VTouchable *  touch,
int  idcal,
int  idin 
)
private

Definition at line 336 of file CaloTrkProcessing.cc.

References detectorLevel(), detectors, detLevels(), detLV(), globals_cff::id1, globals_cff::id2, jets_cff::levels, LogDebug, convertSQLiteXML::ok, and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by fillHits(), and update().

337  {
338  int lastLevel = -1;
339  G4LogicalVolume* lv=nullptr;
340  int id1, id2;
341  if (idcal < 0) {id1 = 0; id2 = static_cast<int>(detectors.size());}
342  else {id1 = idcal; id2 = id1+1;}
343  for (int it1 = id1; it1 < id2; it1++) {
344  if (idin < 0) {
345  for (unsigned int it2 = 0; it2 < detectors[it1].fromDets.size(); it2++) {
346  if (lastLevel != detectors[it1].fromLevels[it2]) {
347  lastLevel = detectors[it1].fromLevels[it2];
348  lv = detLV(touch,lastLevel);
349 #ifdef DebugLog
350  std::string name1 = "Unknown";
351  if (lv != 0) name1 = lv->GetName();
352  LogDebug("CaloSim") << "CaloTrkProcessing: volume " << name1
353  << " at Level " << lastLevel;
354  int levels = detLevels(touch);
355  if (levels > 0) {
356  G4String name2[20]; int copyno2[20];
357  detectorLevel(touch, levels, copyno2, name2);
358  for (int i2=0; i2<levels; i2++)
359  LogDebug("CaloSim") << " " << i2 << " " << name2[i2] << " "
360  << copyno2[i2];
361  }
362 #endif
363  }
364  bool ok = (lv == detectors[it1].fromDetL[it2]);
365  if (ok) return it2;
366  }
367  } else {
368  lastLevel = detectors[it1].fromLevels[idin];
369  lv = detLV(touch,lastLevel);
370 #ifdef DebugLog
371  std::string name1 = "Unknown";
372  if (lv != 0) name1 = lv->GetName();
373  LogDebug("CaloSim") << "CaloTrkProcessing: volume " << name1
374  << " at Level " << lastLevel;
375  int levels = detLevels(touch);
376  if (levels > 0) {
377  G4String name2[20]; int copyno2[20];
378  detectorLevel(touch, levels, copyno2, name2);
379  for (int i2=0; i2<levels; i2++)
380  LogDebug("CaloSim") << " " << i2 << " " << name2[i2] << " "
381  << copyno2[i2];
382  }
383 #endif
384  bool ok = (lv == detectors[it1].fromDetL[idin]);
385  if (ok) return idin;
386  }
387  }
388  return -1;
389 }
#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
bool CaloTrkProcessing::ProcessHits ( G4Step *  ,
G4TouchableHistory *   
)
inlineoverridevirtual

Implements SensitiveDetector.

Definition at line 35 of file CaloTrkProcessing.h.

35  {
36  return true;
37  }
uint32_t CaloTrkProcessing::setDetUnitId ( const G4Step *  step)
inlineoverridevirtual

Implements SensitiveDetector.

Definition at line 38 of file CaloTrkProcessing.h.

38 {return 0;}
void CaloTrkProcessing::update ( const BeginOfEvent )
overrideprivatevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfEvent * >.

Definition at line 166 of file CaloTrkProcessing.cc.

References lastTrackID.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.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(), fillHits(), 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(), MatrixUtil.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().

166  {
167  lastTrackID = -1;
168 }
void CaloTrkProcessing::update ( const G4Step *  )
overrideprivatevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const G4Step * >.

Definition at line 170 of file CaloTrkProcessing.cc.

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

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.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(), MatrixUtil.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().

170  {
171 
172  // define if you are at the surface of CALO
173 
174  G4Track* theTrack = aStep->GetTrack();
175  int id = theTrack->GetTrackID();
176 
177  TrackInformation* trkInfo = dynamic_cast<TrackInformation*>
178  (theTrack->GetUserInformation());
179 
180  if (trkInfo == nullptr) {
181  edm::LogError("CaloSim") << "CaloTrkProcessing: No trk info !!!! abort ";
182  throw cms::Exception("Unknown", "CaloTrkProcessing")
183  << "cannot get trkInfo for Track " << id << "\n";
184  }
185 
186  if (testBeam) {
187  if (trkInfo->getIDonCaloSurface() == 0) {
188 #ifdef DebugLog
189  LogDebug("CaloSim") << "CaloTrkProcessing set IDonCaloSurface to " << id
190  << " at step Number "
191  << theTrack->GetCurrentStepNumber();
192 #endif
193  trkInfo->setIDonCaloSurface(id,0,0,
194  theTrack->GetDefinition()->GetPDGEncoding(),
195  theTrack->GetMomentum().mag());
196  lastTrackID = id;
197  if (theTrack->GetKineticEnergy()/MeV > eMin)
198  trkInfo->putInHistory();
199  }
200  } else {
201  if (putHistory) {
202  trkInfo->putInHistory();
203  // trkInfo->setAncestor();
204  }
205 #ifdef DebugLog
206  LogDebug("CaloSim") << "CaloTrkProcessing Entered for " << id
207  << " at stepNumber "<< theTrack->GetCurrentStepNumber()
208  << " IDonCaloSur.. " << trkInfo->getIDonCaloSurface()
209  << " CaloCheck " << trkInfo->caloIDChecked();
210 #endif
211  if (trkInfo->getIDonCaloSurface() != 0) {
212  if (trkInfo->caloIDChecked() == false) {
213  G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
214  const G4VTouchable* post_touch = postStepPoint->GetTouchable();
215 
216  if (isItInside(post_touch, trkInfo->getIDCaloVolume(),
217  trkInfo->getIDLastVolume()) > 0) {
218  trkInfo->setIDonCaloSurface(0, -1, -1, 0, 0);
219  } else {
220  trkInfo->setCaloIDChecked(true);
221  }
222  }
223  } else {
224  G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
225  const G4VTouchable* post_touch = postStepPoint->GetTouchable();
226  int ical = isItCalo(post_touch);
227  if (ical >= 0) {
228  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
229  const G4VTouchable* pre_touch = preStepPoint->GetTouchable();
230  int inside = isItInside(pre_touch, ical, -1);
231  if (inside >= 0 || (theTrack->GetCurrentStepNumber()==1)) {
232  trkInfo->setIDonCaloSurface(id,ical,inside,
233  theTrack->GetDefinition()->GetPDGEncoding(),
234  theTrack->GetMomentum().mag());
235  trkInfo->setCaloIDChecked(true);
236  lastTrackID = id;
237  if (theTrack->GetKineticEnergy()/MeV > eMin) trkInfo->putInHistory();
238 #ifdef DebugLog
239  LogDebug("CaloSim") <<"CaloTrkProcessing: set ID on Calo " << ical
240  << " surface (Inside " << inside << ") to "
241  << id << " of a Track with Kinetic Energy "
242  << theTrack->GetKineticEnergy()/MeV << " MeV";
243 #endif
244  }
245  }
246  }
247  }
248 }
#define LogDebug(id)
int getIDonCaloSurface() const
int isItCalo(const G4VTouchable *)
const double MeV
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 71 of file CaloTrkProcessing.h.

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

double CaloTrkProcessing::eMin
private

Definition at line 69 of file CaloTrkProcessing.h.

Referenced by CaloTrkProcessing(), and update().

int CaloTrkProcessing::lastTrackID
private

Definition at line 70 of file CaloTrkProcessing.h.

Referenced by update().

const SimTrackManager* CaloTrkProcessing::m_trackManager
private

Definition at line 72 of file CaloTrkProcessing.h.

bool CaloTrkProcessing::putHistory
private

Definition at line 68 of file CaloTrkProcessing.h.

Referenced by CaloTrkProcessing(), and update().

bool CaloTrkProcessing::testBeam
private

Definition at line 68 of file CaloTrkProcessing.h.

Referenced by CaloTrkProcessing(), and update().