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

#include <SG_StepNtuple.h>

Inheritance diagram for G4UA::SG_StepNtuple:
Collaboration diagram for G4UA::SG_StepNtuple:

Public Member Functions

 SG_StepNtuple (const std::vector< int > &)
virtual void BeginOfRunAction (const G4Run *) override
virtual void BeginOfEventAction (const G4Event *) override
virtual void EndOfEventAction (const G4Event *) override
virtual void UserSteppingAction (const G4Step *) override
bool msgLvl (const MSG::Level lvl) const
 Test the output level.
MsgStream & msg () const
 The standard message stream.
MsgStream & msg (const MSG::Level lvl) const
 The standard message stream.
void setLevel (MSG::Level lvl)
 Change the current logging level.

Private Member Functions

void initMessaging () const
 Initialize our message level and MessageSvc.

Private Attributes

NTuple::Item< long > m_nsteps
NTuple::Item< long > m_evtid
NTuple::Array< int > m_pdg
NTuple::Array< int > m_charge
NTuple::Array< int > m_baryon
NTuple::Array< float > m_x1
NTuple::Array< float > m_y1
NTuple::Array< float > m_z1
NTuple::Array< float > m_t1
NTuple::Array< float > m_pt1
NTuple::Array< float > m_x2
NTuple::Array< float > m_y2
NTuple::Array< float > m_z2
NTuple::Array< float > m_t2
NTuple::Array< float > m_pt2
NTuple::Array< float > m_minA
NTuple::Array< float > m_v2
NTuple::Array< float > m_vthresh
NTuple::Array< float > m_vbelowthresh
NTuple::Array< float > m_dep
NTuple::Array< float > m_mass
NTuple::Array< float > m_ke1
NTuple::Array< float > m_ke2
NTuple::Array< int > m_rh
NTuple::Array< int > m_rhid
NTuple::Array< int > m_step
std::set< int > m_rhs
long m_nevents = 0
long m_rhadronIndex = 0
std::string m_nm
 Message source name.
boost::thread_specific_ptr< MsgStream > m_msg_tls
 MsgStream instance (a std::cout like with print-out levels)
std::atomic< IMessageSvc * > m_imsg { nullptr }
 MessageSvc pointer.
std::atomic< MSG::Level > m_lvl { MSG::NIL }
 Current logging level.
std::atomic_flag m_initialized ATLAS_THREAD_SAFE = ATOMIC_FLAG_INIT
 Messaging initialized (initMessaging)

Detailed Description

Definition at line 20 of file SG_StepNtuple.h.

Constructor & Destructor Documentation

◆ SG_StepNtuple()

G4UA::SG_StepNtuple::SG_StepNtuple ( const std::vector< int > & pdgids)

Definition at line 28 of file SG_StepNtuple.cxx.

29 : AthMessaging(Gaudi::svcLocator()->service< IMessageSvc >( "MessageSvc" ),
30 "SG_StepNtuple")
31 {
32 // Load up the PDG IDs
33 for (int i : pdgids ){
34 m_rhs.insert(i);
35 }
36 }
AthMessaging()
Default constructor:
std::set< int > m_rhs

Member Function Documentation

◆ BeginOfEventAction()

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

Definition at line 93 of file SG_StepNtuple.cxx.

94 {
95 m_nsteps=0;
96 m_rhadronIndex=0;//the rhadron index (either the first or second rhadon, usually)
97 m_nevents++;
98 m_evtid=m_nevents;//since it gets cleared out after every fill...
99 }
NTuple::Item< long > m_nsteps
NTuple::Item< long > m_evtid

◆ BeginOfRunAction()

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

Definition at line 38 of file SG_StepNtuple.cxx.

