G4ImportanceProcess Class Reference

#include <G4ImportanceProcess.hh>

Inheritance diagram for G4ImportanceProcess:

G4VProcess G4VTrackTerminator

Public Member Functions

 G4ImportanceProcess (const G4VImportanceAlgorithm &aImportanceAlgorithm, const G4VIStore &aIstore, const G4VTrackTerminator *TrackTerminator, const G4String &aName="ImportanceProcess", G4bool para=false)
virtual ~G4ImportanceProcess ()
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 61 of file G4ImportanceProcess.hh.


Constructor & Destructor Documentation

G4ImportanceProcess::G4ImportanceProcess ( const G4VImportanceAlgorithm aImportanceAlgorithm,
const G4VIStore aIstore,
const G4VTrackTerminator TrackTerminator,
const G4String aName = "ImportanceProcess",
G4bool  para = false 
)

Definition at line 55 of file G4ImportanceProcess.cc.

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

00059  : G4VProcess(aName),
00060    fParticleChange(new G4ParticleChange),
00061    fImportanceAlgorithm(aImportanceAlgorithm),
00062    fIStore(aIstore),
00063    fPostStepAction(0),
00064    fGhostNavigator(0), fNavigatorID(-1), fFieldTrack('0'),
00065    paraflag(para)
00066   
00067 {
00068   if (TrackTerminator)
00069   {
00070     fPostStepAction = new G4SamplingPostStepAction(*TrackTerminator);
00071   }
00072   else
00073   {
00074     fPostStepAction = new G4SamplingPostStepAction(*this);
00075   }
00076   if (!fParticleChange)
00077   {
00078     G4Exception("G4ImportanceProcess::G4ImportanceProcess()",
00079                 "FatalError", FatalException,
00080                 "Failed allocation of G4ParticleChange !");
00081   }
00082   G4VProcess::pParticleChange = fParticleChange;
00083 
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   G4cout << " importance process paraflag is: " << paraflag << G4endl;
00098 
00099 }

G4ImportanceProcess::~G4ImportanceProcess (  )  [virtual]

Definition at line 101 of file G4ImportanceProcess.cc.

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


Member Function Documentation

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

Implements G4VProcess.

Definition at line 442 of file G4ImportanceProcess.cc.

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

00443 {
00444   // Dummy ParticleChange ie: does nothing
00445   // Expecting G4Transportation to move the track
00446   //AH  G4cout << " along step do it " << G4endl;
00447   pParticleChange->Initialize(aTrack);
00448   return pParticleChange; 
00449   //  return 0;
00450 }

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

Implements G4VProcess.

Definition at line 347 of file G4ImportanceProcess.cc.

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

00350 {
00351 
00352   //AH  G4cout << " along step physical interaction length " << G4endl;
00353 
00354   if(paraflag) {
00355     static G4FieldTrack endTrack('0');
00356     static ELimited eLimited;
00357     
00358     *selection = NotCandidateForSelection;
00359     G4double returnedStep = DBL_MAX;
00360     
00361     if (previousStepSize > 0.)
00362       { fGhostSafety -= previousStepSize; }
00363     //  else
00364     //  { fGhostSafety = -1.; }
00365     if (fGhostSafety < 0.) fGhostSafety = 0.0;
00366     
00367     // ------------------------------------------
00368     // Determination of the proposed STEP LENGTH:
00369     // ------------------------------------------
00370     if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
00371       {
00372         // I have no chance to limit
00373         returnedStep = currentMinimumStep;
00374         fOnBoundary = false;
00375         proposedSafety = fGhostSafety - currentMinimumStep;
00376         //AH    G4cout << " step not limited, why? " << G4endl;
00377       }
00378     else // (currentMinimumStep > fGhostSafety: I may limit the Step)
00379       {
00380         G4FieldTrackUpdator::Update(&fFieldTrack,&track);
00381         //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00382         // ComputeStep
00383         //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00384         returnedStep
00385           = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
00386                                      track.GetCurrentStepNumber(),fGhostSafety,eLimited,
00387                                      endTrack,track.GetVolume());
00388         //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00389         if(eLimited == kDoNot)
00390           {
00391             //AH            G4cout << " computestep came back with not-boundary " << G4endl;
00392             // Track is not on the boundary
00393             fOnBoundary = false;
00394             fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition());
00395           }
00396         else
00397           {
00398             // Track is on the boundary
00399             //AH            G4cout << " FOUND A BOUNDARY ! " << track.GetCurrentStepNumber() << G4endl;
00400             fOnBoundary = true;
00401             // proposedSafety = fGhostSafety;
00402           }
00403         proposedSafety = fGhostSafety;
00404         if(eLimited == kUnique || eLimited == kSharedOther) {
00405           *selection = CandidateForSelection;
00406         }else if (eLimited == kSharedTransport) { 
00407           returnedStep *= (1.0 + 1.0e-9);  
00408           // Expand to disable its selection in Step Manager comparison
00409         }
00410         
00411       }
00412 
00413   // ----------------------------------------------
00414   // Returns the fGhostSafety as the proposedSafety
00415   // The SteppingManager will take care of keeping
00416   // the smallest one.
00417   // ----------------------------------------------
00418     return returnedStep;
00419 
00420   } else {
00421 
00422     return DBL_MAX;
00423 
00424   }
00425 
00426 }

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

