Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
extended/optical/OpNovice2/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 //
26 // $Id: SteppingAction.cc 71007 2013-06-09 16:14:59Z maire $
27 //
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 
33 #include "SteppingAction.hh"
34 //#include "EventAction.hh"
35 #include "HistoManager.hh"
36 #include "TrackInformation.hh"
37 #include "Run.hh"
38 
39 #include "G4Cerenkov.hh"
40 #include "G4Scintillation.hh"
41 #include "G4OpBoundaryProcess.hh"
42 
43 #include "G4Step.hh"
44 #include "G4Track.hh"
45 #include "G4OpticalPhoton.hh"
46 #include "G4Event.hh"
47 #include "G4EventManager.hh"
48 #include "G4SteppingManager.hh"
49 #include "G4RunManager.hh"
50 #include "G4ProcessManager.hh"
51 
52 #include "G4SystemOfUnits.hh"
53 
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57 {}
58 
59 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
61 {}
62 
63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65 {
66  static G4ParticleDefinition* opticalphoton =
68  static G4AnalysisManager* analysisMan = G4AnalysisManager::Instance();
69  Run* run = static_cast<Run*>(
71 
72  G4Track* track = step->GetTrack();
73  G4StepPoint* endPoint = step->GetPostStepPoint();
74  G4StepPoint* startPoint = step->GetPreStepPoint();
75 
76  G4String particleName = track->GetDynamicParticle()->
77  GetParticleDefinition()->GetParticleName();
78 
79  TrackInformation* trackInfo =
81 
82  if (particleName == "opticalphoton") {
83  const G4VProcess* pds = endPoint->GetProcessDefinedStep();
84  if (pds->GetProcessName() == "OpAbsorption") {
85  run->AddOpAbsorption();
86  if (trackInfo->GetIsFirstTankX()) {
87  run->AddOpAbsorptionPrior();
88  }
89  }
90  else if (pds->GetProcessName() == "OpRayleigh") {
91  run->AddRayleigh();
92  }
93 
94  // optical process has endpt on bdry,
95  if (endPoint->GetStepStatus() == fGeomBoundary) {
96 
97  const G4DynamicParticle* theParticle = track->GetDynamicParticle();
98 
99  G4ThreeVector oldMomentumDir = theParticle->GetMomentumDirection();
100 
101  G4ThreeVector m0 = startPoint->GetMomentumDirection();
102  G4ThreeVector m1 = endPoint->GetMomentumDirection();
103 
105 
106  G4ProcessManager* OpManager =
108  G4int MAXofPostStepLoops =
109  OpManager->GetPostStepProcessVector()->entries();
110  G4ProcessVector* postStepDoItVector =
112 
113  if (trackInfo->GetIsFirstTankX()) {
114  G4ThreeVector momdir = endPoint->GetMomentumDirection();
115  G4double px1 = momdir.x();
116  G4double py1 = momdir.y();
117  G4double pz1 = momdir.z();
118  if (px1 < 0.) {
119  analysisMan->FillH1(4, px1);
120  analysisMan->FillH1(5, py1);
121  analysisMan->FillH1(6, pz1);
122  } else if (px1 >= 0.) {
123  analysisMan->FillH1(7, px1);
124  analysisMan->FillH1(8, py1);
125  analysisMan->FillH1(9, pz1);
126  }
127 
128  trackInfo->SetIsFirstTankX(false);
129  run->AddTotalSurface();
130 
131  for (G4int i=0; i<MAXofPostStepLoops; ++i) {
132  G4VProcess* currentProcess = (*postStepDoItVector)[i];
133 
134  G4OpBoundaryProcess* opProc =
135  dynamic_cast<G4OpBoundaryProcess*>(currentProcess);
136  if (opProc) {
137  theStatus = opProc->GetStatus();
138  analysisMan->FillH1(3, theStatus);
139  if (theStatus == Transmission) {
140  run->AddTransmission();
141  }
142  else if (theStatus == FresnelRefraction) {
143  run->AddFresnelRefraction();
144  analysisMan->FillH1(10, px1);
145  analysisMan->FillH1(11, py1);
146  analysisMan->FillH1(12, pz1);
147  }
148  else if (theStatus == FresnelReflection) {
149  run->AddFresnelReflection();
150  }
151  else if (theStatus == TotalInternalReflection) {
153  }
154  else if (theStatus == LambertianReflection) {
156  }
157  else if (theStatus == LobeReflection) {
158  run->AddLobeReflection();
159  }
160  else if (theStatus == SpikeReflection) {
161  run->AddSpikeReflection();
162  }
163  else if (theStatus == BackScattering) {
164  run->AddBackScattering();
165  }
166  else if (theStatus == Absorption) {
167  run->AddAbsorption();
168  }
169  else if (theStatus == Detection) {
170  run->AddDetection();
171  }
172  else if (theStatus == NotAtBoundary) {
173  run->AddNotAtBoundary();
174  }
175  else if (theStatus == SameMaterial) {
176  run->AddSameMaterial();
177  }
178  else if (theStatus == StepTooSmall) {
179  run->AddStepTooSmall();
180  }
181  else if (theStatus == NoRINDEX) {
182  run->AddNoRINDEX();
183  }
184  else if (theStatus == PolishedLumirrorAirReflection) {
186  }
187  else if (theStatus == PolishedLumirrorGlueReflection) {
189  }
190  else if (theStatus == PolishedAirReflection) {
192  }
193  else if (theStatus == PolishedTeflonAirReflection) {
195  }
196  else if (theStatus == PolishedTiOAirReflection) {
198  }
199  else if (theStatus == PolishedTyvekAirReflection) {
201  }
202  else if (theStatus == PolishedVM2000AirReflection) {
204  }
205  else if (theStatus == PolishedVM2000GlueReflection) {
207  }
208  else if (theStatus == EtchedLumirrorAirReflection) {
210  }
211  else if (theStatus == EtchedLumirrorGlueReflection) {
213  }
214  else if (theStatus == EtchedAirReflection) {
215  run->AddEtchedAirReflection();
216  }
217  else if (theStatus == EtchedTeflonAirReflection) {
219  }
220  else if (theStatus == EtchedTiOAirReflection) {
222  }
223  else if (theStatus == EtchedTyvekAirReflection) {
225  }
226  else if (theStatus == EtchedVM2000AirReflection) {
228  }
229  else if (theStatus == EtchedVM2000GlueReflection) {
231  }
232  else if (theStatus == GroundLumirrorAirReflection) {
234  }
235  else if (theStatus == GroundLumirrorGlueReflection) {
237  }
238  else if (theStatus == GroundAirReflection) {
239  run->AddGroundAirReflection();
240  }
241  else if (theStatus == GroundTeflonAirReflection) {
243  }
244  else if (theStatus == GroundTiOAirReflection) {
246  }
247  else if (theStatus == GroundTyvekAirReflection) {
249  }
250  else if (theStatus == GroundVM2000AirReflection) {
252  }
253  else if (theStatus == GroundVM2000GlueReflection) {
255  }
256  else if (theStatus == Dichroic) {
257  run->AddDichroic();
258  }
259 
260  else {
261  G4cout << "theStatus: " << theStatus
262  << " was none of the above." << G4endl;
263  }
264 
265  }
266  }
267  }
268  }
269  }
270 
271  else { // particle != opticalphoton
272  const std::vector<const G4Track*>* secondaries =
274 
275  for (auto sec : *secondaries) {
276  if (sec->GetDynamicParticle()->GetParticleDefinition() == opticalphoton){
277  if (sec->GetCreatorProcess()->GetProcessName().compare("Cerenkov")==0){
278  G4double en = sec->GetKineticEnergy();
279  run->AddCerenkovEnergy(en);
280  run->AddCerenkov();
281  G4AnalysisManager::Instance()->FillH1(1, en);
282  }
283  else if (sec->GetCreatorProcess()
284  ->GetProcessName().compare("Scintillation") == 0) {
285  G4double en = sec->GetKineticEnergy();
286  run->AddScintillationEnergy(en);
287  run->AddScintillation();
288  G4AnalysisManager::Instance()->FillH1(2, en);
289  }
290  }
291  }
292  }
293 
294  return;
295 
296 }
297 
298 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
void AddPolishedLumirrorAirReflection(void)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:80
void AddScintillationEnergy(G4double en)
G4int entries() const
G4OpBoundaryProcessStatus
G4StepPoint * GetPreStepPoint() const
G4bool GetIsFirstTankX() const
#define G4endl
Definition: G4ios.hh:61
const G4ThreeVector & GetMomentumDirection() const
const std::vector< const G4Track * > * GetSecondaryInCurrentStep() const
Definition: G4Step.cc:193
double z() const
void AddGroundLumirrorGlueReflection(void)
void AddPolishedTyvekAirReflection(void)
double G4double
Definition: G4Types.hh:76
G4StepStatus GetStepStatus() const
static G4OpticalPhoton * OpticalPhoton()
G4Track * GetTrack() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:411
const G4ThreeVector & GetMomentumDirection() const
static G4OpticalPhoton * OpticalPhotonDefinition()
G4CsvAnalysisManager G4AnalysisManager
Definition: g4csv_defs.hh:77
Definition: G4Step.hh:76
G4StepPoint * GetPostStepPoint() const
void UserSteppingAction(const G4Step *)
G4Run * GetNonConstCurrentRun() const
G4ProcessVector * GetPostStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
void AddCerenkovEnergy(G4double en)
int G4int
Definition: G4Types.hh:78
void AddEtchedLumirrorGlueReflection(void)
void SetIsFirstTankX(G4bool b)
G4ProcessManager * GetProcessManager() const
G4GLOB_DLL std::ostream G4cout
double x() const
G4VUserTrackInformation * GetUserInformation() const
void AddPolishedTeflonAirReflection(void)
double y() const
const G4VProcess * GetProcessDefinedStep() const
void AddPolishedLumirrorGlueReflection(void)
void AddPolishedVM2000AirReflection(void)
const G4DynamicParticle * GetDynamicParticle() const
G4OpBoundaryProcessStatus GetStatus() const
Definition of the TrackInformation class.