Geant4  10.02.p01
PhysicsList.cc
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 // This example is provided by the Geant4-DNA collaboration
27 // Any report or published results obtained using the Geant4-DNA software
28 // shall cite the following Geant4-DNA collaboration publication:
29 // Med. Phys. 37 (2010) 4692-4708
30 // The Geant4-DNA web site is available at http://geant4-dna.org
31 //
32 // $ID$
35 
36 #include "PhysicsList.hh"
37 #include "G4SystemOfUnits.hh"
38 
39 // Geant4-DNA MODELS
40 
41 #include "G4DNAElastic.hh"
44 
45 #include "G4DNAExcitation.hh"
48 
49 #include "G4DNAIonisation.hh"
52 
53 #include "G4DNAChargeDecrease.hh"
55 
56 #include "G4DNAChargeIncrease.hh"
58 
59 #include "G4DNAAttachment.hh"
61 
62 #include "G4DNAVibExcitation.hh"
64 
65 //
66 
67 #include "G4LossTableManager.hh"
68 #include "G4EmConfigurator.hh"
69 #include "G4VEmModel.hh"
70 #include "G4DummyModel.hh"
71 #include "G4eIonisation.hh"
72 #include "G4hIonisation.hh"
73 #include "G4ionIonisation.hh"
74 #include "G4eMultipleScattering.hh"
75 #include "G4hMultipleScattering.hh"
76 #include "G4BraggIonGasModel.hh"
78 #include "G4UrbanMscModel.hh"
79 #include "G4MollerBhabhaModel.hh"
80 #include "G4IonFluctuations.hh"
82 
83 #include "G4ElectronCapture.hh"
84 
85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
86 
88 {
94 
95  SetVerboseLevel(1);
96 }
97 
98 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
99 
101 {}
102 
103 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
104 
106 {
107  ConstructBosons();
110 }
111 
112 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
113 
115 {
116  // gamma
118 }
119 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
120 
122 {
123  // leptons
126 }
127 
128 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
129 
130 //DNA
132 //ENDDNA
133 
135 {
136  // baryons
139 
140  // Geant4 DNA new particles
141  G4DNAGenericIonsManager * genericIonsManager;
142  genericIonsManager=G4DNAGenericIonsManager::Instance();
143  genericIonsManager->GetIon("alpha++");
144  genericIonsManager->GetIon("alpha+");
145  genericIonsManager->GetIon("helium");
146  genericIonsManager->GetIon("hydrogen");
147 }
148 
149 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
150 
152 {
154  ConstructEM();
156 }
157 
158 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
159 
161 {
162 
163  theParticleIterator->reset();
164 
165  while( (*theParticleIterator)() )
166  {
167 
168  G4ParticleDefinition* particle = theParticleIterator->value();
169  G4ProcessManager* pmanager = particle->GetProcessManager();
170  G4String particleName = particle->GetParticleName();
171 
172  // *********************************
173  // 1) Processes for the World region
174  // *********************************
175 
176  if (particleName == "e-") {
177 
178  // STANDARD msc is active in the world
180  pmanager->AddProcess(msc, -1, 1, 1);
181 
182  // STANDARD ionisation is active in the world
183  G4eIonisation* eion = new G4eIonisation();
184  eion->SetEmModel(new G4MollerBhabhaModel(), 1);
185  pmanager->AddProcess(eion, -1, 2, 2);
186 
187  // DNA elastic is not active in the world
188  G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
189  theDNAElasticProcess->SetEmModel(new G4DummyModel(),1);
190  pmanager->AddDiscreteProcess(theDNAElasticProcess);
191 
192  // DNA excitation is not active in the world
193  G4DNAExcitation* dnaex = new G4DNAExcitation("e-_G4DNAExcitation");
194  dnaex->SetEmModel(new G4DummyModel(),1);
195  pmanager->AddDiscreteProcess(dnaex);
196 
197  // DNA ionisation is not active in the world
198  G4DNAIonisation* dnaioni = new G4DNAIonisation("e-_G4DNAIonisation");
199  dnaioni->SetEmModel(new G4DummyModel(),1);
200  pmanager->AddDiscreteProcess(dnaioni);
201 
202  // DNA attachment is not active in the world
203  G4DNAAttachment* dnaatt = new G4DNAAttachment("e-_G4DNAAttachment");
204  dnaatt->SetEmModel(new G4DummyModel(),1);
205  pmanager->AddDiscreteProcess(dnaatt);
206 
207  // DNA vib. excitation is not active in the world
208  G4DNAVibExcitation* dnavib =
209  new G4DNAVibExcitation("e-_G4DNAVibExcitation");
210  dnavib->SetEmModel(new G4DummyModel(),1);
211  pmanager->AddDiscreteProcess(dnavib);
212 
213  // THE FOLLOWING PROCESS WILL KILL ALL ELECTRONS BELOW A
214  // SELECTED ENERY THRESHOLD
215  // Capture of low-energy e-
216  G4ElectronCapture* ecap = new G4ElectronCapture("Target", 5.1*eV);
217  pmanager->AddDiscreteProcess(ecap);
218 
219  } else if ( particleName == "proton" ) {
220 
221  // STANDARD msc is active in the world
223  pmanager->AddProcess(msc, -1, 1, 1);
224 
225  // STANDARD ionisation is active in the world
226  G4hIonisation* hion = new G4hIonisation();
227  hion->SetEmModel(new G4BraggIonGasModel(), 1);
228  hion->SetEmModel(new G4BetheBlochIonGasModel(), 2);
229  pmanager->AddProcess(hion, -1, 2, 2);
230 
231  // DNA excitation is not active in the world
232  G4DNAExcitation* dnaex = new G4DNAExcitation("proton_G4DNAExcitation");
233  dnaex->SetEmModel(new G4DummyModel(),1);
234  dnaex->SetEmModel(new G4DummyModel(),2);
235  pmanager->AddDiscreteProcess(dnaex);
236 
237  // DNA ionisation is not active in the world
238  G4DNAIonisation* dnaioni = new G4DNAIonisation("proton_G4DNAIonisation");
239  dnaioni->SetEmModel(new G4DummyModel(),1);
240  dnaioni->SetEmModel(new G4DummyModel(),2);
241  pmanager->AddDiscreteProcess(dnaioni);
242 
243  // DNA charge decrease is ACTIVE in the world since
244  // no corresponding STANDARD process exist
245  pmanager->AddDiscreteProcess(
246  new G4DNAChargeDecrease("proton_G4DNAChargeDecrease"));
247 
248  } else if ( particleName == "hydrogen" ) {
249 
250  // DNA processes are ACTIVE in the world since
251  // no corresponding STANDARD processes exist
252  pmanager->AddDiscreteProcess(
253  new G4DNAIonisation("hydrogen_G4DNAIonisation"));
254  pmanager->AddDiscreteProcess(
255  new G4DNAExcitation("hydrogen_G4DNAExcitation"));
256  pmanager->AddDiscreteProcess(
257  new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease"));
258 
259  } else if (particleName == "GenericIon")
260  { // THIS IS NEEDED FOR STANDARD ALPHA G4ionIonisation PROCESS
261 
262  // STANDARD msc is active in the world
263  pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
264 
265  // STANDARD ionisation is active in the world
266  G4ionIonisation* hion = new G4ionIonisation();
267  hion->SetEmModel(new G4BraggIonGasModel(),1);
268  hion->SetEmModel(new G4BetheBlochIonGasModel(), 2);
269  pmanager->AddProcess(hion, -1, 2, 2);
270 
271  } else if ( particleName == "alpha" ) {
272 
273  // STANDARD msc is active in the world
275  pmanager->AddProcess(msc, -1, 1, 1);
276 
277  // STANDARD ionisation is active in the world
278  G4ionIonisation* hion = new G4ionIonisation();
279  hion->SetEmModel(new G4BraggIonGasModel(),1);
280  hion->SetEmModel(new G4BetheBlochIonGasModel(), 2);
281  pmanager->AddProcess(hion, -1, 2, 2);
282 
283  // DNA excitation is not active in the world
284  G4DNAExcitation* dnaex = new G4DNAExcitation("alpha_G4DNAExcitation");
285  dnaex->SetEmModel(new G4DummyModel(),1);
286  pmanager->AddDiscreteProcess(dnaex);
287 
288  // DNA ionisation is not active in the world
289  G4DNAIonisation* dnaioni = new G4DNAIonisation("alpha_G4DNAIonisation");
290  dnaioni->SetEmModel(new G4DummyModel(),1);
291  pmanager->AddDiscreteProcess(dnaioni);
292 
293  // DNA charge decrease is ACTIVE in the world since no
294  // corresponding STANDARD process exist
295  pmanager->AddDiscreteProcess(
296  new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease"));
297 
298  } else if ( particleName == "alpha+" ) {
299 
300  // DNA processes are ACTIVE in the world since no
301  // corresponding STANDARD processes exist
302  pmanager->AddDiscreteProcess(
303  new G4DNAExcitation("alpha+_G4DNAExcitation"));
304  pmanager->AddDiscreteProcess(
305  new G4DNAIonisation("alpha+_G4DNAIonisation"));
306  pmanager->AddDiscreteProcess(
307  new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease"));
308  pmanager->AddDiscreteProcess(
309  new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease"));
310 
311  } else if ( particleName == "helium" ) {
312 
313  // DNA processes are ACTIVE in the world since no
314  // corresponding STANDARD processes exist
315  pmanager->AddDiscreteProcess(
316  new G4DNAExcitation("helium_G4DNAExcitation"));
317  pmanager->AddDiscreteProcess(
318  new G4DNAIonisation("helium_G4DNAIonisation"));
319  pmanager->AddDiscreteProcess(
320  new G4DNAChargeIncrease("helium_G4DNAChargeIncrease"));
321 
322  }
323  }
324 
325  // **************************************
326  // 2) Define processes for Target region
327  // **************************************
328 
329  // STANDARD EM processes should be inactivated when
330  // corresponding DNA processes are used
331  // - STANDARD EM e- processes are inactivated below 1 MeV
332  // - STANDARD EM proton & alpha processes are inactivated below
333  // standEnergyLimit
334  G4double standEnergyLimit = 9.9*MeV;
335  //
336 
337  G4double massFactor = 1.0079/4.0026;
340 
341  G4VEmModel* mod;
342 
343  // *** e-
344 
345  // ---> STANDARD EM processes are inactivated below 1 MeV
346 
347  mod = new G4UrbanMscModel();
349  em_config->SetExtraEmModel("e-","msc",mod,"Target");
350 
351  mod = new G4MollerBhabhaModel();
353  em_config->SetExtraEmModel("e-",
354  "eIoni",
355  mod,
356  "Target",
357  0.0,
358  100*TeV,
359  new G4UniversalFluctuation());
360 
361  // ---> DNA processes activated
362 
363  mod = new G4DNAChampionElasticModel();
364  em_config->SetExtraEmModel("e-","e-_G4DNAElastic",
365  mod,"Target",0.0,1*MeV);
366 
367  mod = new G4DNABornIonisationModel();
368  em_config->SetExtraEmModel("e-","e-_G4DNAIonisation",
369  mod,"Target",11*eV,1*MeV);
370 
371  mod = new G4DNABornExcitationModel();
372  em_config->SetExtraEmModel("e-","e-_G4DNAExcitation",
373  mod,"Target",9*eV,1*MeV);
374 
375  mod = new G4DNAMeltonAttachmentModel();
376  em_config->SetExtraEmModel("e-","e-_G4DNAAttachment",
377  mod,"Target",4*eV,13*eV);
378 
379  mod = new G4DNASancheExcitationModel();
380  em_config->SetExtraEmModel("e-","e-_G4DNAVibExcitation",
381  mod,"Target",2*eV,100*eV);
382 
383  // *** proton
384 
385  // ---> STANDARD EM processes inactivated below standEnergyLimit
386 
387  // STANDARD msc is still active
388  // Inactivate following STANDARD processes
389 
390  mod = new G4BraggIonGasModel();
391  mod->SetActivationLowEnergyLimit(standEnergyLimit);
392  em_config->SetExtraEmModel("proton","hIoni",
393  mod,"Target",0.0,2*MeV, new G4IonFluctuations());
394 
395  mod = new G4BetheBlochIonGasModel();
396  mod->SetActivationLowEnergyLimit(standEnergyLimit);
397  em_config->SetExtraEmModel("proton","hIoni",
398  mod,"Target",2*MeV,100*TeV,
399  new G4UniversalFluctuation());
400 
401  // ---> DNA processes activated
402 
403  mod = new G4DNARuddIonisationModel();
404  em_config->SetExtraEmModel("proton","proton_G4DNAIonisation",
405  mod,"Target",0.0,0.5*MeV);
406 
407  mod = new G4DNABornIonisationModel();
408  em_config->SetExtraEmModel("proton","proton_G4DNAIonisation",
409  mod,"Target",0.5*MeV,10*MeV);
410 
412  em_config->SetExtraEmModel("proton","proton_G4DNAExcitation",
413  mod,"Target",10*eV,0.5*MeV);
414 
415  mod = new G4DNABornExcitationModel();
416  em_config->SetExtraEmModel("proton","proton_G4DNAExcitation",
417  mod,"Target",0.5*MeV,10*MeV);
418 
419  // *** alpha
420 
421  // ---> STANDARD EM processes inactivated below standEnergyLimit
422 
423  // STANDARD msc is still active
424  // Inactivate following STANDARD processes
425 
426  mod = new G4BraggIonGasModel();
427  mod->SetActivationLowEnergyLimit(standEnergyLimit);
428  em_config->SetExtraEmModel("alpha","ionIoni",
429  mod,"Target",0.0,2*MeV/massFactor,
430  new G4IonFluctuations());
431 
432  mod = new G4BetheBlochIonGasModel();
433  mod->SetActivationLowEnergyLimit(standEnergyLimit);
434  em_config->SetExtraEmModel("alpha","ionIoni",
435  mod,"Target",2*MeV/massFactor,100*TeV,
436  new G4UniversalFluctuation());
437 
438  // ---> DNA processes activated
439 
440  mod = new G4DNARuddIonisationModel();
441  em_config->SetExtraEmModel("alpha","alpha_G4DNAIonisation",
442  mod,"Target",0.0,10*MeV);
443 
445  em_config->SetExtraEmModel("alpha","alpha_G4DNAExcitation",
446  mod,"Target",1*keV,10*MeV);
447 
448 }
449 
450 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
451 
453 { }
454 
455 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
456 
458 {
459  if (verboseLevel >0)
460  {
461  G4cout << "PhysicsList::SetCuts:";
462  G4cout << "CutLength : " << G4BestUnit(defaultCutValue,"Length") << G4endl;
463  }
464 
465  // set cut values for gamma at first and for e- second and next for e+,
466  // because some processes for e+/e- need cut values for gamma
467  SetCutValue(fCutForGamma, "gamma");
470  SetCutValue(fCutForProton, "proton");
471 
472  if (verboseLevel>0) { DumpCutValuesTable(); }
473 }
G4EmConfigurator * EmConfigurator()
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:88
void ConstructParticle()
Definition: PhysicsList.cc:117
static const double MeV
Definition: G4SIunits.hh:211
#define G4DNABornIonisationModel
static G4LossTableManager * Instance()
void SetCutValue(G4double aCut, const G4String &pname)
G4double fCutForPositron
Definition: PhysicsList.hh:70
void ConstructBosons()
Definition: PhysicsList.cc:64
void ConstructLeptons()
Definition: PhysicsList.cc:71
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4double fCutForElectron
Definition: PhysicsList.hh:69
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4ProcessManager * GetProcessManager() const
G4EmConfigurator em_config
Definition: PhysicsList.hh:65
const G4String & GetParticleName() const
G4DNABornExcitationModel1 G4DNABornExcitationModel
G4double fCutForProton
Definition: PhysicsList.hh:68
void ConstructBarions()
Definition: PhysicsList.cc:80
void SetEmModel(G4VEmModel *, G4int index=1)
void DumpCutValuesTable(G4int flag=1)
G4GLOB_DLL std::ostream G4cout
G4double fCutForGamma
Definition: PhysicsList.hh:64
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
static G4DNAGenericIonsManager * Instance(void)
void SetVerboseLevel(G4int value)
void SetCuts()
Definition: PhysicsList.cc:219
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
static G4Positron * PositronDefinition()
Definition: G4Positron.cc:89
static const double micrometer
Definition: G4SIunits.hh:99
static const double eV
Definition: G4SIunits.hh:212
void SetEmModel(G4VEmModel *, G4int index=1)
#define G4endl
Definition: G4ios.hh:61
static const double TeV
Definition: G4SIunits.hh:215
static const double keV
Definition: G4SIunits.hh:213
double G4double
Definition: G4Types.hh:76
void ConstructProcess()
Definition: PhysicsList.cc:170
void ConstructGeneral()
Definition: PhysicsList.cc:287
void ConstructEM()
Definition: PhysicsList.cc:120
void SetExtraEmModel(const G4String &particleName, const G4String &processName, G4VEmModel *, const G4String &regionName="", G4double emin=0.0, G4double emax=DBL_MAX, G4VEmFluctuationModel *fm=0)
#define theParticleIterator
G4ParticleDefinition * GetIon(const G4String &name)
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81