Geant4  10.01.p02
SteppingAction.cc
Go to the documentation of this file.
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 // Please cite the following paper if you use this software
27 // Nucl.Instrum.Meth.B260:20-27, 2007
28 
29 #include "SteppingAction.hh"
30 #include "Analysis.hh"
31 #include "G4SystemOfUnits.hh"
32 
33 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
34 
36 :fRun(run),fDetector(det)
37 {}
38 
39 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
40 
42 {}
43 
44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
45 
47 
48 {
49 
50 G4AnalysisManager* man = G4AnalysisManager::Instance();
51 
52 
53 if ( (step->GetTrack()->GetDynamicParticle()->GetDefinition() ==
55 
56 /*
57 // for doublet
58 
59  && (step->GetPostStepPoint()->GetPosition().z()/mm>-3230.2)
60  && (step->GetPostStepPoint()->GetPosition().z()/mm<-3229.8)
61 */
62 
63 // for triplet and whole line
64 
65  && (step->GetPostStepPoint()->GetPosition().z()/mm>249.99999)
66  && (step->GetPostStepPoint()->GetPosition().z()/mm<250.00001)
67  && (step->GetPreStepPoint()->GetTouchableHandle()->GetVolume()->
68  GetLogicalVolume()->GetName() == fDetector->GetLogicalVol()->GetName())
70  GetLogicalVolume()->GetName() == fDetector->GetLogicalWorld()->GetName())
71  )
72 
73  {
74  fXIn = step->GetPostStepPoint()->GetPosition().x();
75  fYIn = step->GetPostStepPoint()->GetPosition().y();
76  fZIn = step->GetPostStepPoint()->GetPosition().z();
77  fE = step->GetTrack()->GetKineticEnergy();
78 
79  G4ThreeVector angleIn;
80  angleIn = step->GetTrack()->GetMomentumDirection();
81 
82  fThetaIn = std::asin(angleIn[0]/std::sqrt(angleIn[0]
83  *angleIn[0]+angleIn[1]*angleIn[1]+angleIn[2]*angleIn[2]));
84  fPhiIn = std::asin(angleIn[1]/std::sqrt(angleIn[0]
85  *angleIn[0]+angleIn[1]*angleIn[1]+angleIn[2]*angleIn[2]));
86 
87  G4cout << " =>IMAGE : X(microns)=" << fXIn/micrometer
88  <<" Y(microns)="<< fYIn/micrometer << " THETA(mrad)="
89  << (fThetaIn/mrad) << " PHI(mrad)=" << (fPhiIn/mrad) << G4endl;
90  G4cout << G4endl;
91 
92  if (fDetector->GetCoef()==1)
93  {
94  fRun->AddRow();
95  fRun->AddToXVector(fXIn/um);
96  fRun->AddToYVector(fYIn/um);
99  }
100 
101  //Fill ntuple 3
102  man->FillNtupleDColumn(3,0,fXIn/um);
103  man->FillNtupleDColumn(3,1,fYIn/um);
104  man->FillNtupleDColumn(3,2,fThetaIn/mrad);
105  man->FillNtupleDColumn(3,3,fPhiIn/mrad);
106  man->AddNtupleRow(3);
107 
108  }
109 
110 if (fDetector->GetProfile()==1)
111 {
112 
113  if (
116  && (step->GetPreStepPoint()->GetTouchableHandle()
118  && (step->GetPostStepPoint()->GetTouchableHandle()
120  {
121  fXIn = step->GetPostStepPoint()->GetPosition().x();
122  fYIn = step->GetPostStepPoint()->GetPosition().y();
123  fZIn = step->GetPostStepPoint()->GetPosition().z();
124 
125  //Fill ntuple 1
126  man->FillNtupleDColumn(1,0,fXIn/um);
127  man->FillNtupleDColumn(1,1,fYIn/um);
128  man->FillNtupleDColumn(1,2,fZIn/um);
129  man->AddNtupleRow(1);
130  }
131 }
132 
133 if (fDetector->GetGrid()==1)
134 {
135 
136  if (
139  && (step->GetPreStepPoint()->GetTouchableHandle()
141  && (step->GetPostStepPoint()->GetTouchableHandle()
143  {
144  fXIn = step->GetPostStepPoint()->GetPosition().x();
145  fYIn = step->GetPostStepPoint()->GetPosition().y();
146  fE = step->GetTrack()->GetKineticEnergy();
147 
148  //Fill ntuple 2
149  man->FillNtupleDColumn(2,0,fXIn/um);
150  man->FillNtupleDColumn(2,1,fYIn/um);
151  man->FillNtupleDColumn(2,2,fE/MeV);
152  man->AddNtupleRow(2);
153  }
154  }
155 
156 // end
157 }
G4LogicalVolume * GetLogicalWorld()
void AddToXVector(float v)
Definition: RunAction.hh:55
static const double MeV
Definition: G4SIunits.hh:193
G4String GetName() const
CLHEP::Hep3Vector G4ThreeVector
const G4DynamicParticle * GetDynamicParticle() const
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
G4ParticleDefinition * GetDefinition() const
void UserSteppingAction(const G4Step *)
G4StepPoint * GetPreStepPoint() const
G4LogicalVolume * GetLogicalGrid()
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
const G4ThreeVector & GetPosition() const
ExG4HbookAnalysisManager G4AnalysisManager
Definition: g4hbook_defs.hh:61
Definition: G4Step.hh:76
G4LogicalVolume * GetLogicalVol()
static const double mrad
Definition: G4SIunits.hh:131
DetectorConstruction * fDetector
void AddToPhiVector(float v)
Definition: RunAction.hh:58
static const double micrometer
Definition: G4SIunits.hh:90
void AddRow()
Definition: RunAction.hh:53
const G4ThreeVector & GetMomentumDirection() const
G4LogicalVolume * GetLogicalVolume() const
void AddToThetaVector(float v)
Definition: RunAction.hh:57
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
G4StepPoint * GetPostStepPoint() const
void AddToYVector(float v)
Definition: RunAction.hh:56
Detector construction class to demonstrate various ways of placement.
#define G4endl
Definition: G4ios.hh:61
G4Track * GetTrack() const
static const double mm
Definition: G4SIunits.hh:102
RunAction * fRun
const G4TouchableHandle & GetTouchableHandle() const