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

#include <G4UCNAbsorption.hh>

+ Inheritance diagram for G4UCNAbsorption:

Public Member Functions

 G4UCNAbsorption (const G4String &processName="UCNAbsorption", G4ProcessType type=fUCN)
 
virtual ~G4UCNAbsorption ()
 
G4bool IsApplicable (const G4ParticleDefinition &aParticleType)
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *condition)
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
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 &)
 
virtual G4double GetCrossSection (const G4double, const G4MaterialCutsCouple *)
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool 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 const G4VProcessGetCreatorProcess () const
 
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
 
virtual void ProcessDescription (std::ostream &outfile) 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)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 65 of file G4UCNAbsorption.hh.

Constructor & Destructor Documentation

◆ G4UCNAbsorption()

G4UCNAbsorption::G4UCNAbsorption ( const G4String processName = "UCNAbsorption",
G4ProcessType  type = fUCN 
)

Definition at line 73 of file G4UCNAbsorption.cc.

74 : G4VDiscreteProcess(processName, type)
75{
76 if (verboseLevel>0) G4cout << GetProcessName() << " is created " << G4endl;
77
79}
@ fUCNAbsorption
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4int verboseLevel
Definition: G4VProcess.hh:360
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386

◆ ~G4UCNAbsorption()

G4UCNAbsorption::~G4UCNAbsorption ( )
virtual

Definition at line 89 of file G4UCNAbsorption.cc.

89{}

Member Function Documentation

◆ GetMeanFreePath()

G4double G4UCNAbsorption::GetMeanFreePath ( const G4Track aTrack,
G4double  ,
G4ForceCondition condition 
)
virtual

Implements G4VDiscreteProcess.

Definition at line 119 of file G4UCNAbsorption.cc.

122{
123 G4double AttenuationLength = DBL_MAX;
124
125 const G4Material* aMaterial = aTrack.GetMaterial();
126 G4MaterialPropertiesTable* aMaterialPropertiesTable =
127 aMaterial->GetMaterialPropertiesTable();
128
129 G4double losscs = 0.0;
130 if (aMaterialPropertiesTable) {
131 losscs = aMaterialPropertiesTable->GetConstProperty("ABSCS");
132// if (losscs == 0.0)
133// G4cout << "No UCN Absorption length specified" << G4endl;
134 }
135// else G4cout << "No UCN Absorption length specified" << G4endl;
136
137 if (losscs) {
138
139 // Calculate a UCN absorption length for this cross section
140
141 // *** Thermal boost ***
142
143 // Prepare neutron
144
145 //G4double theA = aMaterial->GetElement(0)->GetN();
146 //G4double theZ = aMaterial->GetElement(0)->GetZ();
147
148 //G4ReactionProduct
149 // theNeutron(const_cast<G4ParticleDefinition *>(aTrack.GetDefinition()));
150 //theNeutron.SetMomentum(aTrack.GetMomentum());
151 //theNeutron.SetKineticEnergy(aTrack.GetKineticEnergy());
152 //G4ThreeVector neuVelo = theNeutron.GetMomentum()/
153 // aTrack.GetDefinition()->GetPDGMass());
154
155 // Prepare properly biased thermal nucleus
156
157 //G4double theA = aMaterial->GetElement(0)->GetN();
158 //G4double theZ = aMaterial->GetElement(0)->GetZ();
159
160 //G4double eps = 0.0001;
161
162 //G4double eleMass =
163 // G4NucleiPropertiesTable::
164 // GetNuclearMass(static_cast<G4int>(theZ+eps),
165 // static_cast<G4int>(theA+eps)))
166 // / G4Neutron::Neutron()->GetPDGMass();
167
168 //G4Nucleus aNuc;
169
170 //G4ReactionProduct aThermalNuc =
171 // aNuc.GetBiasedThermalNucleus(eleMass,
172 // neuVelo,
173 // aMaterial->GetTemperature());
174
175 // Boost to rest system and return
176
177 //G4ReactionProduct boosted;
178 //boosted.Lorentz(theNeutron, aThermalNuc);
179
180 //G4double vel = sqrt(2*boosted.GetKineticEnergy()/
181 // neutron_mass_c2*c_squared);
182
183 G4double density = aMaterial->GetTotNbOfAtomsPerVolume();
184
185 // Calculate cross section for a constant loss
186
187 G4double vel = aTrack.GetVelocity();
188
189 //G4cout << aTrack.GetVelocity()/meter*second << " "
190 // << vel/meter*second << "meters/second" << G4endl;
191
192 // Input data is normally taken from the website:
193 // http://rrdjazz.nist.gov/resources/n-lengths/list.html
194 // and coresponds to 2200 m/s fast neutrons
195
196 G4double crossect = losscs*barn*2200.*meter/second/vel;
197
198 // In principle, if one asks for the MaterialProperty incoherent cross
199 // section, one could put the formula for inelastic up scattering here
200 // and add the cross section to the absorption
201
202 // sigma inelastic = ... ignatovic, p. 174.
203
204 // attenuation length in mm
205 AttenuationLength = 1./density/crossect;
206
207 if (verboseLevel>0) G4cout << "UCNABSORPTION with"
208 << " AttenuationLength: " << AttenuationLength/m << "m"
209 << " CrossSection: " << crossect/barn << "barn" << G4endl;
210 }
211
212 return AttenuationLength;
213}
double G4double
Definition: G4Types.hh:83
G4double GetConstProperty(const G4String &key) const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:251
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
G4double GetVelocity() const
G4Material * GetMaterial() const
#define DBL_MAX
Definition: templates.hh:62

◆ IsApplicable()

G4bool G4UCNAbsorption::IsApplicable ( const G4ParticleDefinition aParticleType)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 113 of file G4UCNAbsorption.hh.

114{
115 return ( &aParticleType == G4Neutron::NeutronDefinition() );
116}
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:98

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 99 of file G4UCNAbsorption.cc.

100{
102
104
105 if ( verboseLevel > 0 ) G4cout << "UCNABSORPTION at: "
106 << aTrack.GetProperTime()/s << "s, "
107 << aTrack.GetGlobalTime()/s << "s. "
108 << ", after track length " << aTrack.GetTrackLength()/cm << "cm, "
109 << "in volume "
111 << G4endl;
112
113 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
114}
@ fStopAndKill
void Initialize(const G4Track &) override
G4VPhysicalVolume * GetPhysicalVolume() const
G4StepPoint * GetPostStepPoint() const
G4double GetTrackLength() const
G4double GetGlobalTime() const
G4double GetProperTime() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void ProposeTrackStatus(G4TrackStatus status)
const G4String & GetName() const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:331

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