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