39 {
40 NTupleFilePtr file1(ntupleSvc(), "/NTUPLES/FILE1");
41
42 SmartDataPtr<NTuple::Directory>
43 ntdir(ntupleSvc(),"/NTUPLES/FILE1/StepNtuple");
44
45 if ( !ntdir ) ntdir = ntupleSvc()->createDirectory(file1,"StepNtuple");
46 //if ( !ntdir ) log << MSG::ERROR << " failed to get ntuple directory" << endmsg;
47 NTuplePtr nt(ntupleSvc(), "/NTUPLES/FILE1/StepNtuple/10");
48 if ( !nt ) { // Check if already booked
49 nt = ntupleSvc()->book (ntdir.ptr(), 10,CLID_ColumnWiseTuple, "GEANT4 Step NTuple");
50 if ( nt ) {
51
52 ATH_MSG_INFO("booked step ntuple ");
53
54 if( nt->addItem ("nstep", m_nsteps,0 ,50000)!=StatusCode::SUCCESS // WARNING!! Force limit to 50k tracks
55 || nt->addItem ("pdg", m_nsteps, m_pdg)!=StatusCode::SUCCESS
56 || nt->addItem ("charge", m_nsteps, m_charge)!=StatusCode::SUCCESS
57 || nt->addItem ("mass", m_nsteps, m_mass)!=StatusCode::SUCCESS
58 || nt->addItem ("baryon", m_nsteps, m_baryon)!=StatusCode::SUCCESS
59 || nt->addItem ("x1", m_nsteps, m_x1)!=StatusCode::SUCCESS
60 || nt->addItem ("y1", m_nsteps, m_y1)!=StatusCode::SUCCESS
61 || nt->addItem ("z1", m_nsteps, m_z1)!=StatusCode::SUCCESS
62 || nt->addItem ("t1", m_nsteps, m_t1)!=StatusCode::SUCCESS
63 || nt->addItem ("x2", m_nsteps, m_x2)!=StatusCode::SUCCESS
64 || nt->addItem ("y2", m_nsteps, m_y2)!=StatusCode::SUCCESS
65 || nt->addItem ("z2", m_nsteps, m_z2)!=StatusCode::SUCCESS
66 || nt->addItem ("t2", m_nsteps, m_t2)!=StatusCode::SUCCESS
67 || nt->addItem ("dep", m_nsteps, m_dep)!=StatusCode::SUCCESS
68 || nt->addItem ("ke1", m_nsteps, m_ke1)!=StatusCode::SUCCESS
69 || nt->addItem ("ke2", m_nsteps, m_ke2)!=StatusCode::SUCCESS
70 || nt->addItem ("rh", m_nsteps, m_rh)!=StatusCode::SUCCESS
71 || nt->addItem ("rhid", m_nsteps, m_rhid)!=StatusCode::SUCCESS
72 || nt->addItem ("step", m_nsteps, m_step)!=StatusCode::SUCCESS
73 || nt->addItem ("pt1", m_nsteps, m_pt1)!=StatusCode::SUCCESS
74 || nt->addItem ("pt2", m_nsteps, m_pt2)!=StatusCode::SUCCESS
75 || nt->addItem ("minA",m_nsteps, m_minA)!=StatusCode::SUCCESS
76 || nt->addItem ("v2",m_nsteps, m_v2)!=StatusCode::SUCCESS
77 || nt->addItem ("vthresh",m_nsteps, m_vthresh)!=StatusCode::SUCCESS
78 || nt->addItem ("vbelowthresh",m_nsteps, m_vbelowthresh)!=StatusCode::SUCCESS
79 || nt->addItem ("evtid", m_evtid)!=StatusCode::SUCCESS)
80
81 ATH_MSG_ERROR("Could not configure branches ");
82
83 }
84
85 else ATH_MSG_ERROR("Could not book step ntuple!! ");
86 }
87
88 //set initial values
89 m_nevents=0;
90
91 }
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
INTupleSvc * ntupleSvc()
NTuple::Array< float > m_x1
NTuple::Array< float > m_y2
NTuple::Array< float > m_pt1
NTuple::Array< float > m_t1
NTuple::Array< float > m_v2
NTuple::Array< int > m_baryon
NTuple::Array< float > m_ke2
NTuple::Array< float > m_dep
NTuple::Array< float > m_z2
NTuple::Array< int > m_step
NTuple::Array< float > m_ke1
NTuple::Array< float > m_vthresh
NTuple::Array< float > m_x2
NTuple::Array< float > m_mass
NTuple::Array< float > m_minA
NTuple::Array< float > m_t2
NTuple::Array< float > m_vbelowthresh
NTuple::Array< int > m_rhid
NTuple::Array< float > m_z1
NTuple::Array< float > m_y1
NTuple::Array< int > m_rh
NTuple::Array< float > m_pt2
NTuple::Array< int > m_pdg
NTuple::Array< int > m_charge

◆ EndOfEventAction()

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

Definition at line 101 of file SG_StepNtuple.cxx.

102 {
103 if(! ntupleSvc()->writeRecord("/NTUPLES/FILE1/StepNtuple/10").isSuccess()) {
104 ATH_MSG_ERROR( " failed to write record for this event" );
105 }
106
107 //this also seems to zero out all the arrays... so beware!
108 }

◆ initMessaging()

void AthMessaging::initMessaging ( ) const
privateinherited

Initialize our message level and MessageSvc.

This method should only be called once.

Definition at line 39 of file AthMessaging.cxx.

