Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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)
 
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 &, 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 G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
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 &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 83 of file G4UCNBoundaryProcess.hh.

Constructor & Destructor Documentation

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 
77  theStatus = Undefined;
78 
79  fMessenger = new G4UCNBoundaryProcessMessenger(this);
80 
81  neV = 1.0e-9*eV;
82 
83  kCarTolerance = G4GeometryTolerance::GetInstance()
85 
86  Material1 = NULL;
87  Material2 = NULL;
88 
89  aMaterialPropertiesTable1 = NULL;
90  aMaterialPropertiesTable2 = NULL;
91 
92  UseMicroRoughnessReflection = false;
93  DoMicroRoughnessReflection = false;
94 
95  nNoMPT = nNoMRT = nNoMRCondition = 0;
96  nAbsorption = nEzero = nFlip = 0;
97  aSpecularReflection = bSpecularReflection = 0;
98  bLambertianReflection = 0;
99  aMRDiffuseReflection = bMRDiffuseReflection = 0;
100  nSnellTransmit = mSnellTransmit = 0;
101  aMRDiffuseTransmit = 0;
102  ftheta_o = fphi_o = 0;
103 }
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double GetSurfaceTolerance() const
G4GLOB_DLL std::ostream G4cout
static constexpr double eV
Definition: G4SIunits.hh:215
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
#define G4endl
Definition: G4ios.hh:61
static G4GeometryTolerance * GetInstance()

Here is the call graph for this function:

G4UCNBoundaryProcess::~G4UCNBoundaryProcess ( )
virtual

Definition at line 105 of file G4UCNBoundaryProcess.cc.

106 {
107  delete fMessenger;
108 }

Member Function Documentation

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: "
959  << aSpecularReflection << G4endl;
960  G4cout << "Sum No. E < V Specular Reflection: "
961  << bSpecularReflection << G4endl;
962  G4cout << "Sum No. E < V Lambertian Reflection: "
963  << bLambertianReflection << G4endl;
964  G4cout << "Sum No. E > V MR DiffuseReflection: "
965  << aMRDiffuseReflection << G4endl;
966  G4cout << "Sum No. E < V MR DiffuseReflection: "
967  << bMRDiffuseReflection << G4endl;
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: "
973  << aMRDiffuseTransmit << G4endl;
974  G4cout << " " << G4endl;
975 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
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
G4bool G4UCNBoundaryProcess::GetMicroRoughness ( )
inline

Definition at line 248 of file G4UCNBoundaryProcess.hh.

249 {
250  return UseMicroRoughnessReflection;
251 }
G4double G4UCNBoundaryProcess::GetPhi_o ( )
inline

Definition at line 213 of file G4UCNBoundaryProcess.hh.

213 {return fphi_o;};
G4UCNBoundaryProcessStatus G4UCNBoundaryProcess::GetStatus ( ) const
inline

Definition at line 228 of file G4UCNBoundaryProcess.hh.

229 {
230  return theStatus;
231 }
G4double G4UCNBoundaryProcess::GetTheta_o ( )
inline

Definition at line 212 of file G4UCNBoundaryProcess.hh.

