Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4EmDNAPhysics_option1 Class Reference

#include <G4EmDNAPhysics_option1.hh>

Inheritance diagram for G4EmDNAPhysics_option1:
Collaboration diagram for G4EmDNAPhysics_option1:

Public Member Functions

 G4EmDNAPhysics_option1 (G4int ver=1, const G4String &name="")
 
virtual ~G4EmDNAPhysics_option1 ()
 
virtual void ConstructParticle ()
 
virtual void ConstructProcess ()
 
- Public Member Functions inherited from G4VPhysicsConstructor
 G4VPhysicsConstructor (const G4String &="")
 
 G4VPhysicsConstructor (const G4String &name, G4int physics_type)
 
virtual ~G4VPhysicsConstructor ()
 
void SetPhysicsName (const G4String &="")
 
const G4StringGetPhysicsName () const
 
void SetPhysicsType (G4int)
 
G4int GetPhysicsType () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
G4int GetInstanceID () const
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VPhysicsConstructor
static const G4VPCManagerGetSubInstanceManager ()
 
- Protected Member Functions inherited from G4VPhysicsConstructor
G4bool RegisterProcess (G4VProcess *process, G4ParticleDefinition *particle)
 
G4ParticleTable::G4PTblDicIteratorGetParticleIterator () const
 
- 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_option1.hh.

Constructor & Destructor Documentation

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

Definition at line 90 of file G4EmDNAPhysics_option1.cc.

91  : G4VPhysicsConstructor("G4EmDNAPhysics_option1"), verbose(ver)
92 {
94  param->SetDefaults();
95  param->SetFluo(true);
96  param->SetAuger(true);
97  param->SetAugerCascade(true);
98  param->SetDeexcitationIgnoreCut(true);
99 
101 }
void SetDeexcitationIgnoreCut(G4bool val)
void SetAuger(G4bool val)
void SetAugerCascade(G4bool val)
static G4EmParameters * Instance()
G4VPhysicsConstructor(const G4String &="")
void SetFluo(G4bool val)

Here is the call graph for this function:

G4EmDNAPhysics_option1::~G4EmDNAPhysics_option1 ( )
virtual

Definition at line 105 of file G4EmDNAPhysics_option1.cc.

106 {}

Member Function Documentation

void G4EmDNAPhysics_option1::ConstructParticle ( void  )
virtual

Implements G4VPhysicsConstructor.

Definition at line 110 of file G4EmDNAPhysics_option1.cc.

111 {
112 // bosons
113  G4Gamma::Gamma();
114 
115 // leptons
118 
119 // baryons
121 
123 
124  G4DNAGenericIonsManager * genericIonsManager;
125  genericIonsManager=G4DNAGenericIonsManager::Instance();
126  genericIonsManager->GetIon("alpha++");
127  genericIonsManager->GetIon("alpha+");
128  genericIonsManager->GetIon("helium");
129  genericIonsManager->GetIon("hydrogen");
130  //genericIonsManager->GetIon("carbon");
131  //genericIonsManager->GetIon("nitrogen");
132  //genericIonsManager->GetIon("oxygen");
133  //genericIonsManager->GetIon("iron");
134 
135 }
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:88
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4DNAGenericIonsManager * Instance(void)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4ParticleDefinition * GetIon(const G4String &name)

Here is the call graph for this function:

void G4EmDNAPhysics_option1::ConstructProcess ( void  )
virtual

Implements G4VPhysicsConstructor.

Definition at line 139 of file G4EmDNAPhysics_option1.cc.