40{
42 // If user did not set an explicit level, set a default
43 if (m_lvl == MSG::NIL) {
44 m_lvl = m_imsg ?
45 static_cast<MSG::Level>( m_imsg.load()->outputLevel(m_nm) ) :
46 MSG::INFO;
47 }
48}
std::string m_nm
Message source name.
std::atomic< IMessageSvc * > m_imsg
MessageSvc pointer.
std::atomic< MSG::Level > m_lvl
Current logging level.
IMessageSvc * getMessageSvc(bool quiet=false)

◆ msg() [1/2]

MsgStream & AthMessaging::msg ( ) const
inlineinherited

The standard message stream.

Returns a reference to the default message stream May not be invoked before sysInitialize() has been invoked.

Definition at line 163 of file AthMessaging.h.

164{
165 MsgStream* ms = m_msg_tls.get();
166 if (!ms) {
167 if (!m_initialized.test_and_set()) initMessaging();
168 ms = new MsgStream(m_imsg,m_nm);
169 m_msg_tls.reset( ms );
170 }
171
172 ms->setLevel (m_lvl);
173 return *ms;
174}
boost::thread_specific_ptr< MsgStream > m_msg_tls
MsgStream instance (a std::cout like with print-out levels)
void initMessaging() const
Initialize our message level and MessageSvc.

◆ msg() [2/2]

MsgStream & AthMessaging::msg ( const MSG::Level lvl) const
inlineinherited

The standard message stream.

Returns a reference to the default message stream May not be invoked before sysInitialize() has been invoked.

Definition at line 178 of file AthMessaging.h.

179{ return msg() << lvl; }
MsgStream & msg() const
The standard message stream.

◆ msgLvl()

bool AthMessaging::msgLvl ( const MSG::Level lvl) const
inlineinherited

Test the output level.

Parameters
lvlThe message level to test against
Returns
boolean Indicating if messages at given level will be printed
Return values
trueMessages at level "lvl" will be printed

Definition at line 151 of file AthMessaging.h.

152{
153 if (m_lvl <= lvl) {
154 msg() << lvl;
155 return true;
156 } else {
157 return false;
158 }
159}

◆ setLevel()

void AthMessaging::setLevel ( MSG::Level lvl)
inherited

Change the current logging level.

Use this rather than msg().setLevel() for proper operation with MT.

Definition at line 28 of file AthMessaging.cxx.

29{
30 m_lvl = lvl;
31}

◆ UserSteppingAction()

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

Definition at line 110 of file SG_StepNtuple.cxx.

