Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
examples/extended/electromagnetic/TestEm1/src/PhysicsList.cc
이 파일의 문서화 페이지로 가기
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 // $Id: PhysicsList.cc 108857 2018-03-12 07:42:27Z gcosmo $
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 
34 #include "PhysicsList.hh"
35 #include "PhysicsListMessenger.hh"
36 
37 #include "PhysListEmStandard.hh"
38 
39 #include "G4EmStandardPhysics.hh"
44 #include "G4EmStandardPhysicsSS.hh"
45 #include "G4EmStandardPhysicsGS.hh"
47 #include "G4EmLivermorePhysics.hh"
48 #include "G4EmPenelopePhysics.hh"
49 #include "G4EmLowEPPhysics.hh"
50 
51 #include "DetectorConstruction.hh"
52 
53 #include "G4LossTableManager.hh"
54 #include "G4UnitsTable.hh"
55 #include "G4SystemOfUnits.hh"
56 
57 // particles
58 
59 #include "G4BosonConstructor.hh"
60 #include "G4LeptonConstructor.hh"
61 #include "G4MesonConstructor.hh"
62 #include "G4BosonConstructor.hh"
63 #include "G4BaryonConstructor.hh"
64 #include "G4IonConstructor.hh"
66 
67 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68 
70  : G4VModularPhysicsList(), fDet(det)
71 {
72  fMessenger = new PhysicsListMessenger(this);
73  SetVerboseLevel(1);
74 
75  // EM physics
76  fEmName = G4String("local");
78 
80  // fix lower limit for cut
83 }
84 
85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
86 
88 {
89  delete fMessenger;
90 }
91 
92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
93 
95 {
96  G4BosonConstructor pBosonConstructor;
97  pBosonConstructor.ConstructParticle();
98 
99  G4LeptonConstructor pLeptonConstructor;
100  pLeptonConstructor.ConstructParticle();
101 
102  G4MesonConstructor pMesonConstructor;
103  pMesonConstructor.ConstructParticle();
104 
105  G4BaryonConstructor pBaryonConstructor;
106  pBaryonConstructor.ConstructParticle();
107 
108  G4IonConstructor pIonConstructor;
109  pIonConstructor.ConstructParticle();
110 
111  G4ShortLivedConstructor pShortLivedConstructor;
112  pShortLivedConstructor.ConstructParticle();
113 }
114 
115 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
116 
118 {
119  // Transportation
120  //
122 
123  // Electromagnetic physics list
124  //
126 
127  // Decay Process
128  //
129  AddDecay();
130 
131  // Decay Process
132  //
134 
135  // step limitation (as a full process)
136  //
137  AddStepMax();
138 
139  // Get process
140  auto process = GetProcess("RadioactiveDecay");
141  if (process != nullptr)
142  G4cout << "\n GetProcess : " << process->GetProcessName() << G4endl;
143 }
144 
145 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
146 
148 {
149  if (verboseLevel>1) {
150  G4cout << "PhysicsList::AddPhysicsList: <" << name << ">" << G4endl;
151  }
152 
153  if (name == fEmName) return;
154 
155  if (name == "local") {
156  fEmName = name;
157  delete fEmPhysicsList;
159 
160  } else if (name == "emstandard_opt0") {
161  fEmName = name;
162  delete fEmPhysicsList;
164 
165  } else if (name == "emstandard_opt1") {
166  fEmName = name;
167  delete fEmPhysicsList;
169 
170  } else if (name == "emstandard_opt2") {
171  fEmName = name;
172  delete fEmPhysicsList;
174 
175  } else if (name == "emstandard_opt3") {
176  fEmName = name;
177  delete fEmPhysicsList;
179 
180  } else if (name == "emstandard_opt4") {
181  fEmName = name;
182  delete fEmPhysicsList;
184 
185  } else if (name == "emstandardSS") {
186  fEmName = name;
187  delete fEmPhysicsList;
189 
190  } else if (name == "emstandardGS") {
191  fEmName = name;
192  delete fEmPhysicsList;
194 
195  } else if (name == "emstandardWVI") {
196  fEmName = name;
197  delete fEmPhysicsList;
199 
200  } else if (name == "emlivermore") {
201  fEmName = name;
202  delete fEmPhysicsList;
204 
205  } else if (name == "empenelope") {
206  fEmName = name;
207  delete fEmPhysicsList;
209 
210  } else if (name == "emlowenergy") {
211  fEmName = name;
212  delete fEmPhysicsList;
214 
215  } else {
216 
217  G4cout << "PhysicsList::AddPhysicsList: <" << name << ">"
218  << " is not defined"
219  << G4endl;
220  }
221 }
222 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
223 
224 #include "G4PhysicsListHelper.hh"
225 #include "G4Decay.hh"
226 
228 {
230 
231  // Decay Process
232  //
233  G4Decay* fDecayProcess = new G4Decay();
234 
236  particleIterator->reset();
237  while( (*particleIterator)() ){
238  G4ParticleDefinition* particle = particleIterator->value();
239  if (fDecayProcess->IsApplicable(*particle))
240  ph->RegisterProcess(fDecayProcess, particle);
241  }
242 }
243 
244 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
245 
246 #include "G4PhysicsListHelper.hh"
247 #include "G4RadioactiveDecay.hh"
248 #include "G4GenericIon.hh"
249 #include "G4NuclideTable.hh"
250 
252 {
253  G4RadioactiveDecay* radioactiveDecay = new G4RadioactiveDecay();
254 
255  radioactiveDecay->SetARM(true); //Atomic Rearangement
256 
258  ph->RegisterProcess(radioactiveDecay, G4GenericIon::GenericIon());
259 
260  // mandatory for G4NuclideTable
261  //
263 }
264 
265 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
266 
267 #include "G4ProcessManager.hh"
268 #include "StepMax.hh"
269 
271 {
272  // Step limitation seen as a process
273  StepMax* stepMaxProcess = new StepMax();
274 
276  particleIterator->reset();
277  while ((*particleIterator)()){
278  G4ParticleDefinition* particle = particleIterator->value();
279  G4ProcessManager* pmanager = particle->GetProcessManager();
280 
281  if (stepMaxProcess->IsApplicable(*particle))
282  {
283  pmanager ->AddDiscreteProcess(stepMaxProcess);
284  }
285  }
286 }
287 
288 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
289 
290 #include "G4Material.hh"
291 
293 {
296  const G4MaterialCutsCouple* couple = lBox->GetMaterialCutsCouple();
297  const G4Material* currMat = lBox->GetMaterial();
298 
300  G4double cut;
301  part = particleTable->FindParticle("e-");
302  cut = G4LossTableManager::Instance()->GetRange(part,val,couple);
303  G4cout << "material : " << currMat->GetName() << G4endl;
304  G4cout << "particle : " << part->GetParticleName() << G4endl;
305  G4cout << "energy : " << G4BestUnit(val,"Energy") << G4endl;
306  G4cout << "range : " << G4BestUnit(cut,"Length") << G4endl;
307 }
308 
309 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
310 
311 #include "G4Gamma.hh"
312 #include "G4Electron.hh"
313 #include "G4Proton.hh"
314 #include "G4GenericIon.hh"
315 
316 G4VProcess* PhysicsList::GetProcess(const G4String& processName) const
317 {
319  G4ProcessVector* procList = particle->GetProcessManager()->GetProcessList();
320  G4int nbProc = particle->GetProcessManager()->GetProcessListLength();
321  for (G4int k=0; k<nbProc; k++) {
322  G4VProcess* process = (*procList)[k];
323  if (process->GetProcessName() == processName) return process;
324  }
325  return nullptr;
326 }
327 
328 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:93
const XML_Char * name
Definition: expat.h:151
static void ConstructParticle()
G4LogicalVolume * GetLogicalVolume() const
static G4ParticleTable * GetParticleTable()
G4int GetProcessListLength() const
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
static constexpr double mm
Definition: G4SIunits.hh:115
#define G4endl
Definition: G4ios.hh:61
static void ConstructParticle()
static void ConstructParticle()
const G4String & GetParticleName() const
virtual G4bool IsApplicable(const G4ParticleDefinition &) override
Definition: G4Decay.cc:89
void SetEnergyRange(G4double lowedge, G4double highedge)
G4ProcessVector * GetProcessList() const
G4Material * GetMaterial() const
const G4String & GetName() const
Definition: G4Material.hh:179
G4VProcess * GetProcess(const G4String &) const
virtual void ConstructProcess()=0
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
double G4double
Definition: G4Types.hh:76
TString part[npart]
Definition: Style.C:32
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
PhysicsList()
Implementation of the PhysicsList class.
const G4String & GetProcessName() const
Definition: G4VProcess.hh:411
static constexpr double eV
Definition: G4SIunits.hh:215
static void ConstructParticle()
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
void SetDefaultCutValue(G4double newCutValue)
static void ConstructParticle()
void SetARM(G4bool arm)
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4ProductionCutsTable * GetProductionCutsTable()
void SetThresholdOfHalfLife(G4double)
static G4PhysicsListHelper * GetPhysicsListHelper()
int G4int
Definition: G4Types.hh:78
virtual G4bool IsApplicable(const G4ParticleDefinition &)
G4ProcessManager * GetProcessManager() const
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
G4GLOB_DLL std::ostream G4cout
static G4LossTableManager * Instance()
static constexpr double picosecond
Definition: G4SIunits.hh:161
static constexpr double GeV
Definition: G4SIunits.hh:217
static G4NuclideTable * GetInstance()
Simple detector construction with a box volume placed in a world.
void SetVerboseLevel(G4int value)