4 // Geant4 implementation of the MonteCarlo interface methods
5 // for access to Geant4 at step level
7 // The public methods that do not implement AliMC methods
8 // are commented as G4 specific
10 #ifndef TG4_STEP_MANAGER_H
11 #define TG4_STEP_MANAGER_H
13 #include <G4ThreeVector.hh>
19 class G4SteppingManager;
23 enum TG4StepStatus { kPreStepPoint, kPostStepPoint};
30 // TG4StepManager(const TG4StepManager& right);
31 virtual ~TG4StepManager();
33 // static access method
34 static TG4StepManager* Instance();
37 void StopTrack(); //new
38 void StopEvent(); //new
39 void Rndm(Float_t* array, const Int_t size) const;
42 void SetStep(G4Step* step, TG4StepStatus status); // G4 specific
43 void SetSteppingManager(G4SteppingManager* manager); // G4 specific
44 void SetMaxStep(Float_t step);
45 void SetMaxNStep(Int_t maxNofSteps); //??
46 void SetUserDecay(Int_t pdg); //NEW
49 G4Step* GetStep() const; // G4 specific
50 TG4StepStatus GetStepStatus() const; // G4 specific
53 Int_t CurrentVolID(Int_t& copyNo) const;
54 Int_t CurrentVolOffID(Int_t off, Int_t& copyNo) const;
55 const char* CurrentVolName() const;
56 const char* CurrentVolOffName(Int_t off) const;
57 Int_t CurrentMaterial(Float_t &a, Float_t &z, Float_t &dens,
58 Float_t &radl, Float_t &absl) const;
59 void Gmtod(Float_t* xm, Float_t* xd, Int_t iflag); //new
60 void Gdtom(Float_t* xd, Float_t* xm, Int_t iflag); //new
61 Float_t MaxStep() const;
62 Int_t GetMaxNStep() const; //??
63 Int_t GetMedium() const; //??
67 void TrackPosition(TLorentzVector& position) const;
68 void TrackMomentum(TLorentzVector& momentum) const;
69 void TrackVertexPosition(TLorentzVector& position) const;
70 void TrackVertexMomentum(TLorentzVector& momentum) const;
71 Float_t TrackStep() const;
72 Float_t TrackLength() const;
73 Float_t TrackTime() const;
76 Int_t TrackPid() const;
77 Float_t TrackCharge() const;
78 Float_t TrackMass() const;
82 Bool_t IsTrackInside() const;
83 Bool_t IsTrackEntering() const;
84 Bool_t IsTrackExiting() const;
85 Bool_t IsTrackOut() const;
86 Bool_t IsTrackDisappeared() const;
87 Bool_t IsTrackStop() const;
88 Bool_t IsTrackAlive() const;
89 Bool_t IsNewTrack() const;
92 Int_t NSecondaries() const;
93 void GetSecondary(Int_t isec, Int_t& particleId,
94 TLorentzVector& position, TLorentzVector& momentum);
95 const char* ProdProcess() const;
98 TG4StepManager(const TG4StepManager& right);
101 TG4StepManager& operator=(const TG4StepManager& right);
105 void SetTLorentzVector(G4ThreeVector xyz, G4double t,
106 TLorentzVector& lv) const;
108 // static data members
109 static TG4StepManager* fgInstance; //this instance
112 G4Step* fStep; //current step
113 TG4StepStatus fStepStatus; //step status that decides whether
114 //track properties will be returned
115 //from PreStepPoint or PostStepPoint
116 G4SteppingManager* fSteppingManager; //G4SteppingManager
121 inline TG4StepManager* TG4StepManager::Instance()
122 { return fgInstance; }
124 inline void TG4StepManager::SetStep(G4Step* step, TG4StepStatus status)
125 { fStep = step; fStepStatus = status; }
127 inline void TG4StepManager::SetSteppingManager(G4SteppingManager* manager)
128 { fSteppingManager = manager; }
130 inline G4Step* TG4StepManager::GetStep() const
133 inline TG4StepStatus TG4StepManager::GetStepStatus() const
134 { return fStepStatus; }
136 #endif //TG4_STEP_MANAGER_H