140 {
141  if(verbose > 1) {
142  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
143  }
145 
146  auto myParticleIterator=GetParticleIterator();
147  myParticleIterator->reset();
148  while( (*myParticleIterator)() )
149  {
150  G4ParticleDefinition* particle = myParticleIterator->value();
151  G4String particleName = particle->GetParticleName();
152 
153  if (particleName == "e-") {
154 
155  ph->RegisterProcess(new G4DNAElectronSolvation("e-_G4DNAElectronSolvation"),
156  particle);
157 
158  // *** Elastic scattering (two alternative models available) ***
159 
160  //G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
161  //theDNAElasticProcess->SetEmModel(new G4DNAChampionElasticModel());
162 
163  // or alternative model
164  //theDNAElasticProcess->SetEmModel(new G4DNAScreenedRutherfordElasticModel());
165 
166  //ph->RegisterProcess(theDNAElasticProcess, particle);
167 
169  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
170  ph->RegisterProcess(msc, particle);
171 
172 
173  // *** Excitation ***
174  ph->RegisterProcess(new G4DNAExcitation("e-_G4DNAExcitation"), particle);
175 
176  // *** Ionisation ***
177  ph->RegisterProcess(new G4DNAIonisation("e-_G4DNAIonisation"), particle);
178 
179  // *** Vibrational excitation ***
180  ph->RegisterProcess(new G4DNAVibExcitation("e-_G4DNAVibExcitation"), particle);
181 
182  // *** Attachment ***
183  ph->RegisterProcess(new G4DNAAttachment("e-_G4DNAAttachment"), particle);
184 
185  } else if ( particleName == "proton" ) {
186 
188  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
189  ph->RegisterProcess(msc, particle);
190 
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 G4DNAExcitation("hydrogen_G4DNAExcitation"), particle);
197  ph->RegisterProcess(new G4DNAIonisation("hydrogen_G4DNAIonisation"), particle);
198  ph->RegisterProcess(new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease"), particle);
199 
200  } else if ( particleName == "alpha" ) {
201 
203  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
204  ph->RegisterProcess(msc, particle);
205 
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 
213  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
214  ph->RegisterProcess(msc, particle);
215 
216  ph->RegisterProcess(new G4DNAExcitation("alpha+_G4DNAExcitation"), particle);
217  ph->RegisterProcess(new G4DNAIonisation("alpha+_G4DNAIonisation"), particle);
218  ph->RegisterProcess(new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease"), particle);
219  ph->RegisterProcess(new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease"), particle);
220 
221  } else if ( particleName == "helium" ) {
222  ph->RegisterProcess(new G4DNAExcitation("helium_G4DNAExcitation"), particle);
223  ph->RegisterProcess(new G4DNAIonisation("helium_G4DNAIonisation"), particle);
224  ph->RegisterProcess(new G4DNAChargeIncrease("helium_G4DNAChargeIncrease"), particle);
225 
226  // Extension to HZE proposed by Z. Francis
227 
228  // S. Incerti
229 
230  } else if ( particleName == "GenericIon" ) {
231 
233  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
234  ph->RegisterProcess(msc, particle);
235 
236  ph->RegisterProcess(new G4DNAIonisation("GenericIon_G4DNAIonisation"), particle);
237 
238  /*
239  } else if ( particleName == "carbon" ) {
240 
241  G4hMultipleScattering* msc = new G4hMultipleScattering();
242  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
243  ph->RegisterProcess(msc, particle);
244 
245  ph->RegisterProcess(new G4DNAIonisation("carbon_G4DNAIonisation"), particle);
246 
247  } else if ( particleName == "nitrogen" ) {
248 
249  G4hMultipleScattering* msc = new G4hMultipleScattering();
250  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
251  ph->RegisterProcess(msc, particle);
252 
253  ph->RegisterProcess(new G4DNAIonisation("nitrogen_G4DNAIonisation"), particle);
254 
255  } else if ( particleName == "oxygen" ) {
256 
257  G4hMultipleScattering* msc = new G4hMultipleScattering();
258  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
259  ph->RegisterProcess(msc, particle);
260 
261  ph->RegisterProcess(new G4DNAIonisation("oxygen_G4DNAIonisation"), particle);
262 
263  } else if ( particleName == "iron" ) {
264 
265  G4hMultipleScattering* msc = new G4hMultipleScattering();
266  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
267  ph->RegisterProcess(msc, particle);
268 
269  ph->RegisterProcess(new G4DNAIonisation("iron_G4DNAIonisation"), particle);
270  */
271  }
272 
273  // Warning : the following particles and processes are needed by EM Physics builders
274  // They are taken from the default Livermore Physics list
275  // These particles are currently not handled by Geant4-DNA
276 
277  // e+
278 
279  else if (particleName == "e+") {
280 
281  // Identical to G4EmStandardPhysics_option3
282 
285  G4eIonisation* eIoni = new G4eIonisation();
286  eIoni->SetStepFunction(0.2, 100*um);
287 
288  ph->RegisterProcess(msc, particle);
289  ph->RegisterProcess(eIoni, particle);
290  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
291  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
292 
293  } else if (particleName == "gamma") {
294 
295  G4double LivermoreHighEnergyLimit = GeV;
296 
297  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
298  G4LivermorePhotoElectricModel* theLivermorePhotoElectricModel =
300  theLivermorePhotoElectricModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
301  thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel);
302  ph->RegisterProcess(thePhotoElectricEffect, particle);
303 
304  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
305  G4LivermoreComptonModel* theLivermoreComptonModel =
307  theLivermoreComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
308  theComptonScattering->AddEmModel(0, theLivermoreComptonModel);
309  ph->RegisterProcess(theComptonScattering, particle);
310 
311  G4GammaConversion* theGammaConversion = new G4GammaConversion();
312  G4LivermoreGammaConversionModel* theLivermoreGammaConversionModel =
314  theLivermoreGammaConversionModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
315  theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel);
316  ph->RegisterProcess(theGammaConversion, particle);
317 
318  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
319  G4LivermoreRayleighModel* theRayleighModel = new G4LivermoreRayleighModel();
320  theRayleighModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
321  theRayleigh->AddEmModel(0, theRayleighModel);
322  ph->RegisterProcess(theRayleigh, particle);
323  }
324 
325  // Warning : end of particles and processes are needed by EM Physics builders
326 
327  }
328 
329  // Deexcitation
330  //
333 }
static G4LossTableManager * Instance()
void SetEmModel(G4VMscModel *, G4int index=1)
void SetStepFunction(G4double v1, G4double v2, G4bool lock=true)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
G4GLOB_DLL std::ostream G4cout
static constexpr double um
Definition: G4SIunits.hh:113
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
const G4String & GetPhysicsName() const
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
static G4PhysicsListHelper * GetPhysicsListHelper()
static constexpr double GeV
Definition: G4SIunits.hh:217
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetAtomDeexcitation(G4VAtomDeexcitation *)
void SetStepLimitType(G4MscStepLimitType val)

Here is the call graph for this function:


The documentation for this class was generated from the following files: