G4ScoreSplittingProcess Class Reference

#include <G4ScoreSplittingProcess.hh>

Inheritance diagram for G4ScoreSplittingProcess:

G4VProcess

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

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.

References G4cout, G4endl, G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4VProcess::GetProcessName(), G4VProcess::pParticleChange, and G4VProcess::verboseLevel.

00054   :G4VProcess(processName,theType),
00055    fOldTouchableH(),  fNewTouchableH(), fInitialTouchableH(), fFinalTouchableH()
00056 {
00057   pParticleChange = &xParticleChange;
00058 
00059   fSplitStep = new G4Step();
00060   fSplitPreStepPoint =  fSplitStep->GetPreStepPoint();
00061   fSplitPostStepPoint = fSplitStep->GetPostStepPoint();
00062 
00063   if (verboseLevel>0)
00064   {
00065     G4cout << GetProcessName() << " is created " << G4endl;
00066   }
00067   fpEnergySplitter = new G4EnergySplitter(); 
00068 }

G4ScoreSplittingProcess::~G4ScoreSplittingProcess (  )  [virtual]

Definition at line 73 of file G4ScoreSplittingProcess.cc.

00074 {
00075   delete fSplitStep;
00076   delete fpEnergySplitter; 
00077 }


Member Function Documentation

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

Implements G4VProcess.

Definition at line 375 of file G4ScoreSplittingProcess.cc.

References G4VParticleChange::Initialize().

00377 {
00378   // Dummy ParticleChange ie: does nothing
00379   // Expecting G4Transportation to move the track
00380   dummyParticleChange.Initialize(track);
00381   return &dummyParticleChange;
00382 }

G4double G4ScoreSplittingProcess::AlongStepGetPhysicalInteractionLength ( const G4Track ,
G4double  ,
G4double  ,
G4double ,
G4GPILSelection  
) [virtual]

Implements G4VProcess.

Definition at line 360 of file G4ScoreSplittingProcess.cc.

References DBL_MAX, and NotCandidateForSelection.

00366 {
00367   *selection = NotCandidateForSelection;
00368   return DBL_MAX;
00369 }

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

Implements G4VProcess.

Definition at line 387 of file G4ScoreSplittingProcess.cc.

References G4VParticleChange::Initialize(), and G4VProcess::pParticleChange.

00390 { 
00391   pParticleChange->Initialize(track);
00392   return pParticleChange;
00393 }

G4double G4ScoreSplittingProcess::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition  
) [virtual]

Implements G4VProcess.

Definition at line 348 of file G4ScoreSplittingProcess.cc.

References DBL_MAX, and NotForced.

00351 {
00352   *condition = NotForced;  // Was Forced
00353   return DBL_MAX;
00354 }

G4VParticleChange * G4ScoreSplittingProcess::PostStepDoIt ( const G4Track ,
const G4Step  
) [virtual]

Implements G4VProcess.

Definition at line 134 of file G4ScoreSplittingProcess.cc.

