Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/electromagnetic/TestEm3/src/SteppingAction.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 electromagnetic/TestEm3/src/SteppingAction.cc
 27 /// \brief Implementation of the SteppingAction class
 28 //
 29 //
 30 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 32 
 33 #include "SteppingAction.hh"
 34 
 35 #include "DetectorConstruction.hh"
 36 #include "EventAction.hh"
 37 #include "HistoManager.hh"
 38 #include "Run.hh"
 39 
 40 #include "G4PhysicalConstants.hh"
 41 #include "G4Positron.hh"
 42 #include "G4RunManager.hh"
 43 
 44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 45 
 46 SteppingAction::SteppingAction(DetectorConstruction* det, EventAction* evt)
 47   : fDetector(det), fEventAct(evt)
 48 {}
 49 
 50 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 51 
 52 void SteppingAction::UserSteppingAction(const G4Step* aStep)
 53 {
 54   // track informations
 55   const G4StepPoint* prePoint = aStep->GetPreStepPoint();
 56 
 57   // if World, return
 58   //
 59   G4VPhysicalVolume* volume = prePoint->GetTouchableHandle()->GetVolume();
 60   // if sum of absorbers do not fill exactly a layer: check material, not volume.
 61   const G4Material* mat = volume->GetLogicalVolume()->GetMaterial();
 62   if (mat == fDetector->GetWorldMaterial()) return;
 63 
 64   const G4StepPoint* endPoint = aStep->GetPostStepPoint();
 65   const G4ParticleDefinition* particle = aStep->GetTrack()->GetDefinition();
 66 
 67   // here we are in an absorber. Locate it
 68   //
 69   G4int absorNum = prePoint->GetTouchableHandle()->GetCopyNumber(0);
 70   G4int layerNum = prePoint->GetTouchableHandle()->GetCopyNumber(1);
 71 
 72   // get Run
 73   Run* run = static_cast<Run*>(G4RunManager::GetRunManager()->GetNonConstCurrentRun());
 74 
 75   // collect energy deposit taking into account track weight
 76   G4double edep = aStep->GetTotalEnergyDeposit() * aStep->GetTrack()->GetWeight();
 77 
 78   // collect step length of charged particles
 79   G4double stepl = 0.;
 80   if (particle->GetPDGCharge() != 0.) {
 81     stepl = aStep->GetStepLength();
 82     run->AddChargedStep();
 83   }
 84   else {
 85     run->AddNeutralStep();
 86   }
 87 
 88   //  G4cout << "Nabs= " << absorNum << "   edep(keV)= " << edep << G4endl;
 89 
 90   // sum up per event
 91   fEventAct->SumEnergy(absorNum, edep, stepl);
 92 
 93   // longitudinal profile of edep per absorber
 94   if (edep > 0.) {
 95     G4AnalysisManager::Instance()->FillH1(kMaxAbsor + absorNum, G4double(layerNum + 1), edep);
 96   }
 97   // energy flow
 98   //
 99   //  unique identificator of layer+absorber
100   G4int Idnow = (fDetector->GetNbOfAbsor()) * layerNum + absorNum;
101   G4int plane;
102   //
103   // leaving the absorber ?
104   if (endPoint->GetStepStatus() == fGeomBoundary) {
105     G4ThreeVector position = endPoint->GetPosition();
106     G4ThreeVector direction = endPoint->GetMomentumDirection();
107     G4double sizeYZ = 0.5 * fDetector->GetCalorSizeYZ();
108     G4double Eflow = endPoint->GetKineticEnergy();
109     if (particle == G4Positron::Positron()) Eflow += 2 * electron_mass_c2;
110     if ((std::abs(position.y()) >= sizeYZ) || (std::abs(position.z()) >= sizeYZ))
111       run->SumLateralEleak(Idnow, Eflow);
112     else if (direction.x() >= 0.)
113       run->SumEnergyFlow(plane = Idnow + 1, Eflow);
114     else
115       run->SumEnergyFlow(plane = Idnow, -Eflow);
116   }
117 
118   ////  example of Birk attenuation
119   /// G4double destep   = aStep->GetTotalEnergyDeposit();
120   /// G4double response = BirksAttenuation(aStep);
121   /// G4cout << " Destep: " << destep/keV << " keV"
122   ///       << " response after Birks: " << response/keV << " keV" << G4endl;
123 }
124 
125 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
126 
127 G4double SteppingAction::BirksAttenuation(const G4Step* aStep)
128 {
129   // Example of Birk attenuation law in organic scintillators.
130   // adapted from Geant3 PHYS337. See MIN 80 (1970) 239-244
131   //
132   const G4Material* material = aStep->GetTrack()->GetMaterial();
133   G4double birk1 = material->GetIonisation()->GetBirksConstant();
134   G4double destep = aStep->GetTotalEnergyDeposit();
135   G4double stepl = aStep->GetStepLength();
136   G4double charge = aStep->GetTrack()->GetDefinition()->GetPDGCharge();
137   //
138   G4double response = destep;
139   if (birk1 * destep * stepl * charge != 0.) {
140     response = destep / (1. + birk1 * destep / stepl);
141   }
142   return response;
143 }
144 
145 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
146