Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmDNAPhysics_option1 Class Reference

#include <G4EmDNAPhysics_option1.hh>

+ Inheritance diagram for G4EmDNAPhysics_option1:

Public Member Functions

 G4EmDNAPhysics_option1 (G4int ver=1, const G4String &name="")
 
virtual ~G4EmDNAPhysics_option1 ()
 
virtual void ConstructParticle ()
 
virtual void ConstructProcess ()
 
- Public Member Functions inherited from G4VPhysicsConstructor
 G4VPhysicsConstructor (const G4String &="")
 
 G4VPhysicsConstructor (const G4String &name, G4int physics_type)
 
virtual ~G4VPhysicsConstructor ()
 
virtual void ConstructParticle ()=0
 
virtual void ConstructProcess ()=0
 
void SetPhysicsName (const G4String &="")
 
const G4StringGetPhysicsName () const
 
void SetPhysicsType (G4int)
 
G4int GetPhysicsType () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
G4int GetInstanceID () const
 
virtual void TerminateWorker ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VPhysicsConstructor
static const G4VPCManagerGetSubInstanceManager ()
 
- Protected Types inherited from G4VPhysicsConstructor
using PhysicsBuilder_V = G4VPCData::PhysicsBuilders_V
 
- Protected Member Functions inherited from G4VPhysicsConstructor
G4bool RegisterProcess (G4VProcess *process, G4ParticleDefinition *particle)
 
G4ParticleTable::G4PTblDicIteratorGetParticleIterator () const
 
PhysicsBuilder_V GetBuilders () const
 
void AddBuilder (G4PhysicsBuilderInterface *bld)
 
- Protected Attributes inherited from G4VPhysicsConstructor
G4int verboseLevel
 
G4String namePhysics
 
G4int typePhysics
 
G4ParticleTabletheParticleTable
 
G4int g4vpcInstanceID
 
- Static Protected Attributes inherited from G4VPhysicsConstructor
static G4RUN_DLL G4VPCManager subInstanceManager
 

Detailed Description

Definition at line 36 of file G4EmDNAPhysics_option1.hh.

Constructor & Destructor Documentation

◆ G4EmDNAPhysics_option1()

G4EmDNAPhysics_option1::G4EmDNAPhysics_option1 ( G4int  ver = 1,
const G4String name = "" 
)
explicit

Definition at line 89 of file G4EmDNAPhysics_option1.cc.

90 : G4VPhysicsConstructor("G4EmDNAPhysics_option1"), verbose(ver)
91{
93 param->SetDefaults();
94 param->SetFluo(true);
95 param->SetAuger(true);
96 param->SetAugerCascade(true);
97 param->SetDeexcitationIgnoreCut(true);
98 param->ActivateDNA();
99
101}
@ bElectromagnetic
static G4EmParameters * Instance()
void SetDeexcitationIgnoreCut(G4bool val)
void SetFluo(G4bool val)
void SetAugerCascade(G4bool val)
void SetAuger(G4bool val)

◆ ~G4EmDNAPhysics_option1()

G4EmDNAPhysics_option1::~G4EmDNAPhysics_option1 ( )
virtual

Definition at line 105 of file G4EmDNAPhysics_option1.cc.

106{}

Member Function Documentation

◆ ConstructParticle()

void G4EmDNAPhysics_option1::ConstructParticle ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 110 of file G4EmDNAPhysics_option1.cc.

111{
112// bosons
114
115// leptons
118
119// baryons
121
123
124 G4DNAGenericIonsManager * genericIonsManager;
125 genericIonsManager=G4DNAGenericIonsManager::Instance();
126 genericIonsManager->GetIon("alpha++");
127 genericIonsManager->GetIon("alpha+");
128 genericIonsManager->GetIon("helium");
129 genericIonsManager->GetIon("hydrogen");
130 //genericIonsManager->GetIon("carbon");
131 //genericIonsManager->GetIon("nitrogen");
132 //genericIonsManager->GetIon("oxygen");
133 //genericIonsManager->GetIon("iron");
134
135}
static G4DNAGenericIonsManager * Instance(void)
G4ParticleDefinition * GetIon(const G4String &name)
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:87
static G4Positron * Positron()
Definition: G4Positron.cc:93
static G4Proton * Proton()
Definition: G4Proton.cc:92

