ATLAS Offline Software
Loading...
Searching...
No Matches
G4UA::SteppingValidation Class Reference

User action to do some basic step-based validation of G4. More...

#include <SteppingValidation.h>

Inheritance diagram for G4UA::SteppingValidation:
Collaboration diagram for G4UA::SteppingValidation:

Public Member Functions

 SteppingValidation ()
user action interface
virtual void BeginOfRunAction (const G4Run *) override
virtual void EndOfEventAction (const G4Event *) override
virtual void BeginOfEventAction (const G4Event *) override
virtual void UserSteppingAction (const G4Step *) override

Protected Attributes

std::string m_path {"/truth/"}
ServiceHandle< ITHistSvc > m_histSvc {"THistSvc", "SimTestHisto"}

Private Attributes

TH1 * m_stepL
TH1 * m_stepProc
TH1 * m_mscAngle
TH1 * m_stepELoss
TH1 * m_secE
TH1 * m_latPhi
TH1 * m_latEta
TH2 * m_EvsR
G4Track * m_prim
G4Track * m_sec
double m_primH
double m_primF
double m_dh
double m_dh2
double m_dp
double m_dp2
double m_nsec

Detailed Description

User action to do some basic step-based validation of G4.

Definition at line 24 of file SteppingValidation.h.

Constructor & Destructor Documentation

◆ SteppingValidation()

G4UA::SteppingValidation::SteppingValidation ( )

Definition at line 23 of file SteppingValidation.cxx.

24 : m_stepL(nullptr), m_stepProc(nullptr), m_mscAngle(nullptr),
25 m_stepELoss(nullptr), m_secE(nullptr),
26 m_latPhi(nullptr), m_latEta(nullptr), m_EvsR(nullptr),
27 m_prim(nullptr), m_sec(nullptr),
28 m_primH(0), m_primF(0), m_dh(0), m_dh2(0), m_dp(0), m_dp2(0), m_nsec(0)
29 {}

Member Function Documentation

◆ BeginOfEventAction()

void G4UA::SteppingValidation::BeginOfEventAction ( const G4Event * )
overridevirtual

Definition at line 67 of file SteppingValidation.cxx.

68 {
69 m_prim = m_sec = 0;
70 m_primH = m_primF = 0;
71 m_dh = m_dh2 = m_dp = m_dp2 = 0;
72 m_nsec = 0;
73 }

◆ BeginOfRunAction()

void G4UA::SteppingValidation::BeginOfRunAction ( const G4Run * )
overridevirtual

Definition at line 32 of file SteppingValidation.cxx.

33 {
34 m_path += "Stepping/";
35
36 // Set up all the histograms...
37 _TH1D_NOCHECK( m_stepL , "stepL" , 100 , -9 , 3 );
38 _TH1D_NOCHECK( m_stepProc , "stepProc" , 15 , -0.5 , 14.5 );
39 _TH1D_NOCHECK( m_mscAngle , "mscAngle" , 100 , 0. , M_PI );
40 _TH1D_NOCHECK( m_stepELoss , "stepELoss" , 100 , -10 , 2 );
41 _TH1D_NOCHECK( m_secE , "secE" , 100 , -3.5 , 5 );
42 _TH2D_NOCHECK( m_EvsR , "EvsR" , 100 , 0 , 1.2 , 100 , 1000. , 20000.);
43 _TH1D_NOCHECK( m_latPhi , "latPhi" , 100 , 0. , 1.2 );
44 _TH1D_NOCHECK( m_latEta , "latEta" , 100 , 0. , 0.6 );
45
46 _SET_TITLE(m_stepL, "Step length distribution", "Step length [log(mm)]", "Steps");
47 _SET_TITLE(m_stepProc, "Step process distribution", "Step process", "Steps");
48 _SET_TITLE(m_mscAngle, "Multiple scattering angle distribution", "Angle", "Steps");
49 _SET_TITLE(m_stepELoss, "Energy loss distribution", "Energy Loss [log(MeV)]", "Steps");
50 _SET_TITLE(m_secE, "Secondary energy distribution", "Energy [log(MeV)]", "Secondaries");
51 _SET_TITLE(m_EvsR, "Primary Energy vs Radius", "Radius [m]", "Primary Energy [MeV]");
52 _SET_TITLE(m_latPhi, "Phi energy distribution", "Energy-weighted #phi RMS", "Primaries");
53 _SET_TITLE(m_latEta, "Eta energy distribution", "Energy-weighted #eta RMS", "Primaries");
54 }
#define M_PI
#define _SET_TITLE(var, title, xaxis, yaxis)
#define _TH1D_NOCHECK(var, name, nbin, xmin, xmax)
#define _TH2D_NOCHECK(var, name, nbinx, xmin, xmax, nbiny, ymin, ymax)
std::string m_path

◆ EndOfEventAction()

void G4UA::SteppingValidation::EndOfEventAction ( const G4Event * )
overridevirtual

Definition at line 57 of file SteppingValidation.cxx.

58 {
59 // Fill lateral energy spread
60 if (m_nsec>0){
61 m_latPhi->Fill( std::sqrt( m_dp2/m_nsec - std::pow(m_dp/m_nsec,2) ) );
62 m_latEta->Fill( std::sqrt( m_dh2/m_nsec - std::pow(m_dh/m_nsec,2) ) );
63 }
64 }

◆ UserSteppingAction()

void G4UA::SteppingValidation::UserSteppingAction ( const G4Step * aStep)
overridevirtual

