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