◆ ConstructProcess()

void G4EmDNAPhysics_option1::ConstructProcess ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 139 of file G4EmDNAPhysics_option1.cc.

140{
141 if(verbose > 1) {
142 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
143 }
145
146 auto myParticleIterator=GetParticleIterator();
147 myParticleIterator->reset();
148 while( (*myParticleIterator)() )
149 {
150 G4ParticleDefinition* particle = myParticleIterator->value();
151 G4String particleName = particle->GetParticleName();
152
153 if (particleName == "e-") {
154
155 ph->RegisterProcess(new G4DNAElectronSolvation("e-_G4DNAElectronSolvation"),
156 particle);
157
158 /*
159 // THESE LINES MUST BE TESTED:
160 G4DNAElectronSolvation* solvation =
161 new G4DNAElectronSolvation("e-_G4DNAElectronSolvation");
162 auto therm = G4DNASolvationModelFactory::GetMacroDefinedModel();
163 therm->SetHighEnergyLimit(7.4*eV); // limit of the Champion's model
164 solvation->SetEmModel(therm);
165 ph->RegisterProcess(solvation, particle);
166 */
167
168 // *** Elastic scattering (two alternative models available) ***
169
170 //G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
171 //theDNAElasticProcess->SetEmModel(new G4DNAChampionElasticModel());
172
173 // or alternative model
174 //theDNAElasticProcess->SetEmModel(new G4DNAScreenedRutherfordElasticModel());
175
176 //ph->RegisterProcess(theDNAElasticProcess, particle);
177
180 ph->RegisterProcess(msc, particle);
181
182
183 // *** Excitation ***
184 ph->RegisterProcess(new G4DNAExcitation("e-_G4DNAExcitation"), particle);
185
186 // *** Ionisation ***
187 ph->RegisterProcess(new G4DNAIonisation("e-_G4DNAIonisation"), particle);
188
189 // *** Vibrational excitation ***
190 ph->RegisterProcess(new G4DNAVibExcitation("e-_G4DNAVibExcitation"), particle);
191
192 // *** Attachment ***
193 ph->RegisterProcess(new G4DNAAttachment("e-_G4DNAAttachment"), particle);
194
195 } else if ( particleName == "proton" ) {
196
199 ph->RegisterProcess(msc, particle);
200
201 ph->RegisterProcess(new G4DNAExcitation("proton_G4DNAExcitation"), particle);
202 ph->RegisterProcess(new G4DNAIonisation("proton_G4DNAIonisation"), particle);
203 ph->RegisterProcess(new G4DNAChargeDecrease("proton_G4DNAChargeDecrease"), particle);
204
205 } else if ( particleName == "hydrogen" ) {
206 ph->RegisterProcess(new G4DNAExcitation("hydrogen_G4DNAExcitation"), particle);
207 ph->RegisterProcess(new G4DNAIonisation("hydrogen_G4DNAIonisation"), particle);
208 ph->RegisterProcess(new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease"), particle);
209
210 } else if ( particleName == "alpha" ) {
211
214 ph->RegisterProcess(msc, particle);
215
216 ph->RegisterProcess(new G4DNAExcitation("alpha_G4DNAExcitation"), particle);
217 ph->RegisterProcess(new G4DNAIonisation("alpha_G4DNAIonisation"), particle);
218 ph->RegisterProcess(new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease"), particle);
219
220 } else if ( particleName == "alpha+" ) {
221
224 ph->RegisterProcess(msc, particle);
225
226 ph->RegisterProcess(new G4DNAExcitation("alpha+_G4DNAExcitation"), particle);
227 ph->RegisterProcess(new G4DNAIonisation("alpha+_G4DNAIonisation"), particle);
228 ph->RegisterProcess(new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease"), particle);
229 ph->RegisterProcess(new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease"), particle);
230
231 } else if ( particleName == "helium" ) {
232 ph->RegisterProcess(new G4DNAExcitation("helium_G4DNAExcitation"), particle);
233 ph->RegisterProcess(new G4DNAIonisation("helium_G4DNAIonisation"), particle);
234 ph->RegisterProcess(new G4DNAChargeIncrease("helium_G4DNAChargeIncrease"), particle);
235
236 // Extension to HZE proposed by Z. Francis
237
238 // S. Incerti
239
240 } else if ( particleName == "GenericIon" ) {
241
243 msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
244 ph->RegisterProcess(msc, particle);
245
246 ph->RegisterProcess(new G4DNAIonisation("GenericIon_G4DNAIonisation"), particle);
247
248 /*
249 } else if ( particleName == "carbon" ) {
250
251 G4hMultipleScattering* msc = new G4hMultipleScattering();
252 msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
253 ph->RegisterProcess(msc, particle);
254
255 ph->RegisterProcess(new G4DNAIonisation("carbon_G4DNAIonisation"), particle);
256
257 } else if ( particleName == "nitrogen" ) {
258
259 G4hMultipleScattering* msc = new G4hMultipleScattering();
260 msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
261 ph->RegisterProcess(msc, particle);
262
263 ph->RegisterProcess(new G4DNAIonisation("nitrogen_G4DNAIonisation"), particle);
264
265 } else if ( particleName == "oxygen" ) {
266
267 G4hMultipleScattering* msc = new G4hMultipleScattering();
268 msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
269 ph->RegisterProcess(msc, particle);
270
271 ph->RegisterProcess(new G4DNAIonisation("oxygen_G4DNAIonisation"), particle);
272
273 } else if ( particleName == "iron" ) {
274
275 G4hMultipleScattering* msc = new G4hMultipleScattering();
276 msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
277 ph->RegisterProcess(msc, particle);
278
279 ph->RegisterProcess(new G4DNAIonisation("iron_G4DNAIonisation"), particle);
280 */
281 }
282
283 // Warning : the following particles and processes are needed by EM Physics builders
284 // They are taken from the default Livermore Physics list
285 // These particles are currently not handled by Geant4-DNA
286
287 // e+
288
289 else if (particleName == "e+") {
290
291 // Identical to G4EmStandardPhysics_option3
292
295 G4eIonisation* eIoni = new G4eIonisation();
296 eIoni->SetStepFunction(0.2, 100*um);
297
298 ph->RegisterProcess(msc, particle);
299 ph->RegisterProcess(eIoni, particle);
300 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
301 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
302
303 } else if (particleName == "gamma") {
304
305 // photoelectric effect - Livermore model only
306 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
307 thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
308 ph->RegisterProcess(thePhotoElectricEffect, particle);
309
310 // Compton scattering - Livermore model only
311 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
312 theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
313 ph->RegisterProcess(theComptonScattering, particle);
314
315 // gamma conversion - Livermore model below 80 GeV
316 G4GammaConversion* theGammaConversion = new G4GammaConversion();
317 theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
318 ph->RegisterProcess(theGammaConversion, particle);
319
320 // default Rayleigh scattering is Livermore
321 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
322 ph->RegisterProcess(theRayleigh, particle);
323 }
324
325 // Warning : end of particles and processes are needed by EM Physics builders
326
327 }
328
329 // Deexcitation
330 //
333}
@ fUseDistanceToBoundary
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void SetAtomDeexcitation(G4VAtomDeexcitation *)
static G4LossTableManager * Instance()
const G4String & GetParticleName() const
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
void SetEmModel(G4VEmModel *, G4int index=0)
void SetStepFunction(G4double v1, G4double v2)
void SetEmModel(G4VMscModel *, size_t index=0)
void SetStepLimitType(G4MscStepLimitType val)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
const G4String & GetPhysicsName() const

The documentation for this class was generated from the following files: