LDMX Software
Public Member Functions | Private Attributes | List of all members
simcore::ScoringPlaneSD Class Reference

Class defining a basic sensitive detector for scoring planes. More...

#include <ScoringPlaneSD.h>

Public Member Functions

 ScoringPlaneSD (const std::string &name, simcore::ConditionsInterface &ci, const framework::config::Parameters &params)
 Constructor.
 
virtual ~ScoringPlaneSD ()=default
 Destructor.
 
bool isSensDet (G4LogicalVolume *volume) const override
 Check if the input logical volume is a scoring plane we should include.
 
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 should process its interaction with us.
 
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.
 
virtual void OnFinishedEvent () override
 Cleanup SD and prepare a new-event state.
 
- Public Member Functions inherited from simcore::SensitiveDetector
 SensitiveDetector (const std::string &name, simcore::ConditionsInterface &ci, const framework::config::Parameters &parameters)
 Constructor.
 
virtual ~SensitiveDetector ()=default
 Destructor.
 
virtual void EndOfEvent (G4HCofThisEvent *) override
 This is Geant4's handle to tell us the event is ending.
 

Private Attributes

std::string match_substr_
 Substring to match to logical volumes.
 
std::string collection_name_
 Name of output collection to add.
 
std::vector< ldmx::SimTrackerHithits_
 The actual output collection.
 

Additional Inherited Members

- Public Types inherited from simcore::SensitiveDetector
using Factory = ::simcore::Factory< SensitiveDetector, SensitiveDetector *, const std::string &, simcore::ConditionsInterface &, const framework::config::Parameters & >
 The SD Factory.
 
- Protected Member Functions inherited from simcore::SensitiveDetector
template<class T >
const T & getCondition (const std::string &condition_name)
 Record the configuration of this detector into the run header.
 
bool isGeantino (const G4Step *step) const
 Check if the passed step is a step of a geantino.
 
const TrackMapgetTrackMap () const
 Get a handle to the current track map.
 

Detailed Description

Class defining a basic sensitive detector for scoring planes.

Definition at line 13 of file ScoringPlaneSD.h.

Constructor & Destructor Documentation

◆ ScoringPlaneSD()

simcore::ScoringPlaneSD::ScoringPlaneSD ( const std::string &  name,
simcore::ConditionsInterface ci,
const framework::config::Parameters params 
)

Constructor.

Parameters
nameThe name of the sensitive detector.
ciConditions interface handle
paramspython configuration parameters

Definition at line 22 of file ScoringPlaneSD.cxx.

25 : SensitiveDetector(name, ci, params) {
26 collection_name_ = params.getParameter<std::string>("collection_name");
27 match_substr_ = params.getParameter<std::string>("match_substr");
28}
T getParameter(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:89
std::string collection_name_
Name of output collection to add.
std::string match_substr_
Substring to match to logical volumes.
SensitiveDetector(const std::string &name, simcore::ConditionsInterface &ci, const framework::config::Parameters &parameters)
Constructor.

References collection_name_, framework::config::Parameters::getParameter(), and match_substr_.

Member Function Documentation

◆ isSensDet()

bool simcore::ScoringPlaneSD::isSensDet ( G4LogicalVolume *  volume) const
inlineoverridevirtual

Check if the input logical volume is a scoring plane we should include.

Implements simcore::SensitiveDetector.

Definition at line 31 of file ScoringPlaneSD.h.

31 {
32 return volume->GetName().contains(match_substr_);
33 }

References match_substr_.

◆ OnFinishedEvent()

virtual void simcore::ScoringPlaneSD::OnFinishedEvent ( )
inlineoverridevirtual

Cleanup SD and prepare a new-event state.

Implements simcore::SensitiveDetector.

Definition at line 53 of file ScoringPlaneSD.h.

53{ hits_.clear(); }
std::vector< ldmx::SimTrackerHit > hits_
The actual output collection.

References hits_.

◆ ProcessHits()

G4bool simcore::ScoringPlaneSD::ProcessHits ( G4Step *  step,
G4TouchableHistory *  hist 
)
overridevirtual

This is Geant4's handle to tell us that a particle has stepped through our sensitive detector and we should process its interaction with us.

Parameters
[in]stepthe step that happened within one of our logical volumes
[in]histthe touchable history of the step

Implements simcore::SensitiveDetector.

Definition at line 30 of file ScoringPlaneSD.cxx.

30 {
31 // Get the edep from the step.
32 G4double edep = step->GetTotalEnergyDeposit();
33
34 // Create a new hit object.
35 ldmx::SimTrackerHit& hit{hits_.emplace_back()};
36
37 // Assign track ID for finding the SimParticle in post event processing.
38 hit.setTrackID(step->GetTrack()->GetTrackID());
39 hit.setPdgID(step->GetTrack()->GetDynamicParticle()->GetPDGcode());
40
41 // Set the edep.
42 hit.setEdep(edep);
43
44 // Set the start position.
45 G4StepPoint* prePoint = step->GetPreStepPoint();
46
47 // Set the end position.
48 G4StepPoint* postPoint = step->GetPostStepPoint();
49
50 G4ThreeVector start = prePoint->GetPosition();
51 G4ThreeVector end = postPoint->GetPosition();
52
53 // Set the mid position.
54 G4ThreeVector mid = 0.5 * (start + end);
55 hit.setPosition(mid.x(), mid.y(), mid.z());
56
57 // Compute path length.
58 G4double pathLength =
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);
62
63 // Set the global time.
64 hit.setTime(step->GetTrack()->GetGlobalTime());
65
66 // Set the momentum
67 G4ThreeVector p = postPoint->GetMomentum();
68 hit.setMomentum(p.x(), p.y(), p.z());
69 hit.setEnergy(postPoint->GetTotalEnergy());
70
71 /*
72 * Set the 32-bit ID on the hit.
73 */
74 int cpNumber = prePoint->GetTouchableHandle()->GetCopyNumber();
76 hit.setID(id.raw());
77
78 return true;
79}
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.

References hits_, ldmx::SimSpecialID::ScoringPlaneID(), and ldmx::SimTrackerHit::setTrackID().

◆ saveHits()

void simcore::ScoringPlaneSD::saveHits ( framework::Event event)
overridevirtual

We are given the event bus here and we must decide now what to persist into the event.

Parameters
[in,out]eventevent bus to add thing(s) to

Implements simcore::SensitiveDetector.

Definition at line 81 of file ScoringPlaneSD.cxx.

81 {
82 event.add(collection_name_, hits_);
83}

References collection_name_, and hits_.

Member Data Documentation

◆ collection_name_

std::string simcore::ScoringPlaneSD::collection_name_
private

Name of output collection to add.

Definition at line 60 of file ScoringPlaneSD.h.

Referenced by saveHits(), and ScoringPlaneSD().

◆ hits_

std::vector<ldmx::SimTrackerHit> simcore::ScoringPlaneSD::hits_
private

The actual output collection.

Definition at line 63 of file ScoringPlaneSD.h.

Referenced by OnFinishedEvent(), ProcessHits(), and saveHits().

◆ match_substr_

std::string simcore::ScoringPlaneSD::match_substr_
private

Substring to match to logical volumes.

Definition at line 57 of file ScoringPlaneSD.h.

Referenced by isSensDet(), and ScoringPlaneSD().


The documentation for this class was generated from the following files: