Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
OpNovicePhysicsList.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 //
28 //
29 //
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 
34 #include "globals.hh"
35 #include "OpNovicePhysicsList.hh"
37 
38 #include "G4ParticleDefinition.hh"
39 #include "G4ParticleTypes.hh"
40 #include "G4ParticleTable.hh"
41 
42 #include "G4BosonConstructor.hh"
43 #include "G4LeptonConstructor.hh"
44 #include "G4MesonConstructor.hh"
45 #include "G4BaryonConstructor.hh"
46 #include "G4IonConstructor.hh"
48 
49 #include "G4ProcessManager.hh"
50 
51 #include "G4Cerenkov.hh"
52 #include "G4Scintillation.hh"
53 #include "G4OpAbsorption.hh"
54 #include "G4OpRayleigh.hh"
55 #include "G4OpMieHG.hh"
56 #include "G4OpBoundaryProcess.hh"
57 
58 #include "G4LossTableManager.hh"
59 #include "G4EmSaturation.hh"
60 
61 G4ThreadLocal G4int OpNovicePhysicsList::fVerboseLevel = 1;
62 G4ThreadLocal G4int OpNovicePhysicsList::fMaxNumPhotonStep = 20;
63 G4ThreadLocal G4Cerenkov* OpNovicePhysicsList::fCerenkovProcess = 0;
64 G4ThreadLocal G4Scintillation* OpNovicePhysicsList::fScintillationProcess = 0;
65 G4ThreadLocal G4OpAbsorption* OpNovicePhysicsList::fAbsorptionProcess = 0;
66 G4ThreadLocal G4OpRayleigh* OpNovicePhysicsList::fRayleighScatteringProcess = 0;
67 G4ThreadLocal G4OpMieHG* OpNovicePhysicsList::fMieHGScatteringProcess = 0;
68 G4ThreadLocal G4OpBoundaryProcess* OpNovicePhysicsList::fBoundaryProcess = 0;
69 
70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71 
74 {
75  fMessenger = new OpNovicePhysicsListMessenger(this);
76 }
77 
78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
79 
81 
82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
83 
85 {
86  // In this method, static member functions should be called
87  // for all particles which you want to use.
88  // This ensures that objects of these particle types will be
89  // created in the program.
90 
91  G4BosonConstructor bConstructor;
92  bConstructor.ConstructParticle();
93 
94  G4LeptonConstructor lConstructor;
95  lConstructor.ConstructParticle();
96 
97  G4MesonConstructor mConstructor;
98  mConstructor.ConstructParticle();
99 
100  G4BaryonConstructor rConstructor;
101  rConstructor.ConstructParticle();
102 
103  G4IonConstructor iConstructor;
104  iConstructor.ConstructParticle();
105 }
106 
107 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
108 
110 {
112  ConstructDecay();
113  ConstructEM();
114  ConstructOp();
115 }
116 
117 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
118 
119 #include "G4Decay.hh"
120 
121 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
122 
124 {
125  // Add Decay Process
126  G4Decay* theDecayProcess = new G4Decay();
128  particleIterator->reset();
129  while( (*particleIterator)() ){
130  G4ParticleDefinition* particle = particleIterator->value();
131  G4ProcessManager* pmanager = particle->GetProcessManager();
132  if (theDecayProcess->IsApplicable(*particle)) {
133  pmanager ->AddProcess(theDecayProcess);
134  // set ordering for PostStepDoIt and AtRestDoIt
135  pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
136  pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
137  }
138  }
139 }
140 
141 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
142 
143 #include "G4ComptonScattering.hh"
144 #include "G4GammaConversion.hh"
145 #include "G4PhotoElectricEffect.hh"
146 
147 #include "G4eMultipleScattering.hh"
148 #include "G4MuMultipleScattering.hh"
149 #include "G4hMultipleScattering.hh"
150 
151 #include "G4eIonisation.hh"
152 #include "G4eBremsstrahlung.hh"
153 #include "G4eplusAnnihilation.hh"
154 
155 #include "G4MuIonisation.hh"
156 #include "G4MuBremsstrahlung.hh"
157 #include "G4MuPairProduction.hh"
158 
159 #include "G4hIonisation.hh"
160 
161 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
162 
164 {
166  particleIterator->reset();
167  while( (*particleIterator)() ){
168  G4ParticleDefinition* particle = particleIterator->value();
169  G4ProcessManager* pmanager = particle->GetProcessManager();
170  G4String particleName = particle->GetParticleName();
171 
172  if (particleName == "gamma") {
173  // gamma
174  // Construct processes for gamma
175  pmanager->AddDiscreteProcess(new G4GammaConversion());
176  pmanager->AddDiscreteProcess(new G4ComptonScattering());
177  pmanager->AddDiscreteProcess(new G4PhotoElectricEffect());
178 
179  } else if (particleName == "e-") {
180  //electron
181  // Construct processes for electron
182  pmanager->AddProcess(new G4eMultipleScattering(),-1, 1, 1);
183  pmanager->AddProcess(new G4eIonisation(), -1, 2, 2);
184  pmanager->AddProcess(new G4eBremsstrahlung(), -1, 3, 3);
185 
186  } else if (particleName == "e+") {
187  //positron
188  // Construct processes for positron
189  pmanager->AddProcess(new G4eMultipleScattering(),-1, 1, 1);
190  pmanager->AddProcess(new G4eIonisation(), -1, 2, 2);
191  pmanager->AddProcess(new G4eBremsstrahlung(), -1, 3, 3);
192  pmanager->AddProcess(new G4eplusAnnihilation(), 0,-1, 4);
193 
194  } else if( particleName == "mu+" ||
195  particleName == "mu-" ) {
196  //muon
197  // Construct processes for muon
198  pmanager->AddProcess(new G4MuMultipleScattering(),-1, 1, 1);
199  pmanager->AddProcess(new G4MuIonisation(), -1, 2, 2);
200  pmanager->AddProcess(new G4MuBremsstrahlung(), -1, 3, 3);
201  pmanager->AddProcess(new G4MuPairProduction(), -1, 4, 4);
202 
203  } else {
204  if ((particle->GetPDGCharge() != 0.0) &&
205  (particle->GetParticleName() != "chargedgeantino") &&
206  !particle->IsShortLived()) {
207  // all others charged particles except geantino
208  pmanager->AddProcess(new G4hMultipleScattering(),-1,1,1);
209  pmanager->AddProcess(new G4hIonisation(), -1,2,2);
210  }
211  }
212  }
213 }
214 
215 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
216 #include "G4Threading.hh"
217 
219 {
220  fCerenkovProcess = new G4Cerenkov("Cerenkov");
221  fCerenkovProcess->SetMaxNumPhotonsPerStep(fMaxNumPhotonStep);
222  fCerenkovProcess->SetMaxBetaChangePerStep(10.0);
223  fCerenkovProcess->SetTrackSecondariesFirst(true);
224  fScintillationProcess = new G4Scintillation("Scintillation");
225  fScintillationProcess->SetScintillationYieldFactor(1.);
226  fScintillationProcess->SetTrackSecondariesFirst(true);
227  fAbsorptionProcess = new G4OpAbsorption();
228  fRayleighScatteringProcess = new G4OpRayleigh();
229  fMieHGScatteringProcess = new G4OpMieHG();
230  fBoundaryProcess = new G4OpBoundaryProcess();
231 
232  fCerenkovProcess->SetVerboseLevel(fVerboseLevel);
233  fScintillationProcess->SetVerboseLevel(fVerboseLevel);
234  fAbsorptionProcess->SetVerboseLevel(fVerboseLevel);
235  fRayleighScatteringProcess->SetVerboseLevel(fVerboseLevel);
236  fMieHGScatteringProcess->SetVerboseLevel(fVerboseLevel);
237  fBoundaryProcess->SetVerboseLevel(fVerboseLevel);
238 
239  // Use Birks Correction in the Scintillation process
241  {
242  G4EmSaturation* emSaturation =
244  fScintillationProcess->AddSaturation(emSaturation);
245  }
246 
248  particleIterator->reset();
249  while( (*particleIterator)() ){
250  G4ParticleDefinition* particle = particleIterator->value();
251  G4ProcessManager* pmanager = particle->GetProcessManager();
252  G4String particleName = particle->GetParticleName();
253  if (fCerenkovProcess->IsApplicable(*particle)) {
254  pmanager->AddProcess(fCerenkovProcess);
255  pmanager->SetProcessOrdering(fCerenkovProcess,idxPostStep);
256  }
257  if (fScintillationProcess->IsApplicable(*particle)) {
258  pmanager->AddProcess(fScintillationProcess);
259  pmanager->SetProcessOrderingToLast(fScintillationProcess, idxAtRest);
260  pmanager->SetProcessOrderingToLast(fScintillationProcess, idxPostStep);
261  }
262  if (particleName == "opticalphoton") {
263  G4cout << " AddDiscreteProcess to OpticalPhoton " << G4endl;
264  pmanager->AddDiscreteProcess(fAbsorptionProcess);
265  pmanager->AddDiscreteProcess(fRayleighScatteringProcess);
266  pmanager->AddDiscreteProcess(fMieHGScatteringProcess);
267  pmanager->AddDiscreteProcess(fBoundaryProcess);
268  }
269  }
270 }
271 
272 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
273 
275 {
276  fVerboseLevel = verbose;
277 
278  fCerenkovProcess->SetVerboseLevel(fVerboseLevel);
279  fScintillationProcess->SetVerboseLevel(fVerboseLevel);
280  fAbsorptionProcess->SetVerboseLevel(fVerboseLevel);
281  fRayleighScatteringProcess->SetVerboseLevel(fVerboseLevel);
282  fMieHGScatteringProcess->SetVerboseLevel(fVerboseLevel);
283  fBoundaryProcess->SetVerboseLevel(fVerboseLevel);
284 }
285 
286 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
287 
289 {
290  fMaxNumPhotonStep = MaxNumber;
291 
292  fCerenkovProcess->SetMaxNumPhotonsPerStep(fMaxNumPhotonStep);
293 }
294 
295 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
296 
298 {
299  // " G4VUserPhysicsList::SetCutsWithDefault" method sets
300  // the default cut value for all particle types
301  //
303 
305 }
306 
307 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
void SetMaxBetaChangePerStep(const G4double d)
Definition: G4Cerenkov.cc:150
static G4LossTableManager * Instance()
void SetNbOfPhotonsCerenkov(G4int)
virtual void ConstructProcess()
void SetTrackSecondariesFirst(const G4bool state)
Definition: G4Cerenkov.cc:145
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4Decay.cc:89
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
Definition of the OpNovicePhysicsListMessenger class.
static void ConstructParticle()
#define G4ThreadLocal
Definition: tls.hh:89
static void ConstructParticle()
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
static void ConstructParticle()
void SetScintillationYieldFactor(const G4double yieldfactor)
void DumpCutValuesTable(G4int flag=1)
G4GLOB_DLL std::ostream G4cout
static void ConstructParticle()
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
G4EmSaturation * EmSaturation()
void SetMaxNumPhotonsPerStep(const G4int NumPhotons)
Definition: G4Cerenkov.cc:155
void SetProcessOrdering(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt, G4int ordDoIt=ordDefault)
static void ConstructParticle()
Definition of the OpNovicePhysicsList class.
void SetTrackSecondariesFirst(const G4bool state)
void SetProcessOrderingToLast(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt)
G4ProcessManager * GetProcessManager() const
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
Definition: G4Cerenkov.cc:133
virtual void ConstructParticle()
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
#define G4endl
Definition: G4ios.hh:61
G4bool IsMasterThread()
Definition: G4Threading.cc:146
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
G4double GetPDGCharge() const
void AddSaturation(G4EmSaturation *sat)
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:437