Geant4 10.7.0
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)
 
virtual ~G4ScoreSplittingProcess ()
 
void StartTracking (G4Track *)
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
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
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
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 void StartTracking (G4Track *)
 
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 72 of file G4ScoreSplittingProcess.hh.

Constructor & Destructor Documentation

◆ G4ScoreSplittingProcess()

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

Definition at line 51 of file G4ScoreSplittingProcess.cc.

53 :G4VProcess(processName,theType),
54 fOldTouchableH(), fNewTouchableH(), fInitialTouchableH(), fFinalTouchableH()
55{
56 pParticleChange = &xParticleChange;
57
58 fSplitStep = new G4Step();
59 fSplitPreStepPoint = fSplitStep->GetPreStepPoint();
60 fSplitPostStepPoint = fSplitStep->GetPostStepPoint();
61
62 if (verboseLevel>0)
63 {
64 G4cout << GetProcessName() << " is created " << G4endl;
65 }
66 fpEnergySplitter = new G4EnergySplitter();
67}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
G4int verboseLevel
Definition: G4VProcess.hh:356
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

◆ ~G4ScoreSplittingProcess()

G4ScoreSplittingProcess::~G4ScoreSplittingProcess ( )
virtual

Definition at line 72 of file G4ScoreSplittingProcess.cc.

73{
74 delete fSplitStep;
75 delete fpEnergySplitter;
76}

Member Function Documentation

◆ AlongStepDoIt()

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

Implements G4VProcess.

Definition at line 374 of file G4ScoreSplittingProcess.cc.

376{
377 // Dummy ParticleChange ie: does nothing
378 // Expecting G4Transportation to move the track
379 dummyParticleChange.Initialize(track);
380 return &dummyParticleChange;
381}
virtual void Initialize(const G4Track &)

◆ AlongStepGetPhysicalInteractionLength()

G4double G4ScoreSplittingProcess::AlongStepGetPhysicalInteractionLength ( const G4Track ,
G4double  ,
G4double  ,
G4double ,
G4GPILSelection selection 
)
virtual

Implements G4VProcess.

Definition at line 359 of file G4ScoreSplittingProcess.cc.

365{
366 *selection = NotCandidateForSelection;
367 return DBL_MAX;
368}
@ NotCandidateForSelection
#define DBL_MAX
Definition: templates.hh:62

◆ AtRestDoIt()

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

Implements G4VProcess.

Definition at line 386 of file G4ScoreSplittingProcess.cc.

389{
391 return pParticleChange;
392}

◆ AtRestGetPhysicalInteractionLength()

G4double G4ScoreSplittingProcess::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 347 of file G4ScoreSplittingProcess.cc.

350{
351 *condition = NotForced; // Was Forced
352 return DBL_MAX;
353}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced

◆ PostStepDoIt()

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

Implements G4VProcess.

Definition at line 133 of file G4ScoreSplittingProcess.cc.

