DD4hep  1.37.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;
44  Position global;
46  if( !sd.useVolumeManager() || !segmentation.isValid() ) {
47  pos = h.avgPosition();
48  global = h.localToGlobal(pos);
49  }
50  else if( !segmentation.cellsSpanVolumes() ) {
51  // Convert the position relative to the local readout volume
52  // to a global position.
53  pos = segmentation.position(cell);
54  global = h.localToGlobal(pos);
55  }
56  else {
57  // The segmentation can gang together multiple volumes.
58  // In this case, we can't use the transformation we get from
59  // the step --- the volume that actually contains the hit
60  // may not be the same volume that the segmentation uses
61  // for the local coordinate system. We need to get the
62  // actual volID used from the segmentation and then look
63  // it up the volume manager to get the proper transformation.
64  VolumeID volID = segmentation.volumeID(cell);
65  VolumeManager vman_glob = VolumeManager::getVolumeManager(sd.detectorDescription());
66  VolumeManager vman = vman_glob.subdetector(sd.id());
67  VolumeManagerContext* vc = vman.lookupContext(volID);
68  // explicit unit conversion; h.localToGlobal does it internally already
69  pos = segmentation.position(cell);
70  global = vc->localToWorld(Position(pos)) / dd4hep::mm;
71  }
72  hit = new Hit(global);
73  hit->cellID = cell;
74  coll.add(cell, hit);
75  Geant4TouchableHandler handler(h.touchable());
76  sd.printM2("%s> CREATE hit with deposit:%e MeV Pos:%8.2f %8.2f %8.2f %s [%s]",
77  sd.c_name(),contrib.deposit,pos.X,pos.Y,pos.Z,handler.path().c_str(),
78  coll.GetName().c_str());
79  if ( 0 == hit->cellID ) { // for debugging only!
80  hit->cellID = cell;
81  sd.except("+++ Invalid CELL ID for hit!");
82  }
83  }
84  hit->truth.emplace_back(contrib);
85  hit->energyDeposit += contrib.deposit;
86  }
87  }
88 
89  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
90  // Geant4SensitiveAction<Geant4VoidSensitive>
91  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
104  m_collectionID = -1;
105  }
106 
108  template <> bool
110  G4TouchableHistory* /* hist */) {
111  return true;
112  }
113 
115  template <> bool
117  G4TouchableHistory* /* hist */) {
118  return true;
119  }
121 
122  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
123  // Geant4SensitiveAction<Geant4Tracker>
124  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
136  m_collectionID = declareReadoutFilteredCollection<Geant4Tracker::Hit>();
137  }
138 
140  template <> bool
141  Geant4SensitiveAction<Geant4Tracker>::process(const G4Step* step,G4TouchableHistory* /* hist */) {
142  typedef Geant4Tracker::Hit Hit;
143  // Note: 1) We store in the hit the hit-direction, which is not the same as the track direction.
144  // 2) The energy deposit is the difference between incoming and outcoming energy.
145  Geant4StepHandler h(step);
146  auto contrib = Hit::extractContribution(step);
147  Direction hit_momentum = 0.5 * (h.preMom() + h.postMom());
148  double hit_deposit = h.deposit();
149  Hit* hit = new Hit(contrib, hit_momentum, hit_deposit);
150 
151  hit->cellID = cellID(step);
152  if ( 0 == hit->cellID ) {
153  hit->cellID = volumeID(step);
154  except("+++ Invalid CELL ID for hit!");
155  }
156  collection(m_collectionID)->add(hit);
157  mark(h.track);
158  print("Hit with deposit:%f Pos:%f %f %f ID=%016X",
159  hit->energyDeposit,hit->position.X(),hit->position.Y(),hit->position.Z(),(void*)hit->cellID);
160  Geant4TouchableHandler handler(step);
161  print(" Geant4 path:%s",handler.path().c_str());
162  return true;
163  }
164 
166  template <> bool
168  G4TouchableHistory* /* hist */)
169  {
170  typedef Geant4Tracker::Hit Hit;
171  Geant4FastSimHandler h(spot);
172  auto contrib = Hit::extractContribution(spot);
173  Hit* hit = new Hit(contrib, h.momentum(), h.deposit());
174 
175  hit->cellID = cellID(h.touchable(), h.avgPositionG4());
176  if ( 0 == hit->cellID ) {
177  hit->cellID = volumeID(h.touchable());
178  except("+++ Invalid CELL ID for hit!");
179  }
180  collection(m_collectionID)->add(hit);
181  mark(h.track);
182  print("Hit with deposit:%f Pos:%f %f %f ID=%016X",
183  hit->energyDeposit,hit->position.X(),hit->position.Y(),hit->position.Z(),(void*)hit->cellID);
184  Geant4TouchableHandler handler(h.touchable());
185  print(" Geant4 path:%s",handler.path().c_str());
186  return true;
187  }
188 
190 
191  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
192  // Geant4SensitiveAction<Geant4PhotonCounter>
193  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
204  struct Geant4OpticalTracker {};
206 
209  m_collectionID = declareReadoutFilteredCollection<Geant4Tracker::Hit>();
210  }
211 
213  template <> bool
214  Geant4SensitiveAction<Geant4OpticalTracker>::process(const G4Step* step,G4TouchableHistory* /* hist */) {
215  // Note: 1) We store in the hit the hit-direction, which is not the same as the track direction.
216  // 2) The energy deposit is the track momentum
217  typedef Geant4Tracker::Hit Hit;
218  Geant4StepHandler h(step);
219  auto contrib = Hit::extractContribution(step);
220  Direction hit_momentum = 0.5 * (h.preMom() + h.postMom());
221  double hit_deposit = contrib.deposit;
222  Hit* hit = new Hit(contrib, hit_momentum, hit_deposit);
223 
224  if (h.trackDef() == G4OpticalPhoton::OpticalPhotonDefinition()) {
225  step->GetTrack()->SetTrackStatus(fStopAndKill);
226  }
227  hit->cellID = cellID(step);
228  if ( 0 == hit->cellID ) {
229  hit->cellID = volumeID( step ) ;
230  except("+++ Invalid CELL ID for hit!");
231  }
232  collection(m_collectionID)->add(hit);
233  mark(h.track);
234  print("Hit with deposit:%f Pos:%f %f %f ID=%016X",
235  hit->energyDeposit,hit->position.X(),hit->position.Y(),hit->position.Z(),(void*)hit->cellID);
236  Geant4TouchableHandler handler(step);
237  print(" Geant4 path:%s",handler.path().c_str());
238  return true;
239  }
240 
242  template <> bool
244  G4TouchableHistory* /* hist */)
245  {
246  typedef Geant4Tracker::Hit Hit;
247  Geant4FastSimHandler h(spot);
248  auto contrib = Hit::extractContribution(spot);
249  Hit* hit = new Hit(contrib, h.momentum(), contrib.deposit);
250 
251  hit->cellID = cellID(h.touchable(), h.avgPositionG4());
252  if ( 0 == hit->cellID ) {
253  hit->cellID = volumeID(h.touchable());
254  except("+++ Invalid CELL ID for hit!");
255  }
256  collection(m_collectionID)->add(hit);
257  mark(h.track);
258  print("Hit with deposit:%f Pos:%f %f %f ID=%016X",
259  hit->energyDeposit,hit->position.X(),hit->position.Y(),hit->position.Z(),(void*)hit->cellID);
260  Geant4TouchableHandler handler(h.touchable());
261  print(" Geant4 path:%s",handler.path().c_str());
262  return true;
263  }
265 
266  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
267  // Geant4SensitiveAction<Calorimeter>
268  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
281  m_collectionID = declareReadoutFilteredCollection<Geant4Calorimeter::Hit>();
282  }
283 
285  template <> bool
286  Geant4SensitiveAction<Geant4Calorimeter>::process(const G4Step* step,G4TouchableHistory*) {
287  typedef Geant4Calorimeter::Hit Hit;
288  Geant4StepHandler h(step);
289  HitContribution contrib = Hit::extractContribution(step);
290  Geant4HitCollection* coll = collection(m_collectionID);
291  VolumeID cell = 0;
292 
293  try {
294  cell = cellID(step);
295  } catch(std::runtime_error &e) {
296  std::stringstream out;
297  out << std::setprecision(20) << std::scientific;
298  out << "ERROR: " << e.what() << std::endl;
299  out << "Position: "
300  << "Pre (" << std::setw(24) << step->GetPreStepPoint()->GetPosition() << ") "
301  << "Post (" << std::setw(24) << step->GetPostStepPoint()->GetPosition() << ") "
302  << std::endl;
303  out << "Momentum: "
304  << " Pre (" <<std::setw(24) << step->GetPreStepPoint() ->GetMomentum() << ") "
305  << " Post (" <<std::setw(24) << step->GetPostStepPoint()->GetMomentum() << ") "
306  << std::endl;
307  std::cout << out.str();
308  return true;
309  }
310 
311  handleCalorimeterHit(cell, contrib, *coll, h, *this, m_segmentation);
312  mark(h.track);
313  return true;
314  }
316  template <> bool
318  G4TouchableHistory* /* hist */)
319  {
320  typedef Geant4Calorimeter::Hit Hit;
321  Geant4FastSimHandler h(spot);
322  HitContribution contrib = Hit::extractContribution(spot);
323  Geant4HitCollection* coll = collection(m_collectionID);
324  VolumeID cell = 0;
325 
326  try {
327  cell = cellID(h.touchable(), h.avgPositionG4());
328  } catch(std::runtime_error &e) {
329  std::stringstream out;
330  out << std::setprecision(20) << std::scientific;
331  out << "ERROR: " << e.what() << std::endl;
332  out << "Position: (" << std::setw(24) << h.avgPositionG4() << ") " << std::endl;
333  out << "Momentum: (" << std::setw(24) << h.momentumG4() << ") " << std::endl;
334  std::cout << out.str();
335  return true;
336  }
337  handleCalorimeterHit(cell, contrib, *coll, h, *this, m_segmentation);
338  mark(h.track);
339  return true;
340  }
341 
343 
344  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
345  // Geant4SensitiveAction<OpticalCalorimeter>
346  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
347 
360  struct Geant4OpticalCalorimeter {};
362 
365  m_collectionID = declareReadoutFilteredCollection<Geant4Calorimeter::Hit>();
366  }
368  template <> bool
369  Geant4SensitiveAction<Geant4OpticalCalorimeter>::process(const G4Step* step,G4TouchableHistory*) {
370  G4Track * track = step->GetTrack();
371  // check that particle is optical photon:
372  if( track->GetDefinition() != G4OpticalPhoton::OpticalPhotonDefinition() ) {
373  return false;
374  }
375  else if ( track->GetCreatorProcess()->G4VProcess::GetProcessName() != "Cerenkov") {
376  track->SetTrackStatus(fStopAndKill);
377  return false;
378  }
379  else {
380  typedef Geant4Calorimeter::Hit Hit;
381  Geant4StepHandler h(step);
382  Geant4HitCollection* coll = collection(m_collectionID);
383  HitContribution contrib = Hit::extractContribution(step);
384  Position pos = h.prePos();
385  Hit* hit = coll->find<Hit>(PositionCompare<Hit,Position>(pos));
386  if ( !hit ) {
387  hit = new Hit(pos);
388  hit->cellID = volumeID(step);
389  coll->add(hit);
390  if ( 0 == hit->cellID ) {
391  hit->cellID = volumeID(step);
392  except("+++ Invalid CELL ID for hit!");
393  }
394  }
395  hit->energyDeposit += contrib.deposit;
396  hit->truth.emplace_back(contrib);
397  track->SetTrackStatus(fStopAndKill); // don't step photon any further
398  mark(h.track);
399  return true;
400  }
401  }
402 
404  template <> bool
406  G4TouchableHistory* /* hist */)
407  {
408  typedef Geant4Calorimeter::Hit Hit;
409  Geant4FastSimHandler h(spot);
410  const G4Track* track = h.track;
411  if( track->GetDefinition() != G4OpticalPhoton::OpticalPhotonDefinition() ) {
412  return false;
413  }
414  else if ( track->GetCreatorProcess()->G4VProcess::GetProcessName() != "Cerenkov") {
415  return false;
416  }
417  else {
418  Geant4HitCollection* coll = collection(m_collectionID);
419  HitContribution contrib = Hit::extractContribution(spot);
420  Position pos = h.avgPosition();
421  Hit* hit = coll->find<Hit>(PositionCompare<Hit,Position>(pos));
422  if ( !hit ) {
423  hit = new Hit(pos);
424  hit->cellID = volumeID(h.touchable());
425  coll->add(hit);
426  if ( 0 == hit->cellID ) {
427  hit->cellID = volumeID(h.touchable());
428  except("+++ Invalid CELL ID for hit!");
429  }
430  }
431  hit->energyDeposit += contrib.deposit;
432  hit->truth.emplace_back(contrib);
433  mark(h.track);
434  return true;
435  }
436  }
438 
439 
440  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
441  // Geant4SensitiveAction<ScintillatorCalorimeter>
442  // For scintillator with Geant4 BirksLaw effect
443  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
444 
445 
446  /* \addtogroup Geant4SDActionPlugin
447  *
448  * @{
449  * \package Geant4ScintillatorCalorimeterAction
450  * \brief Sensitive detector meant for scintillator calorimeters
451  *
452  * This sensitive action will apply Birks' law to the energy deposits
453  *
454  * @}
455  */
458 
461  m_collectionID = declareReadoutFilteredCollection<Geant4Calorimeter::Hit>();
462  }
464  template <> bool
465  Geant4SensitiveAction<Geant4ScintillatorCalorimeter>::process(const G4Step* step,G4TouchableHistory*) {
466  typedef Geant4Calorimeter::Hit Hit;
467  Geant4StepHandler h(step);
468  HitContribution contrib = Hit::extractContribution(step,true);
469  Geant4HitCollection* coll = collection(m_collectionID);
470  VolumeID cell = 0;
471  try {
472  cell = cellID(step);
473  } catch(std::runtime_error &e) {
474  std::stringstream out;
475  out << std::setprecision(20) << std::scientific;
476  out << "ERROR: " << e.what() << std::endl;
477  out << "Position: "
478  << "Pre (" << std::setw(24) << step->GetPreStepPoint()->GetPosition() << ") "
479  << "Post (" << std::setw(24) << step->GetPostStepPoint()->GetPosition() << ") "
480  << std::endl;
481  out << "Momentum: "
482  << " Pre (" <<std::setw(24) << step->GetPreStepPoint() ->GetMomentum() << ") "
483  << " Post (" <<std::setw(24) << step->GetPostStepPoint()->GetMomentum() << ") "
484  << std::endl;
485  std::cout << out.str();
486  return true;
487  }
488  handleCalorimeterHit(cell, contrib, *coll, h, *this, m_segmentation);
489  mark(h.track);
490  return true;
491  }
492 
494  template <> bool
496  G4TouchableHistory* /* hist */)
497  {
498  typedef Geant4Calorimeter::Hit Hit;
499  Geant4FastSimHandler h(spot);
500  HitContribution contrib = Hit::extractContribution(spot);
501  Geant4HitCollection* coll = collection(m_collectionID);
502  VolumeID cell = 0;
503 
504  try {
505  cell = cellID(h.touchable(), h.avgPositionG4());
506  } catch(std::runtime_error &e) {
507  std::stringstream out;
508  out << std::setprecision(20) << std::scientific;
509  out << "ERROR: " << e.what() << std::endl;
510  out << "Position: (" << std::setw(24) << h.avgPositionG4() << ") " << std::endl;
511  out << "Momentum: (" << std::setw(24) << h.momentumG4() << ") " << std::endl;
512  std::cout << out.str();
513  return true;
514  }
515  handleCalorimeterHit(cell, contrib, *coll, h, *this, m_segmentation);
516  mark(h.track);
517  return true;
518  }
519 
521 
535 
544  struct TrackerCombine {
548  G4VPhysicalVolume* firstSpotVolume { nullptr };
549  double mean_time { 0e0 };
550  double e_cut { 0e0 };
551  int current { -1 };
552  int combined { 0 };
553  long long int cell { 0 };
554 
556  }
557 
559  void start_collecting(const G4Track* track) {
560  pre.truth.deposit = 0.0;
562  sensitive->mark(track);
563  mean_pos.SetXYZ(0,0,0);
564  mean_time = 0;
565  post.copyFrom(pre);
566  combined = 0;
567  cell = 0;
568  }
569  void start(const G4Step* step, const G4StepPoint* point) {
570  pre.storePoint(step,point);
571  start_collecting(step->GetTrack());
572  firstSpotVolume = step->GetPreStepPoint()->GetTouchableHandle()->GetVolume();
573  }
574  void start(const Geant4FastSimSpot* spot) {
575  pre.storePoint(spot);
576  start_collecting(spot->primary);
577  firstSpotVolume = spot->volume();
578  }
579 
581  void update_collected_hit(const G4VTouchable* h, G4ThreeVector&& global) {
587  if ( 0 == cell ) {
588  cell = sensitive->cellID(h, global);
589  if ( 0 == cell ) {
590  cell = sensitive->volumeID(h) ;
591  sensitive->except("+++ Invalid CELL ID for hit!");
592  }
593  }
594  ++combined;
595  }
596  void update(const Geant4StepHandler& h) {
597  post.storePoint(h.step, h.post);
598  update_collected_hit(h.preTouchable(), h.avgPositionG4()); // Compute cellID
599  }
600  void update(const Geant4FastSimHandler& h) {
601  post.storePoint(h.spot);
602  update_collected_hit(h.touchable(), h.avgPositionG4()); // Compute cellID
603  }
604 
606  void clear() {
607  mean_pos.SetXYZ(0,0,0);
608  mean_time = 0;
609  post.clear();
610  pre.clear();
611  current = -1;
612  combined = 0;
613  cell = 0;
614  firstSpotVolume = nullptr;
615  }
616 
618  bool mustSaveTrack(const G4Track* tr) const {
619  return current > 0 && current != tr->GetTrackID();
620  }
621 
623  void extractHit(Geant4HitCollection* collection) {
624  if ( current == -1 ) {
625  return;
626  }
627  double depo = pre.truth.deposit;
628  double time = depo != 0 ? mean_time / depo : mean_time;
629  Position pos = depo != 0 ? mean_pos / depo : mean_pos;
630  Momentum mom = 0.5 * (pre.momentum + post.momentum);
631  double path_len = (post.position - pre.position).R();
633  depo, time, path_len, pos, mom);
634  hit->cellID = cell;
635  collection->add(hit);
636  sensitive->printM2("+++ TrackID:%6d [%s] CREATE hit combination with %2d deposit(s):"
637  " %e MeV Pos:%8.2f %8.2f %8.2f",
639  pos.X()/CLHEP::mm,pos.Y()/CLHEP::mm,pos.Z()/CLHEP::mm);
640  clear();
641  }
642 
643 
645  G4bool process(const G4Step* step, G4TouchableHistory* ) {
646  Geant4StepHandler h(step);
647  // std::cout << " process called - pre pos: " << h.prePos() << " post pos " << h.postPos()
648  // << " edep: " << h.deposit() << std::endl ;
649  void *prePV = h.volume(h.pre), *postPV = h.volume(h.post);
650 
653  if ( mustSaveTrack(h.track) ) {
654  extractHit(coll);
655  }
657  if ( current < 0 ) {
658  start(step, h.pre);
659  }
661  update(h);
662 
663  if ( prePV != postPV ) {
664  void* postSD = h.sd(h.post);
665  extractHit(coll);
666  if ( 0 != postSD ) {
667  void* preSD = h.sd(h.pre);
668  if ( preSD == postSD ) {
669  start(step,h.post);
670  }
671  }
672  }
673  else if ( h.track->GetTrackStatus() == fStopAndKill ) {
674  extractHit(coll);
675  }
676  return true;
677  }
678 
680  G4bool process(const Geant4FastSimSpot* spot, G4TouchableHistory* ) {
681  Geant4FastSimHandler h(spot);
682  G4VPhysicalVolume* prePV = firstSpotVolume, *postPV = h.volume();
685  if ( mustSaveTrack(h.track) ) {
686  extractHit(coll);
687  }
689  if ( current < 0 ) {
690  start(spot);
691  }
693  update(h);
694 
695  if ( firstSpotVolume && prePV != postPV ) {
696  void* postSD = h.sd();
697  extractHit(coll);
698  if ( 0 != postSD ) {
699  void* preSD = prePV ? prePV->GetLogicalVolume()->GetSensitiveDetector() : nullptr;
700  if ( preSD == postSD ) {
701  start(spot);
702  }
703  }
704  }
705  else if ( h.track->GetTrackStatus() == fStopAndKill ) {
706  extractHit(coll);
707  }
708  return true;
709  }
710 
712  void endEvent(const G4Event* /* event */) {
713  // We need to add the possibly last added hit to the collection here.
714  // otherwise the last hit would be assigned to the next event and the
715  // MC truth would be screwed.
716  //
717  // Alternatively the 'update' method would become rather CPU consuming,
718  // beacuse the extract action would have to be recalculated over and over.
719  if ( current > 0 ) {
721  extractHit(coll);
722  }
723  }
724  };
725 
728  eventAction().callAtEnd(&m_userData,&TrackerCombine::endEvent);
729  m_userData.e_cut = m_sensitive.energyCutoff();
730  m_userData.sensitive = this;
731  }
732 
735  m_collectionID = declareReadoutFilteredCollection<Geant4Tracker::Hit>();
736  }
737 
739  template <> void Geant4SensitiveAction<TrackerCombine>::clear(G4HCofThisEvent*) {
740  m_userData.clear();
741  }
742 
744  template <> G4bool
745  Geant4SensitiveAction<TrackerCombine>::process(const G4Step* step, G4TouchableHistory* history) {
746  return m_userData.process(step, history);
747  }
748 
750  template <> bool
752  G4TouchableHistory* history) {
753  return m_userData.process(spot, history);
754  }
755 
760  }
761 }
762 
763 using namespace dd4hep::sim;
764 
765 #include <DDG4/Factories.h>
766 // Special void entry point
768 // Standard factories used for simulation
775 
776 // Need these factories for backwards compatibility
dd4hep::sim::TrackerCombine::start
void start(const Geant4FastSimSpot *spot)
Definition: Geant4SDActions.cpp:574
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:437
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:712
dd4hep::sim::Geant4SimpleOpticalCalorimeterAction
Geant4OpticalCalorimeterAction Geant4SimpleOpticalCalorimeterAction
Definition: Geant4SDActions.cpp:759
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:361
dd4hep::sim::Geant4StepHandler::prePos
Position prePos() const
Returns the pre-step position.
Definition: Geant4StepHandler.h:84
dd4hep::sim::Geant4SimpleTrackerAction
Geant4TrackerAction Geant4SimpleTrackerAction
Definition: Geant4SDActions.cpp:757
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:342
dd4hep::sim::Momentum
Position Momentum
Definition: Defs.h:27
dd4hep::sim::Geant4SimpleCalorimeterAction
Geant4CalorimeterAction Geant4SimpleCalorimeterAction
Definition: Geant4SDActions.cpp:758
dd4hep::sim::TrackerCombine::update
void update(const Geant4FastSimHandler &h)
Definition: Geant4SDActions.cpp:600
dd4hep::sim::Geant4StepHandler::pre
G4StepPoint * pre
Definition: Geant4StepHandler.h:49
dd4hep::sim::Geant4TrackerAction
Geant4SensitiveAction< Geant4Tracker > Geant4TrackerAction
Definition: Geant4SDActions.cpp:189
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:596
dd4hep::sim::TrackerCombine::process
G4bool process(const G4Step *step, G4TouchableHistory *)
Method for generating hit(s) using the information of G4Step object.
Definition: Geant4SDActions.cpp:645
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:120
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:581
dd4hep::sim::Geant4HitData::cellID
long long int cellID
cellID
Definition: Geant4Data.h:124
dd4hep::sim::TrackerCombine::sensitive
Geant4Sensitive * sensitive
Definition: Geant4SDActions.cpp:547
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:296
dd4hep::sim::TrackerCombine::cell
long long int cell
Definition: Geant4SDActions.cpp:553
dd4hep::sim::TrackerCombine::TrackerCombine
TrackerCombine()
Definition: Geant4SDActions.cpp:555
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:247
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:569
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:520
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:264
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:284
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:606
dd4hep::sim::TrackerCombine::combined
int combined
Definition: Geant4SDActions.cpp:552
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:680
dd4hep::sim::Geant4Sensitive::id
int id() const
Get the detector identifier (DetElement::id())
Definition: Geant4SensDetAction.cpp:166
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:546
dd4hep::sim::Geant4SensitiveAction
Template class to ease the construction of sensitive detectors using particle template specialization...
Definition: Geant4SensDetAction.h:529
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:549
dd4hep::sim::TrackerCombine::pre
Geant4Tracker::Hit pre
Definition: Geant4SDActions.cpp:545
dd4hep::sim::Geant4TrackerCombineAction
Geant4SensitiveAction< TrackerCombine > Geant4TrackerCombineAction
Definition: Geant4SDActions.cpp:756
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:545
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:80
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:46
dd4hep::sim::Geant4StepHandler::step
const G4Step * step
Definition: Geant4StepHandler.h:48
dd4hep::sim::Geant4Sensitive::useVolumeManager
bool useVolumeManager() const
Access volume manager usage flag. This is a subdetector specific flag.
Definition: Geant4SensDetAction.h:183
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:323
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:548
dd4hep::sim::Geant4ScintillatorCalorimeter
Class to implement the standard sensitive detector for scintillator calorimeters.
Definition: Geant4SDActions.cpp:457
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:205
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:618
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:123
dd4hep::sim::Geant4FastSimHandler::volume
G4VPhysicalVolume * volume() const
Definition: Geant4FastSimHandler.h:93
dd4hep::sim::TrackerCombine::e_cut
double e_cut
Definition: Geant4SDActions.cpp:550
dd4hep::sim::TrackerCombine::current
int current
Definition: Geant4SDActions.cpp:551
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:559
dd4hep::sim::TrackerCombine::extractHit
void extractHit(Geant4HitCollection *collection)
Extract hit information and add the created hit to the collection.
Definition: Geant4SDActions.cpp:623
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 description object.
Definition: Geant4SensDetAction.cpp:237
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:544