Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ScoreSplittingProcess Class Reference

#include <G4ScoreSplittingProcess.hh>

+ Inheritance diagram for G4ScoreSplittingProcess:

Public Member Functions

 G4ScoreSplittingProcess (const G4String &processName="ScoreSplittingProc", G4ProcessType theType=fParameterisation)
 
 ~G4ScoreSplittingProcess () override
 
void StartTracking (G4Track *) override
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *) override
 
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &) override
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *) override
 
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &) override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
void Verbose (const G4Step &) const
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual const G4VProcessGetCreatorProcess () const
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 54 of file G4ScoreSplittingProcess.hh.

Constructor & Destructor Documentation

◆ G4ScoreSplittingProcess()

G4ScoreSplittingProcess::G4ScoreSplittingProcess ( const G4String & processName = "ScoreSplittingProc",
G4ProcessType theType = fParameterisation )

Definition at line 48 of file G4ScoreSplittingProcess.cc.

49 : G4VProcess(processName, theType)
50{
51 pParticleChange = &xParticleChange;
52
53 fSplitStep = new G4Step();
54 fSplitPreStepPoint = fSplitStep->GetPreStepPoint();
55 fSplitPostStepPoint = fSplitStep->GetPostStepPoint();
56
57 if (verboseLevel > 0) {
58 G4cout << GetProcessName() << " is created " << G4endl;
59 }
60 fpEnergySplitter = new G4EnergySplitter();
61}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
G4int verboseLevel
G4VParticleChange * pParticleChange
const G4String & GetProcessName() const

◆ ~G4ScoreSplittingProcess()

G4ScoreSplittingProcess::~G4ScoreSplittingProcess ( )
override

Definition at line 66 of file G4ScoreSplittingProcess.cc.

67{
68 delete fSplitStep;
69 delete fpEnergySplitter;
70}

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4ScoreSplittingProcess::AlongStepDoIt ( const G4Track & track,
const G4Step &  )
overridevirtual

Implements G4VProcess.

Definition at line 360 of file G4ScoreSplittingProcess.cc.

361{
362 // Dummy ParticleChange ie: does nothing
363 // Expecting G4Transportation to move the track
364 dummyParticleChange.Initialize(track);
365 return &dummyParticleChange;
366}
virtual void Initialize(const G4Track &)

◆ AlongStepGetPhysicalInteractionLength()

G4double G4ScoreSplittingProcess::AlongStepGetPhysicalInteractionLength ( const G4Track & ,
G4double ,
G4double ,
G4double & ,
G4GPILSelection * selection )
overridevirtual

Implements G4VProcess.

Definition at line 346 of file G4ScoreSplittingProcess.cc.

351{
352 *selection = NotCandidateForSelection;
353 return DBL_MAX;
354}
@ NotCandidateForSelection
#define DBL_MAX
Definition templates.hh:62

◆ AtRestDoIt()

G4VParticleChange * G4ScoreSplittingProcess::AtRestDoIt ( const G4Track & track,
const G4Step &  )
overridevirtual

Implements G4VProcess.

Definition at line 371 of file G4ScoreSplittingProcess.cc.

372{
374 return pParticleChange;
375}

◆ AtRestGetPhysicalInteractionLength()

G4double G4ScoreSplittingProcess::AtRestGetPhysicalInteractionLength ( const G4Track & ,
G4ForceCondition * condition )
overridevirtual

Implements G4VProcess.

Definition at line 335 of file G4ScoreSplittingProcess.cc.

337{
338 *condition = NotForced; // Was Forced
339 return DBL_MAX;
340}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced

◆ PostStepDoIt()

G4VParticleChange * G4ScoreSplittingProcess::PostStepDoIt ( const G4Track & track,
const G4Step & step )
overridevirtual

Implements G4VProcess.

Definition at line 121 of file G4ScoreSplittingProcess.cc.