136{
137 G4VPhysicalVolume* pCurrentVolume= track.GetVolume();
138 G4LogicalVolume* pLogicalVolume= pCurrentVolume->GetLogicalVolume();
139 G4VSensitiveDetector* ptrSD = pLogicalVolume->GetSensitiveDetector();
140
142 if( ( ! pCurrentVolume->IsRegularStructure() ) || ( !ptrSD )
144 // Set the flag to make sure that Stepping Manager does the scoring
146 } else {
147 G4ThreeVector preStepPosition, postStepPosition, direction, finalPostStepPosition;
149
150 G4double totalEnergyDeposit= step.GetTotalEnergyDeposit();
151 G4StepStatus fullStepStatus= step.GetPostStepPoint()->GetStepStatus();
152
153 CopyStepStart(step);
154 fSplitPreStepPoint->SetSensitiveDetector(ptrSD);
155 fOldTouchableH = fInitialTouchableH;
156 fNewTouchableH= fOldTouchableH;
157 *fSplitPostStepPoint= *(step.GetPreStepPoint());
158
159 // Split the energy
160 // ----------------
161 G4int numberVoxelsInStep= fpEnergySplitter->SplitEnergyInVolumes( &step );
162
163 preStepPosition= step.GetPreStepPoint()->GetPosition();
164 finalPostStepPosition= step.GetPostStepPoint()->GetPosition();
165 direction= (finalPostStepPosition - preStepPosition).unit();
166
167 fFinalTouchableH= track.GetNextTouchableHandle();
168
169 postStepPosition= preStepPosition;
170 // Loop over the sub-parts of this step
171 G4int iStep;
172 for ( iStep=0; iStep < numberVoxelsInStep; iStep++ ){
173 G4int idVoxel= -1; // Voxel ID
174 G4double stepLength=0.0, energyLoss= 0.0;
175
176 *fSplitPreStepPoint = *fSplitPostStepPoint;
177 fOldTouchableH = fNewTouchableH;
178
179 preStepPosition= postStepPosition;
180 fSplitPreStepPoint->SetPosition( preStepPosition );
181 fSplitPreStepPoint->SetTouchableHandle(fOldTouchableH);
182
183 fpEnergySplitter->GetLengthAndEnergyDeposited( iStep, idVoxel, stepLength, energyLoss);
184
185 // Correct the material, so that the track->GetMaterial gives correct answer
186 pLogicalVolume->SetMaterial( fpEnergySplitter->GetVoxelMaterial( iStep) ); // idVoxel) );
187
188 postStepPosition= preStepPosition + stepLength * direction;
189 fSplitPostStepPoint->SetPosition(postStepPosition);
190
191 // Load the Step with the new values
192 fSplitStep->SetStepLength(stepLength);
193 fSplitStep->SetTotalEnergyDeposit(energyLoss);
194 if( iStep < numberVoxelsInStep -1 ){
196 G4int nextVoxelId= -1;
197 fpEnergySplitter->GetVoxelID( iStep+1, nextVoxelId );
198
199 // Create new "next" touchable for each section ??
200 G4VTouchable* fNewTouchablePtr=
201 CreateTouchableForSubStep( nextVoxelId, postStepPosition );
202 fNewTouchableH= G4TouchableHandle(fNewTouchablePtr);
203 fSplitPostStepPoint->SetTouchableHandle( fNewTouchableH );
204 } else {
205 fSplitStep->GetPostStepPoint()->SetStepStatus( fullStepStatus );
206 fSplitPostStepPoint->SetTouchableHandle( fFinalTouchableH );
207 }
208
209
210 // As first approximation, split the NIEL in the same fractions as the energy deposit
211 G4double eLossFraction;
212 eLossFraction= (totalEnergyDeposit>0.0) ? energyLoss / totalEnergyDeposit : 1.0 ;
213 fSplitStep->SetNonIonizingEnergyDeposit(step.GetNonIonizingEnergyDeposit()*eLossFraction);
214
215 fSplitPostStepPoint->SetSensitiveDetector( ptrSD );
216
217 // Call the Sensitive Detector
218 ptrSD->Hit(fSplitStep);
219
220 if (verboseLevel>1) Verbose(step);
221 }
222 }
223
224 // This must change the Stepping Control
225 return pParticleChange;
226}
G4StepStatus
Definition: G4StepStatus.hh:40
@ fGeomBoundary
Definition: G4StepStatus.hh:43
@ 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 
)
virtual

Implements G4VProcess.

Definition at line 109 of file G4ScoreSplittingProcess.cc.

113{
114 // This process must be invoked anyway to score the hit
115 // - to do the scoring if the current volume is a regular structure, or
116 // - else to toggle the flag so that the SteppingManager does the scoring.
118
119 // Future optimisation: check whether in regular structure.
120 // If it is in regular structure, be StronglyForced
121 // If not in regular structure,
122 // ask to be called only if SteppingControl is AvoidHitInvocation
123 // in order to reset it to NormalCondition
124
125 return DBL_MAX;
126}
@ StronglyForced

◆ StartTracking()

void G4ScoreSplittingProcess::StartTracking ( G4Track trk)
virtual

Initialize

Reimplemented from G4VProcess.

Definition at line 83 of file G4ScoreSplittingProcess.cc.

