source: trunk/examples/extended/electromagnetic/TestEm16/src/SteppingAction.cc@ 1036

Last change on this file since 1036 was 807, checked in by garnier, 17 years ago

update

File size: 4.4 KB
RevLine 
[807]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// $Id: SteppingAction.cc,v 1.5 2007/01/18 09:07:20 hbu Exp $
27// GEANT4 tag $Name: $
28//
29//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
30//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
31
32#include "SteppingAction.hh"
33#include "RunAction.hh"
34#include "G4SteppingManager.hh"
35#include "G4VProcess.hh"
36#include "G4ParticleTypes.hh"
37#include "G4UnitsTable.hh"
38#include "HistoManager.hh"
39
40//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
41
42SteppingAction::SteppingAction(RunAction* RuAct, HistoManager* histo)
43:runAction(RuAct), histoManager(histo)
44{ }
45
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
47
48SteppingAction::~SteppingAction()
49{ }
50
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52
53void SteppingAction::UserSteppingAction(const G4Step* aStep)
54{
55 const G4VProcess* process = aStep->GetPostStepPoint()->GetProcessDefinedStep();
56 if (process == 0) return;
57 G4String processName = process->GetProcessName();
58 static G4int iCalled=0;
59 const G4int nprint=0; // set to 10 to get debug print for first 10 calls
60
61 if (processName == "SynRad")
62 {
63 iCalled++;
64 G4StepPoint* PrePoint = aStep->GetPreStepPoint();
65 G4TrackVector* secondary = fpSteppingManager->GetSecondary();
66
67 //(*secondary)[lp-1] points to the last photon generated
68 //
69 size_t lp=(*secondary).size();
70 if (lp)
71 {
72 G4double Egamma = (*secondary)[lp-1]->GetTotalEnergy();
73 runAction->n_gam_sync++;
74 runAction->e_gam_sync += Egamma;
75 runAction->e_gam_sync2 += Egamma*Egamma;
76 if (Egamma > runAction->e_gam_sync_max) runAction->e_gam_sync_max = Egamma;
77 runAction->lam_gam_sync += aStep->GetStepLength();
78 if (iCalled<nprint)
79 {
80 G4double Eelec = PrePoint->GetTotalEnergy();
81 G4ThreeVector Pelec = PrePoint->GetMomentum();
82 G4ThreeVector PGamma = (*secondary)[lp-1]->GetMomentum();
83 G4bool IsGamma =
84 ((*secondary)[lp-1]->GetDefinition() == G4Gamma::GammaDefinition());
85
86 G4cout << "UserSteppingAction processName=" << process->GetProcessName()
87 << " Step Length=" << std::setw(6)
88 << G4BestUnit(aStep->GetStepLength(),"Length")
89 << " Eelec=" << G4BestUnit(Eelec,"Energy")
90 << " Pelec=" << G4BestUnit(Pelec,"Energy")
91 << " IsGamma=" << IsGamma
92 << " Egamma=" << G4BestUnit(Egamma,"Energy")
93 << " PGamma=" << G4BestUnit(PGamma,"Energy")
94 << " #secondaries lp=" << lp
95 << '\n';
96 }
97 histoManager->FillHisto(1,Egamma);
98 histoManager->FillHisto(2,Egamma,Egamma/keV); // power spectrum : gamma weighted with its energy
99 histoManager->FillHisto(3,aStep->GetStepLength());
100 }
101 }
102
103}
104
105//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.