Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
extended/electromagnetic/TestEm15/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 110439 2018-05-23 11:24:51Z gcosmo $
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 
34 #include "SteppingAction.hh"
35 #include "DetectorConstruction.hh"
36 #include "RunAction.hh"
37 #include "HistoManager.hh"
38 #include "G4ParticleTypes.hh"
39 
40 #include "G4RunManager.hh"
41 
42 #include <G4ThreeVector.hh>
43 #include <G4RotationMatrix.hh>
44 
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
46 
48  RunAction* RuAct)
49 :G4UserSteppingAction(),fDetector(det), fRunAction(RuAct)
50 { }
51 
52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53 
55 { }
56 
57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58 
60 {
61  G4StepPoint* prePoint = aStep->GetPreStepPoint();
62 
63  // if World --> return
64  if(prePoint->GetTouchableHandle()->GetVolume()==fDetector->GetWorld()) return;
65 
66  // here we enter in the absorber Box
67  // tag the event to be killed anyway after this step
68  //
70 
71  //count processes and keep only Multiple Scattering or gamma converion
72  //
73  G4StepPoint* endPoint = aStep->GetPostStepPoint();
74  G4String procName = endPoint->GetProcessDefinedStep()->GetProcessName();
75  fRunAction->CountProcesses(procName);
76 
77  if (procName == "msc" || procName == "muMsc" || procName == "stepMax") {
78 
79  //below, only multiple Scattering happens
80  //
81  G4ThreeVector position = endPoint->GetPosition();
82  G4ThreeVector direction = endPoint->GetMomentumDirection();
83 
84  G4double truePathLength = aStep->GetStepLength();
85  G4double geomPathLength = position.x() + 0.5*fDetector->GetBoxSize();
86  G4double ratio = geomPathLength/truePathLength;
87  fRunAction->SumPathLength(truePathLength,geomPathLength);
88  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
89  analysisManager->FillH1(1,truePathLength);
90  analysisManager->FillH1(2,geomPathLength);
91  analysisManager->FillH1(3,ratio);
92 
93  G4double yend = position.y(), zend = position.z();
94  G4double lateralDisplacement = std::sqrt(yend*yend + zend*zend);
95  fRunAction->SumLateralDisplacement(lateralDisplacement);
96  analysisManager->FillH1(4,lateralDisplacement);
97 
98  G4double psi = std::atan(lateralDisplacement/geomPathLength);
99  fRunAction->SumPsi(psi);
100  analysisManager->FillH1(5,psi);
101 
102  G4double xdir = direction.x(), ydir = direction.y(), zdir = direction.z();
103  G4double tetaPlane = std::atan2(ydir, xdir);
104  fRunAction->SumTetaPlane(tetaPlane);
105  analysisManager->FillH1(6,tetaPlane);
106  tetaPlane = std::atan2(zdir, xdir);
107  fRunAction->SumTetaPlane(tetaPlane);
108  analysisManager->FillH1(6,tetaPlane);
109 
110  G4double phiPos = std::atan2(zend, yend);
111  analysisManager->FillH1(7,phiPos);
112  G4double phiDir = std::atan2(zdir, ydir);
113  analysisManager->FillH1(8,phiDir);
114 
115  G4double phiCorrel = 0.;
116  if (lateralDisplacement > 0.)
117  phiCorrel = (yend*ydir + zend*zdir)/lateralDisplacement;
118  fRunAction->SumPhiCorrel(phiCorrel);
119  analysisManager->FillH1(9,phiCorrel);
120  } else if (procName == "conv" ) {
121 
122  // gamma conversion
123 
124  G4StepPoint* PrePoint = aStep->GetPreStepPoint();
125  G4double EGamma = PrePoint->GetTotalEnergy();
126  G4ThreeVector PGamma = PrePoint->GetMomentum();
127  G4ThreeVector PolaGamma = PrePoint->GetPolarization();
128 
129  G4double Eplus=-1;
130  // G4double Eminus=-1;
131  // G4double Erecoil=-1;
132  G4ThreeVector Pplus, Pminus, Precoil;
133  //G4int recPDG;
134 
135  const G4TrackVector* secondary = fpSteppingManager->GetSecondary();
136 
137  for (size_t lp=0; lp< std::min((*secondary).size(),size_t(2) ); lp++) {
138  if ((*secondary)[lp]->GetDefinition()==G4Electron::ElectronDefinition()) {
139  // Eminus = (*secondary)[lp]->GetTotalEnergy();
140  Pminus = (*secondary)[lp]->GetMomentum();
141  } //else {
142  if ((*secondary)[lp]->GetDefinition()==G4Positron::PositronDefinition()) {
143  Eplus = (*secondary)[lp]->GetTotalEnergy();
144  Pplus = (*secondary)[lp]->GetMomentum();
145  }
146  }
147 
148  if ( (*secondary).size() >= 3 ) {
149  // Erecoil = (*secondary)[2]->GetTotalEnergy();
150  Precoil = (*secondary)[2]->GetMomentum();
151  // recPDG = (*secondary)[2]->GetDynamicParticle()->GetPDGcode();
152  } else {
153  // Erecoil = 0.0;
154  Precoil = G4ThreeVector();
155  // recPDG = 0;
156  }
157 
158  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
159 
160  // Fill Histograms
161 
162  G4ThreeVector gammadir = PGamma.unit(); // gamma direction
163 
164  G4ThreeVector z = gammadir;
165  G4ThreeVector x(1.,0.,0.);
166 
167  // pola perpendicular to direction
168  if ( PolaGamma.mag() != 0.0 ) {
169  x = PolaGamma.unit();
170  } else { // Pola = 0 case
171  // Direction (z) is unitary vector
172  // (projection to plane) p_proj = p - (p o d)/(d o d) x d
173  if ( x.howOrthogonal(z) != 0) {
174  x = x - x.dot(z) * z;
175  }
176  if (x.mag() != 0.0 ) {
177  x = x.unit();
178  } else {
179  x.set(0.0,0.0,1.0);
180  }
181  }
182 
183  G4ThreeVector y = z;
184  y = y.cross(x);
185 
186  G4RotationMatrix GtoW(x,y,z); // from gamma ref. sys. to World
187  G4RotationMatrix WtoG = inverseOf(GtoW); // from World to gamma ref. sys.
188 
189 
190  G4double angleE = Pplus.angle(Pminus) * EGamma;
191  analysisManager->FillH1(10,angleE);
192 
193  analysisManager->FillH1(11,std::log10(Precoil.mag()));
194  //analysisManager->FillH1(12,Precoil.rotateUz(gammadir).phi());
195  analysisManager->FillH1(12,Precoil.transform(WtoG).phi());
196 
197  // G4double phiPlus = Pplus.rotateUz(gammadir).phi();
198  // G4double phiMinus = Pminus.rotateUz(gammadir).phi();
199  G4double phiPlus = Pplus.transform(WtoG).phi();
200  G4double phiMinus = Pminus.transform(WtoG).phi();
201  analysisManager->FillH1(13,phiPlus);
202  analysisManager->FillH1(14,std::cos(phiPlus + phiMinus) * -2.0);
203  analysisManager->FillH1(15,Eplus/EGamma);
204 
205  //G4double phiPola = PolaGamma.rotateUz(gammadir).phi();
206  G4double phiPola = PolaGamma.transform(WtoG).phi();
207  analysisManager->FillH1(16, phiPola);
208  }
209 }
210 
211 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Float_t x
Definition: compare.C:6
Hep3Vector & transform(const HepRotation &)
Definition: ThreeVectorR.cc:24
G4double GetTotalEnergy() const
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:80
const DetectorConstruction * fDetector
CLHEP::Hep3Vector G4ThreeVector
G4StepPoint * GetPreStepPoint() const
Float_t y
Definition: compare.C:6
Double_t z
double z() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
void SumPathLength(G4double truepl, G4double geompl)
const G4TouchableHandle & GetTouchableHandle() const
double angle(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76
G4double GetStepLength() const
const G4ThreeVector & GetPosition() const
const G4TrackVector * GetSecondary() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:411
const G4ThreeVector & GetMomentumDirection() const
G4CsvAnalysisManager G4AnalysisManager
Definition: g4csv_defs.hh:77
HepBoost inverseOf(const HepBoost &lt)
Definition: G4Step.hh:76
G4StepPoint * GetPostStepPoint() const
void UserSteppingAction(const G4Step *)
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
Hep3Vector cross(const Hep3Vector &) const
Hep3Vector unit() const
const G4ThreeVector & GetPolarization() const
double mag() const
G4SteppingManager * fpSteppingManager
double phi() const
std::vector< G4Track * > G4TrackVector
double x() const
virtual void AbortEvent()
static G4Positron * PositronDefinition()
Definition: G4Positron.cc:89
G4ThreeVector GetMomentum() const
double y() const
const G4VProcess * GetProcessDefinedStep() const
Simple detector construction with a box volume placed in a world.
T min(const T t1, const T t2)
brief Return the smallest of the two arguments