Geant4 VMC Version 6.6
Loading...
Searching...
No Matches
TG4VSpecialCuts.h
Go to the documentation of this file.
1#ifndef TG4_V_SPECIAL_CUTS_H
2#define TG4_V_SPECIAL_CUTS_H
3
4//------------------------------------------------
5// The Geant4 Virtual Monte Carlo package
6// Copyright (C) 2007 - 2014 Ivana Hrivnacova
7// All rights reserved.
8//
9// For the licensing terms see geant4_vmc/LICENSE.
10// Contact: root-vmc@cern.ch
11//-------------------------------------------------
12
17
18#include <G4VProcess.hh>
19
20class TG4G3CutVector;
21class TG4Limits;
22class TG4TrackManager;
23
24class G4Track;
25class G4LossTableManager;
26
36
38{
39 public:
40 TG4VSpecialCuts(const G4String& processName);
41 virtual ~TG4VSpecialCuts();
42
43 // methods
45 virtual G4double GetMinEkine(
46 const TG4Limits& limits, const G4Track& track) const = 0;
47
48 virtual G4double PostStepGetPhysicalInteractionLength(const G4Track& track,
49 G4double previousStepSize, G4ForceCondition* condition);
50
51 virtual G4VParticleChange* PostStepDoIt(
52 const G4Track& track, const G4Step& step);
53
56 const G4Track&, G4double, G4double, G4double&, G4GPILSelection*)
57 {
58 return -1.0;
59 }
60
62 virtual G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&)
63 {
64 return 0;
65 }
66
69 const G4Track&, G4ForceCondition*)
70 {
71 return -1.0;
72 }
73
75 virtual G4VParticleChange* AtRestDoIt(const G4Track&, const G4Step&)
76 {
77 return 0;
78 }
79
80 private:
87
89 G4LossTableManager* fLossTableManager;
90
93};
94
95#endif // TG4_SPECIAL_CUTS_H
Vector of kinetic energy cut values with convenient set/get methods.
Extended G4UserLimits class.
Definition TG4Limits.h:38
The class for storing G4 tracks in VMC sack.
Abstract base class for a special process that activates kinetic energy cuts.
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
Not implemented.
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
Not implemented.
virtual ~TG4VSpecialCuts()
TG4VSpecialCuts(const TG4VSpecialCuts &right)
Not implemented.
virtual G4double GetMinEkine(const TG4Limits &limits, const G4Track &track) const =0
Return the kinetic energy limit.
TG4TrackManager * fTrackManager
Cached pointer to thread-local track manager.
G4LossTableManager * fLossTableManager
The G4LossTableManager instance.
TG4VSpecialCuts & operator=(const TG4VSpecialCuts &right)
Not implemented.
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
Not implemented.
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
Not implemented.
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &step)
TG4VSpecialCuts()
Not implemented.