Geant4  10.02.p03
G4UCNBoundaryProcess Class Reference

#include <G4UCNBoundaryProcess.hh>

Inheritance diagram for G4UCNBoundaryProcess:
Collaboration diagram for G4UCNBoundaryProcess:

Public Member Functions

 G4UCNBoundaryProcess (const G4String &processName="UCNBoundaryProcess", G4ProcessType type=fUCN)
 
virtual ~G4UCNBoundaryProcess ()
 
G4bool IsApplicable (const G4ParticleDefinition &aParticleType)
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *condition)
 
G4VParticleChange * PostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
G4ThreeVector MRreflect (G4double, G4ThreeVector, G4ThreeVector, G4double, G4double)
 
G4ThreeVector MRreflectHigh (G4double, G4double, G4double, G4ThreeVector, G4ThreeVector, G4double, G4double, G4double &)
 
void SetMicroRoughness (G4bool)
 
G4bool GetMicroRoughness ()
 
G4UCNBoundaryProcessStatus GetStatus () const
 
void BoundaryProcessSummary () const
 
void SetMaterialPropertiesTable1 (G4UCNMaterialPropertiesTable *MPT)
 
void SetMaterialPropertiesTable2 (G4UCNMaterialPropertiesTable *MPT)
 
G4double GetTheta_o ()
 
G4double GetPhi_o ()
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
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 G4VParticleChange * AtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChange * AlongStepDoIt (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 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 &)
 

Private Member Functions

 G4UCNBoundaryProcess (const G4UCNBoundaryProcess &right)
 
G4UCNBoundaryProcessoperator= (const G4UCNBoundaryProcess &right)
 
G4bool High (G4double, G4double)
 
G4bool Loss (G4double, G4double, G4double)
 
G4bool SpinFlip (G4double)
 
G4double Reflectivity (G4double, G4double)
 
G4ThreeVector Reflect (G4double, G4ThreeVector, G4ThreeVector)
 
G4double Transmit (G4double, G4double)
 
G4ThreeVector LDiffRefl (G4ThreeVector)
 
G4ThreeVector MRDiffRefl (G4ThreeVector, G4double, G4double, G4ThreeVector, G4double)
 
G4ThreeVector MRDiffTrans (G4ThreeVector, G4double, G4double, G4ThreeVector, G4double)
 
G4RotationMatrix GetCoordinateTransformMatrix (G4ThreeVector, G4ThreeVector)
 
void BoundaryProcessVerbose () const
 
G4bool InvokeSD (const G4Step *step)
 

Private Attributes

G4UCNBoundaryProcessMessengerfMessenger
 
G4double neV
 
G4double kCarTolerance
 
G4UCNBoundaryProcessStatus theStatus
 
G4MaterialMaterial1
 
G4MaterialMaterial2
 
G4UCNMaterialPropertiesTableaMaterialPropertiesTable1
 
G4UCNMaterialPropertiesTableaMaterialPropertiesTable2
 
G4bool UseMicroRoughnessReflection
 
G4bool DoMicroRoughnessReflection
 
G4int nNoMPT
 
G4int nNoMRT
 
G4int nNoMRCondition
 
G4int nAbsorption
 
G4int nEzero
 
G4int nFlip
 
G4int aSpecularReflection
 
G4int bSpecularReflection
 
G4int bLambertianReflection
 
G4int aMRDiffuseReflection
 
G4int bMRDiffuseReflection
 
G4int nSnellTransmit
 
G4int mSnellTransmit
 
G4int aMRDiffuseTransmit
 
G4double ftheta_o
 
G4double fphi_o
 

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
 
G4VParticleChange * pParticleChange
 
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 83 of file G4UCNBoundaryProcess.hh.

Constructor & Destructor Documentation

◆ G4UCNBoundaryProcess() [1/2]

G4UCNBoundaryProcess::G4UCNBoundaryProcess ( const G4String processName = "UCNBoundaryProcess",
G4ProcessType  type = fUCN 
)

Definition at line 68 of file G4UCNBoundaryProcess.cc.

70  : G4VDiscreteProcess(processName, type)
71 {
72 
73  if (verboseLevel > 0) G4cout << GetProcessName() << " is created " << G4endl;
74 
76 
78 
80 
81  neV = 1.0e-9*eV;
82 
85 
86  Material1 = NULL;
87  Material2 = NULL;
88 
91 
94 
96  nAbsorption = nEzero = nFlip = 0;
101  aMRDiffuseTransmit = 0;
102  ftheta_o = fphi_o = 0;
103 }
G4int verboseLevel
Definition: G4VProcess.hh:368
G4UCNBoundaryProcessMessenger * fMessenger
G4double GetSurfaceTolerance() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4GLOB_DLL std::ostream G4cout
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
G4UCNMaterialPropertiesTable * aMaterialPropertiesTable1
static const double eV
Definition: G4SIunits.hh:212
#define G4endl
Definition: G4ios.hh:61
G4UCNMaterialPropertiesTable * aMaterialPropertiesTable2
G4UCNBoundaryProcessStatus theStatus
static G4GeometryTolerance * GetInstance()
Here is the call graph for this function:

◆ ~G4UCNBoundaryProcess()

G4UCNBoundaryProcess::~G4UCNBoundaryProcess ( )
virtual

Definition at line 105 of file G4UCNBoundaryProcess.cc.

106 {
107  delete fMessenger;
108 }
G4UCNBoundaryProcessMessenger * fMessenger

◆ G4UCNBoundaryProcess() [2/2]

G4UCNBoundaryProcess::G4UCNBoundaryProcess ( const G4UCNBoundaryProcess right)
private

Member Function Documentation

◆ BoundaryProcessSummary()

void G4UCNBoundaryProcess::BoundaryProcessSummary ( void  ) const

Definition at line 944 of file G4UCNBoundaryProcess.cc.

