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

#include <G4UCNBoundaryProcess.hh>

+ Inheritance 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)
 
G4VParticleChangePostStepDoIt (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 &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 82 of file G4UCNBoundaryProcess.hh.

Constructor & Destructor Documentation

◆ G4UCNBoundaryProcess()

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

Definition at line 67 of file G4UCNBoundaryProcess.cc.

69 : G4VDiscreteProcess(processName, type)
70{
71
72 if (verboseLevel > 0) G4cout << GetProcessName() << " is created " << G4endl;
73
75
76 theStatus = Undefined;
77
78 fMessenger = new G4UCNBoundaryProcessMessenger(this);
79
80 neV = 1.0e-9*eV;
81
82 kCarTolerance = G4GeometryTolerance::GetInstance()
84
85 Material1 = NULL;
86 Material2 = NULL;
87
88 aMaterialPropertiesTable1 = NULL;
89 aMaterialPropertiesTable2 = NULL;
90
91 UseMicroRoughnessReflection = false;
92 DoMicroRoughnessReflection = false;
93
94 nNoMPT = nNoMRT = nNoMRCondition = 0;
95 nAbsorption = nEzero = nFlip = 0;
96 aSpecularReflection = bSpecularReflection = 0;
97 bLambertianReflection = 0;
98 aMRDiffuseReflection = bMRDiffuseReflection = 0;
99 nSnellTransmit = mSnellTransmit = 0;
100 aMRDiffuseTransmit = 0;
101 ftheta_o = fphi_o = 0;
102}
@ fUCNBoundary
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
G4int verboseLevel
Definition: G4VProcess.hh:360
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386

◆ ~G4UCNBoundaryProcess()

G4UCNBoundaryProcess::~G4UCNBoundaryProcess ( )
virtual

Definition at line 104 of file G4UCNBoundaryProcess.cc.

105{
106 delete fMessenger;
107}

Member Function Documentation

◆ BoundaryProcessSummary()

void G4UCNBoundaryProcess::BoundaryProcessSummary ( void  ) const

Definition at line 943 of file G4UCNBoundaryProcess.cc.

944{
945 G4cout << "Sum NoMT: "
946 << nNoMPT << G4endl;
947 G4cout << "Sum NoMRT: "
948 << nNoMRT << G4endl;
949 G4cout << "Sum NoMRCondition: "
950 << nNoMRCondition << G4endl;
951 G4cout << "Sum No. E < V Loss: "
952 << nAbsorption << G4endl;
953 G4cout << "Sum No. E > V Ezero: "
954 << nEzero << G4endl;
955 G4cout << "Sum No. E < V SpinFlip: "
956 << nFlip << G4endl;
957 G4cout << "Sum No. E > V Specular Reflection: "
958 << aSpecularReflection << G4endl;
959 G4cout << "Sum No. E < V Specular Reflection: "
960 << bSpecularReflection << G4endl;
961 G4cout << "Sum No. E < V Lambertian Reflection: "
962 << bLambertianReflection << G4endl;
963 G4cout << "Sum No. E > V MR DiffuseReflection: "
964 << aMRDiffuseReflection << G4endl;
965 G4cout << "Sum No. E < V MR DiffuseReflection: "
966 << bMRDiffuseReflection << G4endl;
967 G4cout << "Sum No. E > V SnellTransmit: "
968 << nSnellTransmit << G4endl;
969 G4cout << "Sum No. E > V MR SnellTransmit: "
970 << mSnellTransmit << G4endl;
971 G4cout << "Sum No. E > V DiffuseTransmit: "
972 << aMRDiffuseTransmit << G4endl;
973 G4cout << " " << G4endl;
974}

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 456 of file G4UCNBoundaryProcess.cc.

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

◆ GetMicroRoughness()

G4bool G4UCNBoundaryProcess::GetMicroRoughness ( )
inline

Definition at line 247 of file G4UCNBoundaryProcess.hh.

248{
249 return UseMicroRoughnessReflection;
250}

◆ GetPhi_o()

G4double G4UCNBoundaryProcess::GetPhi_o ( )
inline

Definition at line 212 of file G4UCNBoundaryProcess.hh.

212{return fphi_o;};

◆ GetStatus()

G4UCNBoundaryProcessStatus G4UCNBoundaryProcess::GetStatus ( ) const
inline

Definition at line 227 of file G4UCNBoundaryProcess.hh.

228{
229 return theStatus;
230}

◆ GetTheta_o()

G4double G4UCNBoundaryProcess::GetTheta_o ( )
inline

Definition at line 211 of file G4UCNBoundaryProcess.hh.

211{return ftheta_o;};

◆ IsApplicable()

G4bool G4UCNBoundaryProcess::IsApplicable ( const G4ParticleDefinition aParticleType)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 221 of file G4UCNBoundaryProcess.hh.

222{
223 return ( &aParticleType == G4Neutron::NeutronDefinition() );
224}
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:98

◆ MRreflect()

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

Definition at line 542 of file G4UCNBoundaryProcess.cc.

547{
548 // Only for Enormal <= VFermi
549
550 G4ThreeVector NewMomentum;
551
552 // Checks if the reflection should be non-specular
553
554 if (G4UniformRand() <= pDiffuse) {
555
556 // Reflect diffuse
557
558 // Determines the angles of the non-specular reflection
559
560 NewMomentum =
561 MRDiffRefl(Normal, Energy, FermiPot, OldMomentum, pDiffuse);
562
563 bMRDiffuseReflection++;
564 theStatus = MRDiffuseReflection;
565 if ( verboseLevel > 0 ) BoundaryProcessVerbose();
566
567 return NewMomentum;
568
569 } else {
570
571 // Reflect specular
572
573 G4double PdotN = OldMomentum * Normal;
574
575 NewMomentum = OldMomentum - (2.*PdotN)*Normal;
576 NewMomentum.unit();
577
578 bSpecularReflection++;
579 theStatus = SpecularReflection;
580 if ( verboseLevel > 0 ) BoundaryProcessVerbose();
581
582 return NewMomentum;
583 }
584}
double G4double
Definition: G4Types.hh:83
@ MRDiffuseReflection
@ SpecularReflection
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const

Referenced by PostStepDoIt().

◆ MRreflectHigh()

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

Definition at line 586 of file G4UCNBoundaryProcess.cc.

594{
595 // Only for Enormal > VFermi
596
597 G4double costheta = OldMomentum*Normal;
598
599 G4double Enormal = Energy * (costheta*costheta);
600
601// G4double pSpecular = Reflectivity(Enormal,FermiPot)*
602 G4double pSpecular = Reflectivity(FermiPot,Enormal)*
603 (1.-pDiffuse-pDiffuseTrans-pLoss);
604
605 G4ThreeVector NewMomentum;
606
607 G4double decide = G4UniformRand();
608
609 if (decide < pSpecular) {
610
611 // Reflect specularly
612
613 G4double PdotN = OldMomentum * Normal;
614 NewMomentum = OldMomentum - (2.*PdotN)*Normal;
615 NewMomentum.unit();
616
617 Enew = Energy;
618
619 aSpecularReflection++;
620 theStatus = SpecularReflection;
621 if ( verboseLevel ) BoundaryProcessVerbose();
622
623 return NewMomentum;
624 }
625
626 if (decide < pSpecular + pDiffuse) {
627
628 // Reflect diffusely
629
630 // Determines the angles of the non-specular reflection
631
632 NewMomentum =
633 MRDiffRefl(Normal, Energy, FermiPot, OldMomentum, pDiffuse);
634
635 if (verboseLevel > 0) G4cout << "Diffuse normal " << Normal
636 << ", " << NewMomentum << G4endl;
637 Enew = Energy;
638
639 aMRDiffuseReflection++;
640 theStatus = MRDiffuseReflection;
641 if ( verboseLevel ) BoundaryProcessVerbose();
642
643 return NewMomentum;
644 }
645
646 if (decide < pSpecular + pDiffuse + pDiffuseTrans) {
647
648 // Transmit diffusely
649
650 // Determines the angles of the non-specular transmission
651
652 NewMomentum =
653 MRDiffTrans(Normal, Energy, FermiPot, OldMomentum, pDiffuseTrans);
654
655 Enew = Energy - FermiPot;
656
657 aMRDiffuseTransmit++;
658 theStatus = MRDiffuseTransmit;
659 if ( verboseLevel ) BoundaryProcessVerbose();
660
661 return NewMomentum;
662 }
663
664 if (decide < pSpecular + pDiffuse + pDiffuseTrans + pLoss) {
665
666 // Loss
667
668 Enew = 0.;
669
670 nEzero++;
671 theStatus = Ezero;
672 if ( verboseLevel > 0 ) BoundaryProcessVerbose();
673
674 return NewMomentum;
675 }
676
677 // last case: Refractive transmission
678
679 Enew = Energy - FermiPot;
680
681 G4double palt = std::sqrt(2*neutron_mass_c2/c_squared*Energy);
682 G4double produ = OldMomentum * Normal;
683
684 NewMomentum = palt*OldMomentum-
685 (palt*produ+std::sqrt(palt*palt*produ*produ-2*neutron_mass_c2/
686 c_squared*FermiPot))*Normal;
687
688 mSnellTransmit++;
689 theStatus = SnellTransmit;
690 if ( verboseLevel > 0 ) BoundaryProcessVerbose();
691
692 return NewMomentum.unit();
693}
@ SnellTransmit
@ MRDiffuseTransmit

Referenced by PostStepDoIt().

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 110 of file G4UCNBoundaryProcess.cc.

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

◆ SetMaterialPropertiesTable1()

void G4UCNBoundaryProcess::SetMaterialPropertiesTable1 ( G4UCNMaterialPropertiesTable MPT)
inline

Definition at line 205 of file G4UCNBoundaryProcess.hh.

206 {aMaterialPropertiesTable1 = MPT;}

◆ SetMaterialPropertiesTable2()

void G4UCNBoundaryProcess::SetMaterialPropertiesTable2 ( G4UCNMaterialPropertiesTable MPT)
inline

Definition at line 208 of file G4UCNBoundaryProcess.hh.

209 {aMaterialPropertiesTable2 = MPT;}

◆ SetMicroRoughness()

void G4UCNBoundaryProcess::SetMicroRoughness ( G4bool  active)
inline

Definition at line 241 of file G4UCNBoundaryProcess.hh.

242{
243 UseMicroRoughnessReflection = active;
244}

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