LDMX Software
Public Member Functions | Private Member Functions | Private Attributes | List of all members
biasing::TargetDarkBremFilter Class Reference

This class is meant to filter for events that produce a dark brem matching originating in the target and matching the following parameters. More...

#include <TargetDarkBremFilter.h>

Public Member Functions

 TargetDarkBremFilter (const std::string &name, framework::config::Parameters &parameters)
 Class constructor.
 
 ~TargetDarkBremFilter ()=default
 Class destructor.
 
std::vector< simcore::TYPE > getTypes () override
 Get the types of actions this class can do.
 
void BeginOfEventAction (const G4Event *e) override
 Reset flag signaling finding of A' to false.
 
void stepping (const G4Step *step) override
 Looking for A' while primary is stepping.
 
void EndOfEventAction (const G4Event *event) override
 Check flag signaling finding of A', if false, abort the event so it isn't saved.
 
- Public Member Functions inherited from simcore::UserAction
 UserAction (const std::string &name, framework::config::Parameters &parameters)
 Constructor.
 
virtual ~UserAction ()=default
 Destructor.
 
virtual void BeginOfRunAction (const G4Run *)
 Method called at the beginning of a run.
 
virtual void EndOfRunAction (const G4Run *)
 Method called at the end of a run.
 
virtual void PreUserTrackingAction (const G4Track *)
 Method called before the UserTrackingAction.
 
virtual void PostUserTrackingAction (const G4Track *)
 Method called after the UserTrackingAction.
 
virtual G4ClassificationOfNewTrack ClassifyNewTrack (const G4Track *, const G4ClassificationOfNewTrack &cl)
 Method called when a track is updated.
 
virtual void NewStage ()
 Method called at the beginning of a new stage.
 
virtual void PrepareNewEvent ()
 Method called at the beginning of a new event.
 

Private Member Functions

bool isOutsideTargetRegion (const G4VPhysicalVolume *vol) const
 Check if the volume is outside the target region.
 
bool isOutsideTargetRegion (const G4LogicalVolume *vol) const
 Check if the volume is outside the target region.
 
void AbortEvent (const std::string &reason) const
 Helper to abort an event with a message.
 

Private Attributes

double threshold_
 Minimum energy [MeV] that the A' should have to keep the event.
 
bool found_aprime_
 flag to signal that we saw an A'
 

Additional Inherited Members

- Public Types inherited from simcore::UserAction
using Factory = ::simcore::Factory< UserAction, std::shared_ptr< UserAction >, const std::string &, framework::config::Parameters & >
 factory for user actions
 
- Protected Member Functions inherited from simcore::UserAction
UserEventInformationgetEventInfo () const
 Get a handle to the event information.
 
const std::map< int, ldmx::SimParticle > & getCurrentParticleMap () const
 Get the current particle map.
 
- Protected Attributes inherited from simcore::UserAction
std::string name_ {""}
 Name of the UserAction.
 
framework::config::Parameters parameters_
 The set of parameters used to configure this class.
 

Detailed Description

This class is meant to filter for events that produce a dark brem matching originating in the target and matching the following parameters.

 threshold: minimum energy [MeV] A' needs to have
See also
TargetBremFilter This filter is designed similar to the target brem filter where we check the secondaries of the primary electron if it is stopping within the target or if it is leaving the target region.

Definition at line 44 of file TargetDarkBremFilter.h.

Constructor & Destructor Documentation

◆ TargetDarkBremFilter()

TargetDarkBremFilter::TargetDarkBremFilter ( const std::string &  name,
framework::config::Parameters parameters 
)

Class constructor.

Retrieve the necessary configuration parameters

Definition at line 18 of file TargetDarkBremFilter.cxx.