945 {
946  G4cout << "Sum NoMT: "
947  << nNoMPT << G4endl;
948  G4cout << "Sum NoMRT: "
949  << nNoMRT << G4endl;
950  G4cout << "Sum NoMRCondition: "
951  << nNoMRCondition << G4endl;
952  G4cout << "Sum No. E < V Loss: "
953  << nAbsorption << G4endl;
954  G4cout << "Sum No. E > V Ezero: "
955  << nEzero << G4endl;
956  G4cout << "Sum No. E < V SpinFlip: "
957  << nFlip << G4endl;
958  G4cout << "Sum No. E > V Specular Reflection: "
960  G4cout << "Sum No. E < V Specular Reflection: "
962  G4cout << "Sum No. E < V Lambertian Reflection: "
964  G4cout << "Sum No. E > V MR DiffuseReflection: "
966  G4cout << "Sum No. E < V MR DiffuseReflection: "
968  G4cout << "Sum No. E > V SnellTransmit: "
969  << nSnellTransmit << G4endl;
970  G4cout << "Sum No. E > V MR SnellTransmit: "
971  << mSnellTransmit << G4endl;
972  G4cout << "Sum No. E > V DiffuseTransmit: "
974  G4cout << " " << G4endl;
975 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
Here is the caller graph for this function:

◆ BoundaryProcessVerbose()

void G4UCNBoundaryProcess::BoundaryProcessVerbose ( ) const
private

Definition at line 910 of file G4UCNBoundaryProcess.cc.

911 {
912  if ( theStatus == Undefined )
913  G4cout << " *** Undefined *** " << G4endl;
914  if ( theStatus == NotAtBoundary )
915  G4cout << " *** NotAtBoundary *** " << G4endl;
916  if ( theStatus == SameMaterial )
917  G4cout << " *** SameMaterial *** " << G4endl;
918  if ( theStatus == StepTooSmall )
919  G4cout << " *** StepTooSmall *** " << G4endl;
920  if ( theStatus == NoMPT )
921  G4cout << " *** No G4UCNMaterialPropertiesTable *** " << G4endl;
922  if ( theStatus == NoMRT )
923  G4cout << " *** No MicroRoughness Table *** " << G4endl;
924  if ( theStatus == NoMRCondition )
925  G4cout << " *** MicroRoughness Condition not satisfied *** " << G4endl;
926  if ( theStatus == Absorption )
927  G4cout << " *** Loss on Surface *** " << G4endl;
928  if ( theStatus == Ezero )
929  G4cout << " *** Ezero on Surface *** " << G4endl;
930  if ( theStatus == Flip )
931  G4cout << " *** Spin Flip on Surface *** " << G4endl;
932  if ( theStatus == SpecularReflection )
933  G4cout << " *** Specular Reflection *** " << G4endl;
935  G4cout << " *** LambertianR Reflection *** " << G4endl;
937  G4cout << " *** MR Model Diffuse Reflection *** " << G4endl;
938  if ( theStatus == SnellTransmit )
939  G4cout << " *** Snell Transmission *** " << G4endl;
940  if ( theStatus == MRDiffuseTransmit )
941  G4cout << " *** MR Model Diffuse Transmission *** " << G4endl;
942 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
G4UCNBoundaryProcessStatus theStatus
Here is the caller graph for this function:

◆ GetCoordinateTransformMatrix()

G4RotationMatrix G4UCNBoundaryProcess::GetCoordinateTransformMatrix ( G4ThreeVector  Normal,
G4ThreeVector  direction 
)
private

Definition at line 879 of file G4UCNBoundaryProcess.cc.

881 {
882  // Definition of the local coordinate system (c.s. of the reflection)
883 
884  G4ThreeVector es1, es2, es3;
885 
886  // z-Coordinate is the surface normal, has already length 1
887 
888  es3 = Normal;
889 
890  // Perpendicular part of incident direction w.r.t. normal
891  es1 = direction.perpPart(Normal);
892 
893  // Set to unit length: x-Coordinate
894 
895  es1.setMag(1.);
896  es2 = es1;
897 
898  // y-Coordinate is the pi/2-rotation of x-axis around z-axis
899 
900  es2.rotate(90.*degree, es3);
901 
902  // Transformation matrix consists just of the three coordinate vectors
903  // as the global coordinate system is the 'standard' coordinate system
904 
905  G4RotationMatrix matrix(es1, es2, es3);
906 
907  return matrix;
908 }
Hep3Vector perpPart() const
static const double degree
Definition: G4SIunits.hh:143
void setMag(double)
Definition: ThreeVector.cc:25
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:28
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 457 of file G4UCNBoundaryProcess.cc.

460 {
461  *condition = Forced;
462 
463  return DBL_MAX;
464 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83

◆ GetMicroRoughness()

G4bool G4UCNBoundaryProcess::GetMicroRoughness ( )
inline

Definition at line 248 of file G4UCNBoundaryProcess.hh.

249 {
251 }

◆ GetPhi_o()

G4double G4UCNBoundaryProcess::GetPhi_o ( )
inline

Definition at line 213 of file G4UCNBoundaryProcess.hh.

213 {return fphi_o;};

◆ GetStatus()

G4UCNBoundaryProcessStatus G4UCNBoundaryProcess::GetStatus ( ) const
inline

Definition at line 228 of file G4UCNBoundaryProcess.hh.

229 {
230  return theStatus;
231 }
G4UCNBoundaryProcessStatus theStatus

◆ GetTheta_o()

G4double G4UCNBoundaryProcess::GetTheta_o ( )
inline

Definition at line 212 of file G4UCNBoundaryProcess.hh.

212 {return ftheta_o;};

◆ High()

G4bool G4UCNBoundaryProcess::High ( G4double  Energy,
G4double  FermiPotDiff 
)
inlineprivate

Definition at line 234 of file G4UCNBoundaryProcess.hh.

235 {
236  // Returns true for Energy > Fermi Potential Difference
237 
238  return (Energy > FermiPotDiff);
239 }
Here is the caller graph for this function:

◆ InvokeSD()

G4bool G4UCNBoundaryProcess::InvokeSD ( const G4Step *  step)
private

Definition at line 977 of file G4UCNBoundaryProcess.cc.

978 {
979  G4Step aStep = *pStep;
980 
981  aStep.AddTotalEnergyDeposit(pStep->GetTrack()->GetKineticEnergy());
982 
983  G4VSensitiveDetector* sd = aStep.GetPostStepPoint()->GetSensitiveDetector();
984  if (sd) return sd->Hit(&aStep);
985  else return false;
986 }
G4bool Hit(G4Step *aStep)
Here is the call graph for this function:

◆ IsApplicable()

G4bool G4UCNBoundaryProcess::IsApplicable ( const G4ParticleDefinition aParticleType)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 222 of file G4UCNBoundaryProcess.hh.

223 {
224  return ( &aParticleType == G4Neutron::NeutronDefinition() );
225 }
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:99
Here is the call graph for this function:

◆ LDiffRefl()

G4ThreeVector G4UCNBoundaryProcess::LDiffRefl ( G4ThreeVector  Normal)
private

Definition at line 861 of file G4UCNBoundaryProcess.cc.

862 {
863  G4ThreeVector momentum;
864 
865  // cosine distribution - Lambert's law
866 
867  momentum.setRThetaPhi(1., std::acos(std::sqrt(G4UniformRand())), 2.*pi*G4UniformRand());
868  momentum.rotateUz(Normal);
869 
870  if (momentum*Normal < 0) {
871  momentum*=-1;
872  G4cout << "G4UCNBoundaryProcess::LDiffRefl: !" << G4endl;
873  }
874 
875  return momentum.unit();
876 }
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
Hep3Vector unit() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
void setRThetaPhi(double r, double theta, double phi)
static const double pi
Definition: G4SIunits.hh:74
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Loss()

G4bool G4UCNBoundaryProcess::Loss ( G4double  pUpScatter,
G4double  theVelocityNormal,
G4double  FermiPot 
)
private

Definition at line 466 of file G4UCNBoundaryProcess.cc.

469 {
470  // The surface roughness is not taken into account here.
471  // One could use e.g. ultracold neutrons, R. Golub, p.35,
472  // where mu is increased by roughness parameters sigma and
473  // omega, which are related to the height and the distance of
474  // "disturbances" on the surface
475 
476  G4double vBound = std::sqrt(2.*FermiPot/neutron_mass_c2*c_squared);
477  G4double vRatio = theVelocityNormal/vBound;
478 
479  G4double pLoss = (2*pUpScatter*vRatio)/(std::sqrt(1-(vRatio*vRatio)));
480 
481  // Check, if enhancement for surface roughness should be done
482 
488 
489  // cf. Golub's book p. 35, eq. 2.103
490 
491  pLoss *= std::sqrt(1+2*b*b*vBound*vBound/
492  (hdm*hdm+0.85*hdm*vBound*w+2*vBound*vBound*w*w));
493  }
494  }
495 
496  return (G4UniformRand() <= std::fabs(pLoss));
497 }
float c_squared
Definition: hepunit.py:258
#define G4UniformRand()
Definition: Randomize.hh:97
float neutron_mass_c2
Definition: hepunit.py:276
float hbar_Planck
Definition: hepunit.py:264
G4UCNMaterialPropertiesTable * aMaterialPropertiesTable2
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ MRDiffRefl()

G4ThreeVector G4UCNBoundaryProcess::MRDiffRefl ( G4ThreeVector  Normal,
G4double  Energy,
G4double  FermiPot,
G4ThreeVector  OldMomentum,
G4double  pDiffuse 
)
private

Definition at line 696 of file G4UCNBoundaryProcess.cc.

701 {
702  G4bool accepted = false;
703 
704  G4double theta_o, phi_o;
705 
706  // Polar angle of incidence
707 
708  G4double theta_i = OldMomentum.polarAngle(-Normal);
709 
710  // Azimuthal angle of incidence
711 
712  // G4double phi_i = -OldMomentum.azimAngle(-Normal);
713 
714  // accept-reject method for MR-reflection
715 
716  G4int count = 0;
717  while (!accepted) {
718  theta_o = G4UniformRand()*pi/2;
719  phi_o = G4UniformRand()*pi*2-pi;
720  // Box over distribution is increased by 50% to ensure no value is above
721  if (1.5*G4UniformRand()*
723  GetMRMaxProbability(theta_i, Energy)/pDiffuse <=
725  GetMRProbability(theta_i,Energy,FermiPot,theta_o,phi_o)/pDiffuse)
726 
727  accepted = true;
728 
729  // For the case that the box is nevertheless exceeded
730 
732  GetMRProbability(theta_i, Energy, FermiPot, theta_o, phi_o)/
734  GetMRMaxProbability(theta_i, Energy)) > 1) {
735  G4cout << "MRMax Wahrscheinlichkeitsueberschreitung!" << G4endl;
737  GetMRProbability(theta_i, Energy, FermiPot, theta_o, phi_o)/
739  GetMRMaxProbability(theta_i, Energy)) << G4endl;
741  SetMRMaxProbability(theta_i, Energy,
743  GetMRProbability(theta_i, Energy,
744  FermiPot, theta_o, phi_o));
745  }
746  // Loop checking, 31-Aug-2015, Vladimir Ivanchenko
747  if(++count > 10000) { accepted = true; }
748  }
749 
750  // Creates vector in the local coordinate system of the reflection
751 
752  G4ThreeVector localmomentum;
753  localmomentum.setRThetaPhi(1., theta_o, phi_o);
754 
755  ftheta_o = theta_o;
756  fphi_o = phi_o;
757 
758  // Get coordinate transform matrix
759 
760  G4RotationMatrix TransCoord =
761  GetCoordinateTransformMatrix(Normal, OldMomentum);
762 
763  //transfer to the coordinates of the global system
764 
765  G4ThreeVector momentum = TransCoord*localmomentum;
766 
767  //momentum.rotateUz(Normal);
768 
769  if (momentum * Normal<0) {
770  momentum*=-1;
771  // something has gone wrong...
772  G4cout << "G4UCNBoundaryProcess::MRDiffRefl: !" << G4endl;
773  }
774 
775  return momentum.unit();
776 }
G4RotationMatrix GetCoordinateTransformMatrix(G4ThreeVector, G4ThreeVector)
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
double polarAngle(const Hep3Vector &v2) const
Definition: SpaceVectorD.cc:28
bool G4bool
Definition: G4Types.hh:79
Hep3Vector unit() const
void setRThetaPhi(double r, double theta, double phi)
static const double pi
Definition: G4SIunits.hh:74
#define G4endl
Definition: G4ios.hh:61
G4UCNMaterialPropertiesTable * aMaterialPropertiesTable2
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ MRDiffTrans()

