G4WeightWindowProcess Class Reference

#include <G4WeightWindowProcess.hh>

Inheritance diagram for G4WeightWindowProcess:

G4VProcess G4VTrackTerminator

Public Member Functions

 G4WeightWindowProcess (const G4VWeightWindowAlgorithm &aWeightWindowAlgorithm, const G4VWeightWindowStore &aWWStore, const G4VTrackTerminator *TrackTerminator, G4PlaceOfAction placeOfAction, const G4String &aName="WeightWindowProcess", G4bool para=false)
virtual ~G4WeightWindowProcess ()
void SetParallelWorld (G4String parallelWorldName)
void SetParallelWorld (G4VPhysicalVolume *parallelWorld)
void StartTracking (G4Track *)
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
virtual void KillTrack () const
virtual const G4StringGetName () const
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)

Detailed Description

Definition at line 62 of file G4WeightWindowProcess.hh.


Constructor & Destructor Documentation

G4WeightWindowProcess::G4WeightWindowProcess ( const G4VWeightWindowAlgorithm aWeightWindowAlgorithm,
const G4VWeightWindowStore aWWStore,
const G4VTrackTerminator TrackTerminator,
G4PlaceOfAction  placeOfAction,
const G4String aName = "WeightWindowProcess",
G4bool  para = false 
)

Definition at line 55 of file G4WeightWindowProcess.cc.

