Geant4-11
G4EmDNAPhysics_stationary_option4.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
27
28#include "G4SystemOfUnits.hh"
29
31
32// *** Processes and models for Geant4-DNA
33
34#include "G4DNAElastic.hh"
37
38#include "G4DNAExcitation.hh"
39#include "G4DNAAttachment.hh"
40#include "G4DNAVibExcitation.hh"
41#include "G4DNAIonisation.hh"
44
47
48// particles
49
50#include "G4Electron.hh"
51#include "G4Proton.hh"
52#include "G4GenericIon.hh"
53
54// Warning : the following is needed in order to use EM Physics builders
55// e+
56#include "G4Positron.hh"
58#include "G4eIonisation.hh"
59#include "G4eBremsstrahlung.hh"
61// gamma
62#include "G4Gamma.hh"
67#include "G4GammaConversion.hh"
71
72#include "G4EmParameters.hh"
73// end of warning
74
75#include "G4LossTableManager.hh"
78#include "G4BuilderType.hh"
79
80// factory
82//
84
85
86//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87
89 : G4VPhysicsConstructor("G4EmDNAPhysics_stationary_option4"), verbose(ver)
90{
92 param->SetDefaults();
93 param->SetFluo(true);
94 param->SetAuger(true);
95 param->SetDeexcitationIgnoreCut(true);
96 param->ActivateDNA();
97
99}
100
101//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
102
104const G4String&)
106{}
107
108//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
109
111{}
112
113//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
114
116{
117// bosons
119
120// leptons
123
124// baryons
126
128
129 G4DNAGenericIonsManager * genericIonsManager;
130 genericIonsManager=G4DNAGenericIonsManager::Instance();
131 genericIonsManager->GetIon("alpha++");
132 genericIonsManager->GetIon("alpha+");
133 genericIonsManager->GetIon("helium");
134 genericIonsManager->GetIon("hydrogen");
135
136}
137
138//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
139
141{
142 if(verbose > 1) {
143 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
144 }
146
147 auto myParticleIterator=GetParticleIterator();
148 myParticleIterator->reset();
149 while( (*myParticleIterator)() )
150 {
151 G4ParticleDefinition* particle = myParticleIterator->value();
152 G4String particleName = particle->GetParticleName();
153
154 if (particleName == "e-") {
155
156 // *** Elastic scattering (two alternative models available) ***
157
158 G4DNAElastic* theDNAElasticProcess =
159 new G4DNAElastic("e-_G4DNAElastic");
160 theDNAElasticProcess->SetEmModel(
162 ph->RegisterProcess(theDNAElasticProcess, particle);
163
164 // *** Excitation ***
165
166 G4DNAExcitation* theDNAExcitationProcess =
167 new G4DNAExcitation("e-_G4DNAExcitation");
168 theDNAExcitationProcess->SetEmModel(
171 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
172 ph->RegisterProcess(theDNAExcitationProcess, particle);
173
174 // *** Ionisation ***
175
176 G4DNAIonisation* theDNAIonisationProcess =
177 new G4DNAIonisation("e-_G4DNAIonisation");
178 theDNAIonisationProcess->SetEmModel(
181 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
182 ph->RegisterProcess(theDNAIonisationProcess, particle);
183
184 // *** Vibrational excitation ***
185
186 G4DNAVibExcitation* theDNAVibExcitationProcess =
187 new G4DNAVibExcitation("e-_G4DNAVibExcitation");
188 theDNAVibExcitationProcess->SetEmModel(
191 (theDNAVibExcitationProcess->EmModel()))->SelectStationary(true);
192 ph->RegisterProcess(theDNAVibExcitationProcess, particle);
193
194 // *** Attachment ***
195
196 G4DNAAttachment* theDNAAttachmentProcess =
197 new G4DNAAttachment("e-_G4DNAAttachment");
198 theDNAAttachmentProcess->SetEmModel(
201 (theDNAAttachmentProcess->EmModel()))->SelectStationary(true);
202 ph->RegisterProcess(theDNAAttachmentProcess, particle);
203
204 } else if ( particleName == "proton" ) {
205
206 // *** Elastic ***
207
208 G4DNAElastic* theDNAElasticProcess =
209 new G4DNAElastic("proton_G4DNAElastic");
210 theDNAElasticProcess->SetEmModel(
213 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
214 ph->RegisterProcess(theDNAElasticProcess, particle);
215
216 // *** Excitation ***
217
218 G4DNAExcitation* theDNAExcitationProcess =
219 new G4DNAExcitation("proton_G4DNAExcitation");
220
221 theDNAExcitationProcess->SetEmModel(
223 theDNAExcitationProcess->SetEmModel(
225
227 (theDNAExcitationProcess->EmModel()))->SetLowEnergyLimit(10*eV);
229 (theDNAExcitationProcess->EmModel()))->SetHighEnergyLimit(500*keV);
231 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
232
234 (theDNAExcitationProcess->EmModel(1)))->SetLowEnergyLimit(500*keV);
236 (theDNAExcitationProcess->EmModel(1)))->SetHighEnergyLimit(100*MeV);
238 (theDNAExcitationProcess->EmModel(1)))->SelectStationary(true);
239
240 ph->RegisterProcess(theDNAExcitationProcess, particle);
241
242 // *** Ionisation ***
243
244 G4DNAIonisation* theDNAIonisationProcess =
245 new G4DNAIonisation("proton_G4DNAIonisation");
246
247 theDNAIonisationProcess->SetEmModel(
249 theDNAIonisationProcess->SetEmModel(
251
253 (theDNAIonisationProcess->EmModel()))->SetLowEnergyLimit(0*eV);
255 (theDNAIonisationProcess->EmModel()))->SetHighEnergyLimit(500*keV);
257 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
258
260 (theDNAIonisationProcess->EmModel(1)))->SetLowEnergyLimit(500*keV);
262 (theDNAIonisationProcess->EmModel(1)))->SetHighEnergyLimit(100*MeV);
264 (theDNAIonisationProcess->EmModel(1)))->SelectStationary(true);
265 //
267 (theDNAIonisationProcess->EmModel(1)))->SelectFasterComputation(true);
268 //
269
270 ph->RegisterProcess(theDNAIonisationProcess, particle);
271
272 // *** Charge decrease ***
273
274 G4DNAChargeDecrease* theDNAChargeDecreaseProcess =
275 new G4DNAChargeDecrease("proton_G4DNAChargeDecrease");
276 theDNAChargeDecreaseProcess->SetEmModel(
279 (theDNAChargeDecreaseProcess->EmModel()))->SelectStationary(true);
280 ph->RegisterProcess(theDNAChargeDecreaseProcess, particle);
281
282 } else if ( particleName == "hydrogen" ) {
283
284 // *** Elastic ***
285
286 G4DNAElastic* theDNAElasticProcess =
287 new G4DNAElastic("hydrogen_G4DNAElastic");
288 theDNAElasticProcess->SetEmModel(
291 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
292 ph->RegisterProcess(theDNAElasticProcess, particle);
293
294 // *** Excitation ***
295
296 G4DNAExcitation* theDNAExcitationProcess =
297 new G4DNAExcitation("hydrogen_G4DNAExcitation");
298 theDNAExcitationProcess->SetEmModel(
301 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
302 ph->RegisterProcess(theDNAExcitationProcess, particle);
303
304 // *** Ionisation ***
305
306 G4DNAIonisation* theDNAIonisationProcess =
307 new G4DNAIonisation("hydrogen_G4DNAIonisation");
308 theDNAIonisationProcess->SetEmModel(
311 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
312 ph->RegisterProcess(theDNAIonisationProcess, particle);
313
314 // *** Charge increase ***
315
316 G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
317 new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease");
318 theDNAChargeIncreaseProcess->SetEmModel(
321 (theDNAChargeIncreaseProcess->EmModel()))->SelectStationary(true);
322 ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
323
324 } else if ( particleName == "alpha" ) {
325
326 // *** Elastic ***
327
328 G4DNAElastic* theDNAElasticProcess =
329 new G4DNAElastic("alpha_G4DNAElastic");
330 theDNAElasticProcess->SetEmModel(
333 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
334 ph->RegisterProcess(theDNAElasticProcess, particle);
335
336 // *** Excitation ***
337
338 G4DNAExcitation* theDNAExcitationProcess =
339 new G4DNAExcitation("alpha_G4DNAExcitation");
340 theDNAExcitationProcess->SetEmModel(
343 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
344 ph->RegisterProcess(theDNAExcitationProcess, particle);
345
346 // *** Ionisation ***
347
348 G4DNAIonisation* theDNAIonisationProcess =
349 new G4DNAIonisation("alpha_G4DNAIonisation");
350 theDNAIonisationProcess->SetEmModel(
353 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
354 ph->RegisterProcess(theDNAIonisationProcess, particle);
355
356 // *** Charge decrease ***
357
358 G4DNAChargeDecrease* theDNAChargeDecreaseProcess =
359 new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease");
360 theDNAChargeDecreaseProcess->SetEmModel(
363 (theDNAChargeDecreaseProcess->EmModel()))->SelectStationary(true);
364 ph->RegisterProcess(theDNAChargeDecreaseProcess, particle);
365
366 } else if ( particleName == "alpha+" ) {
367
368 // *** Elastic ***
369
370 G4DNAElastic* theDNAElasticProcess =
371 new G4DNAElastic("alpha+_G4DNAElastic");
372 theDNAElasticProcess->SetEmModel(
375 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
376 ph->RegisterProcess(theDNAElasticProcess, particle);
377
378 // *** Excitation ***
379
380 G4DNAExcitation* theDNAExcitationProcess =
381 new G4DNAExcitation("alpha+_G4DNAExcitation");
382 theDNAExcitationProcess->SetEmModel(
385 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
386 ph->RegisterProcess(theDNAExcitationProcess, particle);
387
388 // *** Ionisation ***
389
390 G4DNAIonisation* theDNAIonisationProcess =
391 new G4DNAIonisation("alpha+_G4DNAIonisation");
392 theDNAIonisationProcess->SetEmModel(
395 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
396 ph->RegisterProcess(theDNAIonisationProcess, particle);
397
398 // *** Charge decrease ***
399
400 G4DNAChargeDecrease* theDNAChargeDecreaseProcess =
401 new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease");
402 theDNAChargeDecreaseProcess->SetEmModel(
405 (theDNAChargeDecreaseProcess->EmModel()))->SelectStationary(true);
406 ph->RegisterProcess(theDNAChargeDecreaseProcess, particle);
407
408 // *** Charge increase ***
409
410 G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
411 new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease");
412 theDNAChargeIncreaseProcess->SetEmModel(
415 (theDNAChargeIncreaseProcess->EmModel()))->SelectStationary(true);
416 ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
417
418 } else if ( particleName == "helium" ) {
419
420 // *** Elastic ***
421
422 G4DNAElastic* theDNAElasticProcess =
423 new G4DNAElastic("helium_G4DNAElastic");
424 theDNAElasticProcess->SetEmModel(
427 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
428 ph->RegisterProcess(theDNAElasticProcess, particle);
429
430 // *** Excitation ***
431
432 G4DNAExcitation* theDNAExcitationProcess =
433 new G4DNAExcitation("helium_G4DNAExcitation");
434 theDNAExcitationProcess->SetEmModel(
437 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
438 ph->RegisterProcess(theDNAExcitationProcess, particle);
439
440 // *** Ionisation ***
441
442 G4DNAIonisation* theDNAIonisationProcess =
443 new G4DNAIonisation("helium_G4DNAIonisation");
444 theDNAIonisationProcess->SetEmModel(
447 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
448 ph->RegisterProcess(theDNAIonisationProcess, particle);
449
450 // *** Charge increase ***
451
452 G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
453 new G4DNAChargeIncrease("helium_G4DNAChargeIncrease");
454 theDNAChargeIncreaseProcess->SetEmModel(
457 (theDNAChargeIncreaseProcess->EmModel()))->SelectStationary(true);
458 ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
459
460 } else if ( particleName == "GenericIon" ) {
461
462 // *** Ionisation ***
463
464 G4DNAIonisation* theDNAIonisationProcess =
465 new G4DNAIonisation("GenericIon_G4DNAIonisation");
466 theDNAIonisationProcess->SetEmModel(
469 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
470 ph->RegisterProcess(theDNAIonisationProcess, particle);
471
472 }
473
474 // Warning : the following particles and processes are needed by EM Physics
475 // builders
476 // They are taken from the default Livermore Physics list
477 // These particles are currently not handled by Geant4-DNA
478
479 // e+
480
481 else if (particleName == "e+") {
482
483 // Identical to G4EmStandardPhysics_stationary
484
487 G4eIonisation* eIoni = new G4eIonisation();
488 eIoni->SetStepFunction(0.2, 100*um);
489
490 ph->RegisterProcess(msc, particle);
491 ph->RegisterProcess(eIoni, particle);
492 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
493 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
494
495 } else if (particleName == "gamma") {
496
497 // photoelectric effect - Livermore model only
498 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
499 thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
500 ph->RegisterProcess(thePhotoElectricEffect, particle);
501
502 // Compton scattering - Livermore model only
503 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
504 theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
505 ph->RegisterProcess(theComptonScattering, particle);
506
507 // gamma conversion - Livermore model below 80 GeV
508 G4GammaConversion* theGammaConversion = new G4GammaConversion();
509 theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
510 ph->RegisterProcess(theGammaConversion, particle);
511
512 // default Rayleigh scattering is Livermore
513 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
514 ph->RegisterProcess(theRayleigh, particle);
515 }
516
517 // Warning : end of particles and processes are needed by EM Physics build.
518
519 }
520
521 // Deexcitation
522 //
525}
526
@ bElectromagnetic
G4DNABornExcitationModel1 G4DNABornExcitationModel
#define G4DNABornIonisationModel
G4_DECLARE_PHYSCONSTR_FACTORY(G4EmDNAPhysics_stationary_option4)
@ fUseDistanceToBoundary
static constexpr double keV
Definition: G4SIunits.hh:202
static constexpr double eV
Definition: G4SIunits.hh:201
static constexpr double um
Definition: G4SIunits.hh:93
static constexpr double MeV
Definition: G4SIunits.hh:200
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4DNAGenericIonsManager * Instance(void)
G4ParticleDefinition * GetIon(const G4String &name)
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4EmParameters * Instance()
void SetDeexcitationIgnoreCut(G4bool val)
void SetFluo(G4bool val)
void SetAuger(G4bool val)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:87
void SetAtomDeexcitation(G4VAtomDeexcitation *)
static G4LossTableManager * Instance()
const G4String & GetParticleName() const
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4Positron * Positron()
Definition: G4Positron.cc:93
static G4Proton * Proton()
Definition: G4Proton.cc:92
G4VEmModel * EmModel(size_t index=0) const
void SetEmModel(G4VEmModel *, G4int index=0)
void SetStepFunction(G4double v1, G4double v2)
void SetStepLimitType(G4MscStepLimitType val)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
const G4String & GetPhysicsName() const