Geant4 Cross Reference |
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 // ClassName: G4EmStandardPhysics_option2 30 // 31 // Author: V.Ivanchenko 09.11.2005 32 // 33 // Modified: 34 // 19.12.2005 V.Ivanchenko rename 71 -> 72 35 // 15.06.2006 V.Ivanchenko use this class as a constructor of fast EM physics 36 // 13.11.2006 V.Ivanchenko use G4hMultipleScattering 37 // 14.11.2006 V.Ivanchenko use sub-cutoff option for all particles 38 // 13.02.2007 V.Ivanchenko use default msc 39 // 15.05.2007 V.Ivanchenko rename to _option2 40 // 13.03.2008 V.Ivanchenko use G4eMultipleScattering 41 // 21.04.2008 V.Ivanchenko add long-lived D and B mesons; use spline 42 // 28.05.2008 V.Ivanchenko linLossLimit=0.01; added hBrem and hPairProd processes 43 // 44 //---------------------------------------------------------------------------- 45 // 46 47 #include "G4EmStandardPhysics_option2.hh" 48 49 #include "G4SystemOfUnits.hh" 50 #include "G4ParticleDefinition.hh" 51 #include "G4LossTableManager.hh" 52 #include "G4EmParameters.hh" 53 #include "G4EmBuilder.hh" 54 55 #include "G4ComptonScattering.hh" 56 #include "G4GammaConversion.hh" 57 #include "G4PhotoElectricEffect.hh" 58 #include "G4LivermorePhotoElectricModel.hh" 59 #include "G4KleinNishinaModel.hh" 60 #include "G4RayleighScattering.hh" 61 62 #include "G4hMultipleScattering.hh" 63 #include "G4CoulombScattering.hh" 64 #include "G4eCoulombScatteringModel.hh" 65 #include "G4UrbanMscModel.hh" 66 #include "G4WentzelVIModel.hh" 67 68 #include "G4eIonisation.hh" 69 #include "G4eBremsstrahlung.hh" 70 #include "G4eBremsstrahlungRelModel.hh" 71 #include "G4eplusAnnihilation.hh" 72 #include "G4Generator2BS.hh" 73 #include "G4SeltzerBergerModel.hh" 74 75 #include "G4hIonisation.hh" 76 #include "G4ionIonisation.hh" 77 78 #include "G4ParticleTable.hh" 79 #include "G4Gamma.hh" 80 #include "G4Electron.hh" 81 #include "G4Positron.hh" 82 #include "G4GenericIon.hh" 83 84 #include "G4PhysicsListHelper.hh" 85 #include "G4BuilderType.hh" 86 #include "G4EmModelActivator.hh" 87 #include "G4GammaGeneralProcess.hh" 88 89 // factory 90 #include "G4PhysicsConstructorFactory.hh" 91 // 92 G4_DECLARE_PHYSCONSTR_FACTORY(G4EmStandardPhysics_option2); 93 94 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 95 96 G4EmStandardPhysics_option2::G4EmStandardPhysics_option2(G4int ver, 97 const G4String&) 98 : G4VPhysicsConstructor("G4EmStandard_opt2") 99 { 100 SetVerboseLevel(ver); 101 G4EmParameters* param = G4EmParameters::Instance(); 102 param->SetDefaults(); 103 param->SetVerbose(ver); 104 param->SetApplyCuts(false); 105 param->SetStepFunction(0.8, 1*CLHEP::mm); 106 param->SetMscRangeFactor(0.2); 107 param->SetLateralDisplacement(false); 108 param->SetMscStepLimitType(fMinimal); 109 SetPhysicsType(bElectromagnetic); 110 } 111 112 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 113 114 G4EmStandardPhysics_option2::~G4EmStandardPhysics_option2() 115 {} 116 117 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 118 119 void G4EmStandardPhysics_option2::ConstructParticle() 120 { 121 // minimal set of particles for EM physics 122 G4EmBuilder::ConstructMinimalEmSet(); 123 } 124 125 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 126 127 void G4EmStandardPhysics_option2::ConstructProcess() 128 { 129 if(verboseLevel > 1) { 130 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl; 131 } 132 G4EmBuilder::PrepareEMPhysics(); 133 134 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper(); 135 136 // processes used by several particles 137 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc"); 138 G4NuclearStopping* pnuc(nullptr); 139 140 // high energy limit for e+- scattering models and bremsstrahlung 141 G4double highEnergyLimit = G4EmParameters::Instance()->MscEnergyLimit(); 142 143 // Add gamma EM Processes 144 G4ParticleDefinition* particle = G4Gamma::Gamma(); 145 146 G4PhotoElectricEffect* pee = new G4PhotoElectricEffect(); 147 pee->SetEmModel(new G4LivermorePhotoElectricModel()); 148 149 if(G4EmParameters::Instance()->GeneralProcessActive()) { 150 G4GammaGeneralProcess* sp = new G4GammaGeneralProcess(); 151 sp->AddEmProcess(pee); 152 sp->AddEmProcess(new G4ComptonScattering()); 153 sp->AddEmProcess(new G4GammaConversion()); 154 G4LossTableManager::Instance()->SetGammaGeneralProcess(sp); 155 ph->RegisterProcess(sp, particle); 156 157 } else { 158 ph->RegisterProcess(pee, particle); 159 ph->RegisterProcess(new G4ComptonScattering(), particle); 160 ph->RegisterProcess(new G4GammaConversion(), particle); 161 } 162 163 // e- 164 particle = G4Electron::Electron(); 165 166 G4eIonisation* eioni = new G4eIonisation(); 167 168 G4UrbanMscModel* msc1 = new G4UrbanMscModel(); 169 G4WentzelVIModel* msc2 = new G4WentzelVIModel(); 170 msc1->SetHighEnergyLimit(highEnergyLimit); 171 msc2->SetLowEnergyLimit(highEnergyLimit); 172 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle); 173 174 G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel(); 175 G4CoulombScattering* ss = new G4CoulombScattering(); 176 ss->SetEmModel(ssm); 177 ss->SetMinKinEnergy(highEnergyLimit); 178 ssm->SetLowEnergyLimit(highEnergyLimit); 179 ssm->SetActivationLowEnergyLimit(highEnergyLimit); 180 181 ph->RegisterProcess(eioni, particle); 182 ph->RegisterProcess(new G4eBremsstrahlung(), particle); 183 ph->RegisterProcess(ss, particle); 184 185 // e+ 186 particle = G4Positron::Positron(); 187 eioni = new G4eIonisation(); 188 189 msc1 = new G4UrbanMscModel(); 190 msc2 = new G4WentzelVIModel(); 191 msc1->SetHighEnergyLimit(highEnergyLimit); 192 msc2->SetLowEnergyLimit(highEnergyLimit); 193 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle); 194 195 ssm = new G4eCoulombScatteringModel(); 196 ss = new G4CoulombScattering(); 197 ss->SetEmModel(ssm); 198 ss->SetMinKinEnergy(highEnergyLimit); 199 ssm->SetLowEnergyLimit(highEnergyLimit); 200 ssm->SetActivationLowEnergyLimit(highEnergyLimit); 201 202 ph->RegisterProcess(eioni, particle); 203 ph->RegisterProcess(new G4eBremsstrahlung(), particle); 204 ph->RegisterProcess(new G4eplusAnnihilation(), particle); 205 ph->RegisterProcess(ss, particle); 206 207 // generic ion 208 particle = G4GenericIon::GenericIon(); 209 G4ionIonisation* ionIoni = new G4ionIonisation(); 210 ph->RegisterProcess(hmsc, particle); 211 ph->RegisterProcess(ionIoni, particle); 212 213 // muons, hadrons ions 214 G4EmBuilder::ConstructCharged(hmsc, pnuc); 215 216 // extra configuration 217 G4EmModelActivator mact(GetPhysicsName()); 218 } 219 220 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 221