Geant4-11
G4EmDNAPhysics_stationary_option2.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"
38
39#include "G4DNAExcitation.hh"
40#include "G4DNAAttachment.hh"
41#include "G4DNAVibExcitation.hh"
42#include "G4DNAIonisation.hh"
45
46// particles
47
48#include "G4Electron.hh"
49#include "G4Proton.hh"
50#include "G4GenericIon.hh"
51
52// Warning : the following is needed in order to use EM Physics builders
53// e+
54#include "G4Positron.hh"
56#include "G4eIonisation.hh"
57#include "G4eBremsstrahlung.hh"
59// gamma
60#include "G4Gamma.hh"
65#include "G4GammaConversion.hh"
69
70#include "G4EmParameters.hh"
71// end of warning
72
73#include "G4LossTableManager.hh"
76#include "G4BuilderType.hh"
77
78// factory
80//
82
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
85
87 : G4VPhysicsConstructor("G4EmDNAPhysics_stationary_option2"), verbose(ver)
88{
90 param->SetDefaults();
91 param->SetFluo(true);
92 param->SetAuger(true);
93 param->SetDeexcitationIgnoreCut(true);
94 param->ActivateDNA();
95
97}
98
99//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
100
102const G4String&)
104{}
105
106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
107
109{}
110
111//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
112
114{
115// bosons
117
118// leptons
121
122// baryons
124
126
127 G4DNAGenericIonsManager * genericIonsManager;
128 genericIonsManager=G4DNAGenericIonsManager::Instance();
129 genericIonsManager->GetIon("alpha++");
130 genericIonsManager->GetIon("alpha+");
131 genericIonsManager->GetIon("helium");
132 genericIonsManager->GetIon("hydrogen");
133
134}
135
136//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
137
139{
140 if(verbose > 1) {
141 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
142 }
144
145 auto myParticleIterator=GetParticleIterator();
146 myParticleIterator->reset();
147 while( (*myParticleIterator)() )
148 {
149 G4ParticleDefinition* particle = myParticleIterator->value();
150 G4String particleName = particle->GetParticleName();
151
152 if (particleName == "e-") {
153
154 // *** Elastic scattering (two alternative models available) ***
155
156 G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
157 theDNAElasticProcess->SetEmModel(new G4DNAChampionElasticModel());
158
159 // or alternative model
160 //theDNAElasticProcess
161 //->SetEmModel(new G4DNAScreenedRutherfordElasticModel());
162
163 ph->RegisterProcess(theDNAElasticProcess, particle);
164
165 // *** Excitation ***
166
167 G4DNAExcitation* theDNAExcitationProcess =
168 new G4DNAExcitation("e-_G4DNAExcitation");
169 theDNAExcitationProcess->SetEmModel(new G4DNABornExcitationModel());
170 ((G4DNABornExcitationModel*)(theDNAExcitationProcess->EmModel()))
171 ->SelectStationary(true);
172 ph->RegisterProcess(theDNAExcitationProcess, particle);
173
174 // *** Ionisation ***
175
176 G4DNAIonisation* theDNAIonisationProcess =
177 new G4DNAIonisation("e-_G4DNAIonisation");
178 theDNAIonisationProcess->SetEmModel(new G4DNABornIonisationModel());
179 ((G4DNABornIonisationModel*)(theDNAIonisationProcess->EmModel()))
180 ->SelectStationary(true);
181 //
182 ((G4DNABornIonisationModel*)(theDNAIonisationProcess->EmModel()))
183 ->SelectFasterComputation(true);
184 //
185 ph->RegisterProcess(theDNAIonisationProcess, particle);
186
187 // *** Vibrational excitation ***
188
189 G4DNAVibExcitation* theDNAVibExcitationProcess =
190 new G4DNAVibExcitation("e-_G4DNAVibExcitation");
191 theDNAVibExcitationProcess->SetEmModel(new G4DNASancheExcitationModel());
192 ((G4DNASancheExcitationModel*)(theDNAVibExcitationProcess->EmModel()))
193 ->SelectStationary(true);
194 ph->RegisterProcess(theDNAVibExcitationProcess, particle);
195
196 // *** Attachment ***
197
198 G4DNAAttachment* theDNAAttachmentProcess =
199 new G4DNAAttachment("e-_G4DNAAttachment");
200 theDNAAttachmentProcess->SetEmModel(new G4DNAMeltonAttachmentModel());
201 ((G4DNAMeltonAttachmentModel*)(theDNAAttachmentProcess->EmModel()))
202 ->SelectStationary(true);
203 ph->RegisterProcess(theDNAAttachmentProcess, particle);
204
205 } else if ( particleName == "proton" ) {
206
207 // *** Elastic ***
208
209 G4DNAElastic* theDNAElasticProcess =
210 new G4DNAElastic("proton_G4DNAElastic");
211 theDNAElasticProcess->SetEmModel(new G4DNAIonElasticModel());
212 ((G4DNAIonElasticModel*)(theDNAElasticProcess->EmModel()))
213 ->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
226 ((G4DNAMillerGreenExcitationModel*)(theDNAExcitationProcess->EmModel()))
227 ->SetLowEnergyLimit(10*eV);
228 ((G4DNAMillerGreenExcitationModel*)(theDNAExcitationProcess->EmModel()))
229 ->SetHighEnergyLimit(500*keV);
230 ((G4DNAMillerGreenExcitationModel*)(theDNAExcitationProcess->EmModel()))
231 ->SelectStationary(true);
232
233 ((G4DNABornExcitationModel*)(theDNAExcitationProcess->EmModel(1)))
234 ->SetLowEnergyLimit(500*keV);
235 ((G4DNABornExcitationModel*)(theDNAExcitationProcess->EmModel(1)))
236 ->SetHighEnergyLimit(100*MeV);
237 ((G4DNABornExcitationModel*)(theDNAExcitationProcess->EmModel(1)))
238 ->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(new G4DNAIonElasticModel());
289 ((G4DNAIonElasticModel*)(theDNAElasticProcess->EmModel()))
290 ->SelectStationary(true);
291 ph->RegisterProcess(theDNAElasticProcess, particle);
292
293 // *** Excitation ***
294
295 G4DNAExcitation* theDNAExcitationProcess =
296 new G4DNAExcitation("hydrogen_G4DNAExcitation");
297 theDNAExcitationProcess->SetEmModel(
299 ((G4DNAMillerGreenExcitationModel*)(theDNAExcitationProcess->EmModel()))
300 ->SelectStationary(true);
301 ph->RegisterProcess(theDNAExcitationProcess, particle);
302
303 // *** Ionisation ***
304
305 G4DNAIonisation* theDNAIonisationProcess =
306 new G4DNAIonisation("hydrogen_G4DNAIonisation");
307 theDNAIonisationProcess->SetEmModel(
310 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
311 ph->RegisterProcess(theDNAIonisationProcess, particle);
312
313 // *** Charge increase ***
314
315 G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
316 new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease");
317 theDNAChargeIncreaseProcess->SetEmModel(
320 (theDNAChargeIncreaseProcess->EmModel()))->SelectStationary(true);
321 ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
322
323 } else if ( particleName == "alpha" ) {
324
325 // *** Elastic ***
326
327 G4DNAElastic* theDNAElasticProcess =
328 new G4DNAElastic("alpha_G4DNAElastic");
329 theDNAElasticProcess->SetEmModel(new G4DNAIonElasticModel());
331 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
332 ph->RegisterProcess(theDNAElasticProcess, particle);
333
334 // *** Excitation ***
335
336 G4DNAExcitation* theDNAExcitationProcess =
337 new G4DNAExcitation("alpha_G4DNAExcitation");
338 theDNAExcitationProcess->SetEmModel(
341 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
342 ph->RegisterProcess(theDNAExcitationProcess, particle);
343
344 // *** Ionisation ***
345
346 G4DNAIonisation* theDNAIonisationProcess =
347 new G4DNAIonisation("alpha_G4DNAIonisation");
348 theDNAIonisationProcess->SetEmModel(
351 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
352 ph->RegisterProcess(theDNAIonisationProcess, particle);
353
354 // *** Charge decrease ***
355
356 G4DNAChargeDecrease* theDNAChargeDecreaseProcess =
357 new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease");
358 theDNAChargeDecreaseProcess->SetEmModel(
361 (theDNAChargeDecreaseProcess->EmModel()))->SelectStationary(true);
362 ph->RegisterProcess(theDNAChargeDecreaseProcess, particle);
363
364 } else if ( particleName == "alpha+" ) {
365
366 // *** Elastic ***
367
368 G4DNAElastic* theDNAElasticProcess =
369 new G4DNAElastic("alpha+_G4DNAElastic");
370 theDNAElasticProcess->SetEmModel(new G4DNAIonElasticModel());
372 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
373 ph->RegisterProcess(theDNAElasticProcess, particle);
374
375 // *** Excitation ***
376
377 G4DNAExcitation* theDNAExcitationProcess =
378 new G4DNAExcitation("alpha+_G4DNAExcitation");
379 theDNAExcitationProcess->SetEmModel(
382 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
383 ph->RegisterProcess(theDNAExcitationProcess, particle);
384
385 // *** Ionisation ***
386
387 G4DNAIonisation* theDNAIonisationProcess =
388 new G4DNAIonisation("alpha+_G4DNAIonisation");
389 theDNAIonisationProcess->SetEmModel(
392 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
393 ph->RegisterProcess(theDNAIonisationProcess, particle);
394
395 // *** Charge decrease ***
396
397 G4DNAChargeDecrease* theDNAChargeDecreaseProcess =
398 new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease");
399 theDNAChargeDecreaseProcess->SetEmModel(
402 (theDNAChargeDecreaseProcess->EmModel()))->SelectStationary(true);
403 ph->RegisterProcess(theDNAChargeDecreaseProcess, particle);
404
405 // *** Charge increase ***
406
407 G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
408 new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease");
409 theDNAChargeIncreaseProcess->SetEmModel(
412 (theDNAChargeIncreaseProcess->EmModel()))->SelectStationary(true);
413 ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
414
415 } else if ( particleName == "helium" ) {
416
417 // *** Elastic ***
418
419 G4DNAElastic* theDNAElasticProcess =
420 new G4DNAElastic("helium_G4DNAElastic");
421 theDNAElasticProcess->SetEmModel(
424 (theDNAElasticProcess->EmModel()))->SelectStationary(true);
425 ph->RegisterProcess(theDNAElasticProcess, particle);
426
427 // *** Excitation ***
428
429 G4DNAExcitation* theDNAExcitationProcess =
430 new G4DNAExcitation("helium_G4DNAExcitation");
431 theDNAExcitationProcess->SetEmModel(
434 (theDNAExcitationProcess->EmModel()))->SelectStationary(true);
435 ph->RegisterProcess(theDNAExcitationProcess, particle);
436
437 // *** Ionisation ***
438
439 G4DNAIonisation* theDNAIonisationProcess =
440 new G4DNAIonisation("helium_G4DNAIonisation");
441 theDNAIonisationProcess->SetEmModel(
444 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
445 ph->RegisterProcess(theDNAIonisationProcess, particle);
446
447 // *** Charge increase ***
448
449 G4DNAChargeIncrease* theDNAChargeIncreaseProcess =
450 new G4DNAChargeIncrease("helium_G4DNAChargeIncrease");
451 theDNAChargeIncreaseProcess->SetEmModel(
454 (theDNAChargeIncreaseProcess->EmModel()))->SelectStationary(true);
455 ph->RegisterProcess(theDNAChargeIncreaseProcess, particle);
456
457 } else if ( particleName == "GenericIon" ) {
458
459 // *** Ionisation ***
460
461 G4DNAIonisation* theDNAIonisationProcess =
462 new G4DNAIonisation("GenericIon_G4DNAIonisation");
463 theDNAIonisationProcess->SetEmModel(
466 (theDNAIonisationProcess->EmModel()))->SelectStationary(true);
467 ph->RegisterProcess(theDNAIonisationProcess, particle);
468
469 }
470
471 // Warning : the following particles and processes are needed by EM Physics
472 // builders
473 // They are taken from the default Livermore Physics list
474 // These particles are currently not handled by Geant4-DNA
475
476 // e+
477
478 else if (particleName == "e+") {
479
480 // Identical to G4EmStandardPhysics_stationary
481
484 G4eIonisation* eIoni = new G4eIonisation();
485 eIoni->SetStepFunction(0.2, 100*um);
486
487 ph->RegisterProcess(msc, particle);
488 ph->RegisterProcess(eIoni, particle);
489 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
490 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
491
492 } else if (particleName == "gamma") {
493
494 // photoelectric effect - Livermore model only
495 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
496 thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
497 ph->RegisterProcess(thePhotoElectricEffect, particle);
498
499 // Compton scattering - Livermore model only
500 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
501 theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
502 ph->RegisterProcess(theComptonScattering, particle);
503
504 // gamma conversion - Livermore model below 80 GeV
505 G4GammaConversion* theGammaConversion = new G4GammaConversion();
506 theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
507 ph->RegisterProcess(theGammaConversion, particle);
508
509 // default Rayleigh scattering is Livermore
510 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
511 ph->RegisterProcess(theRayleigh, particle);
512 }
513
514 // Warning : end of particles and processes are needed by EM Physics build.
515
516 }
517
518 // Deexcitation
519 //
522 de->SetFluo(true);
523}
524
@ bElectromagnetic
G4DNABornExcitationModel1 G4DNABornExcitationModel
#define G4DNABornIonisationModel
G4_DECLARE_PHYSCONSTR_FACTORY(G4EmDNAPhysics_stationary_option2)
@ 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