212 {return ftheta_o;};
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:

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 
564  bMRDiffuseReflection++;
565  theStatus = MRDiffuseReflection;
566  if ( verboseLevel > 0 ) BoundaryProcessVerbose();
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 
579  bSpecularReflection++;
580  theStatus = SpecularReflection;
581  if ( verboseLevel > 0 ) BoundaryProcessVerbose();
582 
583  return NewMomentum;
584  }
585 }
G4int verboseLevel
Definition: G4VProcess.hh:368
#define G4UniformRand()
Definition: Randomize.hh:97
Hep3Vector unit() const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

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 
620  aSpecularReflection++;
621  theStatus = SpecularReflection;
622  if ( verboseLevel ) BoundaryProcessVerbose();
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 
640  aMRDiffuseReflection++;
641  theStatus = MRDiffuseReflection;
642  if ( verboseLevel ) BoundaryProcessVerbose();
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 
658  aMRDiffuseTransmit++;
659  theStatus = MRDiffuseTransmit;
660  if ( verboseLevel ) BoundaryProcessVerbose();
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;
673  if ( verboseLevel > 0 ) BoundaryProcessVerbose();
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++;
690  theStatus = SnellTransmit;
691  if ( verboseLevel > 0 ) BoundaryProcessVerbose();
692 
693  return NewMomentum.unit();
694 }
G4int verboseLevel
Definition: G4VProcess.hh:368
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
static constexpr double neutron_mass_c2
static constexpr double c_squared
Hep3Vector unit() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

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);
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 
123 
124  if (hStep) pStep = hStep;
125 
126  G4bool isOnBoundary =
128 
129  if (isOnBoundary) {
130  Material1 = pStep->GetPreStepPoint()->GetMaterial();
131  Material2 = pStep->GetPostStepPoint()->GetMaterial();
132  } else {
133  theStatus = NotAtBoundary;
134  if ( verboseLevel > 1 ) BoundaryProcessVerbose();
135  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
136  }
137 
138  if (aTrack.GetStepLength()<=kCarTolerance/2) {
139  theStatus = StepTooSmall;
140  if ( verboseLevel > 0 ) BoundaryProcessVerbose();
141  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
142  }
143 
144  aMaterialPropertiesTable1 = (G4UCNMaterialPropertiesTable*)Material1->
145  GetMaterialPropertiesTable();
146  aMaterialPropertiesTable2 = (G4UCNMaterialPropertiesTable*)Material2->
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) {
160  theStatus = SameMaterial;
161  if ( verboseLevel > 0 ) BoundaryProcessVerbose();
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 
235  if (aMaterialPropertiesTable2) {
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.;
243  if (aMaterialPropertiesTable1)
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 
254  DoMicroRoughnessReflection = UseMicroRoughnessReflection;
255 
256  G4double theta_i = 0.;
257 
258  if (!aMaterialPropertiesTable2) {
259 
260  nNoMPT++;
261  theStatus = NoMPT;
262  if ( verboseLevel > 0 ) BoundaryProcessVerbose();
263  DoMicroRoughnessReflection = false;
264 
265  } else {
266 
267  if (!aMaterialPropertiesTable2->GetMicroRoughnessTable()) {
268 
269  nNoMRT++;
270  theStatus = NoMRT;
271  if ( verboseLevel > 0 ) BoundaryProcessVerbose();
272 
273  DoMicroRoughnessReflection = false;
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 
283  if (!aMaterialPropertiesTable2->
284  ConditionsValid(Energy, FermiPotDiff, theta_i)) {
285 
286  nNoMRCondition++;
287  theStatus = NoMRCondition;
288  if ( verboseLevel > 0 ) BoundaryProcessVerbose();
289 
290  DoMicroRoughnessReflection = false;
291  }
292  }
293 
294  G4double MRpDiffuse = 0.;
295  G4double MRpDiffuseTrans = 0.;
296 
297  // If microroughness is available and active for material in question
298 
299  if (DoMicroRoughnessReflection) {
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.
334 
335  nAbsorption++;
336  theStatus = Absorption;
337  if ( verboseLevel > 0 ) BoundaryProcessVerbose();
338 
339  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
340  }
341 
342  // spinflips
343 
344  if (SpinFlip(pSpinFlip)) {
345  nFlip++;
346  theStatus = Flip;
347  if ( verboseLevel > 0 ) BoundaryProcessVerbose();
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 
359  if (DoMicroRoughnessReflection)
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 
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 
384  if (DoMicroRoughnessReflection) {
385 
386  G4double Enew;
387 
388  NewMomentum =
389  MRreflectHigh(MRpDiffuse, MRpDiffuseTrans, 0., OldMomentum,
390  theGlobalNormal, Energy, FermiPotDiff, Enew);
391 
392  if (Enew == 0.) {
395  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
396  } else {
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);
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++;
435  theStatus = SnellTransmit;
436  if ( verboseLevel > 0 ) BoundaryProcessVerbose();
437 
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
double angle(const Hep3Vector &) const
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetVelocity() const
const G4DynamicParticle * GetDynamicParticle() const
G4StepStatus GetStepStatus() const
G4Material * GetMaterial() const
void ProposePolarization(G4double Px, G4double Py, G4double Pz)
int G4int
Definition: G4Types.hh:78
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4StepPoint * GetPreStepPoint() const
G4double GetKineticEnergy() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4VPhysicalVolume * GetPhysicalVolume() const
const G4String & GetName() const
static constexpr double degree
Definition: G4SIunits.hh:144
const G4ThreeVector & GetPosition() const
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
static constexpr double neutron_mass_c2
Definition: G4Step.hh:76
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4TransportationManager * GetTransportationManager()
G4ThreeVector MRreflect(G4double, G4ThreeVector, G4ThreeVector, G4double, G4double)
G4ThreeVector GetMomentum() const
virtual void Initialize(const G4Track &)
const G4ThreeVector & GetMomentumDirection() const
G4ThreeVector MRreflectHigh(G4double, G4double, G4double, G4ThreeVector, G4ThreeVector, G4double, G4double, G4double &)
Hep3Vector unit() const
G4StepPoint * GetPostStepPoint() const
static constexpr double c_light
const G4ThreeVector & GetPolarization() const
void ProposeEnergy(G4double finalEnergy)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
G4double GetEnergy() const
#define G4endl
Definition: G4ios.hh:61
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void ProposeVelocity(G4double finalVelocity)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4double GetStepLength() const
static const G4Step * GetHyperStep()
G4double GetConstProperty(const char *key) const

Here is the call graph for this function:

void G4UCNBoundaryProcess::SetMaterialPropertiesTable1 ( G4UCNMaterialPropertiesTable MPT)
inline

Definition at line 206 of file G4UCNBoundaryProcess.hh.

207  {aMaterialPropertiesTable1 = MPT;}
void G4UCNBoundaryProcess::SetMaterialPropertiesTable2 ( G4UCNMaterialPropertiesTable MPT)
inline

Definition at line 209 of file G4UCNBoundaryProcess.hh.

210  {aMaterialPropertiesTable2 = MPT;}
void G4UCNBoundaryProcess::SetMicroRoughness ( G4bool  active)
inline

Definition at line 242 of file G4UCNBoundaryProcess.hh.

243 {
244  UseMicroRoughnessReflection = active;
245 }

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