References FatalException, G4cout, G4endl, G4Exception(), G4PathFinder::GetInstance(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4VProcess::GetProcessName(), G4TransportationManager::GetTransportationManager(), G4VProcess::pParticleChange, and G4VProcess::verboseLevel.

00061  : G4VProcess(aName),
00062    fParticleChange(new G4ParticleChange),
00063    fWeightWindowAlgorithm(aWeightWindowAlgorithm),
00064    fWeightWindowStore(aWWStore),
00065    fPostStepAction(0),
00066    fPlaceOfAction(placeOfAction),
00067    fGhostNavigator(0), fNavigatorID(-1), fFieldTrack('0')
00068 {
00069   if (TrackTerminator)
00070   {
00071     fPostStepAction = new G4SamplingPostStepAction(*TrackTerminator);
00072   }
00073   else
00074   {
00075     fPostStepAction = new G4SamplingPostStepAction(*this);
00076   }
00077   if (!fParticleChange)
00078   {
00079     G4Exception("G4WeightWindowProcess::G4WeightWindowProcess()",
00080                 "FatalError", FatalException,
00081                 "Failed allocation of G4ParticleChange !");
00082   }
00083   G4VProcess::pParticleChange = fParticleChange;
00084 
00085   fGhostStep = new G4Step();
00086   fGhostPreStepPoint = fGhostStep->GetPreStepPoint();
00087   fGhostPostStepPoint = fGhostStep->GetPostStepPoint();
00088 
00089   fTransportationManager = G4TransportationManager::GetTransportationManager();
00090   fPathFinder = G4PathFinder::GetInstance();
00091 
00092   if (verboseLevel>0)
00093   {
00094     G4cout << GetProcessName() << " is created " << G4endl;
00095   }
00096 
00097   paraflag = para;
00098 
00099 }

G4WeightWindowProcess::~G4WeightWindowProcess (  )  [virtual]

Definition at line 101 of file G4WeightWindowProcess.cc.

00102 {
00103 
00104   delete fPostStepAction;
00105   delete fParticleChange;
00106   delete fGhostStep;
00107 
00108 }


Member Function Documentation

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

Implements G4VProcess.

Definition at line 388 of file G4WeightWindowProcess.cc.

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

00389 {
00390   // Dummy ParticleChange ie: does nothing
00391   // Expecting G4Transportation to move the track
00392   pParticleChange->Initialize(track);
00393   return pParticleChange; 
00394 
00395   //  return 0;
00396 }

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

Implements G4VProcess.

Definition at line 300 of file G4WeightWindowProcess.cc.

References CandidateForSelection, G4Navigator::ComputeSafety(), G4PathFinder::ComputeStep(), DBL_MAX, G4Track::GetCurrentStepNumber(), G4Track::GetVolume(), kDoNot, kSharedOther, kSharedTransport, kUnique, NotCandidateForSelection, and G4FieldTrackUpdator::Update().

00303 {
00304   if(paraflag) {
00305     static G4FieldTrack endTrack('0');
00306     static ELimited eLimited;
00307     
00308     *selection = NotCandidateForSelection;
00309     G4double returnedStep = DBL_MAX;
00310     
00311     if (previousStepSize > 0.)
00312       { fGhostSafety -= previousStepSize; }
00313     //  else
00314     //  { fGhostSafety = -1.; }
00315     if (fGhostSafety < 0.) fGhostSafety = 0.0;
00316     
00317     // ------------------------------------------
00318     // Determination of the proposed STEP LENGTH:
00319     // ------------------------------------------
00320     if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
00321       {
00322         // I have no chance to limit
00323         returnedStep = currentMinimumStep;
00324         fOnBoundary = false;
00325         proposedSafety = fGhostSafety - currentMinimumStep;
00326       }
00327     else // (currentMinimumStep > fGhostSafety: I may limit the Step)
00328       {
00329         G4FieldTrackUpdator::Update(&fFieldTrack,&track);
00330         //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00331         // ComputeStep
00332         //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00333         returnedStep
00334           = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
00335                                      track.GetCurrentStepNumber(),fGhostSafety,eLimited,
00336                                      endTrack,track.GetVolume());
00337         //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00338         if(eLimited == kDoNot)
00339           {
00340             // Track is not on the boundary
00341             fOnBoundary = false;
00342             fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition());
00343           }
00344         else
00345           {
00346             // Track is on the boundary
00347             fOnBoundary = true;
00348             proposedSafety = fGhostSafety;
00349           }
00350         //xbug? proposedSafety = fGhostSafety;
00351         if(eLimited == kUnique || eLimited == kSharedOther) {
00352           *selection = CandidateForSelection;
00353         }else if (eLimited == kSharedTransport) { 
00354           returnedStep *= (1.0 + 1.0e-9);  
00355           // Expand to disable its selection in Step Manager comparison
00356         }
00357         
00358       }
00359 
00360   // ----------------------------------------------
00361   // Returns the fGhostSafety as the proposedSafety
00362   // The SteppingManager will take care of keeping
00363   // the smallest one.
00364   // ----------------------------------------------
00365     return returnedStep;
00366 
00367   } else {
00368     return DBL_MAX;
00369     // not sensible - time goes backwards!    return -1.0;
00370   }
00371 
00372 }

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

Implements G4VProcess.

Definition at line 382 of file G4WeightWindowProcess.cc.

00383 {
00384   return 0;
00385 }

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

Implements G4VProcess.

Definition at line 375 of file G4WeightWindowProcess.cc.

00377 {
00378   return -1.0;
00379 }

const G4String & G4WeightWindowProcess::GetName (  )  const [virtual]

Implements G4VTrackTerminator.

Definition at line 294 of file G4WeightWindowProcess.cc.

References G4VProcess::GetProcessName().

00295 {
00296   return G4VProcess::GetProcessName();
00297 }

void G4WeightWindowProcess::KillTrack (  )  const [virtual]

Implements G4VTrackTerminator.

Definition at line 289 of file G4WeightWindowProcess.cc.

References fStopAndKill, and G4VParticleChange::ProposeTrackStatus().

00290 {
00291   fParticleChange->ProposeTrackStatus(fStopAndKill);
00292 }

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

Implements G4VProcess.

Definition at line 199 of file G4WeightWindowProcess.cc.