References AvoidHitInvocation, fGeomBoundary, G4EnergySplitter::GetLengthAndEnergyDeposited(), G4VPhysicalVolume::GetLogicalVolume(), G4Track::GetNextTouchableHandle(), G4Step::GetNonIonizingEnergyDeposit(), G4StepPoint::GetPosition(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4LogicalVolume::GetSensitiveDetector(), G4RegularNavigationHelper::GetStepLengths(), G4StepPoint::GetStepStatus(), G4Step::GetTotalEnergyDeposit(), G4Track::GetVolume(), G4EnergySplitter::GetVoxelID(), G4EnergySplitter::GetVoxelMaterial(), G4VSensitiveDetector::Hit(), G4VParticleChange::Initialize(), G4RegularNavigationHelper::Instance(), G4VPhysicalVolume::IsRegularStructure(), NormalCondition, G4VProcess::pParticleChange, G4VParticleChange::ProposeSteppingControl(), G4LogicalVolume::SetMaterial(), G4Step::SetNonIonizingEnergyDeposit(), G4StepPoint::SetPosition(), G4StepPoint::SetSensitiveDetector(), G4Step::SetStepLength(), G4StepPoint::SetStepStatus(), G4Step::SetTotalEnergyDeposit(), G4StepPoint::SetTouchableHandle(), G4EnergySplitter::SplitEnergyInVolumes(), Verbose(), and G4VProcess::verboseLevel.

00137 { 
00138   G4VPhysicalVolume*    pCurrentVolume= track.GetVolume(); 
00139   G4LogicalVolume*      pLogicalVolume= pCurrentVolume->GetLogicalVolume(); 
00140   G4VSensitiveDetector* ptrSD = pLogicalVolume->GetSensitiveDetector();
00141 
00142   pParticleChange->Initialize(track); 
00143   if(  ( ! pCurrentVolume->IsRegularStructure() ) || ( !ptrSD ) 
00144     || G4RegularNavigationHelper::Instance()->GetStepLengths().size() <= 1) {
00145      // Set the flag to make sure that Stepping Manager does the scoring
00146      pParticleChange->ProposeSteppingControl( NormalCondition );     
00147   } else { 
00148      G4ThreeVector preStepPosition, postStepPosition, direction, finalPostStepPosition;
00149      pParticleChange->ProposeSteppingControl( AvoidHitInvocation );
00150  
00151      G4double totalEnergyDeposit=  step.GetTotalEnergyDeposit();
00152      G4StepStatus  fullStepStatus= step.GetPostStepPoint()->GetStepStatus(); 
00153 
00154      CopyStepStart(step);
00155      fSplitPreStepPoint->SetSensitiveDetector(ptrSD);
00156      fOldTouchableH = fInitialTouchableH;
00157      fNewTouchableH=  fOldTouchableH; 
00158      *fSplitPostStepPoint= *(step.GetPreStepPoint()); 
00159      
00160      // Split the energy
00161      // ----------------
00162      G4int numberVoxelsInStep= fpEnergySplitter->SplitEnergyInVolumes( &step );
00163 
00164      preStepPosition= step.GetPreStepPoint()->GetPosition();
00165      finalPostStepPosition= step.GetPostStepPoint()->GetPosition(); 
00166      direction= (finalPostStepPosition - preStepPosition).unit(); 
00167 
00168      fFinalTouchableH= track.GetNextTouchableHandle(); 
00169 
00170      postStepPosition= preStepPosition;
00171      // Loop over the sub-parts of this step
00172      G4int iStep;
00173      for ( iStep=0; iStep < numberVoxelsInStep; iStep++ ){ 
00174         G4int     idVoxel=  -1;  // Voxel ID
00175         G4double  stepLength=0.0, energyLoss= 0.0;
00176 
00177         *fSplitPreStepPoint  = *fSplitPostStepPoint;
00178         fOldTouchableH = fNewTouchableH; 
00179 
00180         preStepPosition= postStepPosition;
00181         fSplitPreStepPoint->SetPosition( preStepPosition );
00182         fSplitPreStepPoint->SetTouchableHandle(fOldTouchableH);
00183         
00184         fpEnergySplitter->GetLengthAndEnergyDeposited( iStep, idVoxel, stepLength, energyLoss);
00185 
00186         // Correct the material, so that the track->GetMaterial gives correct answer
00187         pLogicalVolume->SetMaterial( fpEnergySplitter->GetVoxelMaterial( iStep) );  // idVoxel) ); 
00188 
00189         postStepPosition= preStepPosition + stepLength * direction; 
00190         fSplitPostStepPoint->SetPosition(postStepPosition); 
00191 
00192         // Load the Step with the new values
00193         fSplitStep->SetStepLength(stepLength);
00194         fSplitStep->SetTotalEnergyDeposit(energyLoss);
00195         if( iStep < numberVoxelsInStep -1 ){ 
00196           fSplitStep->GetPostStepPoint()->SetStepStatus( fGeomBoundary );
00197           G4int  nextVoxelId= -1;
00198           fpEnergySplitter->GetVoxelID( iStep+1, nextVoxelId );
00199 
00200           // Create new "next" touchable for each section ??
00201           G4VTouchable* fNewTouchablePtr= 
00202               CreateTouchableForSubStep( nextVoxelId, postStepPosition );
00203           fNewTouchableH= G4TouchableHandle(fNewTouchablePtr); 
00204           fSplitPostStepPoint->SetTouchableHandle( fNewTouchableH );
00205         } else { 
00206           fSplitStep->GetPostStepPoint()->SetStepStatus( fullStepStatus );
00207           fSplitPostStepPoint->SetTouchableHandle( fFinalTouchableH );
00208         }
00209 
00210 
00211         // As first approximation, split the NIEL in the same fractions as the energy deposit
00212         G4double eLossFraction; 
00213         eLossFraction= (totalEnergyDeposit>0.0) ? energyLoss / totalEnergyDeposit : 1.0 ; 
00214         fSplitStep->SetNonIonizingEnergyDeposit(step.GetNonIonizingEnergyDeposit()*eLossFraction);
00215         
00216         fSplitPostStepPoint->SetSensitiveDetector( ptrSD ); 
00217 
00218         // Call the Sensitive Detector
00219         ptrSD->Hit(fSplitStep);
00220 
00221         if (verboseLevel>1) Verbose(step);
00222      }
00223   }
00224 
00225   // This must change the Stepping Control
00226   return pParticleChange;
00227 }

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

Implements G4VProcess.

Definition at line 110 of file G4ScoreSplittingProcess.cc.

References DBL_MAX, and StronglyForced.

00114 {
00115   // This process must be invoked anyway to score the hit
00116   //   - to do the scoring if the current volume is a regular structure, or
00117   //   - else to toggle the flag so that the SteppingManager does the scoring.
00118   *condition = StronglyForced;
00119 
00120   // Future optimisation: check whether in regular structure.  
00121   //  If it is in regular structure,  be StronglyForced
00122   //  If  not  in regular structure, 
00123   //         ask to be called only if SteppingControl is AvoidHitInvocation
00124   //         in order to reset it to NormalCondition
00125 
00126   return DBL_MAX;
00127 }

void G4ScoreSplittingProcess::StartTracking ( G4Track  )  [virtual]

Reimplemented from G4VProcess.

Definition at line 84 of file G4ScoreSplittingProcess.cc.

References fUndefined, G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4Track::GetStep(), G4Track::GetTouchableHandle(), G4StepPoint::SetStepStatus(), and G4StepPoint::SetTouchableHandle().

00085 {
00086 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00087 // Setup initial touchables for the first step
00088 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00089   const G4Step* pStep= trk->GetStep();
00090 
00091   fOldTouchableH = trk->GetTouchableHandle();
00092   *fSplitPreStepPoint = *(pStep->GetPreStepPoint()); // Best to copy, so as to initialise
00093   fSplitPreStepPoint->SetTouchableHandle(fOldTouchableH);
00094   fNewTouchableH = fOldTouchableH;
00095   *fSplitPostStepPoint= *(pStep->GetPostStepPoint()); // Best to copy, so as to initialise
00096   fSplitPostStepPoint->SetTouchableHandle(fNewTouchableH);
00097 
00099   fSplitPreStepPoint ->SetStepStatus(fUndefined);
00100   fSplitPostStepPoint->SetStepStatus(fUndefined);
00101 }

void G4ScoreSplittingProcess::Verbose ( const G4Step  )  const

Definition at line 280 of file G4ScoreSplittingProcess.cc.

References G4cout, G4endl, G4Track::GetMomentumDirection(), G4VPhysicalVolume::GetName(), G4StepPoint::GetPhysicalVolume(), G4StepPoint::GetPosition(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4StepPoint::GetProcessDefinedStep(), G4VProcess::GetProcessName(), G4VTouchable::GetReplicaNumber(), G4Step::GetStepLength(), G4Step::GetTotalEnergyDeposit(), G4StepPoint::GetTouchable(), and G4Step::GetTrack().

Referenced by PostStepDoIt().

00281 {
00282   G4cout << "In mass geometry ------------------------------------------------" << G4endl;
00283   G4cout << " StepLength : " << step.GetStepLength()/mm << "      TotalEnergyDeposit : "
00284          << step.GetTotalEnergyDeposit()/MeV << G4endl;
00285   G4cout << " PreStepPoint : "
00286          << step.GetPreStepPoint()->GetPhysicalVolume()->GetName() << " - ";
00287   if(step.GetPreStepPoint()->GetProcessDefinedStep())
00288   { G4cout << step.GetPreStepPoint()->GetProcessDefinedStep()->GetProcessName(); }
00289   else
00290   { G4cout << "NoProcessAssigned"; }
00291   G4cout << G4endl;
00292   G4cout << "                " << step.GetPreStepPoint()->GetPosition() << G4endl;
00293   G4cout << " PostStepPoint : ";
00294   if(step.GetPostStepPoint()->GetPhysicalVolume()) 
00295   { G4cout << step.GetPostStepPoint()->GetPhysicalVolume()->GetName(); }
00296   else
00297   { G4cout << "OutOfWorld"; }
00298   G4cout << " - ";
00299   if(step.GetPostStepPoint()->GetProcessDefinedStep())
00300   { G4cout << step.GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName(); }
00301   else
00302   { G4cout << "NoProcessAssigned"; }
00303   G4cout << G4endl;
00304   G4cout << "                 " << step.GetPostStepPoint()->GetPosition() << G4endl;
00305 
00306   G4cout << "In ghost geometry ------------------------------------------------" << G4endl;
00307   G4cout << " StepLength : " << fSplitStep->GetStepLength()/mm
00308          << "      TotalEnergyDeposit : "
00309          << fSplitStep->GetTotalEnergyDeposit()/MeV << G4endl;
00310   G4cout << " PreStepPoint : "
00311          << fSplitStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() << " ["
00312            << fSplitStep->GetPreStepPoint()->GetTouchable()->GetReplicaNumber()
00313            << " ]" << " - ";
00314   if(fSplitStep->GetPreStepPoint()->GetProcessDefinedStep())
00315   { G4cout << fSplitStep->GetPreStepPoint()->GetProcessDefinedStep()->GetProcessName(); }
00316   else
00317   { G4cout << "NoProcessAssigned"; }
00318   G4cout << G4endl;
00319   G4cout << "                " << fSplitStep->GetPreStepPoint()->GetPosition() << G4endl;
00320   G4cout << " PostStepPoint : ";
00321   if(fSplitStep->GetPostStepPoint()->GetPhysicalVolume()) 
00322   {
00323     G4cout << fSplitStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() << " ["
00324            << fSplitStep->GetPostStepPoint()->GetTouchable()->GetReplicaNumber()
00325            << " ]";
00326   }
00327   else
00328   { G4cout << "OutOfWorld"; }
00329   G4cout << " - ";
00330   if(fSplitStep->GetPostStepPoint()->GetProcessDefinedStep())
00331   { G4cout << fSplitStep->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName(); }
00332   else
00333   { G4cout << "NoProcessAssigned"; }
00334   G4cout << G4endl;
00335   G4cout << "                 " << fSplitStep->GetPostStepPoint()->GetPosition() << " == "
00336          << fSplitStep->GetTrack()->GetMomentumDirection() 
00337          << G4endl;
00338 
00339 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:53:21 2013 for Geant4 by  doxygen 1.4.7