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

#include <G4OpBoundaryProcess.hh>

+ Inheritance diagram for G4OpBoundaryProcess:

Public Member Functions

 G4OpBoundaryProcess (const G4String &processName="OpBoundary", G4ProcessType type=fOptical)
 
 ~G4OpBoundaryProcess ()
 
G4bool IsApplicable (const G4ParticleDefinition &aParticleType)
 
G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *condition)
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
G4OpBoundaryProcessStatus GetStatus () const
 
- 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 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 &)
 
- 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
 
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 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
 

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 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 128 of file G4OpBoundaryProcess.hh.

Constructor & Destructor Documentation

◆ G4OpBoundaryProcess()

G4OpBoundaryProcess::G4OpBoundaryProcess ( const G4String processName = "OpBoundary",
G4ProcessType  type = fOptical 
)

Definition at line 96 of file G4OpBoundaryProcess.cc.

98 : G4VDiscreteProcess(processName, type)
99{
100 if ( verboseLevel > 0) {
101 G4cout << GetProcessName() << " is created " << G4endl;
102 }
103
105
106 theStatus = Undefined;
107 theModel = glisur;
108 theFinish = polished;
109 theReflectivity = 1.;
110 theEfficiency = 0.;
111 theTransmittance = 0.;
112
113 prob_sl = 0.;
114 prob_ss = 0.;
115 prob_bs = 0.;
116
117 PropertyPointer = NULL;
118 PropertyPointer1 = NULL;
119 PropertyPointer2 = NULL;
120
121 Material1 = NULL;
122 Material2 = NULL;
123
124 OpticalSurface = NULL;
125
126 kCarTolerance = G4GeometryTolerance::GetInstance()
128
129 iTE = iTM = 0;
130 thePhotonMomentum = 0.;
131 Rindex1 = Rindex2 = cost1 = cost2 = sint1 = sint2 = 0.;
132}
@ Undefined
@ fOpBoundary
@ glisur
@ polished
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
G4int verboseLevel
Definition: G4VProcess.hh:368
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

◆ ~G4OpBoundaryProcess()

G4OpBoundaryProcess::~G4OpBoundaryProcess ( )

Definition at line 142 of file G4OpBoundaryProcess.cc.

142{}

Member Function Documentation

◆ GetMeanFreePath()

G4double G4OpBoundaryProcess::GetMeanFreePath ( const G4Track ,
G4double  ,
G4ForceCondition condition 
)
virtual

Implements G4VDiscreteProcess.

Definition at line 1051 of file G4OpBoundaryProcess.cc.

1054{
1055 *condition = Forced;
1056
1057 return DBL_MAX;
1058}
G4double condition(const G4ErrorSymMatrix &m)
@ Forced
#define DBL_MAX
Definition: templates.hh:83

◆ GetStatus()

G4OpBoundaryProcessStatus G4OpBoundaryProcess::GetStatus ( ) const
inline

Definition at line 268 of file G4OpBoundaryProcess.hh.

269{
270 return theStatus;
271}

◆ IsApplicable()

G4bool G4OpBoundaryProcess::IsApplicable ( const G4ParticleDefinition aParticleType)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 261 of file G4OpBoundaryProcess.hh.

263{
264 return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
265}
static G4OpticalPhoton * OpticalPhoton()

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 152 of file G4OpBoundaryProcess.cc.

