Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
extended/hadronic/Hadr03/src/SteppingAction.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: SteppingAction.cc 105744 2017-08-16 13:13:16Z gcosmo $
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 
34 #include "SteppingAction.hh"
35 #include "Run.hh"
36 #include "HistoManager.hh"
37 
38 #include "G4ParticleTypes.hh"
39 #include "G4RunManager.hh"
40 #include "G4HadronicProcess.hh"
41 
42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
43 
46 { }
47 
48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
49 
51 { }
52 
53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
54 
56 {
57  Run* run
59 
60  // count processes
61  //
62  const G4StepPoint* endPoint = aStep->GetPostStepPoint();
63  G4VProcess* process =
64  const_cast<G4VProcess*>(endPoint->GetProcessDefinedStep());
65  run->CountProcesses(process);
66 
67  // check that an real interaction occured (eg. not a transportation)
68  G4StepStatus stepStatus = endPoint->GetStepStatus();
69  G4bool transmit = (stepStatus==fGeomBoundary || stepStatus==fWorldBoundary);
70  if (transmit) return;
71 
72  //real processes : sum track length
73  //
75  run->SumTrack(stepLength);
76 
77  //energy-momentum balance initialisation
78  //
79  const G4StepPoint* prePoint = aStep->GetPreStepPoint();
80  G4double Q = - prePoint->GetKineticEnergy();
81  G4ThreeVector Pbalance = - prePoint->GetMomentum();
82 
83  //initialisation of the nuclear channel identification
84  //
85  G4ParticleDefinition* particle = aStep->GetTrack()->GetDefinition();
86  G4String partName = particle->GetParticleName();
87  G4String nuclearChannel = partName;
88  G4HadronicProcess* hproc = dynamic_cast<G4HadronicProcess*>(process);
89  const G4Isotope* target = NULL;
90  if (hproc) target = hproc->GetTargetIsotope();
91  G4String targetName = "XXXX";
92  if (target) targetName = target->GetName();
93  nuclearChannel += " + " + targetName + " --> ";
94  if (targetName == "XXXX") run->SetTargetXXX(true);
95 
96  //scattered primary particle (if any)
97  //
98  G4AnalysisManager* analysis = G4AnalysisManager::Instance();
99  G4int ih = 1;
100  if (aStep->GetTrack()->GetTrackStatus() == fAlive) {
101  G4double energy = endPoint->GetKineticEnergy();
102  analysis->FillH1(ih,energy);
103  //
104  G4ThreeVector momentum = endPoint->GetMomentum();
105  Q += energy;
106  Pbalance += momentum;
107  //
108  nuclearChannel += partName + " + ";
109  }
110 
111  //secondaries
112  //
113  const std::vector<const G4Track*>* secondary
114  = aStep->GetSecondaryInCurrentStep();
115  for (size_t lp=0; lp<(*secondary).size(); lp++) {
116  particle = (*secondary)[lp]->GetDefinition();
117  G4String name = particle->GetParticleName();
118  G4String type = particle->GetParticleType();
119  G4double energy = (*secondary)[lp]->GetKineticEnergy();
120  run->ParticleCount(name,energy);
121  //energy spectrum
122  ih = 0;
123  if (particle == G4Gamma::Gamma()) ih = 2;
124  else if (particle == G4Neutron::Neutron()) ih = 3;
125  else if (particle == G4Proton::Proton()) ih = 4;
126  else if (particle == G4Deuteron::Deuteron()) ih = 5;
127  else if (particle == G4Alpha::Alpha()) ih = 6;
128  else if (type == "nucleus") ih = 7;
129  else if (type == "meson") ih = 8;
130  else if (type == "baryon") ih = 9;
131  if (ih > 0) analysis->FillH1(ih,energy);
132  //atomic mass
133  if (type == "nucleus") {
134  G4int A = particle->GetAtomicMass();
135  analysis->FillH1(12, A);
136  }
137  //energy-momentum balance
138  G4ThreeVector momentum = (*secondary)[lp]->GetMomentum();
139  Q += energy;
140  Pbalance += momentum;
141  //count e- from internal conversion together with gamma
142  if (particle == G4Electron::Electron()) particle = G4Gamma::Gamma();
143  //particle flag
144  fParticleFlag[particle]++;
145  }
146 
147  //energy-momentum balance
148  G4double Pbal = Pbalance.mag();
149  run->Balance(Pbal);
150  ih = 10;
151  analysis->FillH1(ih,Q);
152  ih = 11;
153  analysis->FillH1(ih,Pbal);
154 
155  // nuclear channel
156  const G4int kMax = 16;
157  const G4String conver[] = {"0","","2 ","3 ","4 ","5 ","6 ","7 ","8 ","9 ",
158  "10 ","11 ","12 ","13 ","14 ","15 ","16 "};
159  std::map<G4ParticleDefinition*,G4int>::iterator ip;
160  for (ip = fParticleFlag.begin(); ip != fParticleFlag.end(); ip++) {
161  particle = ip->first;
162  G4String name = particle->GetParticleName();
163  G4int nb = ip->second;
164  if (nb > kMax) nb = kMax;
165  G4String Nb = conver[nb];
166  if (particle == G4Gamma::Gamma()) {
167  run->CountGamma(nb);
168  Nb = "N ";
169  name = "gamma or e-";
170  }
171  if (ip != fParticleFlag.begin()) nuclearChannel += " + ";
172  nuclearChannel += Nb + name;
173  }
174 
176  run->CountNuclearChannel(nuclearChannel, Q);
177 
178  fParticleFlag.clear();
179 
180  // kill event after first interaction
181  //
183 }
184 
185 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
186 
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:80
const XML_Char * name
Definition: expat.h:151
const XML_Char * target
Definition: expat.h:268
const G4String & GetName() const
Definition: G4Isotope.hh:88
G4StepPoint * GetPreStepPoint() const
void Balance(G4double)
const std::vector< const G4Track * > * GetSecondaryInCurrentStep() const
Definition: G4Step.cc:193
const G4String & GetParticleName() const
const G4String & GetParticleType() const
void SetTargetXXX(G4bool)
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4ParticleDefinition * GetDefinition() const
void ParticleCount(G4String, G4double)
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4StepStatus GetStepStatus() const
G4double GetStepLength() const
double energy
Definition: plottest35.C:25
double A(double temperature)
G4Track * GetTrack() const
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
G4StepStatus
Definition: G4StepStatus.hh:49
G4CsvAnalysisManager G4AnalysisManager
Definition: g4csv_defs.hh:77
Definition: G4Step.hh:76
G4StepPoint * GetPostStepPoint() const
void UserSteppingAction(const G4Step *)
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4Run * GetNonConstCurrentRun() const
double mag() const
int G4int
Definition: G4Types.hh:78
void CountProcesses(G4String procName)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4int GetAtomicMass() const
const G4Isotope * GetTargetIsotope()
std::map< G4ParticleDefinition *, G4int > fParticleFlag
virtual void AbortEvent()
G4ThreeVector GetMomentum() const
static double Q[]
const G4VProcess * GetProcessDefinedStep() const
void CountNuclearChannel(G4String, G4double)
G4TrackStatus GetTrackStatus() const
void SumTrack(G4double track)