Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 //
28 //
29 // $Id: SteppingAction.cc 98762 2016-08-09 14:08:07Z gcosmo $
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 
34 #include "SteppingAction.hh"
35 
36 #include "DetectorConstruction.hh"
37 #include "Run.hh"
38 #include "EventAction.hh"
39 #include "HistoManager.hh"
40 
41 #include "G4Positron.hh"
42 #include "G4RunManager.hh"
43 #include "G4PhysicalConstants.hh"
44 
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
46 
48 :G4UserSteppingAction(),fDetector(det),fEventAct(evt)
49 { }
50 
51 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52 
54 { }
55 
56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57 
59 {
60  //track informations
61  const G4StepPoint* prePoint = aStep->GetPreStepPoint();
62  const G4StepPoint* endPoint = aStep->GetPostStepPoint();
63  const G4ParticleDefinition* particle = aStep->GetTrack()->GetDefinition();
64 
65  //if World, return
66  //
67  G4VPhysicalVolume* volume = prePoint->GetTouchableHandle()->GetVolume();
68  //if sum of absorbers do not fill exactly a layer: check material, not volume.
69  G4Material* mat = volume->GetLogicalVolume()->GetMaterial();
70  if (mat == fDetector->GetWorldMaterial()) return;
71 
72  //here we are in an absorber. Locate it
73  //
74  G4int absorNum = prePoint->GetTouchableHandle()->GetCopyNumber(0);
75  G4int layerNum = prePoint->GetTouchableHandle()->GetCopyNumber(1);
76 
77  //get Run
78  Run* run = static_cast<Run*>(
80 
81  // collect energy deposit taking into account track weight
82  G4double edep = aStep->GetTotalEnergyDeposit()*aStep->GetTrack()->GetWeight();
83 
84  // collect step length of charged particles
85  G4double stepl = 0.;
86  if (particle->GetPDGCharge() != 0.) {
87  stepl = aStep->GetStepLength();
88  run->AddChargedStep();
89  } else { run->AddNeutralStep(); }
90 
91  // G4cout << "Nabs= " << absorNum << " edep(keV)= " << edep << G4endl;
92 
93  // sum up per event
94  fEventAct->SumEnergy(absorNum,edep,stepl);
95 
96  //longitudinal profile of edep per absorber
97  if (edep>0.) {
98  G4AnalysisManager::Instance()->FillH1(kMaxAbsor+absorNum,
99  G4double(layerNum+1), edep);
100  }
101  //energy flow
102  //
103  // unique identificator of layer+absorber
104  G4int Idnow = (fDetector->GetNbOfAbsor())*layerNum + absorNum;
105  G4int plane;
106  //
107  //leaving the absorber ?
108  if (endPoint->GetStepStatus() == fGeomBoundary) {
109  G4ThreeVector position = endPoint->GetPosition();
110  G4ThreeVector direction = endPoint->GetMomentumDirection();
111  G4double sizeYZ = 0.5*fDetector->GetCalorSizeYZ();
112  G4double Eflow = endPoint->GetKineticEnergy();
113  if (particle == G4Positron::Positron()) Eflow += 2*electron_mass_c2;
114  if ((std::abs(position.y()) >= sizeYZ) || (std::abs(position.z()) >= sizeYZ))
115  run->SumLateralEleak(Idnow, Eflow);
116  else if (direction.x() >= 0.) run->SumEnergyFlow(plane=Idnow+1, Eflow);
117  else run->SumEnergyFlow(plane=Idnow, -Eflow);
118  }
119 
125 }
126 
127 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
128 
130 {
131  //Example of Birk attenuation law in organic scintillators.
132  //adapted from Geant3 PHYS337. See MIN 80 (1970) 239-244
133  //
134  G4Material* material = aStep->GetTrack()->GetMaterial();
135  G4double birk1 = material->GetIonisation()->GetBirksConstant();
136  G4double destep = aStep->GetTotalEnergyDeposit();
137  G4double stepl = aStep->GetStepLength();
138  G4double charge = aStep->GetTrack()->GetDefinition()->GetPDGCharge();
139  //
140  G4double response = destep;
141  if (birk1*destep*stepl*charge != 0.)
142  {
143  response = destep/(1. + birk1*destep/stepl);
144  }
145  return response;
146 }
147 
148 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
149 
G4ParticleDefinition * GetDefinition() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
G4Material * GetMaterial() const
G4double GetStepLength() const
double x() const
G4StepStatus GetStepStatus() const
void SumEnergyFlow(G4int plane, G4double Eflow)
Definition: Run.cc:105
Event action class.
Definition: EventAction.hh:45
G4int GetCopyNumber(G4int depth=0) const
void UserSteppingAction(const G4Step *)
int G4int
Definition: G4Types.hh:78
void AddChargedStep()
Definition: Run.cc:119
double z() const
G4double GetBirksConstant() const
string material
Definition: eplot.py:19
G4StepPoint * GetPreStepPoint() const
void AddNeutralStep()
Definition: Run.cc:126
const G4ThreeVector & GetMomentumDirection() const
void SumLateralEleak(G4int cell, G4double Eflow)
Definition: Run.cc:112
const G4ThreeVector & GetPosition() const
const G4int kMaxAbsor
float electron_mass_c2
Definition: hepunit.py:274
Definition: G4Step.hh:76
G4double GetTotalEnergyDeposit() const
G4Material * GetMaterial() const
G4double BirksAttenuation(const G4Step *)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
static G4Positron * Positron()
Definition: G4Positron.cc:94
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
G4StepPoint * GetPostStepPoint() const
double y() const
G4double GetWeight() const
void SumEnergy(G4int k, G4double de, G4double dl)
Definition: EventAction.hh:52
Detector construction class to define materials and geometry.
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
Definition: Run.hh:46
G4Track * GetTrack() const
G4double GetPDGCharge() const
G4Material * GetWorldMaterial()
G4Run * GetNonConstCurrentRun() const
const G4TouchableHandle & GetTouchableHandle() const