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

#include <G4EmDNAPhysics_option7.hh>

+ Inheritance diagram for G4EmDNAPhysics_option7:

Public Member Functions

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

Constructor & Destructor Documentation

◆ G4EmDNAPhysics_option7()

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

Definition at line 94 of file G4EmDNAPhysics_option7.cc.

94 :
95 G4VPhysicsConstructor("G4EmDNAPhysics_option7"), 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_option7()

G4EmDNAPhysics_option7::~G4EmDNAPhysics_option7 ( )
virtual

Definition at line 110 of file G4EmDNAPhysics_option7.cc.

111{
112}

Member Function Documentation

◆ ConstructParticle()

void G4EmDNAPhysics_option7::ConstructParticle ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 116 of file G4EmDNAPhysics_option7.cc.

117{
118// bosons
120
121// leptons
124
125// baryons
127
129
130 G4DNAGenericIonsManager * genericIonsManager;
131 genericIonsManager = G4DNAGenericIonsManager::Instance();
132 genericIonsManager->GetIon("alpha++");
133 genericIonsManager->GetIon("alpha+");
134 genericIonsManager->GetIon("helium");
135 genericIonsManager->GetIon("hydrogen");
136
137}
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_option7::ConstructProcess ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 141 of file G4EmDNAPhysics_option7.cc.