111 {
112 if(m_nsteps<50000){
113 const int pdg_id = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
114 bool rhad=false;
115 if (std::find(m_rhs.begin(),m_rhs.end(),std::abs(pdg_id))!=m_rhs.end()) {
116 rhad=true;
117 }
118
119 //
120 if (!rhad && (MC::isSquarkLH(pdg_id) ||
121 pdg_id == 1000021 || // gluino
122 MC::isRHadron(pdg_id))) {
123 ATH_MSG_DEBUG (" TruthUtils classifies "<<pdg_id<<" as an R-Hadron, gluino or LH Squark!");
124 rhad=true;
125 }
126 //
127
128 if (rhad){
129
130 //
131 G4Material * mat = aStep->GetTrack()->GetMaterial();
132 double minA=1500000.;
133 for (unsigned int i=0;i<mat->GetNumberOfElements();++i){
134 if (mat->GetElement(i) && minA>mat->GetElement(i)->GetN()){
135 minA=mat->GetElement(i)->GetN();
136 }
137 }
138 //
139
140 bool firstslow = aStep->GetPostStepPoint()->GetVelocity()<0.15*std::pow(minA,-2./3.)*CLHEP::c_light;
141 //just save the first slow step for the rhadron
142 for (int i=0; i<m_nsteps; ++i){
143 if (m_rhid[i]==m_rhadronIndex && m_vbelowthresh[i]>0) firstslow=false;
144 }
145 if (firstslow || aStep->GetTrack()->GetCurrentStepNumber()<=1 || aStep->GetPostStepPoint()->GetKineticEnergy()==0.){
146
147 //
148 if (MC::isSquarkLH(pdg_id) ||
149 pdg_id == 1000021 || // gluino
150 MC::isRHadron(pdg_id)) {
151 m_rh[m_nsteps] = 1;// TruthUtils classifies this particle as an R-Hadron, gluino or LH squark
152 }
153 else{
154 m_rh[m_nsteps] = 0;// particle only passes the RHadronPDGIDList property check
155 }
156 //
157
158 if (aStep->GetPreStepPoint()->GetGlobalTime()==0) m_rhadronIndex++;
160
161 m_pdg[m_nsteps]=pdg_id;
162 m_charge[m_nsteps]=aStep->GetTrack()->GetDefinition()->GetPDGCharge();
163 m_dep[m_nsteps]=aStep->GetTotalEnergyDeposit();
164 m_mass[m_nsteps]=aStep->GetTrack()->GetDefinition()->GetPDGMass();
165 m_baryon[m_nsteps]=aStep->GetTrack()->GetDefinition()->GetBaryonNumber();
166 m_t1[m_nsteps]=aStep->GetPreStepPoint()->GetGlobalTime();
167 m_t2[m_nsteps]=aStep->GetPostStepPoint()->GetGlobalTime();
168 G4ThreeVector pos1=aStep->GetPreStepPoint()->GetPosition();
169 m_x1[m_nsteps]=pos1.x();
170 m_y1[m_nsteps]=pos1.y();
171 m_z1[m_nsteps]=pos1.z();
172 G4ThreeVector pos2=aStep->GetPostStepPoint()->GetPosition();
173 m_x2[m_nsteps]=pos2.x();
174 m_y2[m_nsteps]=pos2.y();
175 m_z2[m_nsteps]=pos2.z();
176 m_ke1[m_nsteps]=aStep->GetPreStepPoint()->GetKineticEnergy();
177 m_ke2[m_nsteps]=aStep->GetPostStepPoint()->GetKineticEnergy();
178 m_pt1[m_nsteps]=aStep->GetPreStepPoint()->GetMomentum().perp();
179 m_pt2[m_nsteps]=aStep->GetPostStepPoint()->GetMomentum().perp();
180
181 //
182 m_minA[m_nsteps]=minA;
183 m_v2[m_nsteps]=aStep->GetPostStepPoint()->GetVelocity();
184 m_vthresh[m_nsteps]=0.15*std::pow(minA,-2./3.)*CLHEP::c_light;
185 m_vbelowthresh[m_nsteps]=(firstslow?1:0);
186 //
187
189 ++m_nsteps;
190 //std::cout<<"stepping, size is "<<m_nsteps<<std::endl;
191 } // writing the step because it stopped or is the start of the "track"
192 } //rhad true
193 else {
194
195 //KILL the particles here, so we don't waste time in GEANT4 tracking what happens to it!
196 if (MC::isBSM(pdg_id)) { // flag SUSY/BSM particles which are skipped
197 ATH_MSG_DEBUG ("UserSteppingAction(): Killing uninteresting track with pdg_id "<<pdg_id);
198 }
199 aStep->GetTrack()->SetTrackStatus(fKillTrackAndSecondaries);
200 const G4TrackVector *tv = aStep->GetSecondary();
201 //if ((*tv).size()>0) std::cout<<" ... and its "<<(*tv).size()<<" secondaries"<<std::endl;
202 for (unsigned int i=0;i<tv->size();i++){
203 G4Track *t = (*tv)[i];
204 t->SetTrackStatus(fKillTrackAndSecondaries);
205 }
206
207 } // not an rhad
208
209 } //m_nsteps<50000
210 }
#define ATH_MSG_DEBUG(x)
bool isSquarkLH(const T &p)
bool isRHadron(const T &p)
bool isBSM(const T &p)
APID: graviton and all Higgs extensions are BSM.

Member Data Documentation

◆ ATLAS_THREAD_SAFE

std::atomic_flag m_initialized AthMessaging::ATLAS_THREAD_SAFE = ATOMIC_FLAG_INIT
mutableprivateinherited

Messaging initialized (initMessaging)

Definition at line 141 of file AthMessaging.h.

◆ m_baryon

NTuple::Array<int> G4UA::SG_StepNtuple::m_baryon
private

Definition at line 36 of file SG_StepNtuple.h.

◆ m_charge

NTuple::Array<int> G4UA::SG_StepNtuple::m_charge
private

Definition at line 36 of file SG_StepNtuple.h.

◆ m_dep

NTuple::Array<float> G4UA::SG_StepNtuple::m_dep
private

Definition at line 40 of file SG_StepNtuple.h.

◆ m_evtid

NTuple::Item<long> G4UA::SG_StepNtuple::m_evtid
private

Definition at line 35 of file SG_StepNtuple.h.

◆ m_imsg

std::atomic<IMessageSvc*> AthMessaging::m_imsg { nullptr }
mutableprivateinherited

MessageSvc pointer.

Definition at line 135 of file AthMessaging.h.

135{ nullptr };

◆ m_ke1