122{
123 G4VPhysicalVolume* pCurrentVolume = track.GetVolume();
124 G4LogicalVolume* pLogicalVolume = pCurrentVolume->GetLogicalVolume();
125 G4VSensitiveDetector* ptrSD = pLogicalVolume->GetSensitiveDetector();
126
128 if ((!pCurrentVolume->IsRegularStructure()) || (ptrSD == nullptr)
130 {
131 // Set the flag to make sure that Stepping Manager does the scoring
133 }
134 else {
135 G4ThreeVector preStepPosition, postStepPosition, direction, finalPostStepPosition;
137
138 G4double totalEnergyDeposit = step.GetTotalEnergyDeposit();
139 G4StepStatus fullStepStatus = step.GetPostStepPoint()->GetStepStatus();
140
141 CopyStepStart(step);
142 fSplitPreStepPoint->SetSensitiveDetector(ptrSD);
143 fOldTouchableH = fInitialTouchableH;
144 fNewTouchableH = fOldTouchableH;
145 *fSplitPostStepPoint = *(step.GetPreStepPoint());
146
147 // Split the energy
148 // ----------------
149 G4int numberVoxelsInStep = fpEnergySplitter->SplitEnergyInVolumes(&step);
150
151 preStepPosition = step.GetPreStepPoint()->GetPosition();
152 finalPostStepPosition = step.GetPostStepPoint()->GetPosition();
153 direction = (finalPostStepPosition - preStepPosition).unit();
154
155 fFinalTouchableH = track.GetNextTouchableHandle();
156
157 postStepPosition = preStepPosition;
158 // Loop over the sub-parts of this step
159 G4int iStep;
160 for (iStep = 0; iStep < numberVoxelsInStep; iStep++) {
161 G4int idVoxel = -1; // Voxel ID
162 G4double stepLength = 0.0, energyLoss = 0.0;
163
164 *fSplitPreStepPoint = *fSplitPostStepPoint;
165 fOldTouchableH = fNewTouchableH;
166
167 preStepPosition = postStepPosition;
168 fSplitPreStepPoint->SetPosition(preStepPosition);
169 fSplitPreStepPoint->SetTouchableHandle(fOldTouchableH);
170
171 fpEnergySplitter->GetLengthAndEnergyDeposited(iStep, idVoxel, stepLength, energyLoss);
172
173 // Correct the material, so that the track->GetMaterial gives correct answer
174 pLogicalVolume->SetMaterial(fpEnergySplitter->GetVoxelMaterial(iStep)); // idVoxel) );
175
176 postStepPosition = preStepPosition + stepLength * direction;
177 fSplitPostStepPoint->SetPosition(postStepPosition);
178
179 // Load the Step with the new values
180 fSplitStep->SetStepLength(stepLength);
181 fSplitStep->SetTotalEnergyDeposit(energyLoss);
182 if (iStep < numberVoxelsInStep - 1) {
184 G4int nextVoxelId = -1;
185 fpEnergySplitter->GetVoxelID(iStep + 1, nextVoxelId);
186
187 // Create new "next" touchable for each section ??
188 G4VTouchable* fNewTouchablePtr = CreateTouchableForSubStep(nextVoxelId, postStepPosition);
189 fNewTouchableH = G4TouchableHandle(fNewTouchablePtr);
190 fSplitPostStepPoint->SetTouchableHandle(fNewTouchableH);
191 }
192 else {
193 fSplitStep->GetPostStepPoint()->SetStepStatus(fullStepStatus);
194 fSplitPostStepPoint->SetTouchableHandle(fFinalTouchableH);
195 }
196
197 // As first approximation, split the NIEL in the same fractions as the energy deposit
198 G4double eLossFraction;
199 eLossFraction = (totalEnergyDeposit > 0.0) ? energyLoss / totalEnergyDeposit : 1.0;
200 fSplitStep->SetNonIonizingEnergyDeposit(step.GetNonIonizingEnergyDeposit() * eLossFraction);
201
202 fSplitPostStepPoint->SetSensitiveDetector(ptrSD);
203
204 // Call the Sensitive Detector
205 ptrSD->Hit(fSplitStep);
206
207 if (verboseLevel > 1) Verbose(step);
208 }
209 }
210
211 // This must change the Stepping Control
212 return pParticleChange;
213}
G4StepStatus
@ fGeomBoundary
@ AvoidHitInvocation
@ NormalCondition
G4ReferenceCountedHandle< G4VTouchable > G4TouchableHandle
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
void GetLengthAndEnergyDeposited(G4int stepNo, G4int &voxelID, G4double &stepLength, G4double &energyLoss)
G4Material * GetVoxelMaterial(G4int stepNo)
G4int SplitEnergyInVolumes(const G4Step *aStep)
void GetVoxelID(G4int stepNo, G4int &voxelID)
G4VSensitiveDetector * GetSensitiveDetector() const
void SetMaterial(G4Material *pMaterial)
const std::vector< std::pair< G4int, G4double > > & GetStepLengths()
static G4RegularNavigationHelper * Instance()
void Verbose(const G4Step &) const
void SetSensitiveDetector(G4VSensitiveDetector *)
G4StepStatus GetStepStatus() const
void SetStepStatus(const G4StepStatus aValue)
void SetTouchableHandle(const G4TouchableHandle &apValue)
const G4ThreeVector & GetPosition() const
void SetPosition(const G4ThreeVector &aValue)
void SetStepLength(G4double value)
void SetNonIonizingEnergyDeposit(G4double value)
G4double GetNonIonizingEnergyDeposit() const
G4double GetTotalEnergyDeposit() const
void SetTotalEnergyDeposit(G4double value)
G4VPhysicalVolume * GetVolume() const
const G4TouchableHandle & GetNextTouchableHandle() const
void ProposeSteppingControl(G4SteppingControl StepControlFlag)
virtual G4bool IsRegularStructure() const =0
G4LogicalVolume * GetLogicalVolume() const
G4bool Hit(G4Step *aStep)

◆ PostStepGetPhysicalInteractionLength()

G4double G4ScoreSplittingProcess::PostStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overridevirtual

Implements G4VProcess.

Definition at line 99 of file G4ScoreSplittingProcess.cc.

101{
102 // This process must be invoked anyway to score the hit
103 // - to do the scoring if the current volume is a regular structure, or
104 // - else to toggle the flag so that the SteppingManager does the scoring.
106
107 // Future optimisation: check whether in regular structure.
108 // If it is in regular structure, be StronglyForced
109 // If not in regular structure,
110 // ask to be called only if SteppingControl is AvoidHitInvocation
111 // in order to reset it to NormalCondition
112
113 return DBL_MAX;
114}
@ StronglyForced

◆ StartTracking()

void G4ScoreSplittingProcess::StartTracking ( G4Track * trk)
overridevirtual

Initialize

Reimplemented from G4VProcess.

Definition at line 77 of file G4ScoreSplittingProcess.cc.

78{
79 // Setup initial touchables for the first step
80 const G4Step* pStep = trk->GetStep();
81
82 fOldTouchableH = trk->GetTouchableHandle();
83 *fSplitPreStepPoint = *(pStep->GetPreStepPoint()); // Best to copy, so as to initialise
84 fSplitPreStepPoint->SetTouchableHandle(fOldTouchableH);
85 fNewTouchableH = fOldTouchableH;
86 *fSplitPostStepPoint = *(pStep->GetPostStepPoint()); // Best to copy, so as to initialise
87 fSplitPostStepPoint->SetTouchableHandle(fNewTouchableH);
88
89 /// Initialize
90 fSplitPreStepPoint->SetStepStatus(fUndefined);
91 fSplitPostStepPoint->SetStepStatus(fUndefined);
92}
@ fUndefined
const G4TouchableHandle & GetTouchableHandle() const
const G4Step * GetStep() const

◆ Verbose()

void G4ScoreSplittingProcess::Verbose ( const G4Step & step) const

Definition at line 265 of file G4ScoreSplittingProcess.cc.

266{
267 G4cout << "In mass geometry ------------------------------------------------" << G4endl;
268 G4cout << " StepLength : " << step.GetStepLength() / mm
269 << " TotalEnergyDeposit : " << step.GetTotalEnergyDeposit() / MeV << G4endl;
270 G4cout << " PreStepPoint : " << step.GetPreStepPoint()->GetPhysicalVolume()->GetName() << " - ";
271 if (step.GetPreStepPoint()->GetProcessDefinedStep() != nullptr) {
273 }
274 else {
275 G4cout << "NoProcessAssigned";
276 }
277 G4cout << G4endl;
278 G4cout << " " << step.GetPreStepPoint()->GetPosition() << G4endl;
279 G4cout << " PostStepPoint : ";
280 if (step.GetPostStepPoint()->GetPhysicalVolume() != nullptr) {
282 }
283 else {
284 G4cout << "OutOfWorld";
285 }
286 G4cout << " - ";
287 if (step.GetPostStepPoint()->GetProcessDefinedStep() != nullptr) {
289 }
290 else {
291 G4cout << "NoProcessAssigned";
292 }
293 G4cout << G4endl;
294 G4cout << " " << step.GetPostStepPoint()->GetPosition() << G4endl;
295
296 G4cout << "In ghost geometry ------------------------------------------------" << G4endl;
297 G4cout << " StepLength : " << fSplitStep->GetStepLength() / mm
298 << " TotalEnergyDeposit : " << fSplitStep->GetTotalEnergyDeposit() / MeV << G4endl;
299 G4cout << " PreStepPoint : " << fSplitStep->GetPreStepPoint()->GetPhysicalVolume()->GetName()
300 << " [" << fSplitStep->GetPreStepPoint()->GetTouchable()->GetReplicaNumber() << " ]"
301 << " - ";
302 if (fSplitStep->GetPreStepPoint()->GetProcessDefinedStep() != nullptr) {
304 }
305 else {
306 G4cout << "NoProcessAssigned";
307 }
308 G4cout << G4endl;
309 G4cout << " " << fSplitStep->GetPreStepPoint()->GetPosition() << G4endl;
310 G4cout << " PostStepPoint : ";
311 if (fSplitStep->GetPostStepPoint()->GetPhysicalVolume() != nullptr) {
312 G4cout << fSplitStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() << " ["
313 << fSplitStep->GetPostStepPoint()->GetTouchable()->GetReplicaNumber() << " ]";
314 }
315 else {
316 G4cout << "OutOfWorld";
317 }
318 G4cout << " - ";
319 if (fSplitStep->GetPostStepPoint()->GetProcessDefinedStep() != nullptr) {
321 }
322 else {
323 G4cout << "NoProcessAssigned";
324 }
325 G4cout << G4endl;
326 G4cout << " " << fSplitStep->GetPostStepPoint()->GetPosition()
327 << " == " << fSplitStep->GetTrack()->GetMomentumDirection() << G4endl;
328}
const G4VTouchable * GetTouchable() const
const G4VProcess * GetProcessDefinedStep() const
G4VPhysicalVolume * GetPhysicalVolume() const
G4Track * GetTrack() const
G4double GetStepLength() const
virtual G4int GetReplicaNumber(G4int depth=0) const
const G4ThreeVector & GetMomentumDirection() const
const G4String & GetName() const

Referenced by PostStepDoIt().


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