EDepSimDokeBirks.hh
Go to the documentation of this file.
1 #ifndef EDepSim_Quanta_h
2 #define EDepSim_Quanta_h
3 
4 #include "globals.hh"
5 #include "templates.hh"
6 #include "Randomize.hh"
7 #include "G4ThreeVector.hh"
8 #include "G4ParticleMomentum.hh"
9 #include "G4Step.hh"
10 #include "G4VRestDiscreteProcess.hh"
11 #include "G4DynamicParticle.hh"
12 #include "G4Material.hh"
13 #include "G4EmSaturation.hh"
14 
15 /// Implement the Doke-Birk recombination probability found in NEST in an
16 /// optimized form for a LAr. The Doke-Birk parameterization is quite
17 /// accurate for LAr over the important energy ranges for a LArTPC (E more
18 /// than about an MeV).
19 ///
20 /// When this is called for a material other than LAr, then the G4 Birk's Law
21 /// implementation (G4EmSaturation) is used.
22 ///
23 /// This file should always be accompanied by the full NEST implementation
24 /// since this is a direct modification of the NEST code. The NEST code needs
25 /// to be provided since it is the arbitor of the performance of this code.
26 /// Where this code deviates from NEST, NEST is right. Finally, if this is
27 /// used, the NEST authors should be cited. They did all of the physics.
28 /// This is just an adaptation of their work to LAr so it's faster, but with
29 /// much less capability.
30 namespace EDepSim {class DokeBirks;}
31 class EDepSim::DokeBirks : public G4VRestDiscreteProcess //class definition
32 {
33 public:
34 
35  DokeBirks(const G4String& processName = "Doke-Birks",
36  G4ProcessType type = fElectromagnetic);
37  ~DokeBirks();
38 
39  /// Determine which particles this process should be applied too.
40  G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
41 
42  /// Don't limit the step since this process doesn't actually change the
43  /// energy deposit (it returns infinity). The process does set the
44  /// 'StronglyForced' condition so the DoIt is invoked at every step.
45  G4double GetMeanFreePath(const G4Track& aTrack,
46  G4double ,
47  G4ForceCondition* );
48 
49  /// Don't limit the step since this process doesn't actually change the
50  /// energy deposit (it returns infinity). The process does set the
51  /// 'StronglyForced' condition so the DoIt is invoked at every step.
52  G4double GetMeanLifeTime(const G4Track& aTrack,
53  G4ForceCondition* );
54 
55  /// Apply the scintilation process for an in-flight particle.
56  G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
57  const G4Step& aStep);
58 
59  /// Apply the scintilation process to a stopped particle (this just calles
60  /// PostStepDoIt().
61  G4VParticleChange* AtRestDoIt ( const G4Track& aTrack,
62  const G4Step& aStep);
63 
64 private:
65 
66  /// Get the expected electron electron dEdX as a function of energy. This
67  /// is normalized to denstity. The energy must be in MeV.
68  G4double CalculateElectronLET ( G4double E);
69 
70  // The EM Saturation model for when the material is not LAr.
71  G4EmSaturation* fEmSaturation;
72 
73 };
74 #endif
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
Determine which particles this process should be applied too.
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
Apply the scintilation process for an in-flight particle.
Construct a module from components.
Definition: TG4HitSegment.h:10
G4double CalculateElectronLET(G4double E)
G4VParticleChange * AtRestDoIt(const G4Track &aTrack, const G4Step &aStep)
G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *)
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
DokeBirks(const G4String &processName="Doke-Birks", G4ProcessType type=fElectromagnetic)
G4EmSaturation * fEmSaturation