References G4VWeightWindowAlgorithm::Calculate(), G4PathFinder::CreateTouchableHandle(), G4SamplingPostStepAction::DoIt(), fGeomBoundary, G4Track::GetKineticEnergy(), G4VWeightWindowStore::GetLowerWeight(), G4StepPoint::GetPhysicalVolume(), G4Step::GetPostStepPoint(), G4VTouchable::GetReplicaNumber(), G4Step::GetStepLength(), G4StepPoint::GetStepStatus(), G4StepPoint::GetTouchable(), G4StepPoint::GetTouchableHandle(), G4Track::GetWeight(), G4ParticleChange::Initialize(), G4VTrackTerminator::kCarTolerance, onBoundary, onBoundaryAndCollision, onCollision, and G4StepPoint::SetTouchableHandle().

00201 {
00202 
00203   fParticleChange->Initialize(aTrack);
00204 
00205   if(paraflag) {
00206     fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
00207     //xbug?    fOnBoundary = false;
00208     CopyStep(aStep);
00209 
00210     if(fOnBoundary)
00211       {
00212 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00213 // Locate the point and get new touchable
00214 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00215   //??  fPathFinder->Locate(step.GetPostStepPoint()->GetPosition(),
00216   //??                      step.GetPostStepPoint()->GetMomentumDirection());
00217         fNewGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
00218 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00219       }
00220     else
00221       {
00222         // Do I need this ??????????????????????????????????????????????????????????
00223         // fGhostNavigator->LocateGlobalPointWithinVolume(track.GetPosition());
00224         // ?????????????????????????????????????????????????????????????????????????
00225         
00226         // fPathFinder->ReLocate(track.GetPosition());
00227         
00228         // reuse the touchable
00229         fNewGhostTouchable = fOldGhostTouchable;
00230       }
00231 
00232     fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
00233     fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
00234 
00235   }
00236 
00237   if (aStep.GetStepLength() > kCarTolerance)
00238   {
00239 //     if ( ( fPlaceOfAction == onBoundaryAndCollision)
00240 //       || ( (fPlaceOfAction == onBoundary) && 
00241 //            (aStep.GetPostStepPoint()->GetStepStatus() == fGeomBoundary) )
00242 //       || ( (fPlaceOfAction == onCollision) && 
00243 //            (aStep.GetPostStepPoint()->GetStepStatus() != fGeomBoundary) ) )
00244     if(paraflag) {
00245       if ( ( fPlaceOfAction == onBoundaryAndCollision)
00246            || ( (fPlaceOfAction == onBoundary) && 
00247                 (fGhostPostStepPoint->GetStepStatus() == fGeomBoundary) )
00248            || ( (fPlaceOfAction == onCollision) && 
00249                 (fGhostPostStepPoint->GetStepStatus() != fGeomBoundary) ) )
00250         {
00251 
00252 //       G4StepPoint *postpoint = aStep.GetPostStepPoint();
00253 
00254 //       G4GeometryCell postCell(*(postpoint->GetPhysicalVolume()), 
00255 //                              postpoint->GetTouchable()->GetReplicaNumber());
00256 
00257           G4GeometryCell postCell(*(fGhostPostStepPoint->GetPhysicalVolume()), 
00258                                   fGhostPostStepPoint->GetTouchable()->GetReplicaNumber());
00259           G4Nsplit_Weight nw =
00260             fWeightWindowAlgorithm.Calculate(aTrack.GetWeight(),
00261                                              fWeightWindowStore.GetLowerWeight(postCell,
00262                                                     aTrack.GetKineticEnergy()));
00263           fPostStepAction->DoIt(aTrack, fParticleChange, nw);
00264         }
00265     } else {
00266       if ( ( fPlaceOfAction == onBoundaryAndCollision)
00267            || ( (fPlaceOfAction == onBoundary) && 
00268                 (aStep.GetPostStepPoint()->GetStepStatus() == fGeomBoundary) )
00269            || ( (fPlaceOfAction == onCollision) && 
00270                 (aStep.GetPostStepPoint()->GetStepStatus() != fGeomBoundary) ) )
00271         {
00272           
00273           G4StepPoint *postpoint = aStep.GetPostStepPoint();
00274           
00275           G4GeometryCell postCell(*(postpoint->GetPhysicalVolume()), 
00276                                   postpoint->GetTouchable()->GetReplicaNumber());
00277           
00278           G4Nsplit_Weight nw =
00279             fWeightWindowAlgorithm.Calculate(aTrack.GetWeight(),
00280                                              fWeightWindowStore.GetLowerWeight(postCell,
00281                                                                                aTrack.GetKineticEnergy()));
00282           fPostStepAction->DoIt(aTrack, fParticleChange, nw);
00283         }
00284     }
00285   }
00286   return fParticleChange;
00287 }

