]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TGeant4/TG4VSpecialCuts.h
default verbose level change to 1
[u/mrichter/AliRoot.git] / TGeant4 / TG4VSpecialCuts.h
1 // $Id$
2 // Category: physics
3 //
4 // Author: I. Hrivnacova
5 //
6 // Class TG4VSpecialCuts
7 // ---------------------
8 // Abstract base class for a special process that activates 
9 // kinetic energy cuts.
10 // The pure virtual functions GetMinEkine have to be implemented
11 // by derived classes specific for each particle type
12 // (see TG4G3ParticleWSP.h).
13
14 #ifndef TG4_V_SPECIAL_CUTS_H
15 #define TG4_V_SPECIAL_CUTS_H
16
17 #include <G4VProcess.hh>
18
19 class TG4G3CutVector;
20 class TG4Limits;
21
22 class G4Track;
23
24 class TG4VSpecialCuts: public G4VProcess
25 {
26   public:
27     TG4VSpecialCuts(const G4String& processName);
28     // --> protected
29     // TG4VSpecialCuts();                  
30     virtual ~TG4VSpecialCuts();
31
32     // methods
33     virtual G4double GetMinEkine(const TG4Limits& limits,
34                                  const G4Track& track) const = 0;
35     
36     virtual G4double PostStepGetPhysicalInteractionLength(
37                          const G4Track& track, G4double previousStepSize,
38                          G4ForceCondition* condition);
39
40     virtual G4VParticleChange* PostStepDoIt(const G4Track& track, 
41                          const G4Step& step);
42                             
43     virtual G4double AlongStepGetPhysicalInteractionLength(
44                          const G4Track&, G4double, G4double, G4double&,
45                          G4GPILSelection*)
46                          { return -1.0; }
47
48     virtual G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&)
49                          { return 0; }
50
51     virtual G4double AtRestGetPhysicalInteractionLength(const G4Track&,
52                          G4ForceCondition* )
53                          { return -1.0; }
54                             
55     virtual G4VParticleChange* AtRestDoIt(
56                          const G4Track&, const G4Step&)
57                          { return 0; }
58
59   protected:
60     TG4VSpecialCuts();             
61     TG4VSpecialCuts(const TG4VSpecialCuts& right);
62     
63     // operators
64     TG4VSpecialCuts& operator = (const TG4VSpecialCuts& right);
65 };
66
67 #endif //TG4_SPECIAL_CUTS_H
68
69
70