Definition at line 76 of file SteppingValidation.cxx.

77 {
78 // Fill process type
79 m_stepProc->Fill( aStep->GetPostStepPoint()->GetProcessDefinedStep()->
80 GetProcessType() );
81
82 // Fill energy deposit
83 m_stepELoss->Fill( log10(aStep->GetTotalEnergyDeposit()) );
84
85 // Fill step length
86 m_stepL->Fill( log10( aStep->GetStepLength() ) );
87
88 // Fill secondary energies
89 auto* preStepPoint = aStep->GetPreStepPoint();
90 auto* postStepPoint = aStep->GetPostStepPoint();
91 if (m_prim == 0){
92 m_prim = aStep->GetTrack();
93 }
94 // New secondary
95 else if (aStep->GetTrack()!=m_prim && aStep->GetTrack()!=m_sec){
96 m_secE->Fill( log10(preStepPoint->GetTotalEnergy()) );
97 m_sec = aStep->GetTrack();
98 }
99
100 // Defined by transport --> MSC
101 if (postStepPoint->GetProcessDefinedStep()->GetProcessType()==1){
102 if (m_prim==aStep->GetTrack()){
103 double angle = std::acos(
104 preStepPoint->GetMomentumDirection().dot(
105 postStepPoint->GetMomentumDirection() ) );
106 m_mscAngle->Fill(angle);
107 m_EvsR->Fill( postStepPoint->GetPosition().perp()*0.001,
108 aStep->GetTrack()->GetTotalEnergy() );
109 }
110 }
111
112 if (aStep->GetTrack()->GetNextVolume()==0){
113 // Leaving world! Use this for RMS calculations
114 if (m_prim==aStep->GetTrack()){
115 m_primH = postStepPoint->GetPosition().eta();
116 m_primF = postStepPoint->GetPosition().phi();
117 } else {
118 m_dh += (postStepPoint->GetPosition().eta() - m_primH )*aStep->GetTrack()->GetTotalEnergy();
119 m_dh2 += std::pow(postStepPoint->GetPosition().eta() - m_primH,2)*aStep->GetTrack()->GetTotalEnergy();
120 m_dp += (postStepPoint->GetPosition().phi() - m_primF )*aStep->GetTrack()->GetTotalEnergy();
121 m_dp2 += std::pow(postStepPoint->GetPosition().phi() - m_primF,2)*aStep->GetTrack()->GetTotalEnergy();
122 m_nsec += aStep->GetTrack()->GetTotalEnergy();
123 }
124 }
125
126 }
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
double angle(const GeoTrf::Vector2D &a, const GeoTrf::Vector2D &b)

Member Data Documentation

◆ m_dh

double G4UA::SteppingValidation::m_dh
private

Definition at line 48 of file SteppingValidation.h.

◆ m_dh2

double G4UA::SteppingValidation::m_dh2
private

Definition at line 48 of file SteppingValidation.h.

◆ m_dp

double G4UA::SteppingValidation::m_dp
private

Definition at line 48 of file SteppingValidation.h.

◆ m_dp2

double G4UA::SteppingValidation::m_dp2
private

Definition at line 48 of file SteppingValidation.h.

◆ m_EvsR

TH2* G4UA::SteppingValidation::m_EvsR
private

Definition at line 46 of file SteppingValidation.h.

◆ m_histSvc

ServiceHandle<ITHistSvc> SimTestHisto::m_histSvc {"THistSvc", "SimTestHisto"}
protectedinherited

Definition at line 35 of file SimTestHisto.h.

35{"THistSvc", "SimTestHisto"};

◆ m_latEta

TH1 * G4UA::SteppingValidation::m_latEta
private

Definition at line 45 of file SteppingValidation.h.

◆ m_latPhi

TH1 * G4UA::SteppingValidation::m_latPhi
private

Definition at line 45 of file SteppingValidation.h.

◆ m_mscAngle

TH1 * G4UA::SteppingValidation::m_mscAngle
private

Definition at line 45 of file SteppingValidation.h.

◆ m_nsec

double G4UA::SteppingValidation::m_nsec
private

Definition at line 48 of file SteppingValidation.h.

◆ m_path

std::string SimTestHisto::m_path {"/truth/"}
protectedinherited

Definition at line 34 of file SimTestHisto.h.

34{"/truth/"};

◆ m_prim

G4Track* G4UA::SteppingValidation::m_prim
private

Definition at line 47 of file SteppingValidation.h.

◆ m_primF

double G4UA::SteppingValidation::m_primF
private

Definition at line 48 of file SteppingValidation.h.

◆ m_primH

double G4UA::SteppingValidation::m_primH
private

Definition at line 48 of file SteppingValidation.h.

◆ m_sec

G4Track * G4UA::SteppingValidation::m_sec
private

Definition at line 47 of file SteppingValidation.h.

◆ m_secE

TH1 * G4UA::SteppingValidation::m_secE
private

Definition at line 45 of file SteppingValidation.h.

◆ m_stepELoss

TH1 * G4UA::SteppingValidation::m_stepELoss
private

Definition at line 45 of file SteppingValidation.h.

◆ m_stepL

TH1* G4UA::SteppingValidation::m_stepL
private

Definition at line 45 of file SteppingValidation.h.

◆ m_stepProc

TH1 * G4UA::SteppingValidation::m_stepProc
private

Definition at line 45 of file SteppingValidation.h.


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