2#include "SimCore/SDs/ScoringPlaneSD.h"
4#include "DetDescr/SimSpecialID.h"
14#include "G4ChargedGeantino.hh"
15#include "G4Geantino.hh"
16#include "G4SDManager.hh"
18#include "G4StepPoint.hh"
32 G4double edep = step->GetTotalEnergyDeposit();
38 hit.
setTrackID(step->GetTrack()->GetTrackID());
39 hit.setPdgID(step->GetTrack()->GetDynamicParticle()->GetPDGcode());
45 G4StepPoint* prePoint = step->GetPreStepPoint();
48 G4StepPoint* postPoint = step->GetPostStepPoint();
50 G4ThreeVector start = prePoint->GetPosition();
51 G4ThreeVector end = postPoint->GetPosition();
54 G4ThreeVector mid = 0.5 * (start + end);
55 hit.setPosition(mid.x(), mid.y(), mid.z());
59 sqrt(pow(start.x() - end.x(), 2) + pow(start.y() - end.y(), 2) +
60 pow(start.z() - end.z(), 2));
61 hit.setPathLength(pathLength);
64 hit.setTime(step->GetTrack()->GetGlobalTime());
67 G4ThreeVector p = postPoint->GetMomentum();
68 hit.setMomentum(p.x(), p.y(), p.z());
69 hit.setEnergy(postPoint->GetTotalEnergy());
74 int cpNumber = prePoint->GetTouchableHandle()->GetCopyNumber();
Implements an event buffer system for storing event data.
Class encapsulating parameters for configuring a processor.
T getParameter(const std::string &name) const
Retrieve the parameter of the given name.
Implements detector ids for special simulation-derived hits like scoring planes.
static SimSpecialID ScoringPlaneID(int plane)
Create a scoring id from pieces.
Represents a simulated tracker hit in the simulation.
void setTrackID(const int simTrackID)
Set the Sim particle track ID of the hit.
Handle to the conditions system, provided at construction to classes which require it.
Class defining a basic sensitive detector for scoring planes.
std::string collection_name_
Name of output collection to add.
ScoringPlaneSD(const std::string &name, simcore::ConditionsInterface &ci, const framework::config::Parameters ¶ms)
Constructor.
std::vector< ldmx::SimTrackerHit > hits_
The actual output collection.
virtual G4bool ProcessHits(G4Step *step, G4TouchableHistory *hist) override
This is Geant4's handle to tell us that a particle has stepped through our sensitive detector and we ...
virtual void saveHits(framework::Event &event) override
We are given the event bus here and we must decide now what to persist into the event.
std::string match_substr_
Substring to match to logical volumes.
Dynamically loaded Geant4 SensitiveDetector for saving hits in specific volumes within the simulation...