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

#include <G4EmDNAPhysics_option3.hh>

+ Inheritance diagram for G4EmDNAPhysics_option3:

Public Member Functions

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

Constructor & Destructor Documentation

◆ G4EmDNAPhysics_option3()

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

Definition at line 88 of file G4EmDNAPhysics_option3.cc.

89 : G4VPhysicsConstructor("G4EmDNAPhysics_option3"), 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}
@ bElectromagnetic
static G4EmParameters * Instance()
void SetDeexcitationIgnoreCut(G4bool val)
void SetFluo(G4bool val)
void SetAugerCascade(G4bool val)
void SetAuger(G4bool val)

◆ ~G4EmDNAPhysics_option3()

G4EmDNAPhysics_option3::~G4EmDNAPhysics_option3 ( )
virtual

Definition at line 104 of file G4EmDNAPhysics_option3.cc.

105{}

Member Function Documentation

◆ ConstructParticle()

void G4EmDNAPhysics_option3::ConstructParticle ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 109 of file G4EmDNAPhysics_option3.cc.

110{
111// bosons
113
114// leptons
117
118// baryons
120
122
123 G4DNAGenericIonsManager * genericIonsManager;
124 genericIonsManager=G4DNAGenericIonsManager::Instance();
125 genericIonsManager->GetIon("alpha++");
126 genericIonsManager->GetIon("alpha+");
127 genericIonsManager->GetIon("helium");
128 genericIonsManager->GetIon("hydrogen");
129
130}
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_option3::ConstructProcess ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 134 of file G4EmDNAPhysics_option3.cc.

135{
136 if(verbose > 1) {
137 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
138 }
140
141 auto myParticleIterator=GetParticleIterator();
142 myParticleIterator->reset();
143 while( (*myParticleIterator)() )
144 {
145 G4ParticleDefinition* particle = myParticleIterator->value();
146 G4String particleName = particle->GetParticleName();
147
148 if (particleName == "e-") {
149
150 G4DNAElectronSolvation* solvation =
151 new G4DNAElectronSolvation("e-_G4DNAElectronSolvation");
153 therm->SetHighEnergyLimit(7.4*eV); // limit of the Champion's model
154 solvation->SetEmModel(therm);
155 ph->RegisterProcess(solvation, particle);
156
157 // *** Elastic scattering (two alternative models available) ***
158
159 G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
160 theDNAElasticProcess->SetEmModel(new G4DNAChampionElasticModel());
161
162 // or alternative model
163 //theDNAElasticProcess->SetEmModel(new G4DNAScreenedRutherfordElasticModel());
164
165 ph->RegisterProcess(theDNAElasticProcess, particle);
166
167 // *** Excitation ***
168 ph->RegisterProcess(new G4DNAExcitation("e-_G4DNAExcitation"), particle);
169
170 // *** Ionisation ***
171 ph->RegisterProcess(new G4DNAIonisation("e-_G4DNAIonisation"), particle);
172
173 // *** Vibrational excitation ***
174 ph->RegisterProcess(new G4DNAVibExcitation("e-_G4DNAVibExcitation"), particle);
175
176 // *** Attachment ***
177 ph->RegisterProcess(new G4DNAAttachment("e-_G4DNAAttachment"), particle);
178
179 } else if ( particleName == "proton" ) {
180 ph->RegisterProcess(new G4DNAElastic("proton_G4DNAElastic"), particle);
181 ph->RegisterProcess(new G4DNAExcitation("proton_G4DNAExcitation"), particle);
182 ph->RegisterProcess(new G4DNAIonisation("proton_G4DNAIonisation"), particle);
183 ph->RegisterProcess(new G4DNAChargeDecrease("proton_G4DNAChargeDecrease"), particle);
184
185 } else if ( particleName == "hydrogen" ) {
186 ph->RegisterProcess(new G4DNAElastic("hydrogen_G4DNAElastic"), particle);
187 ph->RegisterProcess(new G4DNAExcitation("hydrogen_G4DNAExcitation"), particle);
188 ph->RegisterProcess(new G4DNAIonisation("hydrogen_G4DNAIonisation"), particle);
189 ph->RegisterProcess(new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease"), particle);
190
191 } else if ( particleName == "alpha" ) {
192 ph->RegisterProcess(new G4DNAElastic("alpha_G4DNAElastic"), particle);
193 ph->RegisterProcess(new G4DNAExcitation("alpha_G4DNAExcitation"), particle);
194 ph->RegisterProcess(new G4DNAIonisation("alpha_G4DNAIonisation"), particle);
195 ph->RegisterProcess(new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease"), particle);
196
197 } else if ( particleName == "alpha+" ) {
198 ph->RegisterProcess(new G4DNAElastic("alpha+_G4DNAElastic"), particle);
199 ph->RegisterProcess(new G4DNAExcitation("alpha+_G4DNAExcitation"), particle);
200 ph->RegisterProcess(new G4DNAIonisation("alpha+_G4DNAIonisation"), particle);
201 ph->RegisterProcess(new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease"), particle);
202 ph->RegisterProcess(new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease"), particle);
203
204 } else if ( particleName == "helium" ) {
205 ph->RegisterProcess(new G4DNAElastic("helium_G4DNAElastic"), particle);
206 ph->RegisterProcess(new G4DNAExcitation("helium_G4DNAExcitation"), particle);
207 ph->RegisterProcess(new G4DNAIonisation("helium_G4DNAIonisation"), particle);
208 ph->RegisterProcess(new G4DNAChargeIncrease("helium_G4DNAChargeIncrease"), particle);
209
210 } else if ( particleName == "GenericIon" ) {
211 ph->RegisterProcess(new G4DNAIonisation("GenericIon_G4DNAIonisation"), particle);
212
213 }
214
215 // Warning : the following particles and processes are needed by EM Physics builders
216 // They are taken from the default Livermore Physics list
217 // These particles are currently not handled by Geant4-DNA
218
219 // e+
220
221 else if (particleName == "e+") {
222
223 // Identical to G4EmStandardPhysics_option3
224
227 G4eIonisation* eIoni = new G4eIonisation();
228 eIoni->SetStepFunction(0.2, 100*um);
229
230 ph->RegisterProcess(msc, particle);
231 ph->RegisterProcess(eIoni, particle);
232 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
233 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
234
235 } else if (particleName == "gamma") {
236
237 // photoelectric effect - Livermore model only
238 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
239 thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
240 ph->RegisterProcess(thePhotoElectricEffect, particle);
241
242 // Compton scattering - Livermore model only
243 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
244 theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
245 ph->RegisterProcess(theComptonScattering, particle);
246
247 // gamma conversion - Livermore model below 80 GeV
248 G4GammaConversion* theGammaConversion = new G4GammaConversion();
249 theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
250 ph->RegisterProcess(theGammaConversion, particle);
251
252 // default Rayleigh scattering is Livermore
253 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
254 ph->RegisterProcess(theRayleigh, particle);
255 }
256
257 // Warning : end of particles and processes are needed by EM Physics builders
258
259 }
260
261 // Deexcitation
262 //
265}
@ 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 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: