Belle II Software  release-08-01-10
AWESOMESensitiveDetector.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 /* Own header. */
10 #include <online_book/awesome/simulation/AWESOMESensitiveDetector.h>
11 
12 /* Basf2 headers. */
13 #include <framework/gearbox/Unit.h>
14 #include <framework/geometry/B2Vector3.h>
15 
16 /* Geant4 headers. */
17 #include <G4StepPoint.hh>
18 #include <G4Track.hh>
19 
20 /* CLHEP headers. */
21 #include <CLHEP/Units/SystemOfUnits.h>
22 #include <CLHEP/Vector/ThreeVector.h>
23 
24 using namespace Belle2;
25 using namespace Belle2::AWESOME;
26 
28  Simulation::SensitiveDetectorBase{"AwesomeSensitiveDetector", Const::EDetector::TEST}
29 {
30  /* MCParticles must be optional, not required. */
32  /* Register the simulted hits and all the necessary relations in the datastore. */
33  m_SimHits.registerInDataStore();
35  /*
36  * Register the Relation so that the Geant4 TrackIDs get replaced by the actual
37  * MCParticle indices after simulating the events. This is needed as
38  * secondary particles might not be stored so everything relating to those
39  * particles will be attributed to the last saved mother particle.
40  */
42 }
43 
44 bool AWESOMESensitiveDetector::step(G4Step* step, G4TouchableHistory*)
45 {
46  /* Get some basic information from Geant4. */
47  const G4Track& track = *step->GetTrack();
48  const int trackID = track.GetTrackID();
49  const double energyDep = step->GetTotalEnergyDeposit() * Unit::MeV;
50  G4StepPoint* preStep = step->GetPreStepPoint();
51  G4StepPoint* postStep = step->GetPostStepPoint();
52  /* Ignore everything below 1 eV. */
53  if (energyDep < Unit::eV)
54  return false;
55  /*
56  * Compute other useful quantities to be stored.
57  * We must be sure that the quantities are stored using the proper Belle2 units:
58  * positions in cm, time in ns, etc.
59  */
60  const CLHEP::Hep3Vector position = 0.5 * (preStep->GetPosition() + postStep->GetPosition()) / CLHEP::cm; // Now in cm
61  const double time = 0.5 * (preStep->GetGlobalTime() + postStep->GetGlobalTime()); // Already in ns
62  /* Store the simulated hit. */
63  AWESOMESimHit* simHit = m_SimHits.appendNew();
64  simHit->setEnergyDep(energyDep);
65  simHit->setPosition(B2Vector3{position.x(), position.y(), position.z()});
66  simHit->setTime(time);
67  /*
68  * Add a relation between the current MCParticle and the simulated hit.
69  * Since the MCParticle index is not yet defined we use the trackID from Geant4.
70  */
71  m_MCParticlesToSimHits.add(trackID, simHit->getArrayIndex());
72  return true;
73 }
A Geant4 simulated hit for the AWESOME detector.
Definition: AWESOMESimHit.h:33
void setEnergyDep(float energyDep)
Set the deposited energy.
Definition: AWESOMESimHit.h:71
void setTime(float time)
Set the time.
void setPosition(B2Vector3< float > position)
Set the vector for position.
Definition: AWESOMESimHit.h:89
StoreArray< AWESOMESimHit > m_SimHits
AWESOME simulated hits.
RelationArray m_MCParticlesToSimHits
Relation array between MCParticles and AWESOMESimHits.
StoreArray< MCParticle > m_MCParticles
MC particles.
bool step(G4Step *step, G4TouchableHistory *) override
Step processing method.
A fast and root compatible alternative to TVector3.
Definition: B2Vector3.h:42
DataType x() const
access variable X (= .at(0) without boundary check)
Definition: B2Vector3.h:425
This class provides a set of constants for the framework.
Definition: Const.h:34
void add(index_type from, index_type to, weight_type weight=1.0)
Add a new element to the relation.
int getArrayIndex() const
Returns this object's array index (in StoreArray), or -1 if not found.
static void registerMCParticleRelation(const std::string &name, RelationArray::EConsolidationAction ignoreAction=RelationArray::c_negativeWeight)
Register an relation involving MCParticles.
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
bool registerRelationTo(const StoreArray< TO > &toArray, DataStore::EDurability durability=DataStore::c_Event, DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut, const std::string &namedRelation="") const
Register a relation to the given StoreArray.
Definition: StoreArray.h:140
static const double eV
[electronvolt]
Definition: Unit.h:112
static const double MeV
[megaelectronvolt]
Definition: Unit.h:114
TEST(TestgetDetectorRegion, TestgetDetectorRegion)
Test Constructors.
Abstract base class for different kinds of events.