84{
85//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
86// Setup initial touchables for the first step
87//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
88 const G4Step* pStep= trk->GetStep();
89
90 fOldTouchableH = trk->GetTouchableHandle();
91 *fSplitPreStepPoint = *(pStep->GetPreStepPoint()); // Best to copy, so as to initialise
92 fSplitPreStepPoint->SetTouchableHandle(fOldTouchableH);
93 fNewTouchableH = fOldTouchableH;
94 *fSplitPostStepPoint= *(pStep->GetPostStepPoint()); // Best to copy, so as to initialise
95 fSplitPostStepPoint->SetTouchableHandle(fNewTouchableH);
96
97 /// Initialize
98 fSplitPreStepPoint ->SetStepStatus(fUndefined);
99 fSplitPostStepPoint->SetStepStatus(fUndefined);
100}
@ fUndefined
Definition: G4StepStatus.hh:55
const G4TouchableHandle & GetTouchableHandle() const
const G4Step * GetStep() const

◆ Verbose()

void G4ScoreSplittingProcess::Verbose ( const G4Step step) const

Definition at line 279 of file G4ScoreSplittingProcess.cc.

280{
281 G4cout << "In mass geometry ------------------------------------------------" << G4endl;
282 G4cout << " StepLength : " << step.GetStepLength()/mm << " TotalEnergyDeposit : "
283 << step.GetTotalEnergyDeposit()/MeV << G4endl;
284 G4cout << " PreStepPoint : "
285 << step.GetPreStepPoint()->GetPhysicalVolume()->GetName() << " - ";
288 else
289 { G4cout << "NoProcessAssigned"; }
290 G4cout << G4endl;
291 G4cout << " " << step.GetPreStepPoint()->GetPosition() << G4endl;
292 G4cout << " PostStepPoint : ";
295 else
296 { G4cout << "OutOfWorld"; }
297 G4cout << " - ";
300 else
301 { G4cout << "NoProcessAssigned"; }
302 G4cout << G4endl;
303 G4cout << " " << step.GetPostStepPoint()->GetPosition() << G4endl;
304
305 G4cout << "In ghost geometry ------------------------------------------------" << G4endl;
306 G4cout << " StepLength : " << fSplitStep->GetStepLength()/mm
307 << " TotalEnergyDeposit : "
308 << fSplitStep->GetTotalEnergyDeposit()/MeV << G4endl;
309 G4cout << " PreStepPoint : "
310 << fSplitStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() << " ["
311 << fSplitStep->GetPreStepPoint()->GetTouchable()->GetReplicaNumber()
312 << " ]" << " - ";
313 if(fSplitStep->GetPreStepPoint()->GetProcessDefinedStep())
315 else
316 { G4cout << "NoProcessAssigned"; }
317 G4cout << G4endl;
318 G4cout << " " << fSplitStep->GetPreStepPoint()->GetPosition() << G4endl;
319 G4cout << " PostStepPoint : ";
320 if(fSplitStep->GetPostStepPoint()->GetPhysicalVolume())
321 {
322 G4cout << fSplitStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() << " ["
323 << fSplitStep->GetPostStepPoint()->GetTouchable()->GetReplicaNumber()
324 << " ]";
325 }
326 else
327 { G4cout << "OutOfWorld"; }
328 G4cout << " - ";
329 if(fSplitStep->GetPostStepPoint()->GetProcessDefinedStep())
331 else
332 { G4cout << "NoProcessAssigned"; }
333 G4cout << G4endl;
334 G4cout << " " << fSplitStep->GetPostStepPoint()->GetPosition() << " == "
335 << fSplitStep->GetTrack()->GetMomentumDirection()
336 << G4endl;
337
338}
const G4VTouchable * GetTouchable() const
const G4VProcess * GetProcessDefinedStep() const
G4VPhysicalVolume * GetPhysicalVolume() const
G4Track * GetTrack() const
G4double GetStepLength() const
const G4ThreeVector & GetMomentumDirection() const
const G4String & GetName() const
virtual G4int GetReplicaNumber(G4int depth=0) const
Definition: G4VTouchable.cc:55

Referenced by PostStepDoIt().


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