Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4OpBoundaryProcess Class Reference

#include <G4OpBoundaryProcess.hh>

Inheritance diagram for G4OpBoundaryProcess:
Collaboration 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
 
void SetInvokeSD (G4bool)
 
- 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 131 of file G4OpBoundaryProcess.hh.

Constructor & Destructor Documentation

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

Definition at line 103 of file G4OpBoundaryProcess.cc.

105  : G4VDiscreteProcess(processName, type)
106 {
107  if ( verboseLevel > 0) {
108  G4cout << GetProcessName() << " is created " << G4endl;
109  }
110 
112 
113  theStatus = Undefined;
114  theModel = glisur;
115  theFinish = polished;
116  theReflectivity = 1.;
117  theEfficiency = 0.;
118  theTransmittance = 0.;
119 
120  theSurfaceRoughness = 0.;
121 
122  prob_sl = 0.;
123  prob_ss = 0.;
124  prob_bs = 0.;
125 
126  PropertyPointer = NULL;
127  PropertyPointer1 = NULL;
128  PropertyPointer2 = NULL;
129 
130  Material1 = NULL;
131  Material2 = NULL;
132 
133  OpticalSurface = NULL;
134 
135  kCarTolerance = G4GeometryTolerance::GetInstance()
137 
138  iTE = iTM = 0;
139  thePhotonMomentum = 0.;
140  Rindex1 = Rindex2 = 1.;
141  cost1 = cost2 = sint1 = sint2 = 0.;
142 
143  idx = idy = 0;
144  DichroicVector = NULL;
145 
146  fInvokeSD = true;
147 }
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double GetSurfaceTolerance() const
G4GLOB_DLL std::ostream G4cout
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:

G4OpBoundaryProcess::~G4OpBoundaryProcess ( )

Definition at line 157 of file G4OpBoundaryProcess.cc.

157 {}

Member Function Documentation

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

Implements G4VDiscreteProcess.

Definition at line 1204 of file G4OpBoundaryProcess.cc.

