Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/physics_lists/constructors/electromagnetic/src/G4EmPenelopePhysics.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 #include "G4EmPenelopePhysics.hh"
 28 #include "G4ParticleDefinition.hh"
 29 #include "G4SystemOfUnits.hh"
 30 
 31 // *** Processes and models
 32 
 33 // gamma
 34 #include "G4PhotoElectricEffect.hh"
 35 #include "G4PenelopePhotoElectricModel.hh"
 36 
 37 #include "G4ComptonScattering.hh"
 38 #include "G4PenelopeComptonModel.hh"
 39 
 40 #include "G4GammaConversion.hh"
 41 #include "G4PenelopeGammaConversionModel.hh"
 42 
 43 #include "G4RayleighScattering.hh" 
 44 #include "G4PenelopeRayleighModel.hh"
 45 
 46 #include "G4PEEffectFluoModel.hh"
 47 #include "G4KleinNishinaModel.hh"
 48 
 49 // e- and e+
 50 #include "G4eMultipleScattering.hh"
 51 
 52 #include "G4eIonisation.hh"
 53 #include "G4PenelopeIonisationModel.hh"
 54 
 55 #include "G4eBremsstrahlung.hh"
 56 #include "G4PenelopeBremsstrahlungModel.hh"
 57 
 58 // e+ only
 59 #include "G4eplusAnnihilation.hh"
 60 #include "G4PenelopeAnnihilationModel.hh"
 61 
 62 // hadrons
 63 #include "G4hMultipleScattering.hh"
 64 #include "G4MscStepLimitType.hh"
 65 
 66 #include "G4ePairProduction.hh"
 67 
 68 #include "G4hIonisation.hh"
 69 #include "G4ionIonisation.hh"
 70 #include "G4IonParametrisedLossModel.hh"
 71 #include "G4NuclearStopping.hh"
 72 #include "G4LindhardSorensenIonModel.hh"
 73 
 74 // msc models
 75 #include "G4UrbanMscModel.hh"
 76 #include "G4GoudsmitSaundersonMscModel.hh"
 77 #include "G4WentzelVIModel.hh"
 78 #include "G4CoulombScattering.hh"
 79 #include "G4eCoulombScatteringModel.hh"
 80 
 81 // utilities
 82 #include "G4EmBuilder.hh"
 83 #include "G4EmStandUtil.hh"
 84 
 85 // particles
 86 #include "G4Gamma.hh"
 87 #include "G4Electron.hh"
 88 #include "G4Positron.hh"
 89 #include "G4GenericIon.hh"
 90 
 91 //
 92 #include "G4PhysicsListHelper.hh"
 93 #include "G4BuilderType.hh"
 94 #include "G4EmModelActivator.hh"
 95 
 96 // factory
 97 #include "G4PhysicsConstructorFactory.hh"
 98 //
 99 G4_DECLARE_PHYSCONSTR_FACTORY(G4EmPenelopePhysics);
