Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ScoreSplittingProcess Class Reference

#include <G4ScoreSplittingProcess.hh>

Inheritance diagram for G4ScoreSplittingProcess:
Collaboration 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 ()
 
G4int operator== (const G4VProcess &right) const
 
G4int 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 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
 
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 previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 73 of file G4ScoreSplittingProcess.hh.

Constructor & Destructor Documentation

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

Definition at line 53 of file G4ScoreSplittingProcess.cc.

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

Here is the call graph for this function:

G4ScoreSplittingProcess::~G4ScoreSplittingProcess ( )
virtual

Definition at line 73 of file G4ScoreSplittingProcess.cc.

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

Member Function Documentation

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

Implements G4VProcess.

Definition at line 375 of file G4ScoreSplittingProcess.cc.

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

Here is the call graph for this function:

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

Implements G4VProcess.

Definition at line 360 of file G4ScoreSplittingProcess.cc.

366 {
367  *selection = NotCandidateForSelection;
368  return DBL_MAX;
369 }
#define DBL_MAX
Definition: templates.hh:83
G4VParticleChange * G4ScoreSplittingProcess::AtRestDoIt ( const G4Track track,
const G4Step  
)
virtual

Implements G4VProcess.

Definition at line 387 of file G4ScoreSplittingProcess.cc.

390 {
391  pParticleChange->Initialize(track);
392  return pParticleChange;
393 }
virtual void Initialize(const G4Track &)
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283

Here is the call graph for this function:

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

Implements G4VProcess.

Definition at line 348 of file G4ScoreSplittingProcess.cc.

351 {
352  *condition = NotForced; // Was Forced
353  return DBL_MAX;
354 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83
G4VParticleChange * G4ScoreSplittingProcess::PostStepDoIt ( const G4Track track,
const G4Step step 
)
virtual

Implements G4VProcess.

Definition at line 134 of file G4ScoreSplittingProcess.cc.

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

Here is the call graph for this function:

G4double G4ScoreSplittingProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 110 of file G4ScoreSplittingProcess.cc.

114 {
115  // This process must be invoked anyway to score the hit
116  // - to do the scoring if the current volume is a regular structure, or
117  // - else to toggle the flag so that the SteppingManager does the scoring.
119 
120  // Future optimisation: check whether in regular structure.
121  // If it is in regular structure, be StronglyForced
122  // If not in regular structure,
123  // ask to be called only if SteppingControl is AvoidHitInvocation
124  // in order to reset it to NormalCondition
125 
126  return DBL_MAX;
127 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83
void G4ScoreSplittingProcess::StartTracking ( G4Track trk)
virtual

Initialize

Reimplemented from G4VProcess.

Definition at line 84 of file G4ScoreSplittingProcess.cc.

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

Here is the call graph for this function:

void G4ScoreSplittingProcess::Verbose ( const G4Step step) const

Definition at line 280 of file G4ScoreSplittingProcess.cc.

281 {
282  G4cout << "In mass geometry ------------------------------------------------" << G4endl;
283  G4cout << " StepLength : " << step.GetStepLength()/mm << " TotalEnergyDeposit : "
284  << step.GetTotalEnergyDeposit()/MeV << G4endl;
285  G4cout << " PreStepPoint : "
286  << step.GetPreStepPoint()->GetPhysicalVolume()->GetName() << " - ";
289  else
290  { G4cout << "NoProcessAssigned"; }
291  G4cout << G4endl;
292  G4cout << " " << step.GetPreStepPoint()->GetPosition() << G4endl;
293  G4cout << " PostStepPoint : ";
294  if(step.GetPostStepPoint()->GetPhysicalVolume())
295  { G4cout << step.GetPostStepPoint()->GetPhysicalVolume()->GetName(); }
296  else
297  { G4cout << "OutOfWorld"; }
298  G4cout << " - ";
301  else
302  { G4cout << "NoProcessAssigned"; }
303  G4cout << G4endl;
304  G4cout << " " << step.GetPostStepPoint()->GetPosition() << G4endl;
305 
306  G4cout << "In ghost geometry ------------------------------------------------" << G4endl;
307  G4cout << " StepLength : " << fSplitStep->GetStepLength()/mm
308  << " TotalEnergyDeposit : "
309  << fSplitStep->GetTotalEnergyDeposit()/MeV << G4endl;
310  G4cout << " PreStepPoint : "
311  << fSplitStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() << " ["
312  << fSplitStep->GetPreStepPoint()->GetTouchable()->GetReplicaNumber()
313  << " ]" << " - ";
314  if(fSplitStep->GetPreStepPoint()->GetProcessDefinedStep())
315  { G4cout << fSplitStep->GetPreStepPoint()->GetProcessDefinedStep()->GetProcessName(); }
316  else
317  { G4cout << "NoProcessAssigned"; }
318  G4cout << G4endl;
319  G4cout << " " << fSplitStep->GetPreStepPoint()->GetPosition() << G4endl;
320  G4cout << " PostStepPoint : ";
321  if(fSplitStep->GetPostStepPoint()->GetPhysicalVolume())
322  {
323  G4cout << fSplitStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() << " ["
324  << fSplitStep->GetPostStepPoint()->GetTouchable()->GetReplicaNumber()
325  << " ]";
326  }
327  else
328  { G4cout << "OutOfWorld"; }
329  G4cout << " - ";
330  if(fSplitStep->GetPostStepPoint()->GetProcessDefinedStep())
331  { G4cout << fSplitStep->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName(); }
332  else
333  { G4cout << "NoProcessAssigned"; }
334  G4cout << G4endl;
335  G4cout << " " << fSplitStep->GetPostStepPoint()->GetPosition() << " == "
336  << fSplitStep->GetTrack()->GetMomentumDirection()
337  << G4endl;
338 
339 }
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetStepLength() const
const G4VTouchable * GetTouchable() const
G4StepPoint * GetPreStepPoint() const
G4GLOB_DLL std::ostream G4cout
G4VPhysicalVolume * GetPhysicalVolume() const
const G4String & GetName() const
const G4ThreeVector & GetPosition() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4double GetTotalEnergyDeposit() const
const G4VProcess * GetProcessDefinedStep() const
const G4ThreeVector & GetMomentumDirection() const
G4StepPoint * GetPostStepPoint() const
virtual G4int GetReplicaNumber(G4int depth=0) const
Definition: G4VTouchable.cc:58
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4Track * GetTrack() const

Here is the call graph for this function:

Here is the caller graph for this function:


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