G4ThreeVector G4UCNBoundaryProcess::MRDiffTrans ( G4ThreeVector  Normal,
G4double  Energy,
G4double  FermiPot,
G4ThreeVector  OldMomentum,
G4double  pDiffuseTrans 
)
private

Definition at line 778 of file G4UCNBoundaryProcess.cc.

783 {
784  G4bool accepted = false;
785 
786  G4double theta_o, phi_o;
787 
788  // Polar angle of incidence
789 
790  G4double theta_i = OldMomentum.polarAngle(-Normal);
791 
792  // azimuthal angle of incidence
793 
794  // G4double phi_i = -OldMomentum.azimAngle(-Normal);
795 
796  G4int count = 0;
797  while (!accepted) {
798  theta_o = G4UniformRand()*pi/2;
799  phi_o = G4UniformRand()*pi*2-pi;
800 
801  // Box over distribution is increased by 50% to ensure no value is above
802 
803  if (1.5*G4UniformRand()*
805  GetMRMaxTransProbability(theta_i, Energy)/pDiffuseTrans <=
807  GetMRTransProbability(theta_i,Energy,FermiPot,theta_o,phi_o)/
808  pDiffuseTrans)
809 
810  accepted=true;
811 
812  // For the case that the box is nevertheless exceeded
813 
815  GetMRTransProbability(theta_i, Energy, FermiPot, theta_o, phi_o)/
817  GetMRMaxTransProbability(theta_i, Energy)) > 1) {
818  G4cout << "MRMaxTrans Wahrscheinlichkeitsueberschreitung!" << G4endl;
820  GetMRTransProbability(theta_i, Energy, FermiPot, theta_o, phi_o)/
822  GetMRMaxTransProbability(theta_i, Energy)) << G4endl;
824  SetMRMaxTransProbability(theta_i, Energy,
826  GetMRTransProbability(theta_i, Energy,
827  FermiPot, theta_o, phi_o));
828  }
829  // Loop checking, 31-Aug-2015, Vladimir Ivanchenko
830  if(++count > 10000) { accepted = true; }
831  }
832 
833  // Creates vector in the local coordinate system of the reflection
834 
835  G4ThreeVector localmomentum;
836  localmomentum.setRThetaPhi(1., pi-theta_o, phi_o);
837 
838  // Get coordinate transform matrix
839 
840  G4RotationMatrix TransCoord =
841  GetCoordinateTransformMatrix(Normal, OldMomentum);
842 
843  // Transfer to the coordinates of the global system
844 
845  G4ThreeVector momentum = TransCoord*localmomentum;
846 
847  if (momentum*Normal<0) {
848  // something has gone wrong...
849  momentum*=-1;
850  G4cout << "G4UCNBoundaryProcess::MRDiffTrans: !" << G4endl;
851  }
852 
853  return momentum.unit();
854 }
G4RotationMatrix GetCoordinateTransformMatrix(G4ThreeVector, G4ThreeVector)
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
double polarAngle(const Hep3Vector &v2) const
Definition: SpaceVectorD.cc:28
bool G4bool
Definition: G4Types.hh:79
Hep3Vector unit() const
void setRThetaPhi(double r, double theta, double phi)
static const double pi
Definition: G4SIunits.hh:74
#define G4endl
Definition: G4ios.hh:61
G4UCNMaterialPropertiesTable * aMaterialPropertiesTable2
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ MRreflect()

