Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/electromagnetic/TestEm15/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/TestEm15/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 "HistoManager.hh"
 37 #include "RunAction.hh"
 38 
 39 #include "G4ParticleTypes.hh"
 40 #include "G4RunManager.hh"
 41 
 42 #include <G4RotationMatrix.hh>
 43 #include <G4ThreeVector.hh>
 44 
 45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 46 
 47 SteppingAction::SteppingAction(DetectorConstruction* det, RunAction* RuAct)
 48   : G4UserSteppingAction(), fDetector(det), fRunAction(RuAct)
 49 {}
 50 
 51 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 52 
 53 SteppingAction::~SteppingAction() {}
 54 
 55 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 56 
 57 void SteppingAction::UserSteppingAction(const G4Step* aStep)
 58 {
 59   G4StepPoint* prePoint = aStep->GetPreStepPoint();
 60 
 61   // if World --> return
 62   if (prePoint->GetTouchableHandle()->GetVolume() == fDetector->GetWorld()) return;
 63 
 64   // here we enter in the absorber Box
 65   // tag the event to be killed anyway after this step
 66   //
 67   G4RunManager::GetRunManager()->AbortEvent();
 68 
 69   // count processes and keep only Multiple Scattering or gamma converion
 70   //
 71   G4StepPoint* endPoint = aStep->GetPostStepPoint();
 72   G4String procName = endPoint->GetProcessDefinedStep()->GetProcessName();
 73   fRunAction->CountProcesses(procName);
 74 
 75   if (procName == "msc" || procName == "muMsc" || procName == "stepMax") {
 76     // below, only multiple Scattering happens
 77     //
 78     G4ThreeVector position = endPoint->GetPosition();
 79     G4ThreeVector direction = endPoint->GetMomentumDirection();
 80 
 81     G4double truePathLength = aStep->GetStepLength();
 82     G4double geomPathLength = position.x() + 0.5 * fDetector->GetBoxSize();
 83     G4double ratio = geomPathLength / truePathLength;
 84     fRunAction->SumPathLength(truePathLength, geomPathLength);
 85     G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
 86     analysisManager->FillH1(1, truePathLength);
 87     analysisManager->FillH1(2, geomPathLength);
 88     analysisManager->FillH1(3, ratio);
 89 
 90     G4double yend = position.y(), zend = position.z();
 91     G4double lateralDisplacement = std::sqrt(yend * yend + zend * zend);
 92     fRunAction->SumLateralDisplacement(lateralDisplacement);
 93     analysisManager->FillH1(4, lateralDisplacement);
 94 
 95     G4double psi = std::atan(lateralDisplacement / geomPathLength);
 96     fRunAction->SumPsi(psi);
 97     analysisManager->FillH1(5, psi);
 98 
 99     G4double xdir = direction.x(), ydir = direction.y(), zdir = direction.z();
100     G4double tetaPlane = std::atan2(ydir, xdir);
101     fRunAction->SumTetaPlane(tetaPlane);
102     analysisManager->FillH1(6, tetaPlane);
103     tetaPlane = std::atan2(zdir, xdir);
104     fRunAction->SumTetaPlane(tetaPlane);
105     analysisManager->FillH1(6, tetaPlane);
106 
107     G4double phiPos = std::atan2(zend, yend);
108     analysisManager->FillH1(7, phiPos);
109     G4double phiDir = std::atan2(zdir, ydir);
110     analysisManager->FillH1(8, phiDir);
111 
112     G4double phiCorrel = 0.;
113     if (lateralDisplacement > 0.) phiCorrel = (yend * ydir + zend * zdir) / lateralDisplacement;
114     fRunAction->SumPhiCorrel(phiCorrel);
115     analysisManager->FillH1(9, phiCorrel);
116   }
117   else if (procName == "conv" || procName == "GammaToMuPair") {
118     // gamma conversion
119 
120     G4StepPoint* PrePoint = aStep->GetPreStepPoint();
121     G4double EGamma = PrePoint->GetTotalEnergy();
122     G4ThreeVector PGamma = PrePoint->GetMomentum();
123     G4ThreeVector PolaGamma = PrePoint->GetPolarization();
124 
125     G4double Eplus = -1;
126     G4ThreeVector Pplus, Pminus, Precoil;
127 
128     const G4TrackVector* secondary = fpSteppingManager->GetSecondary();
129 
130     const size_t Nsecondaries = (*secondary).size();
131 
132     // No conversion , E < threshold
133     if (Nsecondaries == 0) return;
134 
135     for (size_t lp = 0; lp < std::min(Nsecondaries, size_t(2)); lp++) {
136       if (((*secondary)[lp]->GetDefinition() == G4Electron::Definition())
137           || ((*secondary)[lp]->GetDefinition() == G4MuonMinus::Definition()))
138       {
139         Pminus = (*secondary)[lp]->GetMomentum();
140       }
141       if (((*secondary)[lp]->GetDefinition() == G4Positron::Definition())
142           || ((*secondary)[lp]->GetDefinition() == G4MuonPlus::Definition()))
143       {
144         Eplus = (*secondary)[lp]->GetTotalEnergy();
145         Pplus = (*secondary)[lp]->GetMomentum();
146       }
147     }
148 
149     if (Nsecondaries >= 3) {
150       Precoil = (*secondary)[2]->GetMomentum();
151     }
152 
153     G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
154 
155     // Fill Histograms
156 
157     G4ThreeVector z = PGamma.unit();  // gamma direction
158     G4ThreeVector x(1., 0., 0.);
159 
160     // pola perpendicular to direction
161 
162     if (PolaGamma.mag() != 0.0) {
163       x = PolaGamma.unit();
164     }
165     else {  // Pola = 0 case
166       x = z.orthogonal().unit();
167     }
168 
169     G4ThreeVector y = z;
170     y = y.cross(x);
171 
172     G4RotationMatrix GtoW(x, y, z);  // from  gamma ref. sys. to World
173     G4RotationMatrix WtoG = inverseOf(GtoW);  // from World to gamma ref. sys.
174 
175     G4double angleE = Pplus.angle(Pminus) * EGamma;
176     analysisManager->FillH1(10, angleE);
177 
178     if (Nsecondaries >= 3) {
179       // recoil returned
180       analysisManager->FillH1(11, std::log10(Precoil.mag()));
181       analysisManager->FillH1(12, Precoil.transform(WtoG).phi());
182     }
183     G4double phiPlus = Pplus.transform(WtoG).phi();
184     G4double phiMinus = Pminus.transform(WtoG).phi();
185     analysisManager->FillH1(13, phiPlus);
186     analysisManager->FillH1(14, std::cos(phiPlus + phiMinus) * -2.0);
187     analysisManager->FillH1(15, Eplus / EGamma);
188 
189     G4double phiPola = PolaGamma.transform(WtoG).phi();
190     analysisManager->FillH1(16, phiPola);
191   }
192 }
193 
194 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
195