Geant4-11
Public Member Functions | Private Member Functions | Private Attributes
G4SamplingPostStepAction Class Reference

#include <G4SamplingPostStepAction.hh>

Public Member Functions

void DoIt (const G4Track &aTrack, G4ParticleChange *aParticleChange, const G4Nsplit_Weight &nw)
 
 G4SamplingPostStepAction (const G4VTrackTerminator &TrackTerminator)
 
 ~G4SamplingPostStepAction ()
 

Private Member Functions

void Split (const G4Track &aTrack, const G4Nsplit_Weight &nw, G4ParticleChange *aParticleChange)
 

Private Attributes

const G4VTrackTerminatorfTrackTerminator
 

Detailed Description

Definition at line 48 of file G4SamplingPostStepAction.hh.

Constructor & Destructor Documentation

◆ G4SamplingPostStepAction()

G4SamplingPostStepAction::G4SamplingPostStepAction ( const G4VTrackTerminator TrackTerminator)
explicit

Definition at line 43 of file G4SamplingPostStepAction.cc.

45 : fTrackTerminator(TrackTerminator)
46{
47}
const G4VTrackTerminator & fTrackTerminator

◆ ~G4SamplingPostStepAction()

G4SamplingPostStepAction::~G4SamplingPostStepAction ( )

Definition at line 49 of file G4SamplingPostStepAction.cc.

50{
51}

Member Function Documentation

◆ DoIt()

void G4SamplingPostStepAction::DoIt ( const G4Track aTrack,
G4ParticleChange aParticleChange,
const G4Nsplit_Weight nw 
)

Definition at line 53 of file G4SamplingPostStepAction.cc.

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}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
void Split(const G4Track &aTrack, const G4Nsplit_Weight &nw, G4ParticleChange *aParticleChange)
void ProposeWeight(G4double finalWeight)
virtual void KillTrack() const =0

References FatalException, G4Nsplit_Weight::fN, fTrackTerminator, G4Nsplit_Weight::fW, G4Exception(), G4VTrackTerminator::KillTrack(), G4VParticleChange::ProposeWeight(), and Split().

Referenced by G4ImportanceProcess::PostStepDoIt(), and G4WeightWindowProcess::PostStepDoIt().

◆ Split()

void G4SamplingPostStepAction::Split ( const G4Track aTrack,
const G4Nsplit_Weight nw,
G4ParticleChange aParticleChange 
)
private

Definition at line 87 of file G4SamplingPostStepAction.cc.

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}
int G4int
Definition: G4Types.hh:85
void AddSecondary(G4Track *aSecondary)
void SetWeight(G4double aValue)
const G4ThreeVector & GetMomentumDirection() const
void SetNumberOfSecondaries(G4int totSecondaries)

References G4ParticleChange::AddSecondary(), FatalException, G4Nsplit_Weight::fN, G4Nsplit_Weight::fW, G4Exception(), G4Track::GetMomentumDirection(), G4VParticleChange::ProposeWeight(), G4VParticleChange::SetNumberOfSecondaries(), and G4Track::SetWeight().

Referenced by DoIt().

Field Documentation

◆ fTrackTerminator

const G4VTrackTerminator& G4SamplingPostStepAction::fTrackTerminator
private

Definition at line 72 of file G4SamplingPostStepAction.hh.

Referenced by DoIt().


The documentation for this class was generated from the following files: