Geant4-11
G4EmDNAPhysics_stationary_option6.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
48
49// particles
50
51#include "G4Electron.hh"
52#include "G4Proton.hh"
53#include "G4GenericIon.hh"
54
55// Warning : the following is needed in order to use EM Physics builders
56// e+
57#include "G4Positron.hh"
59#include "G4eIonisation.hh"
60#include "G4eBremsstrahlung.hh"
62// gamma
63#include "G4Gamma.hh"
68#include "G4GammaConversion.hh"
72
73#include "G4EmParameters.hh"
74// end of warning
75
76#include "G4LossTableManager.hh"
79#include "G4BuilderType.hh"
80
81// factory
83//
85
86
87//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
88
90 : G4VPhysicsConstructor("G4EmDNAPhysics_stationary_option6"), verbose(ver)
91{
93 param->SetDefaults();
94 param->SetFluo(true);
95 param->SetAuger(true);
96 param->SetDeexcitationIgnoreCut(true);
97 param->ActivateDNA();
98
100}
101
102//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
103
105const G4String&)
107{}
108
109//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
110
112{}
113
114//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
115
117{
118// bosons
120
121// leptons
124
125// baryons
127
129
130 G4DNAGenericIonsManager * genericIonsManager;
131 genericIonsManager=G4DNAGenericIonsManager::Instance();
132 genericIonsManager->GetIon("alpha++");
133 genericIonsManager->GetIon("alpha+");
134 genericIonsManager->GetIon("helium");
135 genericIonsManager->GetIon("hydrogen");
136
137}
138
139//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
140
142{
143 if(verbose > 1) {
144 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
145 }
147
148 auto myParticleIterator=GetParticleIterator();
149 myParticleIterator->reset();
150 while( (*myParticleIterator)() )
151 {
152 G4ParticleDefinition* particle = myParticleIterator->value();
153 G4String particleName = particle->GetParticleName();
154
155 if (particleName == "e-") {
156
157 // *** Elastic scattering (two alternative models available) ***
158
159 G4DNAElastic* theDNAElasticProcess =
160 new G4DNAElastic("e-_G4DNAElastic");
161 theDNAElasticProcess->SetEmModel(new G4DNACPA100ElasticModel());
163 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
164 ph->RegisterProcess(theDNAElasticProcess, particle);
165
166 // *** Excitation ***
167
168 G4DNAExcitation* theDNAExcitationProcess =
169 new G4DNAExcitation("e-_G4DNAExcitation");
170 theDNAExcitationProcess->SetEmModel(new G4DNACPA100ExcitationModel());
172 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
173 ph->RegisterProcess(theDNAExcitationProcess, particle);
174
175 // *** Ionisation ***
176
177 G4DNAIonisation* theDNAIonisationProcess =
178 new G4DNAIonisation("e-_G4DNAIonisation");
179 theDNAIonisationProcess->SetEmModel(new G4DNACPA100IonisationModel());
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(new G4DNASancheExcitationModel());
190 (theDNAVibExcitationProcess->EmModel()))->SelectStationary(true);
191 ph->RegisterProcess(theDNAVibExcitationProcess, particle);
192
193 // *** Attachment ***
194
195 G4DNAAttachment* theDNAAttachmentProcess =
196 new G4DNAAttachment("e-_G4DNAAttachment");
197 theDNAAttachmentProcess->SetEmModel(new G4DNAMeltonAttachmentModel());
199 (theDNAAttachmentProcess->EmModel()))->SelectStationary(true);
200 ph->RegisterProcess(theDNAAttachmentProcess, particle);
201
202 } else if ( particleName == "proton" ) {
203
204 // *** Elastic ***
205
206 G4DNAElastic* theDNAElasticProcess =
207 new G4DNAElastic("proton_G4DNAElastic");
208 theDNAElasticProcess->SetEmModel(new G4DNAIonElasticModel());
210 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
211 ph->RegisterProcess(theDNAElasticProcess, particle);
212
213 // *** Excitation ***
214
215 G4DNAExcitation* theDNAExcitationProcess =
216 new G4DNAExcitation("proton_G4DNAExcitation");
217
218 theDNAExcitationProcess->SetEmModel(
220 theDNAExcitationProcess->SetEmModel(
222
224 (theDNAExcitationProcess->EmModel()))->SetLowEnergyLimit(10*eV);
226 (theDNAExcitationProcess->EmModel()))->SetHighEnergyLimit(500*keV);
228 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
229
231 (theDNAExcitationProcess->EmModel(1)))->SetLowEnergyLimit(500*keV);
233 (theDNAExcitationProcess->EmModel(1)))->SetHighEnergyLimit(100*MeV);
235 (theDNAExcitationProcess->EmModel(1)))->SelectStationary(true);
236
237 ph->RegisterProcess(theDNAExcitationProcess, particle);
238
239 // *** Ionisation ***
240
241 G4DNAIonisation* theDNAIonisationProcess =
242 new G4DNAIonisation("proton_G4DNAIonisation");
243
244 theDNAIonisationProcess->SetEmModel(
246 theDNAIonisationProcess->SetEmModel(
248
250 (theDNAIonisationProcess->EmModel()))->SetLowEnergyLimit(0*eV);
252 (theDNAIonisationProcess->EmModel()))->SetHighEnergyLimit(500*keV);
254 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
255
257 (theDNAIonisationProcess->EmModel(1)))->SetLowEnergyLimit(500*keV);
259 (theDNAIonisationProcess->EmModel(1)))->SetHighEnergyLimit(100*MeV);
261 (theDNAIonisationProcess->EmModel(1)))->SelectStationary(true);
262 //
264 (theDNAIonisationProcess->EmModel(1)))->SelectFasterComputation(true);
265 //
266
267 ph->RegisterProcess(theDNAIonisationProcess, particle);
268
269 // *** Charge decrease ***
270
271 G4DNAChargeDecrease* theDNAChargeDecreaseProcess =
272 new G4DNAChargeDecrease("proton_G4DNAChargeDecrease");
273 theDNAChargeDecreaseProcess->SetEmModel(
276 (theDNAChargeDecreaseProcess->EmModel()))->SelectStationary(true);
277 ph->RegisterProcess(theDNAChargeDecreaseProcess, particle);
278
279 } else if ( particleName == "hydrogen" ) {
280
281 // *** Elastic ***
282
283 G4DNAElastic* theDNAElasticProcess =
284 new G4DNAElastic("hydrogen_G4DNAElastic");
285 theDNAElasticProcess->SetEmModel(
288 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
289 ph->RegisterProcess(theDNAElasticProcess, particle);
290
291 // *** Excitation ***
292
293 G4DNAExcitation* theDNAExcitationProcess =
294 new G4DNAExcitation("hydrogen_G4DNAExcitation");
295 theDNAExcitationProcess->SetEmModel(
298 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
299 ph->RegisterProcess(theDNAExcitationProcess, particle);
300
301 // *** Ionisation ***
302
303 G4DNAIonisation* theDNAIonisationProcess =
304 new G4DNAIonisation("hydrogen_G4DNAIonisation");
305 theDNAIonisationProcess->SetEmModel(
308 theDNAIonisationProcess->EmModel()))->SelectStationary(true);
309 ph->RegisterProcess(theDNAIonisationProcess, particle);
310
311 // *** Charge increase ***
312
313 G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
314 new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease");
315 theDNAChargeIncreaseProcess->SetEmModel(
318 theDNAChargeIncreaseProcess->EmModel()))->SelectStationary(true);
319 ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
320
321 } else if ( particleName == "alpha" ) {
322
323 // *** Elastic ***
324
325 G4DNAElastic* theDNAElasticProcess =
326 new G4DNAElastic("alpha_G4DNAElastic");
327 theDNAElasticProcess->SetEmModel(new G4DNAIonElasticModel());
329 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
330 ph->RegisterProcess(theDNAElasticProcess, particle);
331
332 // *** Excitation ***
333
334 G4DNAExcitation* theDNAExcitationProcess =
335 new G4DNAExcitation("alpha_G4DNAExcitation");
336 theDNAExcitationProcess->SetEmModel(
339 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
340 ph->RegisterProcess(theDNAExcitationProcess, particle);
341
342 // *** Ionisation ***
343
344 G4DNAIonisation* theDNAIonisationProcess =
345 new G4DNAIonisation("alpha_G4DNAIonisation");
346 theDNAIonisationProcess->SetEmModel(
349 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
350 ph->RegisterProcess(theDNAIonisationProcess, particle);
351
352 // *** Charge decrease ***
353
354 G4DNAChargeDecrease* theDNAChargeDecreaseProcess =
355 new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease");
356 theDNAChargeDecreaseProcess->SetEmModel(
359 (theDNAChargeDecreaseProcess->EmModel()))->SelectStationary(true);
360 ph->RegisterProcess(theDNAChargeDecreaseProcess, particle);
361
362 } else if ( particleName == "alpha+" ) {
363
364 // *** Elastic ***
365
366 G4DNAElastic* theDNAElasticProcess =
367 new G4DNAElastic("alpha+_G4DNAElastic");
368 theDNAElasticProcess->SetEmModel(
371 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
372 ph->RegisterProcess(theDNAElasticProcess, particle);
373
374 // *** Excitation ***
375
376 G4DNAExcitation* theDNAExcitationProcess =
377 new G4DNAExcitation("alpha+_G4DNAExcitation");
378 theDNAExcitationProcess->SetEmModel(
381 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
382 ph->RegisterProcess(theDNAExcitationProcess, particle);
383
384 // *** Ionisation ***
385
386 G4DNAIonisation* theDNAIonisationProcess =
387 new G4DNAIonisation("alpha+_G4DNAIonisation");
388 theDNAIonisationProcess->SetEmModel(
391 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
392 ph->RegisterProcess(theDNAIonisationProcess, particle);
393
394 // *** Charge decrease ***
395
396 G4DNAChargeDecrease* theDNAChargeDecreaseProcess =
397 new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease");
398 theDNAChargeDecreaseProcess->SetEmModel(
401 (theDNAChargeDecreaseProcess->EmModel()))->SelectStationary(true);
402 ph->RegisterProcess(theDNAChargeDecreaseProcess, particle);
403
404 // *** Charge increase ***
405
406 G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
407 new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease");
408 theDNAChargeIncreaseProcess->SetEmModel(
411 (theDNAChargeIncreaseProcess->EmModel()))->SelectStationary(true);
412 ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
413
414 } else if ( particleName == "helium" ) {
415
416 // *** Elastic ***
417
418 G4DNAElastic* theDNAElasticProcess =
419 new G4DNAElastic("helium_G4DNAElastic");
420 theDNAElasticProcess->SetEmModel(new G4DNAIonElasticModel());
422 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
423 ph->RegisterProcess(theDNAElasticProcess, particle);
424
425 // *** Excitation ***
426
427 G4DNAExcitation* theDNAExcitationProcess =
428 new G4DNAExcitation("helium_G4DNAExcitation");
429 theDNAExcitationProcess->SetEmModel(
432 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
433 ph->RegisterProcess(theDNAExcitationProcess, particle);
434
435 // *** Ionisation ***
436
437 G4DNAIonisation* theDNAIonisationProcess =
438 new G4DNAIonisation("helium_G4DNAIonisation");
439 theDNAIonisationProcess->SetEmModel(
442 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
443 ph->RegisterProcess(theDNAIonisationProcess, particle);
444
445 // *** Charge increase ***
446
447 G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
448 new G4DNAChargeIncrease("helium_G4DNAChargeIncrease");
449 theDNAChargeIncreaseProcess->SetEmModel(
452 (theDNAChargeIncreaseProcess->EmModel()))->SelectStationary(true);
453 ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
454
455 } else if ( particleName == "GenericIon" ) {
456
457 // *** Ionisation ***
458
459 G4DNAIonisation* theDNAIonisationProcess =
460 new G4DNAIonisation("GenericIon_G4DNAIonisation");
461 theDNAIonisationProcess->SetEmModel(
464 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
465 ph->RegisterProcess(theDNAIonisationProcess, particle);
466
467 }
468
469 // Warning : the following particles and processes are needed by EM Physics
470 // builders
471 // They are taken from the default Livermore Physics list
472 // These particles are currently not handled by Geant4-DNA
473
474 // e+
475
476 else if (particleName == "e+") {
477
478 // Identical to G4EmStandardPhysics_stationary
479
482 G4eIonisation* eIoni = new G4eIonisation();
483 eIoni->SetStepFunction(0.2, 100*um);
484
485 ph->RegisterProcess(msc, particle);
486 ph->RegisterProcess(eIoni, particle);
487 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
488 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
489
490 } else if (particleName == "gamma") {
491
492 // photoelectric effect - Livermore model only
493 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
494 thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
495 ph->RegisterProcess(thePhotoElectricEffect, particle);
496
497 // Compton scattering - Livermore model only
498 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
499 theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
500 ph->RegisterProcess(theComptonScattering, particle);
501
502 // gamma conversion - Livermore model below 80 GeV
503 G4GammaConversion* theGammaConversion = new G4GammaConversion();
504 theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
505 ph->RegisterProcess(theGammaConversion, particle);
506
507 // default Rayleigh scattering is Livermore
508 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
509 ph->RegisterProcess(theRayleigh, particle);
510 }
511
512 // Warning : end of particles and processes are needed by EM Physics build.
513
514 }
515
516 // Deexcitation
517 //
520}
521
@ bElectromagnetic
G4DNABornExcitationModel1 G4DNABornExcitationModel
#define G4DNABornIonisationModel
G4_DECLARE_PHYSCONSTR_FACTORY(G4EmDNAPhysics_stationary_option6)
@ 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