G4ThreeVector G4UCNBoundaryProcess::MRreflect ( G4double  pDiffuse,
G4ThreeVector  OldMomentum,
G4ThreeVector  Normal,
G4double  Energy,
G4double  FermiPot 
)

Definition at line 543 of file G4UCNBoundaryProcess.cc.

548 {
549  // Only for Enormal <= VFermi
550 
551  G4ThreeVector NewMomentum;
552 
553  // Checks if the reflection should be non-specular
554 
555  if (G4UniformRand() <= pDiffuse) {
556 
557  // Reflect diffuse
558 
559  // Determines the angles of the non-specular reflection
560 
561  NewMomentum =
562  MRDiffRefl(Normal, Energy, FermiPot, OldMomentum, pDiffuse);
563 
567 
568  return NewMomentum;
569 
570  } else {
571 
572  // Reflect specular
573 
574  G4double PdotN = OldMomentum * Normal;
575 
576  NewMomentum = OldMomentum - (2.*PdotN)*Normal;
577  NewMomentum.unit();
578 
582 
583  return NewMomentum;
584  }
585 }
G4int verboseLevel
Definition: G4VProcess.hh:368
G4ThreeVector MRDiffRefl(G4ThreeVector, G4double, G4double, G4ThreeVector, G4double)
#define G4UniformRand()
Definition: Randomize.hh:97
Hep3Vector unit() const
double G4double
Definition: G4Types.hh:76
G4UCNBoundaryProcessStatus theStatus
Here is the call graph for this function:
Here is the caller graph for this function:

◆ MRreflectHigh()

G4ThreeVector G4UCNBoundaryProcess::MRreflectHigh ( G4double  pDiffuse,
G4double  pDiffuseTrans,
G4double  pLoss,
G4ThreeVector  OldMomentum,
G4ThreeVector  Normal,
G4double  Energy,
G4double  FermiPot,
G4double Enew 
)

Definition at line 587 of file G4UCNBoundaryProcess.cc.

