]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TGeant4/TG4SteppingAction.cxx
added getters for cut and control vectors
[u/mrichter/AliRoot.git] / TGeant4 / TG4SteppingAction.cxx
1 // $Id$
2 // Category: event
3 //
4 // Author: I.Hrivnacova
5 //
6 // Class TG4SteppingAction
7 // -----------------------
8 // See the class description in the header file.
9
10 #include "TG4SteppingAction.h"
11 #include "TG4VSensitiveDetector.h"
12 #include "TG4SDServices.h"
13 #include "TG4Globals.h"
14
15 #include <G4Track.hh>
16 #include <G4SteppingManager.hh>
17
18 //_____________________________________________________________________________
19 TG4SteppingAction::TG4SteppingAction() 
20   : fMaxNofSteps(kMaxNofSteps),
21     fStandardVerboseLevel(0),
22     fLoopVerboseLevel(1),
23     fLoopStepCounter(0)
24  {
25 //
26 }
27
28 //_____________________________________________________________________________
29 TG4SteppingAction::TG4SteppingAction(const TG4SteppingAction& right) {
30 //
31   TG4Globals::Exception("TG4SteppingAction is protected from copying.");
32 }
33
34 //_____________________________________________________________________________
35 TG4SteppingAction::~TG4SteppingAction() {
36 //
37 }
38
39 // operators
40
41 //_____________________________________________________________________________
42 TG4SteppingAction& 
43 TG4SteppingAction::operator=(const TG4SteppingAction &right)
44 {
45   // check assignement to self
46   if (this == &right) return *this;
47   
48   TG4Globals::Exception("TG4SteppingAction is protected from assigning.");
49
50   return *this;
51 }
52
53 // protected methods
54
55 //_____________________________________________________________________________
56 void TG4SteppingAction::PrintTrackInfo(const G4Track* track) const
57 {
58 // Prints the track info
59 // - taken from private G4TrackingManager::Verbose()
60 // and the standard header for verbose tracking
61 // - taken from G4SteppingVerbose::TrackingStarted().
62 // ---
63
64   // print track info
65   G4cout << G4endl;
66   G4cout << "*******************************************************"
67          << "**************************************************"
68          << G4endl;
69   G4cout << "* G4Track Information: " 
70          << "  Particle = " << track->GetDefinition()->GetParticleName() 
71          << "," 
72          << "   Track ID = " << track->GetTrackID() 
73          << "," 
74          << "   Parent ID = " << track->GetParentID() 
75          << G4endl;
76   G4cout << "*******************************************************"
77          << "**************************************************"
78          << G4endl;
79   G4cout << G4endl;
80       
81   // print header    
82 #ifdef G4_USE_G4BESTUNIT_FOR_VERBOSE
83     G4cout << G4std::setw( 5) << "Step#"  << " "
84            << G4std::setw( 8) << "X"      << "     "
85            << G4std::setw( 8) << "Y"      << "     "  
86            << G4std::setw( 8) << "Z"      << "     "
87            << G4std::setw( 9) << "KineE"  << "     "
88            << G4std::setw( 8) << "dE"     << "     "  
89            << G4std::setw(12) << "StepLeng"   << " "  
90            << G4std::setw(12) << "TrackLeng"  << " "
91            << G4std::setw(12) << "NextVolume" << " "
92            << G4std::setw( 8) << "ProcName"   << G4endl;             
93 #else
94     G4cout << G4std::setw( 5) << "Step#"      << " "
95            << G4std::setw( 8) << "X(mm)"      << " "
96            << G4std::setw( 8) << "Y(mm)"      << " "  
97            << G4std::setw( 8) << "Z(mm)"      << " "
98            << G4std::setw( 9) << "KinE(MeV)"  << " "
99            << G4std::setw( 8) << "dE(MeV)"    << " "  
100            << G4std::setw( 8) << "StepLeng"   << " "  
101            << G4std::setw( 9) << "TrackLeng"  << " "
102            << G4std::setw(11) << "NextVolume" << " "
103            << G4std::setw( 8) << "ProcName"   << G4endl;             
104 #endif
105 }
106
107 // public methods
108
109 //_____________________________________________________________________________
110 void TG4SteppingAction::UserSteppingAction(const G4Step* step)
111 {
112 // Called by G4 kernel at the end of each step.
113 // ---
114  
115   G4Track* track = step->GetTrack();  
116
117   // reset parameters at beginning of tracking
118   G4int stepNumber = track->GetCurrentStepNumber();
119   if (stepNumber == 1) { 
120     fStandardVerboseLevel = fpSteppingManager->GetverboseLevel();
121     fLoopStepCounter = 0;
122   }  
123   else if (fLoopStepCounter) {
124       // count steps after detecting looping track
125       fLoopStepCounter++;
126       if (fLoopStepCounter == kMaxNofLoopSteps) {
127
128         // stop the looping track
129         track->SetTrackStatus(fStopAndKill);      
130
131         // reset back parameters
132         fpSteppingManager->SetVerboseLevel(fStandardVerboseLevel);
133         fLoopStepCounter = 0;
134       } 
135     }  
136     else if (stepNumber> fMaxNofSteps) { 
137       
138        // print looping info
139        if (fLoopVerboseLevel > 0) {
140           G4cout << "*** Particle reached max step number ("
141                  << fMaxNofSteps << "). ***" << G4endl;
142           if (fStandardVerboseLevel == 0) PrintTrackInfo(track);
143        }        
144
145        // set loop verbose level 
146        fpSteppingManager->SetVerboseLevel(fLoopVerboseLevel);
147       
148        // start looping counter
149        fLoopStepCounter++;
150     }
151       
152   // call stepping action of derived class
153   SteppingAction(step);
154
155   // let sensitive detector process boundary step
156   // if crossing geometry border
157   // (this ensures compatibility with G3 that
158   // makes boundary step of zero length)
159
160   if (step->GetPostStepPoint()->GetStepStatus() == fGeomBoundary &&
161       step->GetTrack()->GetTrackStatus() == fAlive &&
162       step->GetTrack()->GetNextVolume() != 0) {
163
164 #ifdef TGEANT4_DEBUG
165     TG4VSensitiveDetector* tsd
166       = TG4SDServices::Instance()
167            ->GetSensitiveDetector(
168                 step->GetPostStepPoint()->GetPhysicalVolume()
169                     ->GetLogicalVolume()->GetSensitiveDetector());
170
171     if (tsd) = tsd->ProcessHitsOnBoundary((G4Step*)step);
172 #else
173     TG4VSensitiveDetector* tsd
174       = (TG4VSensitiveDetector*) step->GetPostStepPoint()->GetPhysicalVolume()
175           ->GetLogicalVolume()->GetSensitiveDetector();
176
177     if (tsd) tsd->ProcessHitsOnBoundary((G4Step*)step);
178 #endif     
179   }  
180 }
181