Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
DicomPhysicsList.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: DicomPhysicsList.cc 75679 2013-11-05 09:08:41Z gcosmo $
27 //
28 /// \file medical/DICOM/src/DicomPhysicsList.cc
29 /// \brief Implementation of the DicomPhysicsList class
30 //
31 
32 #include "G4ParticleDefinition.hh"
33 #include "G4ProcessManager.hh"
34 #include "G4ParticleTypes.hh"
35 #include "G4StepLimiter.hh"
36 #include "G4BaryonConstructor.hh"
37 #include "G4IonConstructor.hh"
38 #include "G4MesonConstructor.hh"
39 #include "G4SystemOfUnits.hh"
40 
41 #include "DicomPhysicsList.hh"
42 
43 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
45 {
47  fCutForGamma = defaultCutValue;
48  fCutForElectron = defaultCutValue;
49  fCutForPositron = defaultCutValue;
50  SetVerboseLevel(1);
51 }
52 
53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
55 {}
56 
57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
59 {
63 }
64 
65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67 {
68  // gamma
70 
71  // optical photon
73 }
74 
75 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
77 {
78  // leptons
81 }
82 
83 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85 {
86  // baryons
87  G4BaryonConstructor bConstructor;
88  bConstructor.ConstructParticle();
89 
90  G4IonConstructor iConstructor;
91  iConstructor.ConstructParticle();
92 
93  G4MesonConstructor mConstructor;
94  mConstructor.ConstructParticle();
95 }
96 
97 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
99 {
101  ConstructEM();
102  ConstructHad();
104 }
105 
106 // *** Processes and models
107 
108 // gamma
109 
110 #include "G4PhotoElectricEffect.hh"
112 
113 #include "G4ComptonScattering.hh"
115 
116 #include "G4GammaConversion.hh"
118 
119 #include "G4RayleighScattering.hh"
121 
122 // e-
123 
124 #include "G4eMultipleScattering.hh"
125 #include "G4UniversalFluctuation.hh"
126 
127 #include "G4eIonisation.hh"
129 
130 #include "G4eBremsstrahlung.hh"
132 
133 // e+
134 
135 #include "G4eplusAnnihilation.hh"
136 
137 // mu
138 
139 #include "G4MuIonisation.hh"
140 #include "G4MuBremsstrahlung.hh"
141 #include "G4MuPairProduction.hh"
142 
143 // hadrons
144 
145 #include "G4hMultipleScattering.hh"
146 #include "G4MscStepLimitType.hh"
147 
148 #include "G4hBremsstrahlung.hh"
149 #include "G4hPairProduction.hh"
150 
151 #include "G4hIonisation.hh"
152 #include "G4ionIonisation.hh"
154 
155 //
156 
157 #include "G4LossTableManager.hh"
158 #include "G4EmProcessOptions.hh"
159 
160 
161 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
163 {
164  theParticleIterator->reset();
165 
166  while( (*theParticleIterator)() ){
167 
168  G4ParticleDefinition* particle = theParticleIterator->value();
169  G4ProcessManager* pmanager = particle->GetProcessManager();
170  G4String particleName = particle->GetParticleName();
171 
172  if (particleName == "gamma") {
173 
174 
175  G4PhotoElectricEffect* thePhotoElectricEffect =
176  new G4PhotoElectricEffect();
177  G4LivermorePhotoElectricModel* theLivermorePhotoElectricModel
179  thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel);
180  pmanager->AddDiscreteProcess(thePhotoElectricEffect);
181 
182  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
183  G4LivermoreComptonModel* theLivermoreComptonModel =
185  theComptonScattering->AddEmModel(0, theLivermoreComptonModel);
186  pmanager->AddDiscreteProcess(theComptonScattering);
187 
188  G4GammaConversion* theGammaConversion = new G4GammaConversion();
189  G4LivermoreGammaConversionModel* theLivermoreGammaConversionModel
191  theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel);
192  pmanager->AddDiscreteProcess(theGammaConversion);
193 
194  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
195  G4LivermoreRayleighModel* theRayleighModel
196  = new G4LivermoreRayleighModel();
197  theRayleigh->AddEmModel(0, theRayleighModel);
198  pmanager->AddDiscreteProcess(theRayleigh);
199 
200  pmanager->AddProcess(new G4StepLimiter(), -1, -1, 5);
201 
202  } else if (particleName == "e-") {
203 
206  pmanager->AddProcess(msc, -1, 1, 1);
207 
208  // Ionisation
209  G4eIonisation* eIoni = new G4eIonisation();
210  eIoni->AddEmModel(0, new G4LivermoreIonisationModel(),
211  new G4UniversalFluctuation() );
212  eIoni->SetStepFunction(0.2, 100*um); //
213  pmanager->AddProcess(eIoni, -1, 2, 2);
214 
215  // Bremsstrahlung
216  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
218  pmanager->AddProcess(eBrem, -1,-3, 3);
219 
220  pmanager->AddProcess(new G4StepLimiter(), -1, -1, 4);
221 
222  } else if (particleName == "e+") {
223 
224  // Identical to G4EmStandardPhysics_option3
225 
228  pmanager->AddProcess(msc, -1, 1, 1);
229 
230  G4eIonisation* eIoni = new G4eIonisation();
231  eIoni->SetStepFunction(0.2, 100*um);
232  pmanager->AddProcess(eIoni, -1, 2, 2);
233 
234  pmanager->AddProcess(new G4eBremsstrahlung, -1,-3, 3);
235 
236  pmanager->AddProcess(new G4eplusAnnihilation,0,-1, 4);
237 
238  pmanager->AddProcess(new G4StepLimiter(), -1, -1, 5);
239 
240  } else if (particleName == "GenericIon") {
241 
242  pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
243 
244  G4ionIonisation* ionIoni = new G4ionIonisation();
245  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
246  ionIoni->SetStepFunction(0.1, 20*um);
247  pmanager->AddProcess(ionIoni, -1, 2, 2);
248 
249  pmanager->AddProcess(new G4StepLimiter(), -1, -1, 3);
250 
251  } else if (particleName == "alpha" ||
252  particleName == "He3" ) {
253 
254  // Identical to G4EmStandardPhysics_option3
255 
256  pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
257 
258  G4ionIonisation* ionIoni = new G4ionIonisation();
259  ionIoni->SetStepFunction(0.1, 20*um);
260  pmanager->AddProcess(ionIoni, -1, 2, 2);
261 
262  pmanager->AddProcess(new G4StepLimiter(), -1, -1, 3);
263  }
264 
265  //
266 
267  }
268 }
269 
270 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
271 #include "G4HadronElasticProcess.hh"
272 #include "G4HadronElastic.hh"
273 
276 #include "G4TripathiCrossSection.hh"
277 #include "G4IonsShenCrossSection.hh"
278 
280 {
281 
282  G4HadronElasticProcess * theElasticProcess = new G4HadronElasticProcess;
283  theElasticProcess->RegisterMe( new G4HadronElastic() );
284 
285  theParticleIterator->reset();
286  while( (*theParticleIterator)() )
287  {
288  G4ParticleDefinition* particle = theParticleIterator->value();
289  G4ProcessManager* pManager = particle->GetProcessManager();
290 
291  if (particle->GetParticleName() == "alpha") {
292 
293  // INELASTIC SCATTERING
294  // Binary Cascade
296  // DHW - change lower limit from 80 MeV to zero because
297  // G4LEDeuteronInelastic is deprecated
298  theBC->SetMinEnergy(0.*MeV);
299  theBC->SetMaxEnergy(40.*GeV);
300 
301  // TRIPATHI CROSS SECTION
302  // Implementation of formulas in analogy to NASA technical paper 3621 by
303  // Tripathi, et al. Cross-sections for ion ion scattering
304  G4TripathiCrossSection* TripathiCrossSection = new G4TripathiCrossSection;
305 
306  // IONS SHEN CROSS SECTION
307  // Implementation of formulas
308  // Shen et al. Nuc. Phys. A 491 130 (1989)
309  // Total Reaction Cross Section for Heavy-Ion Collisions
311 
313  theIPalpha->AddDataSet(TripathiCrossSection);
314  theIPalpha->AddDataSet(aShen);
315 
316  // Register the Binary Cascade Model
317  theIPalpha->RegisterMe(theBC);
318 
319  // Activate the alpha inelastic scattering using the binary cascade model
320  pManager -> AddDiscreteProcess(theIPalpha);
321 
322  // Activate the Hadron Elastic Process
323  pManager -> AddDiscreteProcess(theElasticProcess);
324  }
325  }
326 }
327 
328 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
330 { }
331 
332 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
334 {
335  if (verboseLevel >0){
336  G4cout << "DicomPhysicsList::SetCuts:";
337  G4cout << "CutLength : " << G4BestUnit(defaultCutValue,"Length") << G4endl;
338  }
339 
340  // set cut values for gamma at first and for e- second and next for e+,
341  // because some processes for e+/e- need cut values for gamma
342  SetCutValue(fCutForGamma, "gamma");
343  SetCutValue(fCutForElectron, "e-");
344  SetCutValue(fCutForPositron, "e+");
345 
347 
348 }
349 
virtual void SetCuts()
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
Definition of the DicomPhysicsList class.
virtual void ConstructProcess()
void SetCutValue(G4double aCut, const G4String &pname)
void SetStepFunction(G4double v1, G4double v2)
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
static void ConstructParticle()
static void ConstructParticle()
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleName() const
static void ConstructParticle()
void RegisterMe(G4HadronicInteraction *a)
void SetMinEnergy(G4double anEnergy)
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
void DumpCutValuesTable(G4int flag=1)
G4GLOB_DLL std::ostream G4cout
void SetVerboseLevel(G4int value)
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
virtual void ConstructParticle()
static G4Positron * PositronDefinition()
Definition: G4Positron.cc:89
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
void SetEmModel(G4VEmModel *, G4int index=1)
void SetMaxEnergy(const G4double anEnergy)
#define G4endl
Definition: G4ios.hh:61
static G4OpticalPhoton * OpticalPhotonDefinition()
int micrometer
Definition: hepunit.py:34
#define theParticleIterator
void SetStepLimitType(G4MscStepLimitType val)
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81