100 
101 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
102 
103 G4EmPenelopePhysics::G4EmPenelopePhysics(G4int ver, const G4String&)
104   : G4VPhysicsConstructor("G4EmPenelope")
105 {
106   SetVerboseLevel(ver);
107   G4EmParameters* param = G4EmParameters::Instance();
108   param->SetDefaults();
109   param->SetVerbose(ver);
110   param->SetMinEnergy(100*CLHEP::eV);
111   param->SetLowestElectronEnergy(100*CLHEP::eV);
112   param->SetNumberOfBinsPerDecade(20);
113   param->ActivateAngularGeneratorForIonisation(true);
114   param->SetStepFunction(0.2, 10*CLHEP::um);
115   param->SetStepFunctionMuHad(0.1, 50*CLHEP::um);
116   param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
117   param->SetStepFunctionIons(0.1, 1*CLHEP::um);
118   param->SetUseMottCorrection(true);
119   param->SetMscStepLimitType(fUseSafetyPlus);
120   param->SetMscSkin(3);
121   param->SetMscRangeFactor(0.08);
122   param->SetMuHadLateralDisplacement(true);
123   param->SetFluo(true);
124   param->SetUseICRU90Data(true);
125   param->SetFluctuationType(fUrbanFluctuation);
126   param->SetMaxNIELEnergy(1*CLHEP::MeV);
127   param->SetPIXEElectronCrossSectionModel("Penelope");
128   param->SetPositronAtRestModelType(fAllisonPositronium);
129   SetPhysicsType(bElectromagnetic);
130 }
131 
132 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
133 
134 G4EmPenelopePhysics::~G4EmPenelopePhysics() = default;
135 
136 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
137 
138 void G4EmPenelopePhysics::ConstructParticle()
139 {
140   // minimal set of particles for EM physics
141   G4EmBuilder::ConstructMinimalEmSet();
142 }
143 
144 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145 
146 void G4EmPenelopePhysics::ConstructProcess()
147 {
148   if(verboseLevel > 1) {
149     G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
150   }
151   G4EmBuilder::PrepareEMPhysics();
152   G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
153   G4EmParameters* param = G4EmParameters::Instance();
154   
155   // processes used by several particles
156   G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
157 
158   // high energy limit for e+- scattering models
159   G4double highEnergyLimit = param->MscEnergyLimit();
160 
161   // nuclear stopping
162   G4double nielEnergyLimit = param->MaxNIELEnergy();
163   G4NuclearStopping* pnuc = nullptr;
164   if(nielEnergyLimit > 0.0) {
165     pnuc = new G4NuclearStopping();
166     pnuc->SetMaxKinEnergy(nielEnergyLimit);
167   }
168 
169   //Applicability range for Penelope models
170   //for higher energies, the Standard models are used   
171   G4double PenelopeHighEnergyLimit = 1.0*CLHEP::GeV;
172 
173   // Add gamma EM Processes
174   G4ParticleDefinition* particle = G4Gamma::Gamma();
175 
176   //Photo-electric effect
177   G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
178   thePhotoElectricEffect->SetEmModel(new G4PEEffectFluoModel());
179   G4PenelopePhotoElectricModel* thePEPenelopeModel = new G4PenelopePhotoElectricModel();   
180   thePEPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
181   thePhotoElectricEffect->AddEmModel(0, thePEPenelopeModel);
182 
183   //Compton scattering
184   G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
185   G4PenelopeComptonModel* theComptonPenelopeModel = new G4PenelopeComptonModel();
186   theComptonScattering->SetEmModel(new G4KleinNishinaModel());
187   theComptonPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
188   theComptonScattering->AddEmModel(0, theComptonPenelopeModel);
189 
190   //Gamma conversion
191   G4GammaConversion* theGammaConversion = new G4GammaConversion();
192   G4PenelopeGammaConversionModel* theGCPenelopeModel = new G4PenelopeGammaConversionModel();
193   theGCPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
194   theGammaConversion->AddEmModel(0, theGCPenelopeModel);
195 
196   //Rayleigh scattering
197   G4RayleighScattering* theRayleigh = new G4RayleighScattering();
198   G4PenelopeRayleighModel* theRayleighPenelopeModel = new G4PenelopeRayleighModel();
199   theRayleigh->SetEmModel(theRayleighPenelopeModel);
200 
201   ph->RegisterProcess(thePhotoElectricEffect, particle);
202   ph->RegisterProcess(theComptonScattering, particle);
203   ph->RegisterProcess(theGammaConversion, particle);
204   ph->RegisterProcess(theRayleigh, particle);
205 
206   // e-
207   particle = G4Electron::Electron();
208 
209   // multiple scattering
210   G4GoudsmitSaundersonMscModel* msc1 = new G4GoudsmitSaundersonMscModel();
211   G4WentzelVIModel* msc2 = new G4WentzelVIModel();
212   msc1->SetHighEnergyLimit(highEnergyLimit);
213   msc2->SetLowEnergyLimit(highEnergyLimit);
214   G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
215 
216   G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel(); 
217   G4CoulombScattering* ss = new G4CoulombScattering();
218   ss->SetEmModel(ssm); 
219   ss->SetMinKinEnergy(highEnergyLimit);
220   ssm->SetLowEnergyLimit(highEnergyLimit);
221   ssm->SetActivationLowEnergyLimit(highEnergyLimit);
222       
223   //Ionisation
224   G4eIonisation* eioni = new G4eIonisation();
225   eioni->SetFluctModel(G4EmStandUtil::ModelOfFluctuations());
226   G4PenelopeIonisationModel* theIoniPenelope = new G4PenelopeIonisationModel();
227   theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);     
228   eioni->AddEmModel(0, theIoniPenelope);
229       
230   //Bremsstrahlung
231   G4eBremsstrahlung* brem = new G4eBremsstrahlung();
232   G4PenelopeBremsstrahlungModel* theBremPenelope = new G4PenelopeBremsstrahlungModel();
233   theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
234   brem->SetEmModel(theBremPenelope);
235 
236   G4ePairProduction* ee = new G4ePairProduction();
237 
238   // register processes
239   ph->RegisterProcess(eioni, particle);
240   ph->RegisterProcess(brem, particle);
241   ph->RegisterProcess(ee, particle);
242   ph->RegisterProcess(ss, particle);
243 
244   // e+
245   particle = G4Positron::Positron();
246     
247   // multiple scattering
248   msc1 = new G4GoudsmitSaundersonMscModel();
249   msc2 = new G4WentzelVIModel();
250   msc1->SetHighEnergyLimit(highEnergyLimit);
251   msc2->SetLowEnergyLimit(highEnergyLimit);
252   G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
253 
254   ssm = new G4eCoulombScatteringModel(); 
255   ss = new G4CoulombScattering();
256   ss->SetEmModel(ssm); 
257   ss->SetMinKinEnergy(highEnergyLimit);
258   ssm->SetLowEnergyLimit(highEnergyLimit);
259   ssm->SetActivationLowEnergyLimit(highEnergyLimit);
260       
261   //Ionisation
262   eioni = new G4eIonisation();
263   eioni->SetFluctModel(G4EmStandUtil::ModelOfFluctuations());
264   theIoniPenelope = new G4PenelopeIonisationModel();
265   theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
266   eioni->AddEmModel(0,theIoniPenelope);
267 
268   //Bremsstrahlung
269   brem = new G4eBremsstrahlung();
270   theBremPenelope = new G4PenelopeBremsstrahlungModel();
271   theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
272   brem->SetEmModel(theBremPenelope);
273   
274   //Annihilation
275   auto anni = new G4eplusAnnihilation();
276   auto theAnnPenelope = new G4PenelopeAnnihilationModel();
277   theAnnPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
278   anni->AddEmModel(0, theAnnPenelope);
279 
280   // register processes
281   ph->RegisterProcess(eioni, particle);
282   ph->RegisterProcess(brem, particle);
283   ph->RegisterProcess(ee, particle);
284   ph->RegisterProcess(anni, particle);
285   ph->RegisterProcess(ss, particle);
286 
287   // generic ion
288   particle = G4GenericIon::GenericIon();
289   G4ionIonisation* ionIoni = new G4ionIonisation();
290   ionIoni->SetEmModel(new G4LindhardSorensenIonModel());
291   ph->RegisterProcess(hmsc, particle);
292   ph->RegisterProcess(ionIoni, particle);
293   if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
294 
295   // muons, hadrons, ions
296   G4EmBuilder::ConstructCharged(hmsc, pnuc);
297 
298   // extra configuration
299   G4EmModelActivator mact(GetPhysicsName());
300 }
301 
302 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
303