Geant4-11
G4SamplingPostStepAction.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//
26//
27//
28// ----------------------------------------------------------------------
29// GEANT 4 class source file
30//
31// G4SamplingPostStepAction.cc
32//
33// ----------------------------------------------------------------------
34
36#include "G4Track.hh"
37#include "G4ParticleChange.hh"
39#include "G4Nsplit_Weight.hh"
40#include "G4VTrackTerminator.hh"
41#include <sstream>
42
44G4SamplingPostStepAction(const G4VTrackTerminator &TrackTerminator)
45 : fTrackTerminator(TrackTerminator)
46{
47}
48
50{
51}
52
54 G4ParticleChange *aParticleChange,
55 const G4Nsplit_Weight &nw)
56{
57 // evaluate results from sampler
58 if (nw.fN>1)
59 {
60 // split track
61 Split(aTrack, nw, aParticleChange);
62 }
63 else if (nw.fN==1)
64 {
65 // don't split, but weight may be changed !
66 aParticleChange->ProposeWeight(nw.fW);
67 }
68 else if (nw.fN==0)
69 {
70 // kill track
72 }
73 else
74 {
75 // wrong answer
76 std::ostringstream os;
77 os << "Sampler returned nw = "
78 << nw
79 << "\n";
80 G4String msg = os.str();
81
82 G4Exception("G4SamplingPostStepAction::DoIt()",
83 "InvalidCondition", FatalException, msg);
84 }
85}
86
88 const G4Nsplit_Weight &nw,
89 G4ParticleChange *aParticleChange)
90{
91 aParticleChange->ProposeWeight(nw.fW);
92 aParticleChange->SetNumberOfSecondaries(nw.fN-1);
93
94 for (G4int i=1;i<nw.fN;i++)
95 {
96 G4Track *ptrack = new G4Track(aTrack);
97
98 // ptrack->SetCreatorProcess(aTrack.GetCreatorProcess());
99 ptrack->SetWeight(nw.fW);
100
101 if (ptrack->GetMomentumDirection() != aTrack.GetMomentumDirection())
102 {
103 G4Exception("G4SamplingPostStepAction::Split()", "InvalidCondition",
104 FatalException, "Track with same momentum !");
105 }
106 aParticleChange->AddSecondary(ptrack);
107 }
108 return;
109}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
int G4int
Definition: G4Types.hh:85
void AddSecondary(G4Track *aSecondary)
G4SamplingPostStepAction(const G4VTrackTerminator &TrackTerminator)
void DoIt(const G4Track &aTrack, G4ParticleChange *aParticleChange, const G4Nsplit_Weight &nw)
void Split(const G4Track &aTrack, const G4Nsplit_Weight &nw, G4ParticleChange *aParticleChange)
const G4VTrackTerminator & fTrackTerminator
void SetWeight(G4double aValue)
const G4ThreeVector & GetMomentumDirection() const
void ProposeWeight(G4double finalWeight)
void SetNumberOfSecondaries(G4int totSecondaries)
virtual void KillTrack() const =0