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

#include <G4EmDNAPhysics.hh>

+ Inheritance diagram for G4EmDNAPhysics:

Public Member Functions

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

Constructor & Destructor Documentation

◆ G4EmDNAPhysics()

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

Definition at line 88 of file G4EmDNAPhysics.cc.

89 : G4VPhysicsConstructor("G4EmDNAPhysics"), 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()

G4EmDNAPhysics::~G4EmDNAPhysics ( )
virtual

Definition at line 104 of file G4EmDNAPhysics.cc.

105{}

Member Function Documentation

◆ ConstructParticle()

void G4EmDNAPhysics::ConstructParticle ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 109 of file G4EmDNAPhysics.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 //genericIonsManager->GetIon("carbon");
130 //genericIonsManager->GetIon("nitrogen");
131 //genericIonsManager->GetIon("oxygen");
132 //genericIonsManager->GetIon("iron");
133
134}
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::ConstructProcess ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 138 of file G4EmDNAPhysics.cc.

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