Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4UserSpecialCuts Class Reference

#include <G4UserSpecialCuts.hh>

Inheritance diagram for G4UserSpecialCuts:
Collaboration diagram for G4UserSpecialCuts:

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 &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

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.

54  : G4VProcess(aName, fGeneral )
55 {
56  // set Process Sub Type
57  SetProcessSubType(static_cast<int>(USER_SPECIAL_CUTS));
58 
59  if (verboseLevel>0)
60  {
61  G4cout << GetProcessName() << " is created "<< G4endl;
62  }
63  theLossTableManager = G4LossTableManager::Instance();
64 }
static G4LossTableManager * Instance()
G4int verboseLevel
Definition: G4VProcess.hh:368
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition: G4VProcess.cc:52
G4GLOB_DLL std::ostream G4cout
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

G4UserSpecialCuts::~G4UserSpecialCuts ( )
virtual

Definition at line 68 of file G4UserSpecialCuts.cc.

69 {
70 }

Member Function Documentation

virtual G4VParticleChange* G4UserSpecialCuts::AlongStepDoIt ( const G4Track ,
const G4Step  
)
inlinevirtual

Implements G4VProcess.

Definition at line 90 of file G4UserSpecialCuts.hh.

93  {return 0;};
virtual G4double G4UserSpecialCuts::AlongStepGetPhysicalInteractionLength ( const G4Track ,
G4double  ,
G4double  ,
G4double ,
G4GPILSelection  
)
inlinevirtual

Implements G4VProcess.

Definition at line 81 of file G4UserSpecialCuts.hh.

87  { return -1.0; };
virtual G4VParticleChange* G4UserSpecialCuts::AtRestDoIt ( const G4Track ,
const G4Step  
)
inlinevirtual

Implements G4VProcess.

Definition at line 75 of file G4UserSpecialCuts.hh.

78  {return 0;};
virtual G4double G4UserSpecialCuts::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition  
)
inlinevirtual

Implements G4VProcess.

Definition at line 69 of file G4UserSpecialCuts.hh.

72  { return -1.0; };
G4VParticleChange * G4UserSpecialCuts::PostStepDoIt ( const G4Track aTrack,
const G4Step  
)
virtual

Implements G4VProcess.

Definition at line 140 of file G4UserSpecialCuts.cc.

145 {
146  aParticleChange.Initialize(aTrack);
150  return &aParticleChange;
151 }
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double GetKineticEnergy() const
virtual void Initialize(const G4Track &)
void ProposeEnergy(G4double finalEnergy)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
void ProposeTrackStatus(G4TrackStatus status)

Here is the call graph for this function:

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

Implements G4VProcess.

Definition at line 83 of file G4UserSpecialCuts.cc.

86 {
87  // condition is set to "Not Forced"
89 
90  G4double ProposedStep = DBL_MAX;
91  G4UserLimits* pUserLimits =
92  aTrack.GetVolume()->GetLogicalVolume()->GetUserLimits();
93  if (pUserLimits)
94  {
95  // check max kinetic energy first
96  //
97  G4double Ekine = aTrack.GetKineticEnergy();
98  if(Ekine <= pUserLimits->GetUserMinEkine(aTrack)) { return 0.; }
99 
100  // max track length
101  //
102  ProposedStep = (pUserLimits->GetUserMaxTrackLength(aTrack)
103  - aTrack.GetTrackLength());
104  if (ProposedStep < 0.) { return 0.; }
105 
106  // max time limit
107  //
108  G4double tlimit = pUserLimits->GetUserMaxTime(aTrack);
109  if(tlimit < DBL_MAX) {
110  G4double beta = (aTrack.GetDynamicParticle()->GetTotalMomentum())
111  /(aTrack.GetTotalEnergy());
112  G4double dTime = (tlimit - aTrack.GetGlobalTime());
113  G4double temp = beta*c_light*dTime;
114  if (temp < 0.) { return 0.; }
115  if (ProposedStep > temp) { ProposedStep = temp; }
116  }
117 
118  // min remaining range
119  // (only for charged particle except for chargedGeantino)
120  //
121  G4double Rmin = pUserLimits->GetUserMinRange(aTrack);
122  if (Rmin > DBL_MIN) {
123  G4ParticleDefinition* Particle = aTrack.GetDefinition();
124  if ( (Particle->GetPDGCharge() != 0.) && (Particle->GetPDGMass() > 0.0))
125  {
126  const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
127  G4double RangeNow = theLossTableManager->GetRange(Particle,Ekine,couple);
128  G4double temp = RangeNow - Rmin;
129  if (temp < 0.) { return 0.; }
130  if (ProposedStep > temp) { ProposedStep = temp; }
131  }
132  }
133  }
134  return ProposedStep;
135 }
G4double condition(const G4ErrorSymMatrix &m)
virtual G4double GetUserMaxTrackLength(const G4Track &)
virtual G4double GetUserMaxTime(const G4Track &)
G4double GetPDGMass() const
static constexpr double c_light
#define DBL_MIN
Definition: templates.hh:75
virtual G4double GetUserMinRange(const G4Track &)
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)

Here is the call graph for this function:


The documentation for this class was generated from the following files: