Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
extended/medical/fanoCavity/src/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 /// \file medical/fanoCavity/src/SteppingAction.cc
27 /// \brief Implementation of the SteppingAction class
28 //
29 // $Id: SteppingAction.cc 73010 2013-08-15 08:46:58Z 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 "TrackingAction.hh"
38 #include "HistoManager.hh"
39 
40 #include "G4SteppingManager.hh"
41 #include "G4Gamma.hh"
42 #include "G4UnitsTable.hh"
43 
44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
45 
47  TrackingAction* TrAct)
48 :fDetector(det), fRunAction(RuAct), fTrackAction(TrAct),
49  fWall(0), fCavity(0)
50 {
51  first = true;
52  fTrackSegm = 0.;
53  fDirectionIn = G4ThreeVector(0.,0.,0.);
54 }
55 
56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57 
59 { }
60 
61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
62 
64 {
65  //get fDetector pointers
66  if (first) {
67  fWall = fDetector->GetWall();
68  fCavity = fDetector->GetCavity();
69  first = false;
70  }
71 
72  //histograms
73  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
74 
75  //get volume
76  //
77  G4StepPoint* point1 = step->GetPreStepPoint();
78  G4VPhysicalVolume* volume = point1->GetTouchableHandle()->GetVolume();
79 
80  // count processes
81  //
82  G4StepPoint* point2 = step->GetPostStepPoint();
83  const G4VProcess* process = point2->GetProcessDefinedStep();
84  if (process) fRunAction->CountProcesses(process->GetProcessName());
85 
86  //energy deposit in cavity
87  //
88  if (volume == fCavity) {
89  G4double edep = step->GetTotalEnergyDeposit();
90  if (edep > 0.) fTrackAction->AddEdepCavity(edep);
91  }
92 
93  //keep only charged particles
94  //
95  if (step->GetTrack()->GetDefinition() == G4Gamma::Gamma()) return;
96 
97  //step size of charged particles
98  //
99  G4int id;
100  G4double steplen = step->GetStepLength();
101  if (volume == fWall) {fRunAction->StepInWall (steplen); id = 9;}
102  else {fRunAction->StepInCavity(steplen); id = 10;}
103  analysisManager->FillH1(id,steplen);
104 
105  //last step before hitting the cavity
106  //
107  if ((volume == fWall) && (point2->GetStepStatus() == fGeomBoundary)) {
108  fDirectionIn = point1->GetMomentumDirection();
109  }
110 
111  //keep only charged particles within cavity
112  //
113  if (volume == fWall) return;
114 
115  G4double ekin1 = point1->GetKineticEnergy();
116  G4double ekin2 = point2->GetKineticEnergy();
117 
118  //first step in cavity
119  //
120  if (point1->GetStepStatus() == fGeomBoundary) {
121  fTrackSegm = 0.;
122  G4ThreeVector vertex = step->GetTrack()->GetVertexPosition();
123  analysisManager->FillH1(4,vertex.z());
124  fRunAction->FlowInCavity(0,ekin1);
125  analysisManager->FillH1(5,ekin1);
126  if (steplen>0.) {
127  G4ThreeVector directionOut =
128  (point2->GetPosition() - point1->GetPosition()).unit();
129  G4ThreeVector normal = point1->GetTouchableHandle()->GetSolid()
130  ->SurfaceNormal(point1->GetPosition());
131  analysisManager->FillH1(6,std::acos(-fDirectionIn*normal));
132  analysisManager->FillH1(7,std::acos(-directionOut*normal));
133  }
134  }
135 
136  //within cavity
137  //
138  if (step->GetTrack()->GetCurrentStepNumber() == 1) fTrackSegm = 0.;
139  fTrackSegm += steplen;
140  if (ekin2 <= 0.) {
141  fRunAction->AddTrakCavity(fTrackSegm);
142  analysisManager->FillH1(8,fTrackSegm);
143  }
144 
145  //exit cavity
146  //
147  if (point2->GetStepStatus() == fGeomBoundary) {
148  fRunAction->FlowInCavity(1,ekin2);
149  fRunAction->AddTrakCavity(fTrackSegm);
150  analysisManager->FillH1(8,fTrackSegm);
151  }
152 }
153 
154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
155 
156 
G4ParticleDefinition * GetDefinition() const
virtual G4VSolid * GetSolid(G4int depth=0) const
Definition: G4VTouchable.cc:51
CLHEP::Hep3Vector G4ThreeVector
G4double GetStepLength() const
G4StepStatus GetStepStatus() const
void UserSteppingAction(const G4Step *)
int G4int
Definition: G4Types.hh:78
double z() const
G4StepPoint * GetPreStepPoint() const
const G4ThreeVector & GetMomentumDirection() const
G4int GetCurrentStepNumber() const
const G4ThreeVector & GetPosition() const
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0
ExG4HbookAnalysisManager G4AnalysisManager
Definition: g4hbook_defs.hh:46
Definition: G4Step.hh:76
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
const G4ThreeVector & GetVertexPosition() const
G4double GetTotalEnergyDeposit() const
const G4VProcess * GetProcessDefinedStep() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
G4StepPoint * GetPostStepPoint() const
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
G4Track * GetTrack() const
const G4TouchableHandle & GetTouchableHandle() const