20 : simcore::UserAction(name, parameters) {
21 threshold_ = parameters.getParameter<double>("threshold");
22}
double threshold_
Minimum energy [MeV] that the A' should have to keep the event.
T getParameter(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:89
Interface that defines a user action.
Definition UserAction.h:42

References framework::config::Parameters::getParameter(), and threshold_.

Member Function Documentation

◆ AbortEvent()

void TargetDarkBremFilter::AbortEvent ( const std::string &  reason) const
private

Helper to abort an event with a message.

Tells the RunManger to abort the current event after displaying the input message.

Parameters
[in]reasonreason for aborting the event

Definition at line 128 of file TargetDarkBremFilter.cxx.

128 {
129 if (G4RunManager::GetRunManager()->GetVerboseLevel() > 1) {
130 std::cout << "[ TargetDarkBremFilter ]: "
131 << "("
132 << G4EventManager::GetEventManager()
133 ->GetConstCurrentEvent()
134 ->GetEventID()
135 << ") " << reason << " Aborting event." << std::endl;
136 }
137 G4RunManager::GetRunManager()->AbortEvent();
138 return;
139}

Referenced by EndOfEventAction(), and stepping().

◆ BeginOfEventAction()

void TargetDarkBremFilter::BeginOfEventAction ( const G4Event *  e)
overridevirtual

Reset flag signaling finding of A' to false.

Parameters
[in]eevent being started, unused

Reimplemented from simcore::UserAction.

Definition at line 24 of file TargetDarkBremFilter.cxx.

24 {
25 found_aprime_ = false;
26}
bool found_aprime_
flag to signal that we saw an A'

References found_aprime_.

◆ EndOfEventAction()

void TargetDarkBremFilter::EndOfEventAction ( const G4Event *  event)
overridevirtual

Check flag signaling finding of A', if false, abort the event so it isn't saved.

Parameters
[in]eventbeing ended, check if it is already aborted

Reimplemented from simcore::UserAction.

Definition at line 122 of file TargetDarkBremFilter.cxx.

122 {
123 if (not event->IsAborted() and not found_aprime_) {
124 AbortEvent("Did not find an A' after entire event was simulated.");
125 }
126}
void AbortEvent(const std::string &reason) const
Helper to abort an event with a message.

References AbortEvent(), and found_aprime_.

◆ getTypes()

std::vector< simcore::TYPE > biasing::TargetDarkBremFilter::getTypes ( )
inlineoverridevirtual

Get the types of actions this class can do.

Returns
list of action types this class does

Implements simcore::UserAction.

Definition at line 64 of file TargetDarkBremFilter.h.

64 {
65 return {simcore::TYPE::STEPPING, simcore::TYPE::EVENT};
66 }

◆ isOutsideTargetRegion() [1/2]

bool biasing::TargetDarkBremFilter::isOutsideTargetRegion ( const G4LogicalVolume *  vol) const
inlineprivate

Check if the volume is outside the target region.

Note
will return true if vol is nullptr.
Parameters
[in]volG4LogicalVolume to check region
Returns
true if vol is outside target region or nullptr or doesn't have a region

Definition at line 118 of file TargetDarkBremFilter.h.

118 {
119 if (!vol) return true;
120 auto region = vol->GetRegion();
121 return region ? (region->GetName().compareTo("target") != 0) : true;
122 }

◆ isOutsideTargetRegion() [2/2]

bool biasing::TargetDarkBremFilter::isOutsideTargetRegion ( const G4VPhysicalVolume *  vol) const
inlineprivate

Check if the volume is outside the target region.

Note
will return true if vol is nullptr
Parameters
[in]volG4VPhysicalVolume to check region
Returns
true if vol is outside target region or nullptr

Definition at line 105 of file TargetDarkBremFilter.h.

105 {
106 return vol ? isOutsideTargetRegion(vol->GetLogicalVolume()) : true;
107 }
bool isOutsideTargetRegion(const G4VPhysicalVolume *vol) const
Check if the volume is outside the target region.

References isOutsideTargetRegion().

Referenced by isOutsideTargetRegion(), and stepping().

◆ stepping()

void TargetDarkBremFilter::stepping ( const G4Step *  step)
overridevirtual

Looking for A' while primary is stepping.

We make sure that the current track is the primary electron that is within the target region. Then if the track is either stopped or leaving the target region, we look through its secondaries for a good A'.

Parameters
[in]stepcurrent G4Step

check on first step of A' to see if it originated in correct volume this needs to be here because sometimes the electron dark brems in the tagger and misses the target completely std::cout << "[ TargetDarkBremFilter ] : " << "(" << G4EventManager::GetEventManager() ->GetConstCurrentEvent() ->GetEventID() << ") " << "A' originated in " << track->GetLogicalVolumeAtVertex()->GetName() << std::endl;

Reimplemented from simcore::UserAction.

Definition at line 28 of file TargetDarkBremFilter.cxx.

28 {
29 // don't process if event is already aborted
30 if (G4EventManager::GetEventManager()->GetConstCurrentEvent()->IsAborted())
31 return;
32
33 // Get the track associated with this step.
34 auto track{step->GetTrack()};
35
36 // Leave if track is not an electron
37 auto particle_def{track->GetParticleDefinition()};
38 if (particle_def != G4Electron::Electron()) {
39 if (particle_def == G4APrime::APrime() and
40 track->GetCurrentStepNumber() == 1) {
53 if (isOutsideTargetRegion(track->GetLogicalVolumeAtVertex()))
54 AbortEvent("A' was not created within target.");
55 // A' was found and originated in correct region
56 found_aprime_ = true;
57 } // first step of A'
58 return;
59 }
60
61 // Leave if track is not primary
62 if (track->GetParentID() != 0) return;
63
64 // Leave if track is not in target region
65 if (isOutsideTargetRegion(track->GetVolume())) return;
66
67 if (isOutsideTargetRegion(track->GetNextVolume()) // leaving target region
68 or
69 track->GetTrackStatus() == fStopAndKill // stopping within target region
70 or track->GetKineticEnergy() == 0. // stopping within target region
71 ) {
72 // Get the electron secondries
73 const std::vector<G4Track*>* secondaries{step->GetSecondary()};
74 if (not secondaries or secondaries->size() == 0) {
75 AbortEvent("Primary electron did not create secondaries.");
76 return;
77 } else {
78 // check secondaries to see if we made a dark brem
79 for (auto& secondary_track : *secondaries) {
80 if (secondary_track->GetParticleDefinition() == G4APrime::APrime()) {
81 // we found an A', woo-hoo!
82
83 if (secondary_track->GetTotalEnergy() < threshold_) {
85 "A' was not created with total energy above input threshold.");
86 return;
87 }
88
89 if (isOutsideTargetRegion(secondary_track->GetVolume())) {
90 AbortEvent("A' was not created within target.");
91 return;
92 }
93
94 /* debug message
95 std::cout << "[ TargetDarkBremFilter ] : "
96 << "(" <<
97 G4EventManager::GetEventManager()->GetConstCurrentEvent()->GetEventID()
98 << ") "
99 << "Found an A' weighted " << secondary_track->GetWeight()
100 << " originating in " <<
101 secondary_track->GetVolume()->GetLogicalVolume()->GetName()
102 << " with energy " << secondary_track->GetTotalEnergy() << " MeV."
103 << std::endl;
104 */
105
106 // we found a good A', so we can leave
107 return;
108 } // this secondary is an A'
109 } // loop through secondaries
110 } // are there secondaries to loop through
111
112 // got here without finding A'
113 AbortEvent("Primary electron did not create A'.");
114 return;
115 } // should check secondaries of primary
116
117 // get here if we are the primary in the target region,
118 // but shouldn't check the secondaries yet
119 return;
120}

References AbortEvent(), found_aprime_, isOutsideTargetRegion(), and threshold_.

Member Data Documentation

◆ found_aprime_

bool biasing::TargetDarkBremFilter::found_aprime_
private

flag to signal that we saw an A'

Definition at line 148 of file TargetDarkBremFilter.h.

Referenced by BeginOfEventAction(), EndOfEventAction(), and stepping().

◆ threshold_

double biasing::TargetDarkBremFilter::threshold_
private

Minimum energy [MeV] that the A' should have to keep the event.

Also used by PartialEnergySorter to determine which tracks should be processed first.

Parameter Name: 'threshold'

Definition at line 143 of file TargetDarkBremFilter.h.

Referenced by stepping(), and TargetDarkBremFilter().


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