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

Last change on this file since 1357 was 1342, checked in by garnier, 15 years ago

update ti head

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