595 {
596  // Only for Enormal > VFermi
597 
598  G4double costheta = OldMomentum*Normal;
599 
600  G4double Enormal = Energy * (costheta*costheta);
601 
602 // G4double pSpecular = Reflectivity(Enormal,FermiPot)*
603  G4double pSpecular = Reflectivity(FermiPot,Enormal)*
604  (1.-pDiffuse-pDiffuseTrans-pLoss);
605 
606  G4ThreeVector NewMomentum;
607 
608  G4double decide = G4UniformRand();
609 
610  if (decide < pSpecular) {
611 
612  // Reflect specularly
613 
614  G4double PdotN = OldMomentum * Normal;
615  NewMomentum = OldMomentum - (2.*PdotN)*Normal;
616  NewMomentum.unit();
617 
618  Enew = Energy;
619 
623 
624  return NewMomentum;
625  }
626 
627  if (decide < pSpecular + pDiffuse) {
628 
629  // Reflect diffusely
630 
631  // Determines the angles of the non-specular reflection
632 
633  NewMomentum =
634  MRDiffRefl(Normal, Energy, FermiPot, OldMomentum, pDiffuse);
635 
636  if (verboseLevel > 0) G4cout << "Diffuse normal " << Normal
637  << ", " << NewMomentum << G4endl;
638  Enew = Energy;
639 
643 
644  return NewMomentum;
645  }
646 
647  if (decide < pSpecular + pDiffuse + pDiffuseTrans) {
648 
649  // Transmit diffusely
650 
651  // Determines the angles of the non-specular transmission
652 
653  NewMomentum =
654  MRDiffTrans(Normal, Energy, FermiPot, OldMomentum, pDiffuseTrans);
655 
656  Enew = Energy - FermiPot;
657 
661 
662  return NewMomentum;
663  }
664 
665  if (decide < pSpecular + pDiffuse + pDiffuseTrans + pLoss) {
666 
667  // Loss
668 
669  Enew = 0.;
670 
671  nEzero++;
672  theStatus = Ezero;
674 
675  return NewMomentum;
676  }
677 
678  // last case: Refractive transmission
679 
680  Enew = Energy - FermiPot;
681 
682  G4double palt = std::sqrt(2*neutron_mass_c2/c_squared*Energy);
683  G4double produ = OldMomentum * Normal;
684 
685  NewMomentum = palt*OldMomentum-
686  (palt*produ+std::sqrt(palt*palt*produ*produ-2*neutron_mass_c2/
687  c_squared*FermiPot))*Normal;
688 
689  mSnellTransmit++;
692 
693  return NewMomentum.unit();
694 }
G4ThreeVector MRDiffTrans(G4ThreeVector, G4double, G4double, G4ThreeVector, G4double)
G4int verboseLevel
Definition: G4VProcess.hh:368
float c_squared
Definition: hepunit.py:258
G4ThreeVector MRDiffRefl(G4ThreeVector, G4double, G4double, G4ThreeVector, G4double)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
Hep3Vector unit() const
float neutron_mass_c2
Definition: hepunit.py:276
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4UCNBoundaryProcessStatus theStatus
G4double Reflectivity(G4double, G4double)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator=()

G4UCNBoundaryProcess& G4UCNBoundaryProcess::operator= ( const G4UCNBoundaryProcess right)
private

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 111 of file G4UCNBoundaryProcess.cc.

