Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
CCalSteppingAction.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 //
27 // File: CCalSeppingAction.cc
28 // Description: Study profiling during the steps
30 #include <iostream>
31 #include <cmath>
32 
33 #include "CCalSteppingAction.hh"
34 #include "CCalRunAction.hh"
35 #include "CCalAnalysis.hh"
36 
37 #include "G4SystemOfUnits.hh"
38 #include "G4SDManager.hh"
39 #include "G4StepPoint.hh"
40 #include "G4ThreeVector.hh"
41 #include "G4RunManager.hh"
42 
43 #include "CCalAnalysis.hh"
44 
46  runAct(nullptr)
47 {
48  for (int i=0; i<200; i++) {timeDeposit[i] = 0.;}
49  for (int i=0; i<70; i++) {LateralProfile[i] = 0.;}
50 
51 }
52 
53 
55  G4cout << "CCalSteppingAction deleted" << G4endl;
56 }
57 
58 
60 {
61  //thread-local run action
62  if (!runAct)
63  runAct =
64  dynamic_cast<const CCalRunAction*>
66 
67 
69  G4StepPoint* PostStepPoint= aStep->GetPostStepPoint();
70  G4StepPoint* PreStepPoint= aStep->GetPreStepPoint();
71  int TSliceID;
72 
73  if ( PostStepPoint->GetGlobalTime() / nanosecond > 1.0E9 ) TSliceID = 999999999;
74  else TSliceID = static_cast<int>( PostStepPoint->GetGlobalTime() / nanosecond );
75  TSliceID = TSliceID<timeHistoMaxBin ? TSliceID : timeHistoMaxBin-1;
76  timeDeposit[TSliceID] += aStep->GetTotalEnergyDeposit() / GeV;
77 
78  G4ThreeVector HitPoint = 0.5*(PostStepPoint->GetPosition()+
79  PreStepPoint->GetPosition());
80  // Because the beam axis has been defined as the x-axis,
81  // the lateral displacement is given in terms of the y and z positions.
82  double perp = std::sqrt(HitPoint.y()*HitPoint.y()+HitPoint.z()*HitPoint.z());
83  int radialPosition = std::min(69,int(perp/cm));
84  LateralProfile[radialPosition] += aStep->GetTotalEnergyDeposit() / GeV;
85 
86 }
87 
88 
90 
91  G4AnalysisManager* man = G4AnalysisManager::Instance();
92  G4double totalFilledProfileHcal = 0.0;
93 
94  static G4int IDlateralProfile = -1;
95  if (IDlateralProfile < 0)
96  IDlateralProfile = man->GetH1Id("h500");
97 
98  for (G4int i=0; i<70; i++) {
99  man->FillH1(IDlateralProfile+i,LateralProfile[i]);
100 #ifdef debug
101  G4cout << "Fill Profile Hcal histo " << i << " with " << LateralProfile[i] << G4endl;
102 #endif
103  totalFilledProfileHcal += LateralProfile[i];
104  }
105 
106 #ifdef debug
107  G4cout << "CCalAnalysis::InsertLateralProfile: Total filled Profile Hcal"
108  << " histo " << totalFilledProfileHcal << G4endl;
109 #endif
110 
111  static G4int IDTimeHist = -1;
112  if (IDTimeHist < 0)
113  IDTimeHist = man->GetH1Id("h300");
114  G4double totalFilledTimeProfile = 0.0;
115  for (G4int j=0; j<timeHistoMaxBin; j++)
116  {
117  man->FillH1(IDTimeHist+j,timeDeposit[j]);
118 #ifdef debug
119  G4cout << "Fill Time slice histo " << j << " with " << timeDeposit[j] << G4endl;
120 #endif
121  totalFilledTimeProfile += timeDeposit[j];
122 
123  static G4int IDTimeProfile = -1;
124  if (IDTimeProfile < 0)
125  IDTimeProfile = man->GetH1Id("h901");
126  G4double t = j + 0.5;
127  man->FillH1(IDTimeProfile+1,t,timeDeposit[j]);
128 #ifdef debug
129  G4cout << "Fill Time profile histo 1 with " << t << " " << x << G4endl;
130 #endif
131  }
132  #ifdef debug
133  G4cout << "CCalAnalysis::InsertTime: Total filled Time profile histo "
134  << totalFilledTimeProfile << G4endl;
135 #endif
136 
137  int i=0;
138  for (i=0; i<70; i++){LateralProfile[i] = 0.;}
139  for (i=0; i<200; i++){timeDeposit[i] = 0.;}
140 
141  G4cout << " --- End of event --- " << G4endl;
142 
143 }
Float_t x
Definition: compare.C:6
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:80
const G4UserRunAction * GetUserRunAction() const
G4StepPoint * GetPreStepPoint() const
#define G4endl
Definition: G4ios.hh:61
double z() const
G4int GetH1Id(const G4String &name, G4bool warn=true) const
double G4double
Definition: G4Types.hh:76
const G4ThreeVector & GetPosition() const
G4double GetTotalEnergyDeposit() const
G4double GetGlobalTime() const
Definition: G4Step.hh:76
G4StepPoint * GetPostStepPoint() const
static constexpr double nanosecond
Definition: G4SIunits.hh:158
int G4int
Definition: G4Types.hh:78
G4int maxbin() const
static constexpr double cm
Definition: G4SIunits.hh:119
G4GLOB_DLL std::ostream G4cout
const CCalRunAction * runAct
virtual void UserSteppingAction(const G4Step *aStep)
double y() const
G4bool FillH1(G4int id, G4double value, G4double weight=1.0)
static constexpr double GeV
Definition: G4SIunits.hh:217
T min(const T t1, const T t2)
brief Return the smallest of the two arguments