Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4UserSpecialCuts Class Reference

#include <G4UserSpecialCuts.hh>

+ Inheritance 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
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
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
 

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::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}
@ fGeneral
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
static G4LossTableManager * Instance()
G4int verboseLevel
Definition: G4VProcess.hh:368
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

◆ ~G4UserSpecialCuts()

G4UserSpecialCuts::~G4UserSpecialCuts ( )
virtual

Definition at line 68 of file G4UserSpecialCuts.cc.

69{
70}

Member Function Documentation

◆ AlongStepDoIt()

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

Implements G4VProcess.

Definition at line 90 of file G4UserSpecialCuts.hh.

93 {return 0;};

◆ AlongStepGetPhysicalInteractionLength()

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; };

◆ AtRestDoIt()

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

Implements G4VProcess.

Definition at line 75 of file G4UserSpecialCuts.hh.

78 {return 0;};

◆ AtRestGetPhysicalInteractionLength()

virtual G4double G4UserSpecialCuts::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition  
)
inlinevirtual

Implements G4VProcess.

Definition at line 69 of file G4UserSpecialCuts.hh.

72 { return -1.0; };

◆ PostStepDoIt()

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

Implements G4VProcess.

Definition at line 139 of file G4UserSpecialCuts.cc.

144{
149 return &aParticleChange;
150}
@ fStopAndKill
void ProposeEnergy(G4double finalEnergy)
virtual void Initialize(const G4Track &)
G4double GetKineticEnergy() const
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289

◆ PostStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 81 of file G4UserSpecialCuts.cc.

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

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