Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4OpBoundaryProcess Class Reference

#include <G4OpBoundaryProcess.hh>

Inheritance diagram for G4OpBoundaryProcess:
G4VDiscreteProcess G4VProcess

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 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.

References fOpBoundary, G4cout, G4endl, G4GeometryTolerance::GetInstance(), G4VProcess::GetProcessName(), G4GeometryTolerance::GetSurfaceTolerance(), glisur, polished, G4VProcess::SetProcessSubType(), Undefined, and G4VProcess::verboseLevel.

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  prob_sl = 0.;
121  prob_ss = 0.;
122  prob_bs = 0.;
123 
124  PropertyPointer = NULL;
125  PropertyPointer1 = NULL;
126  PropertyPointer2 = NULL;
127 
128  Material1 = NULL;
129  Material2 = NULL;
130 
131  OpticalSurface = NULL;
132 
133  kCarTolerance = G4GeometryTolerance::GetInstance()
135 
136  iTE = iTM = 0;
137  thePhotonMomentum = 0.;
138  Rindex1 = Rindex2 = cost1 = cost2 = sint1 = sint2 = 0.;
139 
140  idx = idy = 0;
141  DichroicVector = NULL;
142 }
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()
G4OpBoundaryProcess::~G4OpBoundaryProcess ( )

Definition at line 152 of file G4OpBoundaryProcess.cc.

152 {}

Member Function Documentation

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

Implements G4VDiscreteProcess.

Definition at line 1121 of file G4OpBoundaryProcess.cc.

References DBL_MAX, and Forced.

1124 {
1125  *condition = Forced;
1126 
1127  return DBL_MAX;
1128 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83
G4OpBoundaryProcessStatus G4OpBoundaryProcess::GetStatus ( ) const
inline

Definition at line 277 of file G4OpBoundaryProcess.hh.

Referenced by WLSSteppingAction::UserSteppingAction().

278 {
279  return theStatus;
280 }
G4bool G4OpBoundaryProcess::IsApplicable ( const G4ParticleDefinition aParticleType)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 270 of file G4OpBoundaryProcess.hh.

References G4OpticalPhoton::OpticalPhoton().

272 {
273  return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
274 }
static G4OpticalPhoton * OpticalPhoton()
G4VParticleChange * G4OpBoundaryProcess::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 163 of file G4OpBoundaryProcess.cc.

References G4VProcess::aParticleChange, Detection, dielectric_dichroic, dielectric_dielectric, dielectric_LUT, dielectric_metal, EventMustBeAborted, fGeomBoundary, FresnelRefraction, fStopAndKill, G4cerr, G4cout, G4endl, G4Exception(), G4Track::GetDynamicParticle(), G4OpticalSurface::GetFinish(), G4ParallelWorldProcess::GetHyperStep(), G4ParallelWorldProcess::GetHypNavigatorID(), G4VPhysicalVolume::GetLogicalVolume(), G4StepPoint::GetMaterial(), G4Material::GetMaterialPropertiesTable(), G4OpticalSurface::GetModel(), G4DynamicParticle::GetMomentumDirection(), G4VPhysicalVolume::GetMotherLogical(), G4VPhysicalVolume::GetName(), G4StepPoint::GetPhysicalVolume(), G4DynamicParticle::GetPolarization(), G4StepPoint::GetPosition(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4MaterialPropertiesTable::GetProperty(), G4Track::GetStepLength(), G4StepPoint::GetStepStatus(), G4LogicalSkinSurface::GetSurface(), G4LogicalBorderSurface::GetSurface(), G4LogicalSurface::GetSurfaceProperty(), G4DynamicParticle::GetTotalMomentum(), G4TransportationManager::GetTransportationManager(), G4SurfaceProperty::GetType(), G4Track::GetVelocity(), glisur, ground, groundbackpainted, groundfrontpainted, G4ParticleChange::Initialize(), LambertianReflection, NoRINDEX, NotAtBoundary, polished, polishedbackpainted, polishedfrontpainted, G4VDiscreteProcess::PostStepDoIt(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChange::ProposeMomentumDirection(), G4ParticleChange::ProposePolarization(), G4VParticleChange::ProposeTrackStatus(), G4ParticleChange::ProposeVelocity(), SameMaterial, StepTooSmall, Undefined, unified, CLHEP::Hep3Vector::unit(), G4PhysicsVector::Value(), and G4VProcess::verboseLevel.

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

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