Loading...
Searching...
No Matches
Public Member Functions | Private Member Functions | Private Attributes | List of all members
RE01StackingAction Class Reference

#include <Doxymodules_runAndEvent.h>

Inheritance diagram for RE01StackingAction:
G4UserStackingAction

Public Member Functions

 RE01StackingAction ()
 
virtual ~RE01StackingAction ()
 
virtual G4ClassificationOfNewTrack ClassifyNewTrack (const G4Track *aTrack)
 
virtual void NewStage ()
 
virtual void PrepareNewEvent ()
 

Private Member Functions

G4VHitsCollectionGetCalCollection ()
 

Private Attributes

G4int fStage
 
G4int fCalorimeterHitsColID
 

Detailed Description

Definition at line 31 of file Doxymodules_runAndEvent.h.

Constructor & Destructor Documentation

◆ RE01StackingAction()

RE01StackingAction::RE01StackingAction ( )

◆ ~RE01StackingAction()

RE01StackingAction::~RE01StackingAction ( )
virtual

Definition at line 53 of file RE01StackingAction.cc.

54{;}

Member Function Documentation

◆ ClassifyNewTrack()

G4ClassificationOfNewTrack RE01StackingAction::ClassifyNewTrack ( const G4Track aTrack)
virtual

Definition at line 57 of file RE01StackingAction.cc.

59{
60 G4ClassificationOfNewTrack classification = fUrgent;
61
62 if(fStage==0)
63 {
64 RE01TrackInformation* trackInfo;
65 if(aTrack->GetTrackStatus()==fSuspend) // Track reached to calorimeter
66 {
67 trackInfo = (RE01TrackInformation*)(aTrack->GetUserInformation());
68 trackInfo->SetTrackingStatus(0);
69 trackInfo->SetSourceTrackInformation(aTrack);
70// G4cout << "Track " << aTrack->GetTrackID()
71// << " (parentID " << aTrack->GetParentID()
72// << ") has reached to calorimeter and has been suspended at "
73// << aTrack->GetPosition() << G4endl;
74 classification = fWaiting;
75 }
76 else if(aTrack->GetParentID()==0) // Primary particle
77 {
78 trackInfo = new RE01TrackInformation(aTrack);
79 trackInfo->SetTrackingStatus(1);
80 G4Track* theTrack = (G4Track*)aTrack;
81 theTrack->SetUserInformation(trackInfo);
82 }
83 }
84 return classification;
85}
void SetSourceTrackInformation(const G4Track *aTrack)

◆ NewStage()

void RE01StackingAction::NewStage ( )
virtual

Definition at line 104 of file RE01StackingAction.cc.

105{
106 // G4cout << "+++++++++++ Stage " << fStage << G4endl;
107 if(fStage==0)
108 {
109 // display trajetory information in the tracking region
110 G4cout << G4endl;
111 G4cout << "Tracks in tracking region have been processed. -- Stage 0 over."
112 << G4endl;
113 G4cout << G4endl;
114 }
115 else
116 {
117 // display calorimeter information caused by a "source" track
118 // in the tracker region
119 // G4cout << G4endl;
120 // G4cout << "Processing one shower originated by a source track "
121 // << "in tracker region is over." << G4endl;
122 // G4cout << G4endl;
125 if(CHC)
126 {
127 int n_hit = CHC->entries();
128 G4double totE = 0;
129 G4int n_hitByATrack = 0;
130 for(int i=0;i<n_hit;i++)
131 {
132 G4double edepByATrack = (*CHC)[i]->GetEdepByATrack();
133 if(edepByATrack>0.)
134 {
135 totE += edepByATrack;
136 if(n_hitByATrack==0)
137 { (*CHC)[i]->GetTrackInformation()->Print(); }
138 n_hitByATrack++;
139 G4cout << "Cell[" << (*CHC)[i]->GetZ() << ","
140 << (*CHC)[i]->GetPhi() << "] "
141 << edepByATrack/GeV << " [GeV]" << G4endl;
142 (*CHC)[i]->ClearEdepByATrack();
143 }
144 }
145 if(n_hitByATrack>0)
146 {
147 G4cout <<
148 "### Total energy deposition in calorimeter by a source track in "
149 << n_hitByATrack << " cells : " << totE / GeV << " (GeV)"
150 << G4endl;
151 G4cout << G4endl;
152 }
153 }
154 }
155
156 if(stackManager->GetNUrgentTrack())
157 {
158 // Transfer all tracks in Urgent stack to Waiting stack, since all tracks
159 // in Waiting stack have already been transfered to Urgent stack before
160 // invokation of this method.
161 stackManager->TransferStackedTracks(fUrgent,fWaiting);
162
163 // Then, transfer only one track to Urgent stack.
164 stackManager->TransferOneStackedTrack(fWaiting,fUrgent);
165
166 fStage++;
167 }
168}
G4THitsCollection< RE01CalorimeterHit > RE01CalorimeterHitsCollection
G4VHitsCollection * GetCalCollection()

◆ PrepareNewEvent()

void RE01StackingAction::PrepareNewEvent ( )
virtual

Definition at line 171 of file RE01StackingAction.cc.

172{
173 fStage = 0;
174}

◆ GetCalCollection()

G4VHitsCollection * RE01StackingAction::GetCalCollection ( )
private

Definition at line 88 of file RE01StackingAction.cc.

89{
90 G4SDManager* SDMan = G4SDManager::GetSDMpointer();
91 G4RunManager* runMan = G4RunManager::GetRunManager();
93 { fCalorimeterHitsColID = SDMan->GetCollectionID("calCollection"); }
95 {
96 const G4Event* currentEvent = runMan->GetCurrentEvent();
97 G4HCofThisEvent* HCE = currentEvent->GetHCofThisEvent();
98 return HCE->GetHC(fCalorimeterHitsColID);
99 }
100 return 0;
101}

Member Data Documentation

◆ fStage

G4int RE01StackingAction::fStage
private

Definition at line 54 of file RE01StackingAction.hh.

◆ fCalorimeterHitsColID

G4int RE01StackingAction::fCalorimeterHitsColID
private

Definition at line 55 of file RE01StackingAction.hh.


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

Applications | User Support | Publications | Collaboration