112 {
113  aParticleChange.Initialize(aTrack);
114  aParticleChange.ProposeVelocity(aTrack.GetVelocity());
115 
116  // Get hyperStep from G4ParallelWorldProcess
117  // NOTE: PostSetpDoIt of this process should be
118  // invoked after G4ParallelWorldProcess!
119 
120  const G4Step* pStep = &aStep;
121 
122  const G4Step* hStep = G4ParallelWorldProcess::GetHyperStep();
123 
124  if (hStep) pStep = hStep;
125 
126  G4bool isOnBoundary =
127  (pStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary);
128 
129  if (isOnBoundary) {
130  Material1 = pStep->GetPreStepPoint()->GetMaterial();
131  Material2 = pStep->GetPostStepPoint()->GetMaterial();
132  } else {
135  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
136  }
137 
138  if (aTrack.GetStepLength()<=kCarTolerance/2) {
141  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
142  }
143 
145  GetMaterialPropertiesTable();
147  GetMaterialPropertiesTable();
148 
149  G4String volnam1 = pStep->GetPreStepPoint() ->GetPhysicalVolume()->GetName();
150  G4String volnam2 = pStep->GetPostStepPoint()->GetPhysicalVolume()->GetName();
151 
152  if (verboseLevel > 0) {
153  G4cout << " UCN at Boundary! " << G4endl;
154  G4cout << " vol1: " << volnam1 << ", vol2: " << volnam2 << G4endl;
155  G4cout << " Ekin: " << aTrack.GetKineticEnergy()/neV <<"neV"<< G4endl;
156  G4cout << " MomDir: " << aTrack.GetMomentumDirection() << G4endl;
157  }
158 
159  if (Material1 == Material2) {
162  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
163  }
164 
165  G4ThreeVector theGlobalPoint = pStep->GetPostStepPoint()->GetPosition();
166 
167  G4bool valid;
168  // Use the new method for Exit Normal in global coordinates,
169  // which provides the normal more reliably.
170 
171  // ID of Navigator which limits step
172 
174  std::vector<G4Navigator*>::iterator iNav =
176  GetActiveNavigatorsIterator();
177 
178  G4ThreeVector theGlobalNormal =
179  (iNav[hNavId])->GetGlobalExitNormal(theGlobalPoint,&valid);
180 
181  if (valid) {
182  theGlobalNormal = -theGlobalNormal;
183  }
184  else
185  {
187  ed << " G4UCNBoundaryProcess/PostStepDoIt(): "
188  << " The Navigator reports that it returned an invalid normal"
189  << G4endl;
190  G4Exception("G4UCNBoundaryProcess::PostStepDoIt", "UCNBoun01",
192  "Invalid Surface Normal - Geometry must return valid surface normal");
193  }
194 
195  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
196 
197  G4ThreeVector OldMomentum = aParticle->GetMomentumDirection();
198 
199  if (OldMomentum * theGlobalNormal > 0.0) {
200 #ifdef G4OPTICAL_DEBUG
202  ed << " G4UCNBoundaryProcess/PostStepDoIt(): "
203  << " theGlobalNormal points in a wrong direction. "
204  << G4endl;
205  ed << " The momentum of the photon arriving at interface (oldMomentum)"
206  << " must exit the volume cross in the step. " << G4endl;
207  ed << " So it MUST have dot < 0 with the normal that Exits the new volume (globalNormal)." << G4endl;
208  ed << " >> The dot product of oldMomentum and global Normal is " << OldMomentum*theGlobalNormal << G4endl;
209  ed << " Old Momentum (during step) = " << OldMomentum << G4endl;
210  ed << " Global Normal (Exiting New Vol) = " << theGlobalNormal << G4endl;
211  ed << G4endl;
212  G4Exception("G4UCNBoundaryProcess::PostStepDoIt", "UCNBoun02",
213  EventMustBeAborted, // Or JustWarning to see if it happens repeatedbly on one ray
214  ed,
215  "Invalid Surface Normal - Geometry must return valid surface normal pointing in the right direction");
216 #else
217  theGlobalNormal = -theGlobalNormal;
218 #endif
219  }
220 
221  G4ThreeVector theNeutronMomentum = aTrack.GetMomentum();
222 
223  G4double theMomentumNormal = theNeutronMomentum*theGlobalNormal;
224  G4double theVelocityNormal = aTrack.GetVelocity() *
225  (OldMomentum * theGlobalNormal);
226 
227  G4double Enormal = theMomentumNormal*theMomentumNormal/2./neutron_mass_c2;
228  G4double Energy = aTrack.GetKineticEnergy();
229 
230  G4double FermiPot2 = 0.;
231  G4double pDiffuse = 0.;
232  G4double pSpinFlip = 0.;
233  G4double pUpScatter = 0.;
234 
236  FermiPot2 = aMaterialPropertiesTable2->GetConstProperty("FERMIPOT")*neV;
237  pDiffuse = aMaterialPropertiesTable2->GetConstProperty("DIFFUSION");
238  pSpinFlip = aMaterialPropertiesTable2->GetConstProperty("SPINFLIP");
239  pUpScatter = aMaterialPropertiesTable2->GetConstProperty("LOSS");
240  }
241 
242  G4double FermiPot1 = 0.;
244  FermiPot1 = aMaterialPropertiesTable1->GetConstProperty("FERMIPOT")*neV;
245 
246  G4double FermiPotDiff = FermiPot2 - FermiPot1;
247 
248  if ( verboseLevel > 1 )
249  G4cout << "UCNBoundaryProcess: new FermiPot: " << FermiPot2/neV
250  << "neV, old FermiPot:" << FermiPot1/neV << "neV" << G4endl;
251 
252  // Use microroughness diffuse reflection behavior if activated
253 
255 
256  G4double theta_i = 0.;
257 
259 
260  nNoMPT++;
261  theStatus = NoMPT;
264 
265  } else {
266 
268 
269  nNoMRT++;
270  theStatus = NoMRT;
272 
274  }
275 
276  // Angle theta_in between surface and momentum direction,
277  // Phi_in is defined to be 0
278 
279  theta_i = OldMomentum.angle(-theGlobalNormal);
280 
281  // Checks the MR-conditions
282 
284  ConditionsValid(Energy, FermiPotDiff, theta_i)) {
285 
286  nNoMRCondition++;
289 
291  }
292  }
293 
294  G4double MRpDiffuse = 0.;
295  G4double MRpDiffuseTrans = 0.;
296 
297  // If microroughness is available and active for material in question
298 
300 
301  // Integral probability for non-specular reflection with microroughness
302 
303  MRpDiffuse = aMaterialPropertiesTable2->
304  GetMRIntProbability(theta_i, Energy);
305 
306  // Integral probability for non-specular transmission with microroughness
307 
308  MRpDiffuseTrans = aMaterialPropertiesTable2->
309  GetMRIntTransProbability(theta_i, Energy);
310 
311  if ( verboseLevel > 1 ) {
312  G4cout << "theta: " << theta_i/degree << "degree" << G4endl;
313  G4cout << "Energy: " << Energy/neV << "neV"
314  << ", Enormal: " << Enormal/neV << "neV" << G4endl;
315  G4cout << "FermiPotDiff: " << FermiPotDiff/neV << "neV " << G4endl;
316  G4cout << "Reflection_prob: " << MRpDiffuse
317  << ", Transmission_prob: " << MRpDiffuseTrans << G4endl;
318  }
319  }
320 
321  if (!High(Enormal, FermiPotDiff)) {
322 
323  // Below critical velocity
324 
325  if (verboseLevel > 0) G4cout << "G4UCNBoundaryProcess -> BELOW critical velocity" << G4endl;
326 
327  // Loss on reflection
328 
329  if (Loss(pUpScatter, theVelocityNormal, FermiPotDiff)) {
330 
331  // kill it.
332  aParticleChange.ProposeTrackStatus(fStopAndKill);
333  aParticleChange.ProposeLocalEnergyDeposit(Energy);
334 
335  nAbsorption++;
338 
339  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
340  }
341 
342  // spinflips
343 
344  if (SpinFlip(pSpinFlip)) {
345  nFlip++;
346  theStatus = Flip;
348 
349  G4ThreeVector NewPolarization = -1. * aParticle->GetPolarization();
350  aParticleChange.ProposePolarization(NewPolarization);
351  }
352 
353  // Reflect from surface
354 
355  G4ThreeVector NewMomentum;
356 
357  // If microroughness is available and active - do non-specular reflection
358 
360  NewMomentum = MRreflect(MRpDiffuse, OldMomentum, theGlobalNormal,
361  Energy, FermiPotDiff);
362  else
363 
364  // Else do it with the Lambert model as implemented by Peter Fierlinger
365 
366  NewMomentum = Reflect(pDiffuse, OldMomentum, theGlobalNormal);
367 
368  aParticleChange.ProposeMomentumDirection(NewMomentum);
369 
370  } else {
371 
372  // Above critical velocity
373 
374  if (verboseLevel > 0) G4cout << "G4UCNBoundaryProcess -> ABOVE critical velocity" << G4endl;
375 
376  // If it is faster than the criticial velocity,
377  // there is a probability to be still reflected.
378  // This formula is (only) valid for low loss materials
379 
380  // If microroughness available and active, do reflection with it
381 
382  G4ThreeVector NewMomentum;
383 
385 
386  G4double Enew;
387 
388  NewMomentum =
389  MRreflectHigh(MRpDiffuse, MRpDiffuseTrans, 0., OldMomentum,
390  theGlobalNormal, Energy, FermiPotDiff, Enew);
391 
392  if (Enew == 0.) {
393  aParticleChange.ProposeTrackStatus(fStopAndKill);
394  aParticleChange.ProposeLocalEnergyDeposit(Energy);
395  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
396  } else {
397  aParticleChange.ProposeEnergy(Enew);
398  aParticleChange.ProposeMomentumDirection(NewMomentum);
399  aParticleChange.ProposeVelocity(sqrt(2*Enew/neutron_mass_c2)*c_light);
400  aParticleChange.ProposeLocalEnergyDeposit(Energy-Enew);
401  }
402 
403  } else {
404 
405  G4double reflectivity = Reflectivity(FermiPotDiff, Enormal);
406 
407  if ( verboseLevel > 1 ) G4cout << "UCNBoundaryProcess: reflectivity "
408  << reflectivity << G4endl;
409 
410  if (G4UniformRand() < reflectivity) {
411 
412  // Reflect from surface
413 
414  NewMomentum = Reflect(pDiffuse, OldMomentum, theGlobalNormal);
415  aParticleChange.ProposeMomentumDirection(NewMomentum);
416 
417  } else {
418 
419  // --- Transmission because it is faster than the critical velocity
420 
421  G4double Enew = Transmit(FermiPotDiff, Energy);
422 
423  // --- Change of the normal momentum component
424  // p = sqrt(2*m*Ekin)
425 
426  G4double mass = -std::sqrt(theMomentumNormal*theMomentumNormal -
427  neutron_mass_c2*2.*FermiPotDiff);
428 
429  // --- Momentum direction in new media
430 
431  NewMomentum =
432  theNeutronMomentum - (theMomentumNormal-mass)*theGlobalNormal;
433 
434  nSnellTransmit++;
437 
438  aParticleChange.ProposeEnergy(Enew);
439  aParticleChange.ProposeMomentumDirection(NewMomentum.unit());
440  aParticleChange.ProposeVelocity(sqrt(2*Enew/neutron_mass_c2)*c_light);
441  aParticleChange.ProposeLocalEnergyDeposit(Energy-Enew);
442 
443  if (verboseLevel > 1) {
444  G4cout << "Energy: " << Energy/neV << "neV, Enormal: "
445  << Enormal/neV << "neV, fpdiff " << FermiPotDiff/neV
446  << "neV, Enew " << Enew/neV << "neV" << G4endl;
447  G4cout << "UCNBoundaryProcess: Transmit and set the new Energy "
448  << aParticleChange.GetEnergy()/neV << "neV" << G4endl;
449  }
450  }
451  }
452  }
453 
454  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
455 }
G4int verboseLevel
Definition: G4VProcess.hh:368
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
G4bool High(G4double, G4double)
G4double Transmit(G4double, G4double)
int G4int
Definition: G4Types.hh:78
double angle(const Hep3Vector &) const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
Hep3Vector unit() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4TransportationManager * GetTransportationManager()
float neutron_mass_c2
Definition: hepunit.py:276
G4bool Loss(G4double, G4double, G4double)
G4UCNMaterialPropertiesTable * aMaterialPropertiesTable1
G4ThreeVector MRreflect(G4double, G4ThreeVector, G4ThreeVector, G4double, G4double)
G4ThreeVector MRreflectHigh(G4double, G4double, G4double, G4ThreeVector, G4ThreeVector, G4double, G4double, G4double &)
const G4ThreeVector & GetMomentumDirection() const
const G4ThreeVector & GetPolarization() const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
static const double degree
Definition: G4SIunits.hh:143
#define G4endl
Definition: G4ios.hh:61
G4ThreeVector Reflect(G4double, G4ThreeVector, G4ThreeVector)
G4UCNMaterialPropertiesTable * aMaterialPropertiesTable2
double G4double
Definition: G4Types.hh:76
G4UCNBoundaryProcessStatus theStatus
G4double GetConstProperty(const char *key)
G4double Reflectivity(G4double, G4double)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
float c_light
Definition: hepunit.py:257
static const G4Step * GetHyperStep()
Here is the call graph for this function:

