Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/muons/include/G4MuBremsstrahlungModel.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 //
 27 // -------------------------------------------------------------------
 28 //
 29 // GEANT4 Class header file
 30 //
 31 //
 32 // File name:     G4MuBremsstrahlungModel
 33 //
 34 // Author:        Vladimir Ivanchenko on base of Laszlo Urban code
 35 // 
 36 // Creation date: 18.05.2002
 37 //
 38 // Modifications:
 39 //
 40 // 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
 41 // 27-01-03 Make models region aware (V.Ivanchenko)
 42 // 13-02-03 Add name (V.Ivanchenko)
 43 // 10-02-04 Add lowestKinEnergy (V.Ivanchenko)
 44 // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
 45 // 13-02-06 Add ComputeCrossSectionPerAtom (mma)
 46 // 11-10-07 Add ignoreCut flag (V.Ivanchenko) 
 47 // 28-02-08 Reorganized protected methods and members (V.Ivanchenko) 
 48 // 06-03-08 Remove obsolete methods and members (V.Ivanchenko) 
 49 // 31-05-13 Use element selectors instead of local data (V.Ivanchenko)
 50 //
 51 // Class Description:
 52 //
 53 // Implementation of bremssrahlung by muons
 54 // A.G. Bogdanov et al., IEEE Trans. Nuc. Sci., Vol.53, No.2, 2006
 55 //
 56 // -------------------------------------------------------------------
 57 //
 58 
 59 #ifndef G4MuBremsstrahlungModel_h
 60 #define G4MuBremsstrahlungModel_h 1
 61 
 62 #include "G4VEmModel.hh"
 63 #include "G4NistManager.hh"
 64 #include "G4Exp.hh"
 65 
 66 class G4Element;
 67 class G4ParticleChangeForLoss;
 68 
 69 class G4MuBremsstrahlungModel : public G4VEmModel
 70 {
 71 
 72 public:
 73 
 74   explicit G4MuBremsstrahlungModel(const G4ParticleDefinition* p = nullptr,
 75                                    const G4String& nam = "MuBrem");
 76 
 77   ~G4MuBremsstrahlungModel() override = default;
 78 
 79   void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
 80 
 81   void InitialiseLocal(const G4ParticleDefinition*, 
 82                        G4VEmModel* masterModel) override;
 83 
 84   G4double MinEnergyCut(const G4ParticleDefinition*,
 85                         const G4MaterialCutsCouple*) override;
 86             
 87   G4double ComputeCrossSectionPerAtom(
 88          const G4ParticleDefinition*,
 89          G4double kineticEnergy,
 90          G4double Z, G4double A,
 91          G4double cutEnergy,
 92          G4double maxEnergy) override;
 93                      
 94   G4double ComputeDEDXPerVolume(const G4Material*,
 95                                 const G4ParticleDefinition*,
 96                                 G4double kineticEnergy,
 97                                 G4double cutEnergy) override;
 98             
 99   void SampleSecondaries(std::vector<G4DynamicParticle*>*,
100                          const G4MaterialCutsCouple*,
101        const G4DynamicParticle*,
102        G4double tmin,
103        G4double maxEnergy) override;
104 
105   inline void SetLowestKineticEnergy(G4double e);
106 
107   G4double MinPrimaryEnergy(const G4Material*, const G4ParticleDefinition*, 
108                             G4double) override;
109 
110   // hide assignment operator
111   G4MuBremsstrahlungModel & 
112     operator=(const  G4MuBremsstrahlungModel &right) = delete;
113   G4MuBremsstrahlungModel(const  G4MuBremsstrahlungModel&) = delete;
114 
115 protected:
116 
117   G4double ComputMuBremLoss(G4double Z, G4double tkin, G4double cut);
118   
119   G4double ComputeMicroscopicCrossSection(G4double tkin,
120             G4double Z,
121             G4double cut);
122 
123   virtual G4double ComputeDMicroscopicCrossSection(G4double tkin,
124                G4double Z,
125                G4double gammaEnergy);
126 
127   void SetParticle(const G4ParticleDefinition*);
128 
129 protected:
130 
131   const G4ParticleDefinition* particle = nullptr;
132   G4ParticleDefinition* theGamma = nullptr;
133   G4ParticleChangeForLoss* fParticleChange = nullptr;
134   G4NistManager* nist = nullptr;
135 
136   G4double mass = 1.0;
137   G4double rmass = 1.0;
138   G4double cc = 1.0;
139   G4double coeff = 1.0;
140   G4double sqrte;
141   G4double bh = 202.4;
142   G4double bh1 = 446.;
143   G4double btf = 183.;
144   G4double btf1 = 1429.;
145 
146   G4double lowestKinEnergy;
147   G4double minThreshold;
148 
149   static const G4double xgi[6];
150   static const G4double wgi[6];
151   static G4double fDN[93];
152 };
153 
154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
155 
156 inline void G4MuBremsstrahlungModel::SetLowestKineticEnergy(G4double e) 
157 {
158   lowestKinEnergy = e;
159 }
160 
161 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
162 
163 #endif
164