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

#include <G4EmDNAPhysics_option2.hh>

+ Inheritance diagram for G4EmDNAPhysics_option2:

Public Member Functions

 G4EmDNAPhysics_option2 (G4int ver=1, const G4String &name="")
 
virtual ~G4EmDNAPhysics_option2 ()
 
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_option2.hh.

Constructor & Destructor Documentation

◆ G4EmDNAPhysics_option2()

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

Definition at line 89 of file G4EmDNAPhysics_option2.cc.

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

◆ ~G4EmDNAPhysics_option2()

G4EmDNAPhysics_option2::~G4EmDNAPhysics_option2 ( )
virtual

Definition at line 106 of file G4EmDNAPhysics_option2.cc.

107{}

Member Function Documentation

◆ ConstructParticle()

void G4EmDNAPhysics_option2::ConstructParticle ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 111 of file G4EmDNAPhysics_option2.cc.

112{
113// bosons
115
116// leptons
119
120// baryons
122
124
125 G4DNAGenericIonsManager * genericIonsManager;
126 genericIonsManager=G4DNAGenericIonsManager::Instance();
127 genericIonsManager->GetIon("alpha++");
128 genericIonsManager->GetIon("alpha+");
129 genericIonsManager->GetIon("helium");
130 genericIonsManager->GetIon("hydrogen");
131
132}
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_option2::ConstructProcess ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 136 of file G4EmDNAPhysics_option2.cc.