G4double G4WeightWindowProcess::PostStepGetPhysicalInteractionLength ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
) [virtual]

Implements G4VProcess.

Definition at line 186 of file G4WeightWindowProcess.cc.

References DBL_MAX, and Forced.

00189 {
00190 //   *condition = Forced;
00191 //   return kInfinity;
00192 
00193 //  *condition = StronglyForced;
00194   *condition = Forced;
00195   return DBL_MAX;
00196 }

void G4WeightWindowProcess::SetParallelWorld ( G4VPhysicalVolume parallelWorld  ) 

Definition at line 129 of file G4WeightWindowProcess.cc.

References G4VPhysicalVolume::GetName(), and G4TransportationManager::GetNavigator().

00130 {
00131 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00132 // Get pointer of navigator
00133 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00134   fGhostWorldName = parallelWorld->GetName();
00135   fGhostWorld = parallelWorld;
00136   fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
00137 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00138 }

void G4WeightWindowProcess::SetParallelWorld ( G4String  parallelWorldName  ) 

Definition at line 117 of file G4WeightWindowProcess.cc.

References G4TransportationManager::GetNavigator(), and G4TransportationManager::GetParallelWorld().

Referenced by G4WeightWindowConfigurator::Configure().

00118 {
00119 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00120 // Get pointers of the parallel world and its navigator
00121 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00122   fGhostWorldName = parallelWorldName;
00123   fGhostWorld = fTransportationManager->GetParallelWorld(fGhostWorldName);
00124   fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
00125 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00126 }

void G4WeightWindowProcess::StartTracking ( G4Track  )  [virtual]

Reimplemented from G4VProcess.

Definition at line 145 of file G4WeightWindowProcess.cc.

References G4TransportationManager::ActivateNavigator(), G4PathFinder::CreateTouchableHandle(), FatalException, G4Exception(), G4Track::GetMomentumDirection(), G4Track::GetPosition(), G4PathFinder::PrepareNewTrack(), and G4StepPoint::SetTouchableHandle().

00146 {
00147 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00148 // Activate navigator and get the navigator ID
00149 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00150 // G4cout << " G4ParallelWorldScoringProcess::StartTracking" << G4endl;
00151 
00152   if(paraflag) {
00153     if(fGhostNavigator)
00154       { fNavigatorID = fTransportationManager->ActivateNavigator(fGhostNavigator); }
00155     else
00156       {
00157         G4Exception("G4WeightWindowProcess::StartTracking",
00158                     "ProcParaWorld000",FatalException,
00159                     "G4WeightWindowProcess is used for tracking without having a parallel world assigned");
00160       }
00161 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00162 
00163 // G4cout << "G4ParallelWorldScoringProcess::StartTracking <<<<<<<<<<<<<<<<<< " << G4endl;
00164 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00165 // Let PathFinder initialize
00166 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00167     fPathFinder->PrepareNewTrack(trk->GetPosition(),trk->GetMomentumDirection());
00168 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00169 
00170 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00171 // Setup initial touchables for the first step
00172 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00173     fOldGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
00174     fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
00175     fNewGhostTouchable = fOldGhostTouchable;
00176     fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
00177 
00178   // Initialize 
00179     fGhostSafety = -1.;
00180     fOnBoundary = false;
00181   }
00182 }


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