NTuple::Array<float> G4UA::SG_StepNtuple::m_ke1
private

Definition at line 41 of file SG_StepNtuple.h.

◆ m_ke2

NTuple::Array<float> G4UA::SG_StepNtuple::m_ke2
private

Definition at line 41 of file SG_StepNtuple.h.

◆ m_lvl

std::atomic<MSG::Level> AthMessaging::m_lvl { MSG::NIL }
mutableprivateinherited

Current logging level.

Definition at line 138 of file AthMessaging.h.

138{ MSG::NIL };

◆ m_mass

NTuple::Array<float> G4UA::SG_StepNtuple::m_mass
private

Definition at line 40 of file SG_StepNtuple.h.

◆ m_minA

NTuple::Array<float> G4UA::SG_StepNtuple::m_minA
private

Definition at line 39 of file SG_StepNtuple.h.

◆ m_msg_tls

boost::thread_specific_ptr<MsgStream> AthMessaging::m_msg_tls
mutableprivateinherited

MsgStream instance (a std::cout like with print-out levels)

Definition at line 132 of file AthMessaging.h.

◆ m_nevents

long G4UA::SG_StepNtuple::m_nevents = 0
private

Definition at line 44 of file SG_StepNtuple.h.

◆ m_nm

std::string AthMessaging::m_nm
privateinherited

Message source name.

Definition at line 129 of file AthMessaging.h.

◆ m_nsteps

NTuple::Item<long> G4UA::SG_StepNtuple::m_nsteps
private

Definition at line 35 of file SG_StepNtuple.h.

◆ m_pdg

NTuple::Array<int> G4UA::SG_StepNtuple::m_pdg
private

Definition at line 36 of file SG_StepNtuple.h.

◆ m_pt1

NTuple::Array<float> G4UA::SG_StepNtuple::m_pt1
private

Definition at line 37 of file SG_StepNtuple.h.

◆ m_pt2

NTuple::Array<float> G4UA::SG_StepNtuple::m_pt2
private

Definition at line 38 of file SG_StepNtuple.h.

◆ m_rh

NTuple::Array<int> G4UA::SG_StepNtuple::m_rh
private

Definition at line 42 of file SG_StepNtuple.h.

◆ m_rhadronIndex

long G4UA::SG_StepNtuple::m_rhadronIndex = 0
private

Definition at line 45 of file SG_StepNtuple.h.

◆ m_rhid

NTuple::Array<int> G4UA::SG_StepNtuple::m_rhid
private

Definition at line 42 of file SG_StepNtuple.h.

◆ m_rhs

std::set<int> G4UA::SG_StepNtuple::m_rhs
private

Definition at line 43 of file SG_StepNtuple.h.

◆ m_step

NTuple::Array<int> G4UA::SG_StepNtuple::m_step
private

Definition at line 42 of file SG_StepNtuple.h.

◆ m_t1

NTuple::Array<float> G4UA::SG_StepNtuple::m_t1
private

Definition at line 37 of file SG_StepNtuple.h.

◆ m_t2

NTuple::Array<float> G4UA::SG_StepNtuple::m_t2
private

Definition at line 38 of file SG_StepNtuple.h.

◆ m_v2

NTuple::Array<float> G4UA::SG_StepNtuple::m_v2
private

Definition at line 39 of file SG_StepNtuple.h.

◆ m_vbelowthresh

NTuple::Array<float> G4UA::SG_StepNtuple::m_vbelowthresh
private

Definition at line 39 of file SG_StepNtuple.h.

◆ m_vthresh

NTuple::Array<float> G4UA::SG_StepNtuple::m_vthresh
private

Definition at line 39 of file SG_StepNtuple.h.

◆ m_x1

NTuple::Array<float> G4UA::SG_StepNtuple::m_x1
private

Definition at line 37 of file SG_StepNtuple.h.

◆ m_x2

NTuple::Array<float> G4UA::SG_StepNtuple::m_x2
private

Definition at line 38 of file SG_StepNtuple.h.

◆ m_y1

NTuple::Array<float> G4UA::SG_StepNtuple::m_y1
private

Definition at line 37 of file SG_StepNtuple.h.

◆ m_y2

NTuple::Array<float> G4UA::SG_StepNtuple::m_y2
private

Definition at line 38 of file SG_StepNtuple.h.

◆ m_z1

NTuple::Array<float> G4UA::SG_StepNtuple::m_z1
private

Definition at line 37 of file SG_StepNtuple.h.

◆ m_z2

NTuple::Array<float> G4UA::SG_StepNtuple::m_z2
private

Definition at line 38 of file SG_StepNtuple.h.


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