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

#include <G4EmDNAPhysics_option4.hh>

+ Inheritance diagram for G4EmDNAPhysics_option4:

Public Member Functions

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

Constructor & Destructor Documentation

◆ G4EmDNAPhysics_option4()

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

Definition at line 94 of file G4EmDNAPhysics_option4.cc.

95 : G4VPhysicsConstructor("G4EmDNAPhysics_option4"), verbose(ver)
96{
98 param->SetDefaults();
99 param->SetFluo(true);
100 param->SetAuger(true);
101 param->SetAugerCascade(true);
102 param->SetDeexcitationIgnoreCut(true);
103 param->ActivateDNA();
104
106}
@ bElectromagnetic
static G4EmParameters * Instance()
void SetDeexcitationIgnoreCut(G4bool val)
void SetFluo(G4bool val)
void SetAugerCascade(G4bool val)
void SetAuger(G4bool val)

◆ ~G4EmDNAPhysics_option4()

G4EmDNAPhysics_option4::~G4EmDNAPhysics_option4 ( )
virtual

Definition at line 110 of file G4EmDNAPhysics_option4.cc.

111{}

Member Function Documentation

◆ ConstructParticle()

void G4EmDNAPhysics_option4::ConstructParticle ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 115 of file G4EmDNAPhysics_option4.cc.

116{
117// bosons
119
120// leptons
123
124// baryons
126
128
129 G4DNAGenericIonsManager * genericIonsManager;
130 genericIonsManager=G4DNAGenericIonsManager::Instance();
131 genericIonsManager->GetIon("alpha++");
132 genericIonsManager->GetIon("alpha+");
133 genericIonsManager->GetIon("helium");
134 genericIonsManager->GetIon("hydrogen");
135
136}
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_option4::ConstructProcess ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 140 of file G4EmDNAPhysics_option4.cc.

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