Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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->SetAugerCascade(true);
96 param->SetDeexcitationIgnoreCut(true);
97 param->ActivateDNA();
98
100}
101
102//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
103
105const G4String&)
106 : G4VPhysicsConstructor("G4EmDNAPhysics_stationary_option4"), verbose(ver)
107{
109 param->SetDefaults();
110 param->SetFluo(true);
111 param->SetAuger(true);
112 param->SetAugerCascade(true);
113 param->SetDeexcitationIgnoreCut(true);
114
116}
117
118//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
119
121{}
122
123//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
124
126{
127// bosons
129
130// leptons
133
134// baryons
136
138
139 G4DNAGenericIonsManager * genericIonsManager;
140 genericIonsManager=G4DNAGenericIonsManager::Instance();
141 genericIonsManager->GetIon("alpha++");
142 genericIonsManager->GetIon("alpha+");
143 genericIonsManager->GetIon("helium");
144 genericIonsManager->GetIon("hydrogen");
145
146}
147
148//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
149
151{
152 if(verbose > 1) {
153 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
154 }
156
157 auto myParticleIterator=GetParticleIterator();
158 myParticleIterator->reset();
159 while( (*myParticleIterator)() )
160 {
161 G4ParticleDefinition* particle = myParticleIterator->value();
162 G4String particleName = particle->GetParticleName();
163
164 if (particleName == "e-") {
165
166 // *** Elastic scattering (two alternative models available) ***
167
168 G4DNAElastic* theDNAElasticProcess =
169 new G4DNAElastic("e-_G4DNAElastic");
170 theDNAElasticProcess->SetEmModel(
172 ph->RegisterProcess(theDNAElasticProcess, particle);
173
174 // *** Excitation ***
175
176 G4DNAExcitation* theDNAExcitationProcess =
177 new G4DNAExcitation("e-_G4DNAExcitation");
178 theDNAExcitationProcess->SetEmModel(
181 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
182 ph->RegisterProcess(theDNAExcitationProcess, particle);
183
184 // *** Ionisation ***
185
186 G4DNAIonisation* theDNAIonisationProcess =
187 new G4DNAIonisation("e-_G4DNAIonisation");
188 theDNAIonisationProcess->SetEmModel(
191 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
192 ph->RegisterProcess(theDNAIonisationProcess, particle);
193
194 // *** Vibrational excitation ***
195
196 G4DNAVibExcitation* theDNAVibExcitationProcess =
197 new G4DNAVibExcitation("e-_G4DNAVibExcitation");
198 theDNAVibExcitationProcess->SetEmModel(
201 (theDNAVibExcitationProcess->EmModel()))->SelectStationary(true);
202 ph->RegisterProcess(theDNAVibExcitationProcess, particle);
203
204 // *** Attachment ***
205
206 G4DNAAttachment* theDNAAttachmentProcess =
207 new G4DNAAttachment("e-_G4DNAAttachment");
208 theDNAAttachmentProcess->SetEmModel(
211 (theDNAAttachmentProcess->EmModel()))->SelectStationary(true);
212 ph->RegisterProcess(theDNAAttachmentProcess, particle);
213
214 } else if ( particleName == "proton" ) {
215
216 // *** Elastic ***
217
218 G4DNAElastic* theDNAElasticProcess =
219 new G4DNAElastic("proton_G4DNAElastic");
220 theDNAElasticProcess->SetEmModel(
223 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
224 ph->RegisterProcess(theDNAElasticProcess, particle);
225
226 // *** Excitation ***
227
228 G4DNAExcitation* theDNAExcitationProcess =
229 new G4DNAExcitation("proton_G4DNAExcitation");
230
231 theDNAExcitationProcess->SetEmModel(
233 theDNAExcitationProcess->SetEmModel(
235
237 (theDNAExcitationProcess->EmModel()))->SetLowEnergyLimit(10*eV);
239 (theDNAExcitationProcess->EmModel()))->SetHighEnergyLimit(500*keV);
241 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
242
244 (theDNAExcitationProcess->EmModel(1)))->SetLowEnergyLimit(500*keV);
246 (theDNAExcitationProcess->EmModel(1)))->SetHighEnergyLimit(100*MeV);
248 (theDNAExcitationProcess->EmModel(1)))->SelectStationary(true);
249
250 ph->RegisterProcess(theDNAExcitationProcess, particle);
251
252 // *** Ionisation ***
253
254 G4DNAIonisation* theDNAIonisationProcess =
255 new G4DNAIonisation("proton_G4DNAIonisation");
256
257 theDNAIonisationProcess->SetEmModel(
259 theDNAIonisationProcess->SetEmModel(
261
263 (theDNAIonisationProcess->EmModel()))->SetLowEnergyLimit(0*eV);
265 (theDNAIonisationProcess->EmModel()))->SetHighEnergyLimit(500*keV);
267 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
268
270 (theDNAIonisationProcess->EmModel(1)))->SetLowEnergyLimit(500*keV);
272 (theDNAIonisationProcess->EmModel(1)))->SetHighEnergyLimit(100*MeV);
274 (theDNAIonisationProcess->EmModel(1)))->SelectStationary(true);
275 //
277 (theDNAIonisationProcess->EmModel(1)))->SelectFasterComputation(true);
278 //
279
280 ph->RegisterProcess(theDNAIonisationProcess, particle);
281
282 // *** Charge decrease ***
283
284 G4DNAChargeDecrease* theDNAChargeDecreaseProcess =
285 new G4DNAChargeDecrease("proton_G4DNAChargeDecrease");
286 theDNAChargeDecreaseProcess->SetEmModel(
289 (theDNAChargeDecreaseProcess->EmModel()))->SelectStationary(true);
290 ph->RegisterProcess(theDNAChargeDecreaseProcess, particle);
291
292 } else if ( particleName == "hydrogen" ) {
293
294 // *** Elastic ***
295
296 G4DNAElastic* theDNAElasticProcess =
297 new G4DNAElastic("hydrogen_G4DNAElastic");
298 theDNAElasticProcess->SetEmModel(
301 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
302 ph->RegisterProcess(theDNAElasticProcess, particle);
303
304 // *** Excitation ***
305
306 G4DNAExcitation* theDNAExcitationProcess =
307 new G4DNAExcitation("hydrogen_G4DNAExcitation");
308 theDNAExcitationProcess->SetEmModel(
311 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
312 ph->RegisterProcess(theDNAExcitationProcess, particle);
313
314 // *** Ionisation ***
315
316 G4DNAIonisation* theDNAIonisationProcess =
317 new G4DNAIonisation("hydrogen_G4DNAIonisation");
318 theDNAIonisationProcess->SetEmModel(
321 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
322 ph->RegisterProcess(theDNAIonisationProcess, particle);
323
324 // *** Charge increase ***
325
326 G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
327 new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease");
328 theDNAChargeIncreaseProcess->SetEmModel(
331 (theDNAChargeIncreaseProcess->EmModel()))->SelectStationary(true);
332 ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
333
334 } else if ( particleName == "alpha" ) {
335
336 // *** Elastic ***
337
338 G4DNAElastic* theDNAElasticProcess =
339 new G4DNAElastic("alpha_G4DNAElastic");
340 theDNAElasticProcess->SetEmModel(
343 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
344 ph->RegisterProcess(theDNAElasticProcess, particle);
345
346 // *** Excitation ***
347
348 G4DNAExcitation* theDNAExcitationProcess =
349 new G4DNAExcitation("alpha_G4DNAExcitation");
350 theDNAExcitationProcess->SetEmModel(
353 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
354 ph->RegisterProcess(theDNAExcitationProcess, particle);
355
356 // *** Ionisation ***
357
358 G4DNAIonisation* theDNAIonisationProcess =
359 new G4DNAIonisation("alpha_G4DNAIonisation");
360 theDNAIonisationProcess->SetEmModel(
363 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
364 ph->RegisterProcess(theDNAIonisationProcess, particle);
365
366 // *** Charge decrease ***
367
368 G4DNAChargeDecrease* theDNAChargeDecreaseProcess =
369 new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease");
370 theDNAChargeDecreaseProcess->SetEmModel(
373 (theDNAChargeDecreaseProcess->EmModel()))->SelectStationary(true);
374 ph->RegisterProcess(theDNAChargeDecreaseProcess, particle);
375
376 } else if ( particleName == "alpha+" ) {
377
378 // *** Elastic ***
379
380 G4DNAElastic* theDNAElasticProcess =
381 new G4DNAElastic("alpha+_G4DNAElastic");
382 theDNAElasticProcess->SetEmModel(
385 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
386 ph->RegisterProcess(theDNAElasticProcess, particle);
387
388 // *** Excitation ***
389
390 G4DNAExcitation* theDNAExcitationProcess =
391 new G4DNAExcitation("alpha+_G4DNAExcitation");
392 theDNAExcitationProcess->SetEmModel(
395 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
396 ph->RegisterProcess(theDNAExcitationProcess, particle);
397
398 // *** Ionisation ***
399
400 G4DNAIonisation* theDNAIonisationProcess =
401 new G4DNAIonisation("alpha+_G4DNAIonisation");
402 theDNAIonisationProcess->SetEmModel(
405 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
406 ph->RegisterProcess(theDNAIonisationProcess, particle);
407
408 // *** Charge decrease ***
409
410 G4DNAChargeDecrease* theDNAChargeDecreaseProcess =
411 new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease");
412 theDNAChargeDecreaseProcess->SetEmModel(
415 (theDNAChargeDecreaseProcess->EmModel()))->SelectStationary(true);
416 ph->RegisterProcess(theDNAChargeDecreaseProcess, particle);
417
418 // *** Charge increase ***
419
420 G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
421 new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease");
422 theDNAChargeIncreaseProcess->SetEmModel(
425 (theDNAChargeIncreaseProcess->EmModel()))->SelectStationary(true);
426 ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
427
428 } else if ( particleName == "helium" ) {
429
430 // *** Elastic ***
431
432 G4DNAElastic* theDNAElasticProcess =
433 new G4DNAElastic("helium_G4DNAElastic");
434 theDNAElasticProcess->SetEmModel(
437 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
438 ph->RegisterProcess(theDNAElasticProcess, particle);
439
440 // *** Excitation ***
441
442 G4DNAExcitation* theDNAExcitationProcess =
443 new G4DNAExcitation("helium_G4DNAExcitation");
444 theDNAExcitationProcess->SetEmModel(
447 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
448 ph->RegisterProcess(theDNAExcitationProcess, particle);
449
450 // *** Ionisation ***
451
452 G4DNAIonisation* theDNAIonisationProcess =
453 new G4DNAIonisation("helium_G4DNAIonisation");
454 theDNAIonisationProcess->SetEmModel(
457 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
458 ph->RegisterProcess(theDNAIonisationProcess, particle);
459
460 // *** Charge increase ***
461
462 G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
463 new G4DNAChargeIncrease("helium_G4DNAChargeIncrease");
464 theDNAChargeIncreaseProcess->SetEmModel(
467 (theDNAChargeIncreaseProcess->EmModel()))->SelectStationary(true);
468 ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
469
470 } else if ( particleName == "GenericIon" ) {
471
472 // *** Ionisation ***
473
474 G4DNAIonisation* theDNAIonisationProcess =
475 new G4DNAIonisation("GenericIon_G4DNAIonisation");
476 theDNAIonisationProcess->SetEmModel(
479 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
480 ph->RegisterProcess(theDNAIonisationProcess, particle);
481
482 }
483
484 // Warning : the following particles and processes are needed by EM Physics
485 // builders
486 // They are taken from the default Livermore Physics list
487 // These particles are currently not handled by Geant4-DNA
488
489 // e+
490
491 else if (particleName == "e+") {
492
493 // Identical to G4EmStandardPhysics_stationary
494
497 G4eIonisation* eIoni = new G4eIonisation();
498 eIoni->SetStepFunction(0.2, 100*um);
499
500 ph->RegisterProcess(msc, particle);
501 ph->RegisterProcess(eIoni, particle);
502 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
503 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
504
505 } else if (particleName == "gamma") {
506
507 // photoelectric effect - Livermore model only
508 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
509 thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
510 ph->RegisterProcess(thePhotoElectricEffect, particle);
511
512 // Compton scattering - Livermore model only
513 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
514 theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
515 ph->RegisterProcess(theComptonScattering, particle);
516
517 // gamma conversion - Livermore model below 80 GeV
518 G4GammaConversion* theGammaConversion = new G4GammaConversion();
519 theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
520 ph->RegisterProcess(theGammaConversion, particle);
521
522 // default Rayleigh scattering is Livermore
523 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
524 ph->RegisterProcess(theRayleigh, particle);
525 }
526
527 // Warning : end of particles and processes are needed by EM Physics build.
528
529 }
530
531 // Deexcitation
532 //
535}
536
@ bElectromagnetic
G4DNABornExcitationModel1 G4DNABornExcitationModel
#define G4DNABornIonisationModel
@ fUseDistanceToBoundary
#define G4_DECLARE_PHYSCONSTR_FACTORY(physics_constructor)
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 SetAugerCascade(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