◆ Reflect()

G4ThreeVector G4UCNBoundaryProcess::Reflect ( G4double  pDiffuse,
G4ThreeVector  OldMomentum,
G4ThreeVector  Normal 
)
private

Definition at line 512 of file G4UCNBoundaryProcess.cc.

515 {
516  G4double PdotN = OldMomentum * Normal;
517 
518  G4ThreeVector NewMomentum = OldMomentum - (2.*PdotN)*Normal;
519  NewMomentum.unit();
520 
521  // Reflect diffuse
522 
523  if (NewMomentum == OldMomentum || G4UniformRand() < pDiffuse) {
524 
525  NewMomentum = LDiffRefl(Normal);
526 
530 
531  return NewMomentum;
532  }
533 
534  // Reflect specular
535 
539 
540  return NewMomentum;
541 }
G4int verboseLevel
Definition: G4VProcess.hh:368
#define G4UniformRand()
Definition: Randomize.hh:97
Hep3Vector unit() const
G4ThreeVector LDiffRefl(G4ThreeVector)
double G4double
Definition: G4Types.hh:76
G4UCNBoundaryProcessStatus theStatus
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Reflectivity()

G4double G4UCNBoundaryProcess::Reflectivity ( G4double  FermiPot,
G4double  Enormal 
)
private

Definition at line 504 of file G4UCNBoundaryProcess.cc.

505 {
506  G4double r = (std::sqrt(Enormal) - std::sqrt(Enormal - FermiPot)) /
507  (std::sqrt(Enormal) + std::sqrt(Enormal - FermiPot));
508 
509  return r*r;
510 }
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

