DD4hep  1.30.0
Detector Description Toolkit for High Energy Physics
Geant4SDActions.cpp
Go to the documentation of this file.
1 //==========================================================================
2 // AIDA Detector description implementation
3 //--------------------------------------------------------------------------
4 // Copyright (C) Organisation europeenne pour la Recherche nucleaire (CERN)
5 // All rights reserved.
6 //
7 // For the licensing terms see $DD4hepINSTALL/LICENSE.
8 // For the list of contributors see $DD4hepINSTALL/doc/CREDITS.
9 //
10 // Author : M.Frank
11 //
12 //==========================================================================
13 
14 // Framework include files
15 #include <DDG4/Geant4SensDetAction.inl>
17 #include <DDG4/Geant4EventAction.h>
18 #include <G4OpticalPhoton.hh>
19 #include <G4VProcess.hh>
20 
21 
23 namespace dd4hep {
24 
26  namespace sim {
27 
28  namespace {
29  struct Geant4VoidSensitive {};
30 
32  template <class HANDLER>
33  void handleCalorimeterHit (VolumeID cell,
34  const HitContribution& contrib,
35  Geant4HitCollection& coll,
36  const HANDLER& h,
37  const Geant4Sensitive& sd,
38  const Segmentation& segmentation)
39  {
40  typedef Geant4Calorimeter::Hit Hit;
41  Hit* hit = coll.findByKey<Hit>(cell);
42  if ( !hit ) {
43  DDSegmentation::Vector3D pos = segmentation.position(cell);
44  Position global;
45 
46  // Convert the position relative to the local readout volume
47  // to a global position.
48  if (!segmentation.cellsSpanVolumes()) {
49  global = h.localToGlobal(pos);
50  }
51  else {
52  // The segmentation can gang together multiple volumes.
53  // In this case, we can't use the transformation we get from
54  // the step --- the volume that actually contains the hit
55  // may not be the same volume that the segmentation uses
56  // for the local coordinate system. We need to get the
57  // actual volID used from the segmentation and then look
58  // it up the volume manager to get the proper transformation.
59  VolumeID volID = segmentation.volumeID(cell);
60  VolumeManager vman = VolumeManager::getVolumeManager(sd.detectorDescription());
61  VolumeManagerContext* vc = vman.lookupContext(volID);
62  // explicit unit conversion; h.localToGlobal does it internally already
63  global = vc->localToWorld(Position(pos)) / dd4hep::mm;
64  }
65  hit = new Hit(global);
66  hit->cellID = cell;
67  coll.add(cell, hit);
68  Geant4TouchableHandler handler(h.touchable());
69  sd.printM2("%s> CREATE hit with deposit:%e MeV Pos:%8.2f %8.2f %8.2f %s [%s]",
70  sd.c_name(),contrib.deposit,pos.X,pos.Y,pos.Z,handler.path().c_str(),
71  coll.GetName().c_str());
72  if ( 0 == hit->cellID ) { // for debugging only!
73  hit->cellID = cell;
74  sd.except("+++ Invalid CELL ID for hit!");
75  }
76  }
77  hit->truth.emplace_back(contrib);
78  hit->energyDeposit += contrib.deposit;
79  }
80  }
81 
82  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
83  // Geant4SensitiveAction<Geant4VoidSensitive>
84  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
97  m_collectionID = -1;
98  }
99 
101  template <> bool
103  G4TouchableHistory* /* hist */) {
104  return true;
105  }
106 
108  template <> bool
110  G4TouchableHistory* /* hist */) {
111  return true;
112  }
114 
115  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
116  // Geant4SensitiveAction<Geant4Tracker>
117  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
129  m_collectionID = declareReadoutFilteredCollection<Geant4Tracker::Hit>();
130  }
131 
133  template <> bool
134  Geant4SensitiveAction<Geant4Tracker>::process(const G4Step* step,G4TouchableHistory* /* hist */) {
135  typedef Geant4Tracker::Hit Hit;
136  // Note: 1) We store in the hit the hit-direction, which is not the same as the track direction.
137  // 2) The energy deposit is the difference between incoming and outcoming energy.
138  Geant4StepHandler h(step);
139  auto contrib = Hit::extractContribution(step);
140  Direction hit_momentum = 0.5 * (h.preMom() + h.postMom());
141  double hit_deposit = h.deposit();
142  Hit* hit = new Hit(contrib, hit_momentum, hit_deposit);
143 
144  hit->cellID = cellID(step);
145  if ( 0 == hit->cellID ) {
146  hit->cellID = volumeID(step);
147  except("+++ Invalid CELL ID for hit!");
148  }
149  collection(m_collectionID)->add(hit);
150  mark(h.track);
151  print("Hit with deposit:%f Pos:%f %f %f ID=%016X",
152  hit->energyDeposit,hit->position.X(),hit->position.Y(),hit->position.Z(),(void*)hit->cellID);
153  Geant4TouchableHandler handler(step);
154  print(" Geant4 path:%s",handler.path().c_str());
155  return true;
156  }
157 
159  template <> bool
161  G4TouchableHistory* /* hist */)
162  {
163  typedef Geant4Tracker::Hit Hit;
164  Geant4FastSimHandler h(spot);
165  auto contrib = Hit::extractContribution(spot);
166  Hit* hit = new Hit(contrib, h.momentum(), h.deposit());
167 
168  hit->cellID = cellID(h.touchable(), h.avgPositionG4());
169  if ( 0 == hit->cellID ) {
170  hit->cellID = volumeID(h.touchable());
171  except("+++ Invalid CELL ID for hit!");
172  }
173  collection(m_collectionID)->add(hit);
174  mark(h.track);
175  print("Hit with deposit:%f Pos:%f %f %f ID=%016X",
176  hit->energyDeposit,hit->position.X(),hit->position.Y(),hit->position.Z(),(void*)hit->cellID);
177  Geant4TouchableHandler handler(h.touchable());
178  print(" Geant4 path:%s",handler.path().c_str());
179  return true;
180  }
181 
183 
184  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
185  // Geant4SensitiveAction<Geant4PhotonCounter>
186  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
197  struct Geant4OpticalTracker {};
199 
202  m_collectionID = declareReadoutFilteredCollection<Geant4Tracker::Hit>();
203  }
204 
206  template <> bool
207  Geant4SensitiveAction<Geant4OpticalTracker>::process(const G4Step* step,G4TouchableHistory* /* hist */) {
208  // Note: 1) We store in the hit the hit-direction, which is not the same as the track direction.
209  // 2) The energy deposit is the track momentum
210  typedef Geant4Tracker::Hit Hit;
211  Geant4StepHandler h(step);
212  auto contrib = Hit::extractContribution(step);
213  Direction hit_momentum = 0.5 * (h.preMom() + h.postMom());
214  double hit_deposit = contrib.deposit;
215  Hit* hit = new Hit(contrib, hit_momentum, hit_deposit);
216 
217  if (h.trackDef() == G4OpticalPhoton::OpticalPhotonDefinition()) {
218  step->GetTrack()->SetTrackStatus(fStopAndKill);
219  }
220  hit->cellID = cellID(step);
221  if ( 0 == hit->cellID ) {
222  hit->cellID = volumeID( step ) ;
223  except("+++ Invalid CELL ID for hit!");
224  }
225  collection(m_collectionID)->add(hit);
226  mark(h.track);
227  print("Hit with deposit:%f Pos:%f %f %f ID=%016X",
228  hit->energyDeposit,hit->position.X(),hit->position.Y(),hit->position.Z(),(void*)hit->cellID);
229  Geant4TouchableHandler handler(step);
230  print(" Geant4 path:%s",handler.path().c_str());
231  return true;
232  }
233 
235  template <> bool
237  G4TouchableHistory* /* hist */)
238  {
239  typedef Geant4Tracker::Hit Hit;
240  Geant4FastSimHandler h(spot);
241  auto contrib = Hit::extractContribution(spot);
242  Hit* hit = new Hit(contrib, h.momentum(), contrib.deposit);
243 
244  hit->cellID = cellID(h.touchable(), h.avgPositionG4());
245  if ( 0 == hit->cellID ) {
246  hit->cellID = volumeID(h.touchable());
247  except("+++ Invalid CELL ID for hit!");
248  }
249  collection(m_collectionID)->add(hit);
250  mark(h.track);
251  print("Hit with deposit:%f Pos:%f %f %f ID=%016X",
252  hit->energyDeposit,hit->position.X(),hit->position.Y(),hit->position.Z(),(void*)hit->cellID);
253  Geant4TouchableHandler handler(h.touchable());
254  print(" Geant4 path:%s",handler.path().c_str());
255  return true;
256  }
258 
259  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
260  // Geant4SensitiveAction<Calorimeter>
261  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
274  m_collectionID = declareReadoutFilteredCollection<Geant4Calorimeter::Hit>();
275  }
276 
278  template <> bool
279  Geant4SensitiveAction<Geant4Calorimeter>::process(const G4Step* step,G4TouchableHistory*) {
280  typedef Geant4Calorimeter::Hit Hit;
281  Geant4StepHandler h(step);
282  HitContribution contrib = Hit::extractContribution(step);
283  Geant4HitCollection* coll = collection(m_collectionID);
284  VolumeID cell = 0;
285 
286  try {
287  cell = cellID(step);
288  } catch(std::runtime_error &e) {
289  std::stringstream out;
290  out << std::setprecision(20) << std::scientific;
291  out << "ERROR: " << e.what() << std::endl;
292  out << "Position: "
293  << "Pre (" << std::setw(24) << step->GetPreStepPoint()->GetPosition() << ") "
294  << "Post (" << std::setw(24) << step->GetPostStepPoint()->GetPosition() << ") "
295  << std::endl;
296  out << "Momentum: "
297  << " Pre (" <<std::setw(24) << step->GetPreStepPoint() ->GetMomentum() << ") "
298  << " Post (" <<std::setw(24) << step->GetPostStepPoint()->GetMomentum() << ") "
299  << std::endl;
300  std::cout << out.str();
301  return true;
302  }
303 
304  handleCalorimeterHit(cell, contrib, *coll, h, *this, m_segmentation);
305  mark(h.track);
306  return true;
307  }
309  template <> bool
311  G4TouchableHistory* /* hist */)
312  {
313  typedef Geant4Calorimeter::Hit Hit;
314  Geant4FastSimHandler h(spot);
315  HitContribution contrib = Hit::extractContribution(spot);
316  Geant4HitCollection* coll = collection(m_collectionID);
317  VolumeID cell = 0;
318 
319  try {
320  cell = cellID(h.touchable(), h.avgPositionG4());
321  } catch(std::runtime_error &e) {
322  std::stringstream out;
323  out << std::setprecision(20) << std::scientific;
324  out << "ERROR: " << e.what() << std::endl;
325  out << "Position: (" << std::setw(24) << h.avgPositionG4() << ") " << std::endl;
326  out << "Momentum: (" << std::setw(24) << h.momentumG4() << ") " << std::endl;
327  std::cout << out.str();
328  return true;
329  }
330  handleCalorimeterHit(cell, contrib, *coll, h, *this, m_segmentation);
331  mark(h.track);
332  return true;
333  }
334 
336 
337  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
338  // Geant4SensitiveAction<OpticalCalorimeter>
339  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
340 
353  struct Geant4OpticalCalorimeter {};
355 
358  m_collectionID = declareReadoutFilteredCollection<Geant4Calorimeter::Hit>();
359  }
361  template <> bool
362  Geant4SensitiveAction<Geant4OpticalCalorimeter>::process(const G4Step* step,G4TouchableHistory*) {
363  G4Track * track = step->GetTrack();
364  // check that particle is optical photon:
365  if( track->GetDefinition() != G4OpticalPhoton::OpticalPhotonDefinition() ) {
366  return false;
367  }
368  else if ( track->GetCreatorProcess()->G4VProcess::GetProcessName() != "Cerenkov") {
369  track->SetTrackStatus(fStopAndKill);
370  return false;
371  }
372  else {
373  typedef Geant4Calorimeter::Hit Hit;
374  Geant4StepHandler h(step);
375  Geant4HitCollection* coll = collection(m_collectionID);
376  HitContribution contrib = Hit::extractContribution(step);
377  Position pos = h.prePos();
378  Hit* hit = coll->find<Hit>(PositionCompare<Hit,Position>(pos));
379  if ( !hit ) {
380  hit = new Hit(pos);
381  hit->cellID = volumeID(step);
382  coll->add(hit);
383  if ( 0 == hit->cellID ) {
384  hit->cellID = volumeID(step);
385  except("+++ Invalid CELL ID for hit!");
386  }
387  }
388  hit->energyDeposit += contrib.deposit;
389  hit->truth.emplace_back(contrib);
390  track->SetTrackStatus(fStopAndKill); // don't step photon any further
391  mark(h.track);
392  return true;
393  }
394  }
395 
397  template <> bool
399  G4TouchableHistory* /* hist */)
400  {
401  typedef Geant4Calorimeter::Hit Hit;
402  Geant4FastSimHandler h(spot);
403  const G4Track* track = h.track;
404  if( track->GetDefinition() != G4OpticalPhoton::OpticalPhotonDefinition() ) {
405  return false;
406  }
407  else if ( track->GetCreatorProcess()->G4VProcess::GetProcessName() != "Cerenkov") {
408  return false;
409  }
410  else {
411  Geant4HitCollection* coll = collection(m_collectionID);
412  HitContribution contrib = Hit::extractContribution(spot);
413  Position pos = h.avgPosition();
414  Hit* hit = coll->find<Hit>(PositionCompare<Hit,Position>(pos));
415  if ( !hit ) {
416  hit = new Hit(pos);
417  hit->cellID = volumeID(h.touchable());
418  coll->add(hit);
419  if ( 0 == hit->cellID ) {
420  hit->cellID = volumeID(h.touchable());
421  except("+++ Invalid CELL ID for hit!");
422  }
423  }
424  hit->energyDeposit += contrib.deposit;
425  hit->truth.emplace_back(contrib);
426  mark(h.track);
427  return true;
428  }
429  }
431 
432 
433  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
434  // Geant4SensitiveAction<ScintillatorCalorimeter>
435  // For scintillator with Geant4 BirksLaw effect
436  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
437 
438 
439  /* \addtogroup Geant4SDActionPlugin
440  *
441  * @{
442  * \package Geant4ScintillatorCalorimeterAction
443  * \brief Sensitive detector meant for scintillator calorimeters
444  *
445  * This sensitive action will apply Birks' law to the energy deposits
446  *
447  * @}
448  */
451 
454  m_collectionID = declareReadoutFilteredCollection<Geant4Calorimeter::Hit>();
455  }
457  template <> bool
458  Geant4SensitiveAction<Geant4ScintillatorCalorimeter>::process(const G4Step* step,G4TouchableHistory*) {
459  typedef Geant4Calorimeter::Hit Hit;
460  Geant4StepHandler h(step);
461  HitContribution contrib = Hit::extractContribution(step,true);
462  Geant4HitCollection* coll = collection(m_collectionID);
463  VolumeID cell = 0;
464  try {
465  cell = cellID(step);
466  } catch(std::runtime_error &e) {
467  std::stringstream out;
468  out << std::setprecision(20) << std::scientific;
469  out << "ERROR: " << e.what() << std::endl;
470  out << "Position: "
471  << "Pre (" << std::setw(24) << step->GetPreStepPoint()->GetPosition() << ") "
472  << "Post (" << std::setw(24) << step->GetPostStepPoint()->GetPosition() << ") "
473  << std::endl;
474  out << "Momentum: "
475  << " Pre (" <<std::setw(24) << step->GetPreStepPoint() ->GetMomentum() << ") "
476  << " Post (" <<std::setw(24) << step->GetPostStepPoint()->GetMomentum() << ") "
477  << std::endl;
478  std::cout << out.str();
479  return true;
480  }
481  handleCalorimeterHit(cell, contrib, *coll, h, *this, m_segmentation);
482  mark(h.track);
483  return true;
484  }
485 
487  template <> bool
489  G4TouchableHistory* /* hist */)
490  {
491  typedef Geant4Calorimeter::Hit Hit;
492  Geant4FastSimHandler h(spot);
493  HitContribution contrib = Hit::extractContribution(spot);
494  Geant4HitCollection* coll = collection(m_collectionID);
495  VolumeID cell = 0;
496 
497  try {
498  cell = cellID(h.touchable(), h.avgPositionG4());
499  } catch(std::runtime_error &e) {
500  std::stringstream out;
501  out << std::setprecision(20) << std::scientific;
502  out << "ERROR: " << e.what() << std::endl;
503  out << "Position: (" << std::setw(24) << h.avgPositionG4() << ") " << std::endl;
504  out << "Momentum: (" << std::setw(24) << h.momentumG4() << ") " << std::endl;
505  std::cout << out.str();
506  return true;
507  }
508  handleCalorimeterHit(cell, contrib, *coll, h, *this, m_segmentation);
509  mark(h.track);
510  return true;
511  }
512 
514 
528 
537  struct TrackerCombine {
541  G4VPhysicalVolume* firstSpotVolume { nullptr };
542  double mean_time { 0e0 };
543  double e_cut { 0e0 };
544  int current { -1 };
545  int combined { 0 };
546  long long int cell { 0 };
547 
549  }
550 
552  void start_collecting(const G4Track* track) {
553  pre.truth.deposit = 0.0;
555  sensitive->mark(track);
556  mean_pos.SetXYZ(0,0,0);
557  mean_time = 0;
558  post.copyFrom(pre);
559  combined = 0;
560  cell = 0;
561  }
562  void start(const G4Step* step, const G4StepPoint* point) {
563  pre.storePoint(step,point);
564  start_collecting(step->GetTrack());
565  firstSpotVolume = step->GetPreStepPoint()->GetTouchableHandle()->GetVolume();
566  }
567  void start(const Geant4FastSimSpot* spot) {
568  pre.storePoint(spot);
569  start_collecting(spot->primary);
570  firstSpotVolume = spot->volume();
571  }
572 
574  void update_collected_hit(const G4VTouchable* h, G4ThreeVector&& global) {
580  if ( 0 == cell ) {
581  cell = sensitive->cellID(h, global);
582  if ( 0 == cell ) {
583  cell = sensitive->volumeID(h) ;
584  sensitive->except("+++ Invalid CELL ID for hit!");
585  }
586  }
587  ++combined;
588  }
589  void update(const Geant4StepHandler& h) {
590  post.storePoint(h.step, h.post);
591  update_collected_hit(h.preTouchable(), h.avgPositionG4()); // Compute cellID
592  }
593  void update(const Geant4FastSimHandler& h) {
594  post.storePoint(h.spot);
595  update_collected_hit(h.touchable(), h.avgPositionG4()); // Compute cellID
596  }
597 
599  void clear() {
600  mean_pos.SetXYZ(0,0,0);
601  mean_time = 0;
602  post.clear();
603  pre.clear();
604  current = -1;
605  combined = 0;
606  cell = 0;
607  firstSpotVolume = nullptr;
608  }
609 
611  bool mustSaveTrack(const G4Track* tr) const {
612  return current > 0 && current != tr->GetTrackID();
613  }
614 
616  void extractHit(Geant4HitCollection* collection) {
617  if ( current == -1 ) {
618  return;
619  }
620  double depo = pre.truth.deposit;
621  double time = depo != 0 ? mean_time / depo : mean_time;
622  Position pos = depo != 0 ? mean_pos / depo : mean_pos;
623  Momentum mom = 0.5 * (pre.momentum + post.momentum);
624  double path_len = (post.position - pre.position).R();
626  depo, time, path_len, pos, mom);
627  hit->cellID = cell;
628  collection->add(hit);
629  sensitive->printM2("+++ TrackID:%6d [%s] CREATE hit combination with %2d deposit(s):"
630  " %e MeV Pos:%8.2f %8.2f %8.2f",
632  pos.X()/CLHEP::mm,pos.Y()/CLHEP::mm,pos.Z()/CLHEP::mm);
633  clear();
634  }
635 
636 
638  G4bool process(const G4Step* step, G4TouchableHistory* ) {
639  Geant4StepHandler h(step);
640  // std::cout << " process called - pre pos: " << h.prePos() << " post pos " << h.postPos()
641  // << " edep: " << h.deposit() << std::endl ;
642  void *prePV = h.volume(h.pre), *postPV = h.volume(h.post);
643 
646  if ( mustSaveTrack(h.track) ) {
647  extractHit(coll);
648  }
650  if ( current < 0 ) {
651  start(step, h.pre);
652  }
654  update(h);
655 
656  if ( prePV != postPV ) {
657  void* postSD = h.sd(h.post);
658  extractHit(coll);
659  if ( 0 != postSD ) {
660  void* preSD = h.sd(h.pre);
661  if ( preSD == postSD ) {
662  start(step,h.post);
663  }
664  }
665  }
666  else if ( h.track->GetTrackStatus() == fStopAndKill ) {
667  extractHit(coll);
668  }
669  return true;
670  }
671 
673  G4bool process(const Geant4FastSimSpot* spot, G4TouchableHistory* ) {
674  Geant4FastSimHandler h(spot);
675  G4VPhysicalVolume* prePV = firstSpotVolume, *postPV = h.volume();
678  if ( mustSaveTrack(h.track) ) {
679  extractHit(coll);
680  }
682  if ( current < 0 ) {
683  start(spot);
684  }
686  update(h);
687 
688  if ( firstSpotVolume && prePV != postPV ) {
689  void* postSD = h.sd();
690  extractHit(coll);
691  if ( 0 != postSD ) {
692  void* preSD = prePV ? prePV->GetLogicalVolume()->GetSensitiveDetector() : nullptr;
693  if ( preSD == postSD ) {
694  start(spot);
695  }
696  }
697  }
698  else if ( h.track->GetTrackStatus() == fStopAndKill ) {
699  extractHit(coll);
700  }
701  return true;
702  }
703 
705  void endEvent(const G4Event* /* event */) {
706  // We need to add the possibly last added hit to the collection here.
707  // otherwise the last hit would be assigned to the next event and the
708  // MC truth would be screwed.
709  //
710  // Alternatively the 'update' method would become rather CPU consuming,
711  // beacuse the extract action would have to be recalculated over and over.
712  if ( current > 0 ) {
714  extractHit(coll);
715  }
716  }
717  };
718 
721  eventAction().callAtEnd(&m_userData,&TrackerCombine::endEvent);
722  m_userData.e_cut = m_sensitive.energyCutoff();
723  m_userData.sensitive = this;
724  }
725 
728  m_collectionID = declareReadoutFilteredCollection<Geant4Tracker::Hit>();
729  }
730 
732  template <> void Geant4SensitiveAction<TrackerCombine>::clear(G4HCofThisEvent*) {
733  m_userData.clear();
734  }
735 
737  template <> G4bool
738  Geant4SensitiveAction<TrackerCombine>::process(const G4Step* step, G4TouchableHistory* history) {
739  return m_userData.process(step, history);
740  }
741 
743  template <> bool
745  G4TouchableHistory* history) {
746  return m_userData.process(spot, history);
747  }
748 
753  }
754 }
755 
756 using namespace dd4hep::sim;
757 
758 #include <DDG4/Factories.h>
759 // Special void entry point
761 // Standard factories used for simulation
768 
769 // Need these factories for backwards compatibility
dd4hep::sim::TrackerCombine::start
void start(const Geant4FastSimSpot *spot)
Definition: Geant4SDActions.cpp:567
DECLARE_GEANT4SENSITIVE
#define DECLARE_GEANT4SENSITIVE(name)
Definition: Factories.h:204
dd4hep::sim::Geant4Tracker::Hit::momentum
Direction momentum
Hit direction.
Definition: Geant4Data.h:270
dd4hep::sim::Geant4OpticalCalorimeterAction
Geant4SensitiveAction< Geant4OpticalCalorimeter > Geant4OpticalCalorimeterAction
Definition: Geant4SDActions.cpp:430
dd4hep::sim::Geant4Tracker::Hit::position
Position position
Hit position.
Definition: Geant4Data.h:268
dd4hep::sim::TrackerCombine::endEvent
void endEvent(const G4Event *)
Post-event action callback.
Definition: Geant4SDActions.cpp:705
dd4hep::sim::Geant4SimpleOpticalCalorimeterAction
Geant4OpticalCalorimeterAction Geant4SimpleOpticalCalorimeterAction
Definition: Geant4SDActions.cpp:752
dd4hep::sim::Geant4TouchableHandler
Helper class to ease the extraction of information from a G4Touchable object.
Definition: Geant4TouchableHandler.h:44
dd4hep::sim::Geant4SensitiveAction::defineCollections
virtual void defineCollections()
Define collections created by this sensitivie action object.
dd4hep::sim::Geant4OpticalCalorimeter
Helper class to define properties of optical calorimeters. UNTESTED.
Definition: Geant4SDActions.cpp:354
dd4hep::sim::Geant4StepHandler::prePos
Position prePos() const
Returns the pre-step position.
Definition: Geant4StepHandler.h:84
dd4hep::sim::Geant4SimpleTrackerAction
Geant4TrackerAction Geant4SimpleTrackerAction
Definition: Geant4SDActions.cpp:750
dd4hep::sim::Geant4HitData::MonteCarloContrib::pdgID
int pdgID
Particle ID from the PDG table.
Definition: Geant4Data.h:143
dd4hep::sim::Geant4FastSimHandler::momentumG4
G4ThreeVector momentumG4() const
Returns the track momentum as a G4ThreeVector.
Definition: Geant4FastSimHandler.h:86
dd4hep::sim::Geant4CalorimeterAction
Geant4SensitiveAction< Geant4Calorimeter > Geant4CalorimeterAction
Definition: Geant4SDActions.cpp:335
dd4hep::sim::Momentum
Position Momentum
Definition: Defs.h:27
dd4hep::sim::Geant4SimpleCalorimeterAction
Geant4CalorimeterAction Geant4SimpleCalorimeterAction
Definition: Geant4SDActions.cpp:751
dd4hep::sim::TrackerCombine::update
void update(const Geant4FastSimHandler &h)
Definition: Geant4SDActions.cpp:593
dd4hep::sim::Geant4StepHandler::pre
G4StepPoint * pre
Definition: Geant4StepHandler.h:49
dd4hep::sim::Geant4TrackerAction
Geant4SensitiveAction< Geant4Tracker > Geant4TrackerAction
Definition: Geant4SDActions.cpp:182
dd4hep::sim::Geant4StepHandler::preMom
Momentum preMom() const
Definition: Geant4StepHandler.h:127
dd4hep::sim::Geant4SensitiveAction::initialize
virtual void initialize() final
Initialization overload for specialization.
dd4hep::sim::Geant4FastSimHandler
Helper class to ease the extraction of information from a G4FastSimSpot object.
Definition: Geant4FastSimHandler.h:43
dd4hep::sim::Geant4HitCollection
Generic hit container class using Geant4HitWrapper objects.
Definition: Geant4HitCollection.h:201
Geant4EventAction.h
Geant4TrackerCombineAction
Sensitive detector meant for tracking detectors will combine multiple steps of the same track in the ...
dd4hep::sim::TrackerCombine::update
void update(const Geant4StepHandler &h)
Definition: Geant4SDActions.cpp:589
dd4hep::sim::TrackerCombine::process
G4bool process(const G4Step *step, G4TouchableHistory *)
Method for generating hit(s) using the information of G4Step object.
Definition: Geant4SDActions.cpp:638
dd4hep::sim::Geant4StepHandler::postMom
Momentum postMom() const
Definition: Geant4StepHandler.h:131
dd4hep::sim::Geant4HitCollection::find
TYPE * find(const Compare &cmp)
Find hits in a collection by comparison of attributes.
Definition: Geant4HitCollection.h:350
dd4hep::sim::Geant4Tracker::Hit::copyFrom
void copyFrom(const Hit &c)
Explicit assignment operation.
Definition: Geant4Data.cpp:133
dd4hep::sim::Geant4VoidSensitiveAction
Geant4SensitiveAction< Geant4VoidSensitive > Geant4VoidSensitiveAction
Definition: Geant4SDActions.cpp:113
Geant4TrackerAction
Sensitive detector meant for tracking detectors, will produce one hit per step.
dd4hep::sim::Direction
Position Direction
Definition: Defs.h:26
dd4hep::sim::TrackerCombine::update_collected_hit
void update_collected_hit(const G4VTouchable *h, G4ThreeVector &&global)
Update energy and track information during hit info accumulation.
Definition: Geant4SDActions.cpp:574
dd4hep::sim::Geant4HitData::cellID
long long int cellID
cellID
Definition: Geant4Data.h:124
dd4hep::sim::TrackerCombine::sensitive
Geant4Sensitive * sensitive
Definition: Geant4SDActions.cpp:540
dd4hep::sim::Geant4StepHandler::volume
G4VPhysicalVolume * volume(const G4StepPoint *p) const
Definition: Geant4StepHandler.h:151
dd4hep::sim::Geant4HitHandler::trackDef
G4ParticleDefinition * trackDef() const
Definition: Geant4HitHandler.h:66
dd4hep::sim::Geant4FastSimSpot
Spot definition for fast simulation and GFlash.
Definition: Geant4FastSimSpot.h:71
dd4hep::sim::Geant4Sensitive::volumeID
long long int volumeID(const G4Step *step)
Returns the volumeID of the sensitive volume corresponding to the step.
Definition: Geant4SensDetAction.cpp:245
dd4hep::sim::TrackerCombine::cell
long long int cell
Definition: Geant4SDActions.cpp:546
dd4hep::sim::TrackerCombine::TrackerCombine
TrackerCombine()
Definition: Geant4SDActions.cpp:548
dd4hep::sim::Geant4Sensitive::collection
Geant4HitCollection * collection(std::size_t which)
Retrieve the hits collection associated with this detector by its serial number.
Definition: Geant4SensDetAction.cpp:196
dd4hep::sim::Geant4SensitiveAction::process
virtual bool process(const G4Step *step, G4TouchableHistory *history) final
G4VSensitiveDetector interface: Method for generating hit(s) using the G4Step object.
dd4hep::sim::Geant4StepHandler
Helper class to ease the extraction of information from a G4Step object.
Definition: Geant4StepHandler.h:46
dd4hep::sim::TrackerCombine::start
void start(const G4Step *step, const G4StepPoint *point)
Definition: Geant4SDActions.cpp:562
dd4hep::sim::Geant4SensitiveAction::clear
virtual void clear(G4HCofThisEvent *hce) final
G4VSensitiveDetector interface: Method invoked if the event was aborted.
dd4hep::sim::PositionCompare
Specialized hit selector based on the hit's position.
Definition: Geant4HitCollection.h:395
dd4hep::sim::Geant4Action::except
void except(const char *fmt,...) const
Support of exceptions: Print fatal message and throw runtime_error.
Definition: Geant4Action.cpp:256
dd4hep::sim::Geant4ScintillatorCalorimeterAction
Geant4SensitiveAction< Geant4ScintillatorCalorimeter > Geant4ScintillatorCalorimeterAction
Definition: Geant4SDActions.cpp:513
dd4hep::sim::Geant4FastSimHandler::deposit
double deposit() const
Access the enery deposit of the energy spot.
Definition: Geant4FastSimHandler.h:90
dd4hep::sim::Geant4StepHandler::avgPositionG4
G4ThreeVector avgPositionG4() const
Average position vector of the step.
Definition: Geant4StepHandler.h:102
dd4hep::sim::Geant4StepHandler::post
G4StepPoint * post
Definition: Geant4StepHandler.h:50
dd4hep::sim::Geant4OpticalTrackerAction
Geant4SensitiveAction< Geant4OpticalTracker > Geant4OpticalTrackerAction
Definition: Geant4SDActions.cpp:257
dd4hep::sim::Geant4FastSimSpot::primary
const G4Track * primary
Definition: Geant4FastSimSpot.h:120
dd4hep::sim::Geant4Sensitive::mark
void mark(const G4Track *track) const
Mark the track to be kept for MC truth propagation during hit processing.
Definition: Geant4SensDetAction.cpp:233
dd4hep::sim::Geant4Tracker::Hit
DDG4 tracker hit class used by the generic DDG4 tracker sensitive detector.
Definition: Geant4Data.h:263
dd4hep::sim::Geant4HitHandler::touchable
const G4VTouchable * touchable() const
Definition: Geant4HitHandler.h:63
dd4hep::VolumeManager::getVolumeManager
static VolumeManager getVolumeManager(const Detector &description)
static accessor calling DD4hepVolumeManagerPlugin if necessary
Definition: VolumeManager.cpp:420
dd4hep::sim::TrackerCombine::clear
void clear()
Clear collected information and restart for new hit.
Definition: Geant4SDActions.cpp:599
dd4hep::sim::TrackerCombine::combined
int combined
Definition: Geant4SDActions.cpp:545
dd4hep::sim::Geant4HitData::MonteCarloContrib::deposit
double deposit
Total energy deposit in this hit.
Definition: Geant4Data.h:145
dd4hep::sim::TrackerCombine::process
G4bool process(const Geant4FastSimSpot *spot, G4TouchableHistory *)
Method for generating hit(s) using the information of fast simulation spot object.
Definition: Geant4SDActions.cpp:673
dd4hep::sim::Geant4HitCollection::findByKey
TYPE * findByKey(VolumeID key)
Find hits in a collection by comparison of key value.
Definition: Geant4HitCollection.h:354
Geant4FastSimHandler.h
dd4hep::sim::Geant4FastSimHandler::avgPositionG4
G4ThreeVector avgPositionG4() const
Returns the post-step position as a G4ThreeVector.
Definition: Geant4FastSimHandler.h:77
dd4hep::sim::TrackerCombine::mean_pos
Position mean_pos
Definition: Geant4SDActions.cpp:539
dd4hep::sim::Geant4SensitiveAction
Template class to ease the construction of sensitive detectors using particle template specialization...
Definition: Geant4SensDetAction.h:514
dd4hep::sim::Geant4Tracker::Hit::truth
Contribution truth
Monte Carlo / Geant4 information.
Definition: Geant4Data.h:276
Geant4CalorimeterAction
Sensitive detector meant for calorimeters.
dd4hep::sim::TrackerCombine::mean_time
double mean_time
Definition: Geant4SDActions.cpp:542
dd4hep::sim::TrackerCombine::pre
Geant4Tracker::Hit pre
Definition: Geant4SDActions.cpp:538
dd4hep::sim::Geant4TrackerCombineAction
Geant4SensitiveAction< TrackerCombine > Geant4TrackerCombineAction
Definition: Geant4SDActions.cpp:749
dd4hep::sim::Geant4Calorimeter::Hit
DDG4 calorimeter hit class used by the generic DDG4 calorimeter sensitive detector.
Definition: Geant4Data.h:323
dd4hep::sim::Geant4HitData::MonteCarloContrib::trackID
int trackID
Geant 4 Track identifier.
Definition: Geant4Data.h:141
dd4hep::sim::TrackerCombine::post
Geant4Tracker::Hit post
Definition: Geant4SDActions.cpp:538
dd4hep::sim::Geant4FastSimHandler::spot
const Geant4FastSimSpot * spot
Definition: Geant4FastSimHandler.h:45
dd4hep::sim::Geant4StepHandler::deposit
double deposit() const
Definition: Geant4StepHandler.h:135
dd4hep::Position
ROOT::Math::XYZVector Position
Definition: Objects.h:81
Factories.h
dd4hep::sim::Geant4FastSimSpot::volume
G4VPhysicalVolume * volume() const
Access the physical volume of the hit.
Definition: Geant4FastSimSpot.h:87
dd4hep::sim::Geant4FastSimHandler::sd
G4VSensitiveDetector * sd() const
Definition: Geant4FastSimHandler.h:102
VolumeID
dd4hep::DDSegmentation::VolumeID VolumeID
Definition: SegmentationDictionary.h:50
dd4hep::sim
Namespace for the Geant4 based simulation part of the AIDA detector description toolkit.
Definition: EDM4hepFileReader.cpp:41
dd4hep::sim::Geant4StepHandler::step
const G4Step * step
Definition: Geant4StepHandler.h:48
dd4hep::sim::Geant4Action::c_name
const char * c_name() const
Access name of the action.
Definition: Geant4Action.h:284
dd4hep::sim::Geant4Sensitive::cellID
long long int cellID(const G4Step *step)
Returns the cellID of the sensitive volume corresponding to the step.
Definition: Geant4SensDetAction.cpp:260
dd4hep::sim::Geant4TouchableHandler::path
std::string path() const
Helper: Access the placement path of a Geant4 touchable object as a string.
Definition: Geant4TouchableHandler.cpp:74
dd4hep::sim::TrackerCombine::firstSpotVolume
G4VPhysicalVolume * firstSpotVolume
Definition: Geant4SDActions.cpp:541
dd4hep::sim::Geant4ScintillatorCalorimeter
Class to implement the standard sensitive detector for scintillator calorimeters.
Definition: Geant4SDActions.cpp:450
dd4hep
Namespace for the AIDA detector description toolkit.
Definition: AlignmentsCalib.h:28
dd4hep::sim::Geant4FastSimHandler::avgPosition
Position avgPosition() const
Returns the pre-gflashspot position.
Definition: Geant4FastSimHandler.h:72
dd4hep::sim::Geant4Tracker::Hit::storePoint
Hit & storePoint(const G4Step *step, const G4StepPoint *point)
Store Geant4 point and step information into tracker hit structure.
Definition: Geant4Data.cpp:154
dd4hep::sim::Geant4OpticalTracker
Helper class to define properties of optical trackers. UNTESTED.
Definition: Geant4SDActions.cpp:198
dd4hep::sim::Geant4FastSimHandler::momentum
Momentum momentum() const
Return track momentum in DD4hep notation.
Definition: Geant4FastSimHandler.h:67
dd4hep::sim::TrackerCombine::mustSaveTrack
bool mustSaveTrack(const G4Track *tr) const
Helper function to decide if the hit has to be extracted and saved in the collection.
Definition: Geant4SDActions.cpp:611
dd4hep::sim::Geant4HitCollection::add
void add(TYPE *hit_pointer)
Add a new hit with a check, that the hit is of the same type.
Definition: Geant4HitCollection.h:333
dd4hep::sim::Geant4Sensitive
The base class for Geant4 sensitive detector actions implemented by users.
Definition: Geant4SensDetAction.h:121
dd4hep::sim::Geant4FastSimHandler::volume
G4VPhysicalVolume * volume() const
Definition: Geant4FastSimHandler.h:93
dd4hep::sim::TrackerCombine::e_cut
double e_cut
Definition: Geant4SDActions.cpp:543
dd4hep::sim::TrackerCombine::current
int current
Definition: Geant4SDActions.cpp:544
dd4hep::sim::Geant4HitHandler::track
const G4Track * track
Definition: Geant4HitHandler.h:45
dd4hep::sim::Geant4StepHandler::preTouchable
const G4VTouchable * preTouchable() const
Definition: Geant4StepHandler.h:141
dd4hep::sim::Geant4StepHandler::sd
G4VSensitiveDetector * sd(const G4StepPoint *p) const
Definition: Geant4StepHandler.h:164
Geant4OpticalCalorimeterAction
Sensitive detector meant for optical calorimeters.
dd4hep::sim::TrackerCombine::start_collecting
void start_collecting(const G4Track *track)
Start a new hit.
Definition: Geant4SDActions.cpp:552
dd4hep::sim::TrackerCombine::extractHit
void extractHit(Geant4HitCollection *collection)
Extract hit information and add the created hit to the collection.
Definition: Geant4SDActions.cpp:616
dd4hep::sim::Geant4Tracker::Hit::clear
Hit & clear()
Clear hit content.
Definition: Geant4Data.cpp:144
dd4hep::sim::Geant4SensitiveAction::processFastSim
virtual bool processFastSim(const Geant4FastSimSpot *spot, G4TouchableHistory *history) final
GFLASH/FastSim interface: Method for generating hit(s) using the information of the fast simulation s...
dd4hep::sim::Geant4HitData::MonteCarloContrib::time
double time
Timestamp when this energy was deposited.
Definition: Geant4Data.h:147
dd4hep::sim::Geant4Action::printM2
void printM2(const char *fmt,...) const
Support for messages with variable output level using output level-2.
Definition: Geant4Action.cpp:166
Geant4OpticalTrackerAction
Sensitive detector meant for photon detectors, will produce one hit per step for regular particles,...
dd4hep::sim::Geant4Sensitive::detectorDescription
Detector & detectorDescription() const
Access the detector desciption object.
Definition: Geant4SensDetAction.cpp:186
Geant4VoidSensitiveAction
Void Sensitive detector action to skip the processing of a detector without changing the entire DDG4 ...
dd4hep::sim::TrackerCombine
Geant4 sensitive detector combining all deposits of one G4Track within one sensitive element.
Definition: Geant4SDActions.cpp:537