1207 {
1208  *condition = Forced;
1209 
1210  return DBL_MAX;
1211 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83
G4OpBoundaryProcessStatus G4OpBoundaryProcess::GetStatus ( ) const
inline

Definition at line 285 of file G4OpBoundaryProcess.hh.

286 {
287  return theStatus;
288 }

Here is the caller graph for this function:

G4bool G4OpBoundaryProcess::IsApplicable ( const G4ParticleDefinition aParticleType)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 278 of file G4OpBoundaryProcess.hh.

280 {
281  return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
282 }
static G4OpticalPhoton * OpticalPhoton()

Here is the call graph for this function:

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

Reimplemented from G4VDiscreteProcess.

Definition at line 168 of file G4OpBoundaryProcess.cc.

169 {
170  theStatus = Undefined;
171 
172  aParticleChange.Initialize(aTrack);
174 
175  // Get hyperStep from G4ParallelWorldProcess
176  // NOTE: PostSetpDoIt of this process should be
177  // invoked after G4ParallelWorldProcess!
178 
179  const G4Step* pStep = &aStep;
180 
182 
183  if (hStep) pStep = hStep;
184 
185  G4bool isOnBoundary =
187 
188  if (isOnBoundary) {
189  Material1 = pStep->GetPreStepPoint()->GetMaterial();
190  Material2 = pStep->GetPostStepPoint()->GetMaterial();
191  } else {
192  theStatus = NotAtBoundary;
193  if ( verboseLevel > 0) BoundaryProcessVerbose();
194  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
195  }
196 
197  G4VPhysicalVolume* thePrePV =
198  pStep->GetPreStepPoint() ->GetPhysicalVolume();
199  G4VPhysicalVolume* thePostPV =
201 
202  if ( verboseLevel > 0 ) {
203  G4cout << " Photon at Boundary! " << G4endl;
204  if (thePrePV) G4cout << " thePrePV: " << thePrePV->GetName() << G4endl;
205  if (thePostPV) G4cout << " thePostPV: " << thePostPV->GetName() << G4endl;
206  }
207 
208  if (aTrack.GetStepLength()<=kCarTolerance/2){
209  theStatus = StepTooSmall;
210  if ( verboseLevel > 0) BoundaryProcessVerbose();
211  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
212  }
213 
214  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
215 
216  thePhotonMomentum = aParticle->GetTotalMomentum();
217  OldMomentum = aParticle->GetMomentumDirection();
218  OldPolarization = aParticle->GetPolarization();
219 
220  if ( verboseLevel > 0 ) {
221  G4cout << " Old Momentum Direction: " << OldMomentum << G4endl;
222  G4cout << " Old Polarization: " << OldPolarization << G4endl;
223  }
224 
225  G4ThreeVector theGlobalPoint = pStep->GetPostStepPoint()->GetPosition();
226 
227  G4bool valid;
228  // Use the new method for Exit Normal in global coordinates,
229  // which provides the normal more reliably.
230 
231  // ID of Navigator which limits step
232 
234  std::vector<G4Navigator*>::iterator iNav =
236  GetActiveNavigatorsIterator();
237  theGlobalNormal =
238  (iNav[hNavId])->GetGlobalExitNormal(theGlobalPoint,&valid);
239 
240  if (valid) {
241  theGlobalNormal = -theGlobalNormal;
242  }
243  else
244  {
246  ed << " G4OpBoundaryProcess/PostStepDoIt(): "
247  << " The Navigator reports that it returned an invalid normal"
248  << G4endl;
249  G4Exception("G4OpBoundaryProcess::PostStepDoIt", "OpBoun01",
251  "Invalid Surface Normal - Geometry must return valid surface normal");
252  }
253 
254  if (OldMomentum * theGlobalNormal > 0.0) {
255 #ifdef G4OPTICAL_DEBUG
257  ed << " G4OpBoundaryProcess/PostStepDoIt(): "
258  << " theGlobalNormal points in a wrong direction. "
259  << G4endl;
260  ed << " The momentum of the photon arriving at interface (oldMomentum)"
261  << " must exit the volume cross in the step. " << G4endl;
262  ed << " So it MUST have dot < 0 with the normal that Exits the new volume (globalNormal)." << G4endl;
263  ed << " >> The dot product of oldMomentum and global Normal is " << OldMomentum*theGlobalNormal << G4endl;
264  ed << " Old Momentum (during step) = " << OldMomentum << G4endl;
265  ed << " Global Normal (Exiting New Vol) = " << theGlobalNormal << G4endl;
266  ed << G4endl;
267  G4Exception("G4OpBoundaryProcess::PostStepDoIt", "OpBoun02",
268  EventMustBeAborted, // Or JustWarning to see if it happens repeatedbly on one ray
269  ed,
270  "Invalid Surface Normal - Geometry must return valid surface normal pointing in the right direction");
271 #else
272  theGlobalNormal = -theGlobalNormal;
273 #endif
274  }
275 
276  G4MaterialPropertiesTable* aMaterialPropertiesTable;
277  G4MaterialPropertyVector* Rindex;
278 
279  aMaterialPropertiesTable = Material1->GetMaterialPropertiesTable();
280  if (aMaterialPropertiesTable) {
281  Rindex = aMaterialPropertiesTable->GetProperty("RINDEX");
282  }
283  else {
284  theStatus = NoRINDEX;
285  if ( verboseLevel > 0) BoundaryProcessVerbose();
286  aParticleChange.ProposeLocalEnergyDeposit(thePhotonMomentum);
288  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
289  }
290 
291  if (Rindex) {
292  Rindex1 = Rindex->Value(thePhotonMomentum);
293  }
294  else {
295  theStatus = NoRINDEX;
296  if ( verboseLevel > 0) BoundaryProcessVerbose();
297  aParticleChange.ProposeLocalEnergyDeposit(thePhotonMomentum);
299  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
300  }
301 
302  theReflectivity = 1.;
303  theEfficiency = 0.;
304  theTransmittance = 0.;
305 
306  theSurfaceRoughness = 0.;
307 
308  theModel = glisur;
309  theFinish = polished;
310 
312 
313  Rindex = NULL;
314  OpticalSurface = NULL;
315 
316  G4LogicalSurface* Surface = NULL;
317 
318  Surface = G4LogicalBorderSurface::GetSurface(thePrePV, thePostPV);
319 
320  if (Surface == NULL){
321  G4bool enteredDaughter= (thePostPV->GetMotherLogical() ==
322  thePrePV ->GetLogicalVolume());
323  if(enteredDaughter){
324  Surface =
326  if(Surface == NULL)
327  Surface =
329  }
330  else {
331  Surface =
333  if(Surface == NULL)
334  Surface =
336  }
337  }
338 
339  if (Surface) OpticalSurface =
340  dynamic_cast <G4OpticalSurface*> (Surface->GetSurfaceProperty());
341 
342  if (OpticalSurface) {
343 
344  type = OpticalSurface->GetType();
345  theModel = OpticalSurface->GetModel();
346  theFinish = OpticalSurface->GetFinish();
347 
348  aMaterialPropertiesTable = OpticalSurface->
349  GetMaterialPropertiesTable();
350 
351  if (aMaterialPropertiesTable) {
352 
353  if (theFinish == polishedbackpainted ||
354  theFinish == groundbackpainted ) {
355  Rindex = aMaterialPropertiesTable->GetProperty("RINDEX");
356  if (Rindex) {
357  Rindex2 = Rindex->Value(thePhotonMomentum);
358  }
359  else {
360  theStatus = NoRINDEX;
361  if ( verboseLevel > 0) BoundaryProcessVerbose();
362  aParticleChange.ProposeLocalEnergyDeposit(thePhotonMomentum);
364  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
365  }
366  }
367 
368  PropertyPointer =
369  aMaterialPropertiesTable->GetProperty("REFLECTIVITY");
370  PropertyPointer1 =
371  aMaterialPropertiesTable->GetProperty("REALRINDEX");
372  PropertyPointer2 =
373  aMaterialPropertiesTable->GetProperty("IMAGINARYRINDEX");
374 
375  iTE = 1;
376  iTM = 1;
377 
378  if (PropertyPointer) {
379 
380  theReflectivity =
381  PropertyPointer->Value(thePhotonMomentum);
382 
383  } else if (PropertyPointer1 && PropertyPointer2) {
384 
385  CalculateReflectivity();
386 
387  }
388 
389  PropertyPointer =
390  aMaterialPropertiesTable->GetProperty("EFFICIENCY");
391  if (PropertyPointer) {
392  theEfficiency =
393  PropertyPointer->Value(thePhotonMomentum);
394  }
395 
396  PropertyPointer =
397  aMaterialPropertiesTable->GetProperty("TRANSMITTANCE");
398  if (PropertyPointer) {
399  theTransmittance =
400  PropertyPointer->Value(thePhotonMomentum);
401  }
402 
403  if (aMaterialPropertiesTable->
404  ConstPropertyExists("SURFACEROUGHNESS"))
405  theSurfaceRoughness = aMaterialPropertiesTable->
406  GetConstProperty("SURFACEROUGHNESS");
407 
408  if ( theModel == unified ) {
409  PropertyPointer =
410  aMaterialPropertiesTable->GetProperty("SPECULARLOBECONSTANT");
411  if (PropertyPointer) {
412  prob_sl =
413  PropertyPointer->Value(thePhotonMomentum);
414  } else {
415  prob_sl = 0.0;
416  }
417 
418  PropertyPointer =
419  aMaterialPropertiesTable->GetProperty("SPECULARSPIKECONSTANT");
420  if (PropertyPointer) {
421  prob_ss =
422  PropertyPointer->Value(thePhotonMomentum);
423  } else {
424  prob_ss = 0.0;
425  }
426 
427  PropertyPointer =
428  aMaterialPropertiesTable->GetProperty("BACKSCATTERCONSTANT");
429  if (PropertyPointer) {
430  prob_bs =
431  PropertyPointer->Value(thePhotonMomentum);
432  } else {
433  prob_bs = 0.0;
434  }
435  }
436  }
437  else if (theFinish == polishedbackpainted ||
438  theFinish == groundbackpainted ) {
439  aParticleChange.ProposeLocalEnergyDeposit(thePhotonMomentum);
441  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
442  }
443  }
444 
445  if (type == dielectric_dielectric ) {
446  if (theFinish == polished || theFinish == ground ) {
447 
448  if (Material1 == Material2){
449  theStatus = SameMaterial;
450  if ( verboseLevel > 0) BoundaryProcessVerbose();
451  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
452  }
453  aMaterialPropertiesTable =
454  Material2->GetMaterialPropertiesTable();
455  if (aMaterialPropertiesTable)
456  Rindex = aMaterialPropertiesTable->GetProperty("RINDEX");
457  if (Rindex) {
458  Rindex2 = Rindex->Value(thePhotonMomentum);
459  }
460  else {
461  theStatus = NoRINDEX;
462  if ( verboseLevel > 0) BoundaryProcessVerbose();
463  aParticleChange.ProposeLocalEnergyDeposit(thePhotonMomentum);
465  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
466  }
467  }
468  }
469 
470  if (type == dielectric_metal) {
471 
472  DielectricMetal();
473 
474  }
475  else if (type == dielectric_LUT) {
476 
477  DielectricLUT();
478 
479  }
480  else if (type == dielectric_dichroic) {
481 
482  DielectricDichroic();
483 
484  }
485  else if (type == dielectric_dielectric) {
486 
487  if ( theFinish == polishedbackpainted ||
488  theFinish == groundbackpainted ) {
489  DielectricDielectric();
490  }
491  else {
492  G4double rand = G4UniformRand();
493  if ( rand > theReflectivity ) {
494  if (rand > theReflectivity + theTransmittance) {
495  DoAbsorption();
496  } else {
497  theStatus = Transmission;
498  NewMomentum = OldMomentum;
499  NewPolarization = OldPolarization;
500  }
501  }
502  else {
503  if ( theFinish == polishedfrontpainted ) {
504  DoReflection();
505  }
506  else if ( theFinish == groundfrontpainted ) {
507  theStatus = LambertianReflection;
508  DoReflection();
509  }
510  else {
511  DielectricDielectric();
512  }
513  }
514  }
515  }
516  else {
517 
518  G4cerr << " Error: G4BoundaryProcess: illegal boundary type " << G4endl;
519  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
520 
521  }
522 
523  NewMomentum = NewMomentum.unit();
524  NewPolarization = NewPolarization.unit();
525 
526  if ( verboseLevel > 0) {
527  G4cout << " New Momentum Direction: " << NewMomentum << G4endl;
528  G4cout << " New Polarization: " << NewPolarization << G4endl;
529  BoundaryProcessVerbose();
530  }
531 
533  aParticleChange.ProposePolarization(NewPolarization);
534 
535  if ( theStatus == FresnelRefraction || theStatus == Transmission ) {
536  G4MaterialPropertyVector* groupvel =
537  Material2->GetMaterialPropertiesTable()->GetProperty("GROUPVEL");
538  G4double finalVelocity = groupvel->Value(thePhotonMomentum);
539  aParticleChange.ProposeVelocity(finalVelocity);
540  }
541 
542  if ( theStatus == Detection && fInvokeSD ) InvokeSD(pStep);
543 
544  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
545 }
G4OpticalSurfaceModel GetModel() const
G4int verboseLevel
Definition: G4VProcess.hh:368
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)
const G4SurfaceType & GetType() const
G4SurfaceType
int G4int
Definition: G4Types.hh:78
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
static G4LogicalBorderSurface * GetSurface(const G4VPhysicalVolume *vol1, const G4VPhysicalVolume *vol2)
G4double GetTotalMomentum() const
G4StepPoint * GetPreStepPoint() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4VPhysicalVolume * GetPhysicalVolume() const
const G4String & GetName() const
const G4ThreeVector & GetPosition() const
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
Definition: G4Step.hh:76
G4LogicalVolume * GetMotherLogical() const
G4double Value(G4double theEnergy, size_t &lastidx) const
G4OpticalSurfaceFinish GetFinish() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4TransportationManager * GetTransportationManager()
virtual void Initialize(const G4Track &)
G4LogicalVolume * GetLogicalVolume() const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
Hep3Vector unit() const
G4StepPoint * GetPostStepPoint() const
const G4ThreeVector & GetPolarization() const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
#define G4endl
Definition: G4ios.hh:61
G4SurfaceProperty * GetSurfaceProperty() const
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 &)
G4MaterialPropertyVector * GetProperty(const char *key)
G4double GetStepLength() const
static const G4Step * GetHyperStep()
G4GLOB_DLL std::ostream G4cerr
static G4LogicalSkinSurface * GetSurface(const G4LogicalVolume *vol)

Here is the call graph for this function:

void G4OpBoundaryProcess::SetInvokeSD ( G4bool  flag)
inline

Definition at line 291 of file G4OpBoundaryProcess.hh.

292 {
293  fInvokeSD = flag;
294 }

Here is the caller graph for this function:


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