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

#include <G4UCNMultiScattering.hh>

+ Inheritance diagram for G4UCNMultiScattering:

Public Member Functions

 G4UCNMultiScattering (const G4String &processName="UCNMultiScattering", G4ProcessType type=fUCN)
 
virtual ~G4UCNMultiScattering ()
 
G4bool IsApplicable (const G4ParticleDefinition &aParticleType)
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *condition)
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
G4ThreeVector Scatter ()
 
- 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 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
 
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 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)
 
- Protected Member Functions inherited from G4VDiscreteProcess
- 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 G4UCNMultiScattering.hh.

Constructor & Destructor Documentation

◆ G4UCNMultiScattering()

G4UCNMultiScattering::G4UCNMultiScattering ( const G4String & processName = "UCNMultiScattering",
G4ProcessType type = fUCN )

Definition at line 68 of file G4UCNMultiScattering.cc.

70 : G4VDiscreteProcess(processName, type)
71{
72 if (verboseLevel>0) G4cout << GetProcessName() << " is created " << G4endl;
73
75}
@ fUCNMultiScattering
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4int verboseLevel
void SetProcessSubType(G4int)
const G4String & GetProcessName() const

◆ ~G4UCNMultiScattering()

G4UCNMultiScattering::~G4UCNMultiScattering ( )
virtual

Definition at line 85 of file G4UCNMultiScattering.cc.

85{}

Member Function Documentation

◆ GetMeanFreePath()

G4double G4UCNMultiScattering::GetMeanFreePath ( const G4Track & aTrack,
G4double ,
G4ForceCondition * condition )
virtual

Implements G4VDiscreteProcess.

Definition at line 117 of file G4UCNMultiScattering.cc.

120{
121 G4double AttenuationLength = DBL_MAX;
122
123 const G4Material* aMaterial = aTrack.GetMaterial();
124 G4MaterialPropertiesTable* aMaterialPropertiesTable =
125 aMaterial->GetMaterialPropertiesTable();
126
127 G4double crossect = 0.0;
128 if (aMaterialPropertiesTable) {
129 crossect = aMaterialPropertiesTable->GetConstProperty("SCATCS");
130// if(crossect == 0.0)
131// G4cout << "No UCN MultiScattering length specified" << G4endl;
132 }
133// else G4cout << "No UCN MultiScattering length specified" << G4endl;
134
135 if (crossect) {
136
137 // Calculate a UCN MultiScattering length for this cross section
138
139 G4double density = aMaterial->GetTotNbOfAtomsPerVolume();
140
141 crossect *= barn;
142
143 // attenuation length in mm
144 AttenuationLength = 1./density/crossect;
145 }
146
147 return AttenuationLength;
148}
double G4double
Definition G4Types.hh:83
G4double GetConstProperty(const G4String &key) const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
G4double GetTotNbOfAtomsPerVolume() const
G4Material * GetMaterial() const
#define DBL_MAX
Definition templates.hh:62

◆ IsApplicable()

G4bool G4UCNMultiScattering::IsApplicable ( const G4ParticleDefinition & aParticleType)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 115 of file G4UCNMultiScattering.hh.

116{
117 return ( &aParticleType == G4Neutron::NeutronDefinition() );
118}
static G4Neutron * NeutronDefinition()
Definition G4Neutron.cc:96

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 95 of file G4UCNMultiScattering.cc.

96{
98
99 if ( verboseLevel > 0 ) G4cout << "UCNMULTISCATTER at: "
100 << aTrack.GetProperTime()/s << "s, "
101 << aTrack.GetGlobalTime()/s << "s. "
102 << ", after track length " << aTrack.GetTrackLength()/cm << "cm, "
103 << "in volume "
105 << G4endl;
106
107 G4ThreeVector scattered = Scatter();
108
110
111 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
112}
void Initialize(const G4Track &) override
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4VPhysicalVolume * GetPhysicalVolume() const
G4StepPoint * GetPostStepPoint() const
G4double GetTrackLength() const
G4double GetGlobalTime() const
G4double GetProperTime() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
const G4String & GetName() const
G4ParticleChange aParticleChange

◆ Scatter()

G4ThreeVector G4UCNMultiScattering::Scatter ( )

Definition at line 150 of file G4UCNMultiScattering.cc.

151{
152
153 G4ThreeVector final(0.,0.,1.);
154
155 // Make a simple uniform distribution in 4 pi
156 // apply scattering, calculate angle phi, theta
157
158 G4double theta = std::acos(2*G4UniformRand()-1);
159 G4double phi = G4UniformRand() * 2 * pi;
160
161 final.rotateY(theta);
162 final.rotateZ(phi);
163 final = final.unit();
164
165 return final;
166}
#define G4UniformRand()
Definition Randomize.hh:52

Referenced by PostStepDoIt().


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