Implements G4VProcess.

Definition at line 436 of file G4ImportanceProcess.cc.

00437 {
00438   return 0;
00439 }

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

Implements G4VProcess.

Definition at line 429 of file G4ImportanceProcess.cc.

00431 {
00432   return -1.0;
00433 }

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

Implements G4VTrackTerminator.

Definition at line 341 of file G4ImportanceProcess.cc.

References G4VProcess::theProcessName.

00342 {
00343   return theProcessName;
00344 }

void G4ImportanceProcess::KillTrack (  )  const [virtual]

Implements G4VTrackTerminator.

Definition at line 336 of file G4ImportanceProcess.cc.

References fStopAndKill, and G4VParticleChange::ProposeTrackStatus().

00337 {
00338   fParticleChange->ProposeTrackStatus(fStopAndKill);
00339 }

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

Implements G4VProcess.

Definition at line 201 of file G4ImportanceProcess.cc.

References G4PathFinder::CreateTouchableHandle(), G4SamplingPostStepAction::DoIt(), fGeomBoundary, fStopAndKill, G4cout, G4endl, G4VIStore::GetImportance(), G4StepPoint::GetPhysicalVolume(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4VTouchable::GetReplicaNumber(), G4Step::GetStepLength(), G4StepPoint::GetStepStatus(), G4StepPoint::GetTouchable(), G4StepPoint::GetTouchableHandle(), G4Track::GetTrackStatus(), G4Track::GetWeight(), G4ParticleChange::Initialize(), G4VTrackTerminator::kCarTolerance, and G4StepPoint::SetTouchableHandle().

00203 {
00204   fParticleChange->Initialize(aTrack);
00205 
00206   if(paraflag) {
00207     fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
00208     //xbug?    fOnBoundary = false;
00209     CopyStep(aStep);
00210     
00211     if(fOnBoundary)
00212       {
00213 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00214 // Locate the point and get new touchable
00215 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00216   //??  fPathFinder->Locate(step.GetPostStepPoint()->GetPosition(),
00217   //??                      step.GetPostStepPoint()->GetMomentumDirection());
00218         fNewGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
00219         //AH    G4cout << " on boundary " << G4endl;
00220 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00221       }
00222     else
00223       {
00224         // Do I need this ??????????????????????????????????????????????????????????
00225         // fGhostNavigator->LocateGlobalPointWithinVolume(track.GetPosition());
00226         // ?????????????????????????????????????????????????????????????????????????
00227         
00228         // fPathFinder->ReLocate(track.GetPosition());
00229         
00230         // reuse the touchable
00231         fNewGhostTouchable = fOldGhostTouchable;
00232         //AH    G4cout << " NOT on boundary " << G4endl;
00233       }
00234     
00235     fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
00236     fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
00237     
00238 //   if ( (aStep.GetPostStepPoint()->GetStepStatus() == fGeomBoundary)
00239 //     && (aStep.GetStepLength() > kCarTolerance) )
00240 //   {
00241 //     if (aTrack.GetTrackStatus()==fStopAndKill)
00242 //     {
00243 //       G4cout << "WARNING - G4ImportanceProcess::PostStepDoIt()"
00244 //              << "          StopAndKill track." << G4endl;
00245 //     }
00246 
00247 //     G4StepPoint *prepoint  = aStep.GetPreStepPoint();
00248 //     G4StepPoint *postpoint = aStep.GetPostStepPoint();
00249 //     G4GeometryCell prekey(*(prepoint->GetPhysicalVolume()), 
00250 //                          prepoint->GetTouchable()->GetReplicaNumber());
00251 //     G4GeometryCell postkey(*(postpoint->GetPhysicalVolume()), 
00252 //                           postpoint->GetTouchable()->GetReplicaNumber());
00253 
00254 
00255     if ( (fGhostPostStepPoint->GetStepStatus() == fGeomBoundary)
00256          && (aStep.GetStepLength() > kCarTolerance) )
00257       {
00258         if (aTrack.GetTrackStatus()==fStopAndKill)
00259           {
00260             G4cout << "WARNING - G4ImportanceProcess::PostStepDoIt()"
00261                    << "          StopAndKill track. on boundary" << G4endl;
00262           }
00263         
00264         G4GeometryCell prekey(*(fGhostPreStepPoint->GetPhysicalVolume()), 
00265                               fGhostPreStepPoint->GetTouchable()->GetReplicaNumber());
00266         G4GeometryCell postkey(*(fGhostPostStepPoint->GetPhysicalVolume()), 
00267                                fGhostPostStepPoint->GetTouchable()->GetReplicaNumber());
00268         
00269         //AH
00270         /*
00271         G4cout << G4endl;
00272         G4cout << G4endl;
00273         G4cout << " inside parallel importance process " << aTrack.GetCurrentStepNumber() << G4endl;
00274         G4cout << G4endl;
00275         G4cout << G4endl;
00276         G4cout << " prekey: " << fGhostPreStepPoint->GetPhysicalVolume()->GetName() << " replica:" 
00277                << fGhostPreStepPoint->GetTouchable()->GetReplicaNumber() << G4endl;
00278         G4cout << " prekey ISTORE: " << fIStore.GetImportance(prekey) << G4endl;
00279         G4cout << " postkey: " << G4endl;
00280         G4cout << " postkey ISTORE: " << fIStore.GetImportance(postkey) << G4endl;
00281         */
00282         //AH
00283         G4Nsplit_Weight nw = fImportanceAlgorithm.
00284           Calculate(fIStore.GetImportance(prekey),
00285                     fIStore.GetImportance(postkey), 
00286                     aTrack.GetWeight());
00287         //AH
00288         /*
00289         G4cout << " prekey weight: " << fIStore.GetImportance(prekey) 
00290                << " postkey weight: " << fIStore.GetImportance(postkey) 
00291                << " split weight: " << nw << G4endl;
00292         G4cout << " before poststepaction " << G4endl;
00293         */
00294         //AH
00295         fPostStepAction->DoIt(aTrack, fParticleChange, nw);
00296         //AH    G4cout << " after post step do it " << G4endl;
00297       }
00298   } else {
00299     if ( (aStep.GetPostStepPoint()->GetStepStatus() == fGeomBoundary)
00300          && (aStep.GetStepLength() > kCarTolerance) )
00301       {
00302         //AH    G4cout << " inside non-parallel importance process " << G4endl;
00303         if (aTrack.GetTrackStatus()==fStopAndKill)
00304           {
00305             G4cout << "WARNING - G4ImportanceProcess::PostStepDoIt()"
00306                    << "          StopAndKill track. on boundary non-parallel" << G4endl;
00307           }
00308         
00309         G4StepPoint *prepoint  = aStep.GetPreStepPoint();
00310         G4StepPoint *postpoint = aStep.GetPostStepPoint();
00311         
00312         G4GeometryCell prekey(*(prepoint->GetPhysicalVolume()), 
00313                               prepoint->GetTouchable()->GetReplicaNumber());
00314         G4GeometryCell postkey(*(postpoint->GetPhysicalVolume()), 
00315                                postpoint->GetTouchable()->GetReplicaNumber());
00316         
00317         G4Nsplit_Weight nw = fImportanceAlgorithm.
00318           Calculate(fIStore.GetImportance(prekey),
00319                     fIStore.GetImportance(postkey), 
00320                     aTrack.GetWeight());
00321         //AH
00322         /*
00323         G4cout << " prekey weight: " << fIStore.GetImportance(prekey) 
00324                << " postkey weight: " << fIStore.GetImportance(postkey) 
00325                << " split weight: " << nw << G4endl;
00326         G4cout << " before poststepaction 2 " << G4endl;
00327         */
00328         //AH
00329         fPostStepAction->DoIt(aTrack, fParticleChange, nw);
00330         //AH    G4cout << " after poststepaction 2 " << G4endl;
00331       }
00332   }
00333   return fParticleChange;
00334 }

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

Implements G4VProcess.

Definition at line 188 of file G4ImportanceProcess.cc.

References DBL_MAX, and Forced.

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

void G4ImportanceProcess::SetParallelWorld ( G4VPhysicalVolume parallelWorld  ) 

Definition at line 130 of file G4ImportanceProcess.cc.

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

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

void G4ImportanceProcess::SetParallelWorld ( G4String  parallelWorldName  ) 

Definition at line 118 of file G4ImportanceProcess.cc.

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

Referenced by G4ImportanceConfigurator::Configure().

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

void G4ImportanceProcess::StartTracking ( G4Track  )  [virtual]

Reimplemented from G4VProcess.

Definition at line 146 of file G4ImportanceProcess.cc.

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

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


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