MyTrackerSDAction.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 "MyTrackerHit.h"
16 #include "DDG4/Geant4SensDetAction.inl"
17 #include "DDG4/Factories.h"
18 
19 //using namespace CLHEP;
20 
21 namespace SomeExperiment {
22 
23  class MyTrackerSD {
24  public:
25  typedef MyTrackerHit Hit;
26  // If we need special data to personalize the action, be put it here
27  int mumDeposits = 0;
28  double integratedDeposit = 0;
29  };
30 }
31 
33 namespace dd4hep {
34 
36  namespace sim {
37 
38  using namespace SomeExperiment;
39 
40  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
41  // Geant4SensitiveAction<MyTrackerSD>
42  // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
52  template <> void Geant4SensitiveAction<MyTrackerSD>::defineCollections() {
54  m_collectionID = declareReadoutFilteredCollection<MyTrackerSD::Hit>();
55  }
56 
58  template <> bool Geant4SensitiveAction<MyTrackerSD>::process(G4Step* step,G4TouchableHistory* /*hist*/ ) {
59  Geant4StepHandler h(step);
60  Position prePos = h.prePos();
61  Position postPos = h.postPos();
62  Position direction = postPos - prePos;
63  Position position = mean_direction(prePos,postPos);
64  double hit_len = direction.R();
65 
66  // Somehow extract here the physics you want
67  MyTrackerSD::Hit* hit = new MyTrackerSD::Hit(h.trkID(), h.trkPdgID(), h.deposit(), h.track->GetGlobalTime());
68  Geant4HitData::MonteCarloContrib contrib = Geant4HitData::extractContribution(step);
69  hit->cellID = cellID(step);
70  hit->energyDeposit = contrib.deposit;
71  hit->position = position;
72  hit->momentum = 0.5*(h. preMom() + h.postMom());
73  hit->length = hit_len;
74  collection(m_collectionID)->add(hit);
75  mark(h.track);
76  if ( 0 == hit->cellID ) {
77  hit->cellID = volumeID(step);
78  except("+++ Invalid CELL ID for hit!");
79  }
80  printP1("Hit with deposit:%f Pos:%f %f %f ID=%016X",
81  step->GetTotalEnergyDeposit(),position.X(),position.Y(),position.Z(),
82  (void*)hit->cellID);
83  Geant4TouchableHandler handler(step);
84  print(" Geant4 path:%s",handler.path().c_str());
85 
86  // Do something with my personal data (can be also something more clever ;-):
87  m_userData.integratedDeposit += contrib.deposit;
88  ++m_userData.mumDeposits;
89  return true;
90  }
91 
92  }
93 }
94 
95 //--- Factory declaration
96 namespace dd4hep { namespace sim {
97  typedef Geant4SensitiveAction<MyTrackerSD> MyTrackerSDAction;
98  }}
99 DECLARE_GEANT4SENSITIVE(MyTrackerSDAction)
double energyDeposit
Energy deposit in the tracker hit.
Definition: MyTrackerHit.h:59
This is the hit definition.
Definition: MyTrackerHit.h:39
double length
Length of the track segment contributing to this hit.
Definition: MyTrackerHit.h:55
long long int cellID
dd4hep::sim::Geant4HitData: cellID
Definition: MyTrackerHit.h:43
void print(OStream &os, const parsing_result &result)
prints parsing result
Definition: clipp.h:6161
dd4hep::Direction momentum
Hit direction.
Definition: MyTrackerHit.h:53
Framework include files.
Definition: MyTrackerHit.h:19
Sensitive detector meant for tracking detectors, will produce one hit per step.
Namespace for the AIDA detector description toolkit.
Geant4SensitiveAction< MyTrackerSD > MyTrackerSDAction
dd4hep::Position position
Hit position.
Definition: MyTrackerHit.h:51