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

#include <G4EmDNAPhysics_option8.hh>

+ Inheritance diagram for G4EmDNAPhysics_option8:

Public Member Functions

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

Constructor & Destructor Documentation

◆ G4EmDNAPhysics_option8()

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

Definition at line 90 of file G4EmDNAPhysics_option8.cc.

91 : G4VPhysicsConstructor("G4EmDNAPhysics_option8"), verbose(ver)
92{
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_option8()

G4EmDNAPhysics_option8::~G4EmDNAPhysics_option8 ( )
virtual

Definition at line 106 of file G4EmDNAPhysics_option8.cc.

107{}

Member Function Documentation

◆ ConstructParticle()

void G4EmDNAPhysics_option8::ConstructParticle ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 111 of file G4EmDNAPhysics_option8.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_option8::ConstructProcess ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 136 of file G4EmDNAPhysics_option8.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 // *** Solvation ***
153
154 G4DNAElectronSolvation* solvation =
155 new G4DNAElectronSolvation("e-_G4DNAElectronSolvation");
157 therm->SetHighEnergyLimit(11.*eV); // limit of the CPA100 model
158 solvation->SetEmModel(therm);
159 ph->RegisterProcess(solvation, particle);
160
161 // *** Elastic scattering using CPA100 model ***
162
163 G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
164 G4VEmModel* model1 = new G4DNACPA100ElasticModel();
166
167 model1->SetActivationLowEnergyLimit(11*eV);
168 model1->SetActivationHighEnergyLimit(255.955*keV);
169 theDNAElasticProcess->SetEmModel(model1);
170
171 model2->SetActivationLowEnergyLimit(255.955*keV);
172 model2->SetActivationHighEnergyLimit(1*MeV);
173 theDNAElasticProcess->AddEmModel(2,model2);
174
175 ph->RegisterProcess(theDNAElasticProcess, particle);
176
177 // *** Excitation ***
178 ph->RegisterProcess(new G4DNAExcitation("e-_G4DNAExcitation"), particle);
179
180 // *** Ionisation ***
181 ph->RegisterProcess(new G4DNAIonisation("e-_G4DNAIonisation"), particle);
182
183 // *** Vibrational excitation ***
184 ph->RegisterProcess(new G4DNAVibExcitation("e-_G4DNAVibExcitation"), particle);
185
186 // *** Attachment ***
187 ph->RegisterProcess(new G4DNAAttachment("e-_G4DNAAttachment"), particle);
188
189 } else if ( particleName == "proton" ) {
190 ph->RegisterProcess(new G4DNAElastic("proton_G4DNAElastic"), particle);
191 ph->RegisterProcess(new G4DNAExcitation("proton_G4DNAExcitation"), particle);
192 ph->RegisterProcess(new G4DNAIonisation("proton_G4DNAIonisation"), particle);
193 ph->RegisterProcess(new G4DNAChargeDecrease("proton_G4DNAChargeDecrease"), particle);
194
195 } else if ( particleName == "hydrogen" ) {
196 ph->RegisterProcess(new G4DNAElastic("hydrogen_G4DNAElastic"), particle);
197 ph->RegisterProcess(new G4DNAExcitation("hydrogen_G4DNAExcitation"), particle);
198 ph->RegisterProcess(new G4DNAIonisation("hydrogen_G4DNAIonisation"), particle);
199 ph->RegisterProcess(new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease"), particle);
200
201 } else if ( particleName == "alpha" ) {
202 ph->RegisterProcess(new G4DNAElastic("alpha_G4DNAElastic"), particle);
203 ph->RegisterProcess(new G4DNAExcitation("alpha_G4DNAExcitation"), particle);
204 ph->RegisterProcess(new G4DNAIonisation("alpha_G4DNAIonisation"), particle);
205 ph->RegisterProcess(new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease"), particle);
206
207 } else if ( particleName == "alpha+" ) {
208 ph->RegisterProcess(new G4DNAElastic("alpha+_G4DNAElastic"), particle);
209 ph->RegisterProcess(new G4DNAExcitation("alpha+_G4DNAExcitation"), particle);
210 ph->RegisterProcess(new G4DNAIonisation("alpha+_G4DNAIonisation"), particle);
211 ph->RegisterProcess(new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease"), particle);
212 ph->RegisterProcess(new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease"), particle);
213
214 } else if ( particleName == "helium" ) {
215 ph->RegisterProcess(new G4DNAElastic("helium_G4DNAElastic"), particle);
216 ph->RegisterProcess(new G4DNAExcitation("helium_G4DNAExcitation"), particle);
217 ph->RegisterProcess(new G4DNAIonisation("helium_G4DNAIonisation"), particle);
218 ph->RegisterProcess(new G4DNAChargeIncrease("helium_G4DNAChargeIncrease"), particle);
219
220 // Extension to HZE proposed by Z. Francis
221
222 } else if ( particleName == "GenericIon" ) {
223 ph->RegisterProcess(new G4DNAIonisation("GenericIon_G4DNAIonisation"), particle);
224 }
225
226 // Warning : the following particles and processes are needed by EM Physics builders
227 // They are taken from the default Livermore Physics list
228 // These particles are currently not handled by Geant4-DNA
229
230 // e+
231
232 else if (particleName == "e+") {
233
234 // Identical to G4EmStandardPhysics_option3
235
238 G4eIonisation* eIoni = new G4eIonisation();
239 eIoni->SetStepFunction(0.2, 100*um);
240
241 ph->RegisterProcess(msc, particle);
242 ph->RegisterProcess(eIoni, particle);
243 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
244 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
245
246 } else if (particleName == "gamma") {
247
248 // photoelectric effect - Livermore model only
249 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
250 thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
251 ph->RegisterProcess(thePhotoElectricEffect, particle);
252
253 // Compton scattering - Livermore model only
254 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
255 theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
256 ph->RegisterProcess(theComptonScattering, particle);
257
258 // gamma conversion - Livermore model below 80 GeV
259 G4GammaConversion* theGammaConversion = new G4GammaConversion();
260 theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
261 ph->RegisterProcess(theGammaConversion, particle);
262
263 // default Rayleigh scattering is Livermore
264 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
265 ph->RegisterProcess(theRayleigh, particle);
266 }
267
268 // Warning : end of particles and processes are needed by EM Physics builders
269
270 }
271
272 // Deexcitation
273 //
276}
@ 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 SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:778
void SetActivationHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:771
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
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: