G4UserSpecialCuts Class Reference

#include <G4UserSpecialCuts.hh>

Inheritance diagram for G4UserSpecialCuts:

G4VProcess

Public Member Functions

 G4UserSpecialCuts (const G4String &processName="UserSpecialCut")
virtual ~G4UserSpecialCuts ()
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)

Detailed Description

Definition at line 47 of file G4UserSpecialCuts.hh.


Constructor & Destructor Documentation

G4UserSpecialCuts::G4UserSpecialCuts ( const G4String processName = "UserSpecialCut"  ) 

Definition at line 53 of file G4UserSpecialCuts.cc.

References G4cout, G4endl, G4VProcess::GetProcessName(), G4LossTableManager::Instance(), G4VProcess::SetProcessSubType(), USER_SPECIAL_CUTS, and G4VProcess::verboseLevel.

00054   : G4VProcess(aName, fGeneral  )
00055 {
00056   // set Process Sub Type
00057   SetProcessSubType(static_cast<int>(USER_SPECIAL_CUTS));
00058 
00059   if (verboseLevel>0)
00060     {
00061       G4cout << GetProcessName() << " is created "<< G4endl;
00062     }
00063   theLossTableManager = G4LossTableManager::Instance();
00064 }

G4UserSpecialCuts::~G4UserSpecialCuts (  )  [virtual]

Definition at line 68 of file G4UserSpecialCuts.cc.

00069 {
00070 }


Member Function Documentation

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

Implements G4VProcess.

Definition at line 90 of file G4UserSpecialCuts.hh.

00093                               {return 0;};

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

Implements G4VProcess.

Definition at line 81 of file G4UserSpecialCuts.hh.

00087                              { return -1.0; };

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

Implements G4VProcess.

Definition at line 75 of file G4UserSpecialCuts.hh.

00078                              {return 0;};

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

Implements G4VProcess.

Definition at line 69 of file G4UserSpecialCuts.hh.

00072                              { return -1.0; };

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

Implements G4VProcess.

Definition at line 139 of file G4UserSpecialCuts.cc.

References G4VProcess::aParticleChange, fStopAndKill, G4Track::GetKineticEnergy(), G4ParticleChange::Initialize(), G4ParticleChange::ProposeEnergy(), G4VParticleChange::ProposeLocalEnergyDeposit(), and G4VParticleChange::ProposeTrackStatus().

00144 {
00145    aParticleChange.Initialize(aTrack);
00146    aParticleChange.ProposeEnergy(0.) ;
00147    aParticleChange.ProposeLocalEnergyDeposit(aTrack.GetKineticEnergy()) ;
00148    aParticleChange.ProposeTrackStatus(fStopAndKill);
00149    return &aParticleChange;
00150 }

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

Implements G4VProcess.

Definition at line 82 of file G4UserSpecialCuts.cc.

References DBL_MAX, DBL_MIN, G4Track::GetKineticEnergy(), G4VPhysicalVolume::GetLogicalVolume(), G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), G4LossTableManager::GetRange(), G4LogicalVolume::GetUserLimits(), G4UserLimits::GetUserMaxTime(), G4UserLimits::GetUserMaxTrackLength(), G4UserLimits::GetUserMinRange(), G4Track::GetVolume(), and NotForced.

00085 {
00086   // condition is set to "Not Forced"
00087   *condition = NotForced;
00088 
00089    G4double ProposedStep = DBL_MAX;
00090    G4UserLimits* pUserLimits =
00091                  aTrack.GetVolume()->GetLogicalVolume()->GetUserLimits();
00092    if (pUserLimits)
00093    {
00094      // check max kinetic energy first
00095      //
00096      G4double Ekine = aTrack.GetKineticEnergy();
00097      if(Ekine <= pUserLimits->GetUserMinEkine(aTrack)) { return 0.; }
00098 
00099      // max track length
00100      //
00101      ProposedStep = (pUserLimits->GetUserMaxTrackLength(aTrack)
00102                    - aTrack.GetTrackLength());
00103      if (ProposedStep < 0.) { return 0.; }
00104 
00105      // max time limit
00106      //
00107      G4double tlimit = pUserLimits->GetUserMaxTime(aTrack);
00108      if(tlimit < DBL_MAX) {
00109        G4double beta  = (aTrack.GetDynamicParticle()->GetTotalMomentum())
00110          /(aTrack.GetTotalEnergy());
00111        G4double dTime = (tlimit - aTrack.GetGlobalTime());
00112        G4double temp  = beta*c_light*dTime;
00113        if (temp < 0.) { return 0.; }
00114        if (ProposedStep > temp) { ProposedStep = temp; }
00115      }
00116                  
00117      // min remaining range 
00118      // (only for charged particle except for chargedGeantino)
00119      //
00120      G4double Rmin = pUserLimits->GetUserMinRange(aTrack);
00121      if (Rmin > DBL_MIN) {
00122        G4ParticleDefinition* Particle = aTrack.GetDefinition();
00123        if ( (Particle->GetPDGCharge() != 0.) && (Particle->GetPDGMass() > 0.0))
00124        {
00125          const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
00126          G4double RangeNow = theLossTableManager->GetRange(Particle,Ekine,couple);
00127          G4double temp = RangeNow - Rmin;
00128          if (temp < 0.) { return 0.; }
00129          if (ProposedStep > temp) { ProposedStep = temp; }
00130        }         
00131      }
00132    }
00133    return ProposedStep;
00134 }


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