142{
143 if(verbose > 1)
144 {
145 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
146 }
148
149 auto myParticleIterator=GetParticleIterator();
150 myParticleIterator->reset();
151 while( (*myParticleIterator)() )
152 {
153 G4ParticleDefinition* particle = myParticleIterator->value();
154 G4String particleName = particle->GetParticleName();
155
156 if (particleName == "e-")
157 {
158 // *** Electron solvation ***
159 // Either kills the electron track or transform the electron to
160 // a solvated electron, depending on whether the chemistry is
161 // activated
162
163 G4DNAElectronSolvation* solvation =
164 new G4DNAElectronSolvation("e-_G4DNAElectronSolvation");
166 therm->SetHighEnergyLimit(10.*eV); // limit of the Uehara's model
167 solvation->SetEmModel(therm);
168 ph->RegisterProcess(solvation, particle);
169
170 // *** Elastic scattering ***
171
172 G4DNAElastic* theDNAElasticProcess =
173 new G4DNAElastic("e-_G4DNAElastic");
176 emElast->SetHighEnergyLimit(1*MeV);
177 theDNAElasticProcess->SetEmModel(emElast);
178 ph->RegisterProcess(theDNAElasticProcess, particle);
179
180 // *** Excitation ***
181
182 G4DNAExcitation* theDNAExcitationProcess =
183 new G4DNAExcitation("e-_G4DNAExcitation");
184
185 // 1-st model should be Set
187 emExc->SetActivationLowEnergyLimit(8*eV);
188 emExc->SetActivationHighEnergyLimit(10.*keV);
189 theDNAExcitationProcess->SetEmModel(emExc);
190
191 // 2-nd model should be Add
193 bornExc->SetActivationLowEnergyLimit(10*keV);
194 bornExc->SetActivationHighEnergyLimit(1.*MeV);
195 theDNAExcitationProcess->AddEmModel(2, bornExc);
196
197 ph->RegisterProcess(theDNAExcitationProcess, particle);
198
199 // *** Ionisation ***
200
201 G4DNAIonisation* theDNAIonisationProcess =
202 new G4DNAIonisation("e-_G4DNAIonisation");
203
204 // 1-st model should be Set
206 emIonModel->SetActivationLowEnergyLimit(10.*eV);
207 emIonModel->SetActivationHighEnergyLimit(10.*keV);
208 theDNAIonisationProcess->SetEmModel(emIonModel);
209
210 // 2-nd model should be Add
212 bornIon->SetActivationLowEnergyLimit(10*keV);
213 bornIon->SetActivationHighEnergyLimit(1.*MeV);
214 theDNAIonisationProcess->AddEmModel(2,bornIon);
215
216 ph->RegisterProcess(theDNAIonisationProcess, particle);
217
218 // *** Vibrational excitation ***
219 ph->RegisterProcess(new G4DNAVibExcitation("e-_G4DNAVibExcitation"), particle);
220
221 // *** Attachment ***
222 ph->RegisterProcess(new G4DNAAttachment("e-_G4DNAAttachment"), particle);
223
224 } else if ( particleName == "proton" ) {
225 ph->RegisterProcess(new G4DNAElastic("proton_G4DNAElastic"), particle);
226 ph->RegisterProcess(new G4DNAExcitation("proton_G4DNAExcitation"), particle);
227 ph->RegisterProcess(new G4DNAIonisation("proton_G4DNAIonisation"), particle);
228 ph->RegisterProcess(new G4DNAChargeDecrease("proton_G4DNAChargeDecrease"), particle);
229
230 } else if ( particleName == "hydrogen" ) {
231 ph->RegisterProcess(new G4DNAElastic("hydrogen_G4DNAElastic"), particle);
232 ph->RegisterProcess(new G4DNAExcitation("hydrogen_G4DNAExcitation"), particle);
233 ph->RegisterProcess(new G4DNAIonisation("hydrogen_G4DNAIonisation"), particle);
234 ph->RegisterProcess(new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease"), particle);
235
236 } else if ( particleName == "alpha" ) {
237 ph->RegisterProcess(new G4DNAElastic("alpha_G4DNAElastic"), particle);
238 ph->RegisterProcess(new G4DNAExcitation("alpha_G4DNAExcitation"), particle);
239 ph->RegisterProcess(new G4DNAIonisation("alpha_G4DNAIonisation"), particle);
240 ph->RegisterProcess(new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease"), particle);
241
242 } else if ( particleName == "alpha+" ) {
243 ph->RegisterProcess(new G4DNAElastic("alpha+_G4DNAElastic"), particle);
244 ph->RegisterProcess(new G4DNAExcitation("alpha+_G4DNAExcitation"), particle);
245 ph->RegisterProcess(new G4DNAIonisation("alpha+_G4DNAIonisation"), particle);
246 ph->RegisterProcess(new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease"), particle);
247 ph->RegisterProcess(new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease"), particle);
248
249 } else if ( particleName == "helium" ) {
250 ph->RegisterProcess(new G4DNAElastic("helium_G4DNAElastic"), particle);
251 ph->RegisterProcess(new G4DNAExcitation("helium_G4DNAExcitation"), particle);
252 ph->RegisterProcess(new G4DNAIonisation("helium_G4DNAIonisation"), particle);
253 ph->RegisterProcess(new G4DNAChargeIncrease("helium_G4DNAChargeIncrease"), particle);
254 }
255 // Extension to HZE proposed by Z. Francis
256 else if ( particleName == "GenericIon" ) {
257 ph->RegisterProcess(new G4DNAIonisation("GenericIon_G4DNAIonisation"), particle);
258 }
259
260 // Warning : the following particles and processes are needed by EM Physics builders
261 // They are taken from the default Livermore Physics list
262 // These particles are currently not handled by Geant4-DNA
263
264 // e+
265
266 else if (particleName == "e+") {
267
268 // Identical to G4EmStandardPhysics_option3
269
272 G4eIonisation* eIoni = new G4eIonisation();
273 eIoni->SetStepFunction(0.2, 100*um);
274
275 ph->RegisterProcess(msc, particle);
276 ph->RegisterProcess(eIoni, particle);
277 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
278 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
279
280 } else if (particleName == "gamma") {
281
282 // photoelectric effect - Livermore model only
283 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
284 thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
285 ph->RegisterProcess(thePhotoElectricEffect, particle);
286
287 // Compton scattering - Livermore model only
288 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
289 theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
290 ph->RegisterProcess(theComptonScattering, particle);
291
292 // gamma conversion - Livermore model below 80 GeV
293 G4GammaConversion* theGammaConversion = new G4GammaConversion();
294 theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
295 ph->RegisterProcess(theGammaConversion, particle);
296
297 // default Rayleigh scattering is Livermore
298 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
299 ph->RegisterProcess(theRayleigh, particle);
300 }
301
302 // Warning : end of particles and processes are needed by EM Physics builders
303
304 }
305
306 // Deexcitation
307 //
310}
G4DNABornExcitationModel1 G4DNABornExcitationModel
#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 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: