G4ErrorTrackLengthTarget Class Reference

#include <G4ErrorTrackLengthTarget.hh>

Inheritance diagram for G4ErrorTrackLengthTarget:

G4VDiscreteProcess G4ErrorTarget G4VProcess

Public Member Functions

 G4ErrorTrackLengthTarget (const G4double maxTrkLength)
virtual ~G4ErrorTrackLengthTarget ()
virtual G4double GetDistanceFromPoint (const G4ThreeVector &, const G4ThreeVector &) const
virtual G4double GetDistanceFromPoint (const G4ThreeVector &) const
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual G4double GetMeanFreePath (const class G4Track &track, G4double, G4ForceCondition *)
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
virtual void Dump (const G4String &msg) const

Detailed Description

Definition at line 54 of file G4ErrorTrackLengthTarget.hh.


Constructor & Destructor Documentation

G4ErrorTrackLengthTarget::G4ErrorTrackLengthTarget ( const G4double  maxTrkLength  ) 

Definition at line 47 of file G4ErrorTrackLengthTarget.cc.

References G4ProcessManager::AddDiscreteProcess(), G4ErrorTarget_TrkL, G4Exception(), G4ParticleTable::GetIterator(), G4ParticleDefinition::GetParticleName(), G4ParticleTable::GetParticleTable(), G4ProcessManager::GetProcessList(), G4ParticleDefinition::GetProcessManager(), G4VProcess::GetProcessName(), G4ParticleDefinition::IsShortLived(), G4ProcessManager::RemoveProcess(), G4ParticleTableIterator< K, V >::reset(), RunMustBeAborted, G4ProcessVector::size(), G4ErrorTarget::theType, and G4ParticleTableIterator< K, V >::value().

00048   : G4VDiscreteProcess ("G4ErrorTrackLengthTarget"),
00049     theMaximumTrackLength( maxTrkLength )
00050 {
00051   theType = G4ErrorTarget_TrkL;
00052 
00053    G4ParticleTable::G4PTblDicIterator* theParticleIterator
00054      = G4ParticleTable::GetParticleTable()->GetIterator();
00055 
00056   // loop over all particles in G4ParticleTable
00057 
00058   theParticleIterator->reset();
00059   while( (*theParticleIterator)() )
00060   {
00061     G4ParticleDefinition* particle = theParticleIterator->value();
00062     G4ProcessManager* pmanager = particle->GetProcessManager();
00063     if (!particle->IsShortLived())
00064     {
00065       // Add transportation process for all particles other than  "shortlived"
00066       if ( pmanager == 0)
00067       {
00068         // Error !! no process manager
00069         G4String particleName = particle->GetParticleName();
00070         G4Exception("G4ErrorTrackLengthTarget::G4ErrorTrackLengthTarget",
00071                     "No process manager", RunMustBeAborted, particleName );
00072       }
00073       else
00074       {
00075         G4ProcessVector* procvec = pmanager->GetProcessList();
00076         size_t isiz = procvec->size();
00077 
00078         for( size_t ii=0; ii < isiz; ii++ )
00079         {
00080           if( ((*procvec)[ii])->GetProcessName() == "G4ErrorTrackLengthTarget")
00081           {
00082             pmanager->RemoveProcess( (*procvec)[ii] );
00083           }
00084         }
00085         pmanager ->AddDiscreteProcess(this,4);
00086         isiz = procvec->size();
00087       }
00088     }
00089     else
00090     {
00091       // shortlived particle case
00092     }
00093   }
00094 }

virtual G4ErrorTrackLengthTarget::~G4ErrorTrackLengthTarget (  )  [inline, virtual]

Definition at line 61 of file G4ErrorTrackLengthTarget.hh.

00061 {}


Member Function Documentation

void G4ErrorTrackLengthTarget::Dump ( const G4String msg  )  const [virtual]

Implements G4ErrorTarget.

Definition at line 132 of file G4ErrorTrackLengthTarget.cc.

References G4cout, and G4endl.

00133 {
00134   G4cout << msg << "G4ErrorTrackLengthTarget: max track length = "
00135          << theMaximumTrackLength << G4endl;
00136 }

virtual G4double G4ErrorTrackLengthTarget::GetDistanceFromPoint ( const G4ThreeVector  )  const [inline, virtual]

Reimplemented from G4ErrorTarget.

Definition at line 69 of file G4ErrorTrackLengthTarget.hh.

References DBL_MAX.

00070     { return DBL_MAX; } 

virtual G4double G4ErrorTrackLengthTarget::GetDistanceFromPoint ( const G4ThreeVector ,
const G4ThreeVector  
) const [inline, virtual]

Reimplemented from G4ErrorTarget.

Definition at line 66 of file G4ErrorTrackLengthTarget.hh.

References DBL_MAX.

00068     { return DBL_MAX; }   

G4double G4ErrorTrackLengthTarget::GetMeanFreePath ( const class G4Track track,
G4double  ,
G4ForceCondition  
) [virtual]

Definition at line 109 of file G4ErrorTrackLengthTarget.cc.

References G4cout, G4endl, and G4ErrorPropagatorData::verbose().

Referenced by PostStepGetPhysicalInteractionLength().

00110 {
00111 #ifdef G4VERBOSE
00112   if(G4ErrorPropagatorData::verbose() >= 3 )
00113   { 
00114     G4cout << " G4ErrorTrackLengthTarget::GetMeanFreePath "
00115            << theMaximumTrackLength - track.GetTrackLength() << G4endl;
00116   }
00117 #endif
00118 
00119   return theMaximumTrackLength - track.GetTrackLength();
00120 }

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

Reimplemented from G4VDiscreteProcess.

Definition at line 124 of file G4ErrorTrackLengthTarget.cc.

References G4VParticleChange::Initialize().

00125 {
00126   theParticleChange.Initialize(aTrack);
00127   return &theParticleChange;
00128 }

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

Reimplemented from G4VDiscreteProcess.

Definition at line 99 of file G4ErrorTrackLengthTarget.cc.

References GetMeanFreePath(), and NotForced.

00101 {
00102   *condition = NotForced;
00103   return GetMeanFreePath( track, 0., condition );
00104 }


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