Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/standard/src/G4eBremsstrahlung.cc

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 file
 30 //
 31 //
 32 // File name:     G4eBremsstrahlung
 33 //
 34 // Author:        Michel Maire
 35 //
 36 // Creation date: 26.06.1996
 37 //
 38 // Modified by Michel Maire, Vladimir Ivanchenko, Andreas Schaelicke  
 39 //
 40 // -------------------------------------------------------------------
 41 //
 42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 43 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 44 
 45 #include "G4eBremsstrahlung.hh"
 46 #include "G4SystemOfUnits.hh"
 47 #include "G4Gamma.hh"
 48 #include "G4SeltzerBergerModel.hh"
 49 #include "G4eBremsstrahlungRelModel.hh"
 50 #include "G4UnitsTable.hh"
 51 #include "G4LossTableManager.hh"
 52 
 53 #include "G4ProductionCutsTable.hh"
 54 #include "G4MaterialCutsCouple.hh"
 55 #include "G4EmParameters.hh"
 56 
 57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 58 
 59 using namespace std;
 60  
 61 G4eBremsstrahlung::G4eBremsstrahlung(const G4String& name):
 62   G4VEnergyLossProcess(name)
 63 {
 64   SetProcessSubType(fBremsstrahlung);
 65   SetSecondaryParticle(G4Gamma::Gamma());
 66   SetIonisation(false);
 67   SetCrossSectionType(fEmTwoPeaks);
 68 }
 69 
 70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 71 
 72 G4eBremsstrahlung::~G4eBremsstrahlung() = default;
 73 
 74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 75 
 76 G4bool G4eBremsstrahlung::IsApplicable(const G4ParticleDefinition& p)
 77 {
 78   return (&p == G4Electron::Electron() || &p == G4Positron::Positron());
 79 }
 80 
 81 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 82 
 83 void 
 84 G4eBremsstrahlung::InitialiseEnergyLossProcess(const G4ParticleDefinition*,
 85                  const G4ParticleDefinition*)
 86 {
 87   if(!isInitialised) {
 88     G4EmParameters* param = G4EmParameters::Instance();
 89 
 90     G4double emax = param->MaxKinEnergy();
 91     G4VEmFluctuationModel* fm = nullptr;
 92 
 93     if (nullptr == EmModel(0)) { SetEmModel(new G4SeltzerBergerModel()); }
 94     G4double energyLimit = std::min(EmModel(0)->HighEnergyLimit(), CLHEP::GeV);
 95     EmModel(0)->SetHighEnergyLimit(energyLimit);
 96     EmModel(0)->SetSecondaryThreshold(param->BremsstrahlungTh());
 97     AddEmModel(1, EmModel(0), fm);
 98 
 99     if(emax > energyLimit) {
100       if (nullptr == EmModel(1)) { 
101   SetEmModel(new G4eBremsstrahlungRelModel());
102       }
103       EmModel(1)->SetLowEnergyLimit(energyLimit);
104       EmModel(1)->SetHighEnergyLimit(emax); 
105       EmModel(1)->SetSecondaryThreshold(param->BremsstrahlungTh());
106       AddEmModel(1, EmModel(1), fm);
107     }
108     isInitialised = true;
109   }
110 }
111 
112 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
113 
114 void G4eBremsstrahlung::StreamProcessInfo(std::ostream& out) const
115 {
116   if(nullptr != EmModel(0)) {
117     G4EmParameters* param = G4EmParameters::Instance();
118     G4double eth = param->BremsstrahlungTh(); 
119     out << "      LPM flag: " << param->LPM() << " for E > " 
120   << EmModel(0)->HighEnergyLimit()/GeV << " GeV";
121     if(eth < DBL_MAX) { 
122       out << ",  VertexHighEnergyTh(GeV)= " << eth/GeV; 
123     }
124     out << G4endl;
125   }
126 }
127 
128 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 
129 
130 void G4eBremsstrahlung::ProcessDescription(std::ostream& out) const
131 {
132   out << "  Bremsstrahlung";
133   G4VEnergyLossProcess::ProcessDescription(out);
134 }
135 
136 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 
137