Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4EmCaptureCascade.hh
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4EmCaptureCascade.hh 101422 2016-11-17 10:41:23Z gcosmo $
27 //
28 //-----------------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 // File name: G4EmCaptureCascade
33 //
34 // Author: V.Ivanchenko (Vladimir.Ivantchenko@cern.ch)
35 //
36 // Creation date: 22 April 2012 on base of G4MuMinusCaptureCascade
37 //
38 // Class Description:
39 //
40 // Simulation of electromagnetic cascade from capture level to K-shell
41 // of the mesonic atom
42 //
43 // Probabilities of gamma and Auger transitions from
44 // N.C.Mukhopadhyay Phys. Rep. 30 (1977) 1.
45 //
46 //-----------------------------------------------------------------------------
47 //
48 // Modifications:
49 //
50 //-----------------------------------------------------------------------------
51 
52 #ifndef G4EmCaptureCascade_h
53 #define G4EmCaptureCascade_h 1
54 
55 #include "globals.hh"
56 #include "G4Nucleus.hh"
57 #include "G4Track.hh"
58 #include "G4HadProjectile.hh"
59 #include "G4HadSecondary.hh"
60 #include "G4HadFinalState.hh"
61 #include "G4HadronicInteraction.hh"
62 #include "G4RandomDirection.hh"
63 #include "G4ParticleDefinition.hh"
64 #include "G4DynamicParticle.hh"
65 #include "G4ThreeVector.hh"
66 
68 {
69 public:
70 
71  explicit G4EmCaptureCascade();
72 
73  virtual ~G4EmCaptureCascade();
74 
75  virtual G4HadFinalState* ApplyYourself(const G4HadProjectile &aTrack,
76  G4Nucleus & targetNucleus );
77 
78  virtual void ModelDescription(std::ostream& outFile) const;
79 
80 private:
81 
82  inline void AddNewParticle(G4ParticleDefinition* aParticle,
83  G4double kinEnergy);
84 
85  // hide assignment operator as private
86  G4EmCaptureCascade& operator=(const G4EmCaptureCascade &right) = delete;
87  G4EmCaptureCascade(const G4EmCaptureCascade& ) = delete;
88 
89  G4HadFinalState result;
90  G4ParticleDefinition* theElectron;
91  G4ParticleDefinition* theGamma;
92  G4double fMuMass;
93  G4double fTime;
94  G4double fLevelEnergy[14];
95  G4double fKLevelEnergy[93];
96 };
97 
98 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
99 
100 inline void
101 G4EmCaptureCascade::AddNewParticle(G4ParticleDefinition* aParticle,
102  G4double kinEnergy)
103 {
104  G4DynamicParticle* dp = new G4DynamicParticle(aParticle,
106  kinEnergy);
107  G4HadSecondary hs(dp);
108  hs.SetTime(fTime);
109  result.AddSecondary(hs);
110 }
111 
112 #endif
113 
114 
virtual void ModelDescription(std::ostream &outFile) const
G4ThreeVector G4RandomDirection()
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
double G4double
Definition: G4Types.hh:76