137{
138 if(verbose > 1) {
139 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
140 }
142
143 auto myParticleIterator=GetParticleIterator();
144 myParticleIterator->reset();
145 while( (*myParticleIterator)() )
146 {
147 G4ParticleDefinition* particle = myParticleIterator->value();
148 G4String particleName = particle->GetParticleName();
149
150 if (particleName == "e-") {
151
152 G4DNAElectronSolvation* solvation =
153 new G4DNAElectronSolvation("e-_G4DNAElectronSolvation");
154
156 therm->SetHighEnergyLimit(7.4*eV); // limit of the Champion's model
157 //therm->SetHighEnergyLimit(10*eV); // limit of the ELSEPA model
158 solvation->SetEmModel(therm);
159 ph->RegisterProcess(solvation, particle);
160
161 // *** Elastic scattering (two alternative models available) ***
162
163 G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
164 theDNAElasticProcess->SetEmModel(new G4DNAChampionElasticModel());
165
166 // or alternative model
167 //theDNAElasticProcess->SetEmModel(new G4DNAELSEPAElasticModel());
168 //theDNAElasticProcess->SetEmModel(new G4DNAScreenedRutherfordElasticModel());
169
170 ph->RegisterProcess(theDNAElasticProcess, particle);
171
172 // *** Excitation ***
173 ph->RegisterProcess(new G4DNAExcitation("e-_G4DNAExcitation"), particle);
174
175 // *** Ionisation ***
176 G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("e-_G4DNAIonisation");
178 mod->SelectFasterComputation(true);
179 theDNAIonisationProcess->SetEmModel(mod);
180 ph->RegisterProcess(theDNAIonisationProcess, particle);
181
182 // *** Vibrational excitation ***
183 ph->RegisterProcess(new G4DNAVibExcitation("e-_G4DNAVibExcitation"), particle);
184
185 // *** Attachment ***
186 ph->RegisterProcess(new G4DNAAttachment("e-_G4DNAAttachment"), particle);
187
188 } else if ( particleName == "proton" ) {
189
190 ph->RegisterProcess(new G4DNAElastic("proton_G4DNAElastic"), particle);
191
192 ph->RegisterProcess(new G4DNAExcitation("proton_G4DNAExcitation"), particle);
193
194 G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("proton_G4DNAIonisation");
195
197 mod1->SetLowEnergyLimit(0*eV);
198 mod1->SetHighEnergyLimit(500*keV);
199
201 mod2->SetLowEnergyLimit(500*keV);
202 mod2->SetHighEnergyLimit(100*MeV);
203 mod2->SelectFasterComputation(true);
204
205 theDNAIonisationProcess->SetEmModel(mod1);
206 theDNAIonisationProcess->SetEmModel(mod2);
207
208 ph->RegisterProcess(theDNAIonisationProcess, particle);
209
210 ph->RegisterProcess(new G4DNAChargeDecrease("proton_G4DNAChargeDecrease"), particle);
211
212 } else if ( particleName == "hydrogen" ) {
213
214 ph->RegisterProcess(new G4DNAElastic("hydrogen_G4DNAElastic"), particle);
215
216 ph->RegisterProcess(new G4DNAExcitation("hydrogen_G4DNAExcitation"), particle);
217
218 //ph->RegisterProcess(new G4DNAIonisation("hydrogen_G4DNAIonisation"), particle);
219 G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("hydrogen_G4DNAIonisation");
220 theDNAIonisationProcess->SetEmModel(new G4DNARuddIonisationExtendedModel());
221 ph->RegisterProcess(theDNAIonisationProcess, particle);
222
223 ph->RegisterProcess(new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease"), particle);
224
225 } else if ( particleName == "alpha" ) {
226
227 ph->RegisterProcess(new G4DNAElastic("alpha_G4DNAElastic"), particle);
228
229 ph->RegisterProcess(new G4DNAExcitation("alpha_G4DNAExcitation"), particle);
230
231 G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("alpha_G4DNAIonisation");
232 theDNAIonisationProcess->SetEmModel(new G4DNARuddIonisationExtendedModel());
233 ph->RegisterProcess(theDNAIonisationProcess, particle);
234
235 ph->RegisterProcess(new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease"), particle);
236
237 } else if ( particleName == "alpha+" ) {
238
239 ph->RegisterProcess(new G4DNAElastic("alpha+_G4DNAElastic"), particle);
240
241 ph->RegisterProcess(new G4DNAExcitation("alpha+_G4DNAExcitation"), particle);
242
243 G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("alpha+_G4DNAIonisation");
244 theDNAIonisationProcess->SetEmModel(new G4DNARuddIonisationExtendedModel());
245 ph->RegisterProcess(theDNAIonisationProcess, particle);
246
247 ph->RegisterProcess(new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease"), particle);
248 ph->RegisterProcess(new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease"), particle);
249
250 } else if ( particleName == "helium" ) {
251
252 ph->RegisterProcess(new G4DNAElastic("helium_G4DNAElastic"), particle);
253
254 ph->RegisterProcess(new G4DNAExcitation("helium_G4DNAExcitation"), particle);
255
256 G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("helium_G4DNAIonisation");
257 theDNAIonisationProcess->SetEmModel(new G4DNARuddIonisationExtendedModel());
258 ph->RegisterProcess(theDNAIonisationProcess, particle);
259
260 ph->RegisterProcess(new G4DNAChargeIncrease("helium_G4DNAChargeIncrease"), particle);
261
262 } else if ( particleName == "GenericIon" ) {
263 ph->RegisterProcess(new G4DNAIonisation("GenericIon_G4DNAIonisation"), particle);
264 }
265
266 // Warning : the following particles and processes are needed by EM Physics builders
267 // They are taken from the default Livermore Physics list
268 // These particles are currently not handled by Geant4-DNA
269
270 // e+
271
272 else if (particleName == "e+") {
273
276 G4eIonisation* eIoni = new G4eIonisation();
277 eIoni->SetStepFunction(0.2, 100*um);
278
279 ph->RegisterProcess(msc, particle);
280 ph->RegisterProcess(eIoni, particle);
281 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
282 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
283
284 } else if (particleName == "gamma") {
285
286 // photoelectric effect - Livermore model only
287 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
288 thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
289 ph->RegisterProcess(thePhotoElectricEffect, particle);
290
291 // Compton scattering - Livermore model only
292 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
293 theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
294 ph->RegisterProcess(theComptonScattering, particle);
295
296 // gamma conversion - Livermore model below 80 GeV
297 G4GammaConversion* theGammaConversion = new G4GammaConversion();
298 theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
299 ph->RegisterProcess(theGammaConversion, particle);
300
301 // default Rayleigh scattering is Livermore
302 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
303 ph->RegisterProcess(theRayleigh, particle);
304 }
305
306 // Warning : end of particles and processes are needed by EM Physics builders
307
308 }
309
310 // Deexcitation
311 //
314}
#define G4DNABornIonisationModel
@ fUseDistanceToBoundary
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4VEmModel * GetMacroDefinedModel()
One step thermalization model can be chosen via macro using /process/dna/e-SolvationSubType Ritchie19...
void SetAtomDeexcitation(G4VAtomDeexcitation *)
static G4LossTableManager * Instance()
const G4String & GetParticleName() const
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764
void SetEmModel(G4VEmModel *, G4int index=0)
void SetStepFunction(G4double v1, G4double v2)
void SetStepLimitType(G4MscStepLimitType val)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
const G4String & GetPhysicsName() const

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