◆ SetMaterialPropertiesTable1()

void G4UCNBoundaryProcess::SetMaterialPropertiesTable1 ( G4UCNMaterialPropertiesTable MPT)
inline

Definition at line 206 of file G4UCNBoundaryProcess.hh.

G4UCNMaterialPropertiesTable * aMaterialPropertiesTable1

◆ SetMaterialPropertiesTable2()

void G4UCNBoundaryProcess::SetMaterialPropertiesTable2 ( G4UCNMaterialPropertiesTable MPT)
inline

Definition at line 209 of file G4UCNBoundaryProcess.hh.

G4UCNMaterialPropertiesTable * aMaterialPropertiesTable2

◆ SetMicroRoughness()

void G4UCNBoundaryProcess::SetMicroRoughness ( G4bool  active)
inline

Definition at line 242 of file G4UCNBoundaryProcess.hh.

Here is the caller graph for this function:

◆ SpinFlip()

G4bool G4UCNBoundaryProcess::SpinFlip ( G4double  pSpinFlip)
private

Definition at line 499 of file G4UCNBoundaryProcess.cc.

500 {
501  return (G4UniformRand() <= pSpinFlip);
502 }
#define G4UniformRand()
Definition: Randomize.hh:97
Here is the caller graph for this function:

◆ Transmit()

G4double G4UCNBoundaryProcess::Transmit ( G4double  FermiPot,
G4double  Energy 
)
private

Definition at line 856 of file G4UCNBoundaryProcess.cc.

857 {
858  return (Energy - FermiPot);
859 }
Here is the caller graph for this function:

Member Data Documentation

◆ aMaterialPropertiesTable1

G4UCNMaterialPropertiesTable* G4UCNBoundaryProcess::aMaterialPropertiesTable1
private

Definition at line 136 of file G4UCNBoundaryProcess.hh.

◆ aMaterialPropertiesTable2

G4UCNMaterialPropertiesTable* G4UCNBoundaryProcess::aMaterialPropertiesTable2
private

Definition at line 138 of file G4UCNBoundaryProcess.hh.

◆ aMRDiffuseReflection

G4int G4UCNBoundaryProcess::aMRDiffuseReflection
private

Definition at line 191 of file G4UCNBoundaryProcess.hh.

◆ aMRDiffuseTransmit

G4int G4UCNBoundaryProcess::aMRDiffuseTransmit
private

Definition at line 193 of file G4UCNBoundaryProcess.hh.

◆ aSpecularReflection

G4int G4UCNBoundaryProcess::aSpecularReflection
private

Definition at line 189 of file G4UCNBoundaryProcess.hh.

◆ bLambertianReflection

G4int G4UCNBoundaryProcess::bLambertianReflection
private

Definition at line 190 of file G4UCNBoundaryProcess.hh.

◆ bMRDiffuseReflection

G4int G4UCNBoundaryProcess::bMRDiffuseReflection
private

Definition at line 191 of file G4UCNBoundaryProcess.hh.

◆ bSpecularReflection

G4int G4UCNBoundaryProcess::bSpecularReflection
private

Definition at line 189 of file G4UCNBoundaryProcess.hh.

◆ DoMicroRoughnessReflection

G4bool G4UCNBoundaryProcess::DoMicroRoughnessReflection
private

Definition at line 141 of file G4UCNBoundaryProcess.hh.

◆ fMessenger

G4UCNBoundaryProcessMessenger* G4UCNBoundaryProcess::fMessenger
private

Definition at line 124 of file G4UCNBoundaryProcess.hh.

◆ fphi_o

G4double G4UCNBoundaryProcess::fphi_o
private

Definition at line 195 of file G4UCNBoundaryProcess.hh.

◆ ftheta_o

G4double G4UCNBoundaryProcess::ftheta_o
private

Definition at line 195 of file G4UCNBoundaryProcess.hh.

◆ kCarTolerance

G4double G4UCNBoundaryProcess::kCarTolerance
private

Definition at line 128 of file G4UCNBoundaryProcess.hh.

◆ Material1

G4Material* G4UCNBoundaryProcess::Material1
private

Definition at line 132 of file G4UCNBoundaryProcess.hh.

◆ Material2

G4Material* G4UCNBoundaryProcess::Material2
private

Definition at line 133 of file G4UCNBoundaryProcess.hh.

◆ mSnellTransmit

G4int G4UCNBoundaryProcess::mSnellTransmit
private

Definition at line 192 of file G4UCNBoundaryProcess.hh.

◆ nAbsorption

G4int G4UCNBoundaryProcess::nAbsorption
private

Definition at line 188 of file G4UCNBoundaryProcess.hh.

◆ neV

G4double G4UCNBoundaryProcess::neV
private

Definition at line 126 of file G4UCNBoundaryProcess.hh.

◆ nEzero

G4int G4UCNBoundaryProcess::nEzero
private

Definition at line 188 of file G4UCNBoundaryProcess.hh.

◆ nFlip

G4int G4UCNBoundaryProcess::nFlip
private

Definition at line 188 of file G4UCNBoundaryProcess.hh.

◆ nNoMPT

G4int G4UCNBoundaryProcess::nNoMPT
private

Definition at line 187 of file G4UCNBoundaryProcess.hh.

◆ nNoMRCondition

G4int G4UCNBoundaryProcess::nNoMRCondition
private

Definition at line 187 of file G4UCNBoundaryProcess.hh.

◆ nNoMRT

G4int G4UCNBoundaryProcess::nNoMRT
private

Definition at line 187 of file G4UCNBoundaryProcess.hh.

◆ nSnellTransmit

G4int G4UCNBoundaryProcess::nSnellTransmit
private

Definition at line 192 of file G4UCNBoundaryProcess.hh.

◆ theStatus

G4UCNBoundaryProcessStatus G4UCNBoundaryProcess::theStatus
private

Definition at line 130 of file G4UCNBoundaryProcess.hh.

◆ UseMicroRoughnessReflection

G4bool G4UCNBoundaryProcess::UseMicroRoughnessReflection
private

Definition at line 140 of file G4UCNBoundaryProcess.hh.


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