153{
154 theStatus = Undefined;
155
158
159 G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
160 G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
161
162 if ( verboseLevel > 0 ) {
163 G4cout << " Photon at Boundary! " << G4endl;
164 G4VPhysicalVolume* thePrePV = pPreStepPoint->GetPhysicalVolume();
165 G4VPhysicalVolume* thePostPV = pPostStepPoint->GetPhysicalVolume();
166 if (thePrePV) G4cout << " thePrePV: " << thePrePV->GetName() << G4endl;
167 if (thePostPV) G4cout << " thePostPV: " << thePostPV->GetName() << G4endl;
168 }
169
170 if (pPostStepPoint->GetStepStatus() != fGeomBoundary){
171 theStatus = NotAtBoundary;
172 if ( verboseLevel > 0) BoundaryProcessVerbose();
173 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
174 }
175 if (aTrack.GetStepLength()<=kCarTolerance/2){
176 theStatus = StepTooSmall;
177 if ( verboseLevel > 0) BoundaryProcessVerbose();
178 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
179 }
180
181 Material1 = pPreStepPoint -> GetMaterial();
182 Material2 = pPostStepPoint -> GetMaterial();
183
184 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
185
186 thePhotonMomentum = aParticle->GetTotalMomentum();
187 OldMomentum = aParticle->GetMomentumDirection();
188 OldPolarization = aParticle->GetPolarization();
189
190 if ( verboseLevel > 0 ) {
191 G4cout << " Old Momentum Direction: " << OldMomentum << G4endl;
192 G4cout << " Old Polarization: " << OldPolarization << G4endl;
193 }
194
195 G4ThreeVector theGlobalPoint = pPostStepPoint->GetPosition();
196
197 G4Navigator* theNavigator =
199 GetNavigatorForTracking();
200
201 G4bool valid;
202 // Use the new method for Exit Normal in global coordinates,
203 // which provides the normal more reliably.
204 theGlobalNormal =
205 theNavigator->GetGlobalExitNormal(theGlobalPoint,&valid);
206
207 if (valid) {
208 theGlobalNormal = -theGlobalNormal;
209 }
210 else
211 {
213 ed << " G4OpBoundaryProcess/PostStepDoIt(): "
214 << " The Navigator reports that it returned an invalid normal"
215 << G4endl;
216 G4Exception("G4OpBoundaryProcess::PostStepDoIt", "OpBoun01",
218 "Invalid Surface Normal - Geometry must return valid surface normal");
219 }
220
221 if (OldMomentum * theGlobalNormal > 0.0) {
222#ifdef G4OPTICAL_DEBUG
224 ed << " G4OpBoundaryProcess/PostStepDoIt(): "
225 << " theGlobalNormal points in a wrong direction. "
226 << G4endl;
227 ed << " The momentum of the photon arriving at interface (oldMomentum)"
228 << " must exit the volume cross in the step. " << G4endl;
229 ed << " So it MUST have dot < 0 with the normal that Exits the new volume (globalNormal)." << G4endl;
230 ed << " >> The dot product of oldMomentum and global Normal is " << OldMomentum*theGlobalNormal << G4endl;
231 ed << " Old Momentum (during step) = " << OldMomentum << G4endl;
232 ed << " Global Normal (Exiting New Vol) = " << theGlobalNormal << G4endl;
233 ed << G4endl;
234 G4Exception("G4OpBoundaryProcess::PostStepDoIt", "OpBoun02",
235 EventMustBeAborted, // Or JustWarning to see if it happens repeatedbly on one ray
236 ed,
237 "Invalid Surface Normal - Geometry must return valid surface normal pointing in the right direction");
238#else
239 theGlobalNormal = -theGlobalNormal;
240#endif
241 }
242
243 G4MaterialPropertiesTable* aMaterialPropertiesTable;
245
246 aMaterialPropertiesTable = Material1->GetMaterialPropertiesTable();
247 if (aMaterialPropertiesTable) {
248 Rindex = aMaterialPropertiesTable->GetProperty("RINDEX");
249 }
250 else {
251 theStatus = NoRINDEX;
252 if ( verboseLevel > 0) BoundaryProcessVerbose();
255 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
256 }
257
258 if (Rindex) {
259 Rindex1 = Rindex->Value(thePhotonMomentum);
260 }
261 else {
262 theStatus = NoRINDEX;
263 if ( verboseLevel > 0) BoundaryProcessVerbose();
266 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
267 }
268
269 theReflectivity = 1.;
270 theEfficiency = 0.;
271 theTransmittance = 0.;
272
273 theModel = glisur;
274 theFinish = polished;
275
277
278 Rindex = NULL;
279 OpticalSurface = NULL;
280
281 G4LogicalSurface* Surface = NULL;
282
284 (pPreStepPoint ->GetPhysicalVolume(),
285 pPostStepPoint->GetPhysicalVolume());
286
287 if (Surface == NULL){
288 G4bool enteredDaughter=(pPostStepPoint->GetPhysicalVolume()
289 ->GetMotherLogical() ==
290 pPreStepPoint->GetPhysicalVolume()
291 ->GetLogicalVolume());
292 if(enteredDaughter){
294 (pPostStepPoint->GetPhysicalVolume()->
295 GetLogicalVolume());
296 if(Surface == NULL)
298 (pPreStepPoint->GetPhysicalVolume()->
299 GetLogicalVolume());
300 }
301 else {
303 (pPreStepPoint->GetPhysicalVolume()->
304 GetLogicalVolume());
305 if(Surface == NULL)
307 (pPostStepPoint->GetPhysicalVolume()->
308 GetLogicalVolume());
309 }
310 }
311
312 if (Surface) OpticalSurface =
313 dynamic_cast <G4OpticalSurface*> (Surface->GetSurfaceProperty());
314
315 if (OpticalSurface) {
316
317 type = OpticalSurface->GetType();
318 theModel = OpticalSurface->GetModel();
319 theFinish = OpticalSurface->GetFinish();
320
321 aMaterialPropertiesTable = OpticalSurface->
322 GetMaterialPropertiesTable();
323
324 if (aMaterialPropertiesTable) {
325
326 if (theFinish == polishedbackpainted ||
327 theFinish == groundbackpainted ) {
328 Rindex = aMaterialPropertiesTable->GetProperty("RINDEX");
329 if (Rindex) {
330 Rindex2 = Rindex->Value(thePhotonMomentum);
331 }
332 else {
333 theStatus = NoRINDEX;
334 if ( verboseLevel > 0) BoundaryProcessVerbose();
337 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
338 }
339 }
340
341 PropertyPointer =
342 aMaterialPropertiesTable->GetProperty("REFLECTIVITY");
343 PropertyPointer1 =
344 aMaterialPropertiesTable->GetProperty("REALRINDEX");
345 PropertyPointer2 =
346 aMaterialPropertiesTable->GetProperty("IMAGINARYRINDEX");
347
348 iTE = 1;
349 iTM = 1;
350
351 if (PropertyPointer) {
352
353 theReflectivity =
354 PropertyPointer->Value(thePhotonMomentum);
355
356 } else if (PropertyPointer1 && PropertyPointer2) {
357
358 CalculateReflectivity();
359
360 }
361
362 PropertyPointer =
363 aMaterialPropertiesTable->GetProperty("EFFICIENCY");
364 if (PropertyPointer) {
365 theEfficiency =
366 PropertyPointer->Value(thePhotonMomentum);
367 }
368
369 PropertyPointer =
370 aMaterialPropertiesTable->GetProperty("TRANSMITTANCE");
371 if (PropertyPointer) {
372 theTransmittance =
373 PropertyPointer->Value(thePhotonMomentum);
374 }
375
376 if ( theModel == unified ) {
377 PropertyPointer =
378 aMaterialPropertiesTable->GetProperty("SPECULARLOBECONSTANT");
379 if (PropertyPointer) {
380 prob_sl =
381 PropertyPointer->Value(thePhotonMomentum);
382 } else {
383 prob_sl = 0.0;
384 }
385
386 PropertyPointer =
387 aMaterialPropertiesTable->GetProperty("SPECULARSPIKECONSTANT");
388 if (PropertyPointer) {
389 prob_ss =
390 PropertyPointer->Value(thePhotonMomentum);
391 } else {
392 prob_ss = 0.0;
393 }
394
395 PropertyPointer =
396 aMaterialPropertiesTable->GetProperty("BACKSCATTERCONSTANT");
397 if (PropertyPointer) {
398 prob_bs =
399 PropertyPointer->Value(thePhotonMomentum);
400 } else {
401 prob_bs = 0.0;
402 }
403 }
404 }
405 else if (theFinish == polishedbackpainted ||
406 theFinish == groundbackpainted ) {
409 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
410 }
411 }
412
413 if (type == dielectric_dielectric ) {
414 if (theFinish == polished || theFinish == ground ) {
415
416 if (Material1 == Material2){
417 theStatus = SameMaterial;
418 if ( verboseLevel > 0) BoundaryProcessVerbose();
419 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
420 }
421 aMaterialPropertiesTable =
422 Material2->GetMaterialPropertiesTable();
423 if (aMaterialPropertiesTable)
424 Rindex = aMaterialPropertiesTable->GetProperty("RINDEX");
425 if (Rindex) {
426 Rindex2 = Rindex->Value(thePhotonMomentum);
427 }
428 else {
429 theStatus = NoRINDEX;
430 if ( verboseLevel > 0) BoundaryProcessVerbose();
433 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
434 }
435 }
436 }
437
438 if (type == dielectric_metal) {
439
440 DielectricMetal();
441
442 // Uncomment the following lines if you wish to have
443 // Transmission instead of Absorption
444 // if (theStatus == Absorption) {
445 // return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
446 // }
447
448 }
449 else if (type == dielectric_LUT) {
450
451 DielectricLUT();
452
453 }
454 else if (type == dielectric_dielectric) {
455
456 if ( theFinish == polishedbackpainted ||
457 theFinish == groundbackpainted ) {
458 DielectricDielectric();
459 }
460 else {
461 if ( !G4BooleanRand(theReflectivity) ) {
462 DoAbsorption();
463 }
464 else {
465 if ( theFinish == polishedfrontpainted ) {
466 DoReflection();
467 }
468 else if ( theFinish == groundfrontpainted ) {
469 theStatus = LambertianReflection;
470 DoReflection();
471 }
472 else {
473 DielectricDielectric();
474 }
475 }
476 }
477 }
478 else {
479
480 G4cerr << " Error: G4BoundaryProcess: illegal boundary type " << G4endl;
481 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
482
483 }
484
485 NewMomentum = NewMomentum.unit();
486 NewPolarization = NewPolarization.unit();
487
488 if ( verboseLevel > 0) {
489 G4cout << " New Momentum Direction: " << NewMomentum << G4endl;
490 G4cout << " New Polarization: " << NewPolarization << G4endl;
491 BoundaryProcessVerbose();
492 }
493
495 aParticleChange.ProposePolarization(NewPolarization);
496
497 if ( theStatus == FresnelRefraction ) {
498 G4MaterialPropertyVector* groupvel =
499 Material2->GetMaterialPropertiesTable()->GetProperty("GROUPVEL");
500 G4double finalVelocity = groupvel->Value(thePhotonMomentum);
501 aParticleChange.ProposeVelocity(finalVelocity);
502 }
503
504 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
505}
@ EventMustBeAborted
@ NotAtBoundary
@ SameMaterial
@ StepTooSmall
@ LambertianReflection
@ FresnelRefraction
@ unified
@ groundfrontpainted
@ polishedbackpainted
@ ground
@ polishedfrontpainted
@ groundbackpainted
@ fGeomBoundary
Definition: G4StepStatus.hh:54
G4SurfaceType
@ dielectric_metal
@ dielectric_LUT
@ dielectric_dielectric
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
G4DLLIMPORT std::ostream G4cerr
Hep3Vector unit() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetTotalMomentum() const
const G4ThreeVector & GetPolarization() const
static G4LogicalBorderSurface * GetSurface(const G4VPhysicalVolume *vol1, const G4VPhysicalVolume *vol2)
static G4LogicalSkinSurface * GetSurface(const G4LogicalVolume *vol)
G4SurfaceProperty * GetSurfaceProperty() const
G4MaterialPropertyVector * GetProperty(const char *key)
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:251
virtual G4ThreeVector GetGlobalExitNormal(const G4ThreeVector &point, G4bool *valid)
G4OpticalSurfaceModel GetModel() const
G4OpticalSurfaceFinish GetFinish() const
void ProposePolarization(G4double Px, G4double Py, G4double Pz)
void ProposeVelocity(G4double finalVelocity)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
virtual void Initialize(const G4Track &)
G4double Value(G4double theEnergy)
G4StepStatus GetStepStatus() const
const G4ThreeVector & GetPosition() const
G4VPhysicalVolume * GetPhysicalVolume() const
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
const G4SurfaceType & GetType() const
G4double GetVelocity() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetStepLength() const
static G4TransportationManager * GetTransportationManager()
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4LogicalVolume * GetMotherLogical() const
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

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