Geant4 10.7.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)
 
virtual ~G4OpBoundaryProcess ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &aParticleType) override
 
virtual G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *condition) override
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
virtual G4OpBoundaryProcessStatus GetStatus () const
 
virtual void SetInvokeSD (G4bool)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &) override
 
virtual void Initialise ()
 
- 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 &)
 
- 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 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 116 of file G4OpBoundaryProcess.hh.

Constructor & Destructor Documentation

◆ G4OpBoundaryProcess()

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

Definition at line 91 of file G4OpBoundaryProcess.cc.

93 : G4VDiscreteProcess(processName, type)
94{
95 Initialise();
96
97 if(verboseLevel > 0)
98 {
99 G4cout << GetProcessName() << " is created " << G4endl;
100 }
102
103 theStatus = Undefined;
104 theModel = glisur;
105 theFinish = polished;
106 theReflectivity = 1.;
107 theEfficiency = 0.;
108 theTransmittance = 0.;
109 theSurfaceRoughness = 0.;
110 prob_sl = 0.;
111 prob_ss = 0.;
112 prob_bs = 0.;
113
114 fRealRIndexMPV = nullptr;
115 fImagRIndexMPV = nullptr;
116 Material1 = nullptr;
117 Material2 = nullptr;
118 OpticalSurface = nullptr;
120
121 iTE = iTM = 0;
122 thePhotonMomentum = 0.;
123 Rindex1 = Rindex2 = 1.;
124 cost1 = cost2 = sint1 = sint2 = 0.;
125 idx = idy = 0;
126 DichroicVector = nullptr;
127}
@ Undefined
@ fOpBoundary
@ glisur
@ polished
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
G4int verboseLevel
Definition: G4VProcess.hh:356
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

◆ ~G4OpBoundaryProcess()

G4OpBoundaryProcess::~G4OpBoundaryProcess ( )
virtual

Definition at line 130 of file G4OpBoundaryProcess.cc.

130{}

Member Function Documentation

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 1295 of file G4OpBoundaryProcess.cc.

1297{
1298 *condition = Forced;
1299 return DBL_MAX;
1300}
G4double condition(const G4ErrorSymMatrix &m)
@ Forced
#define DBL_MAX
Definition: templates.hh:62

◆ GetStatus()

G4OpBoundaryProcessStatus G4OpBoundaryProcess::GetStatus ( ) const
inlinevirtual

Definition at line 260 of file G4OpBoundaryProcess.hh.

261{
262 return theStatus;
263}

◆ Initialise()

void G4OpBoundaryProcess::Initialise ( )
virtual

Definition at line 139 of file G4OpBoundaryProcess.cc.

140{
144}
virtual void SetInvokeSD(G4bool)
G4bool GetBoundaryInvokeSD() const
G4int GetBoundaryVerboseLevel() const
static G4OpticalParameters * Instance()
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:412

Referenced by G4OpBoundaryProcess(), and PreparePhysicsTable().

◆ IsApplicable()

G4bool G4OpBoundaryProcess::IsApplicable ( const G4ParticleDefinition aParticleType)
inlineoverridevirtual

Reimplemented from G4VProcess.

Definition at line 254 of file G4OpBoundaryProcess.hh.

256{
257 return (&aParticleType == G4OpticalPhoton::OpticalPhoton());
258}
static G4OpticalPhoton * OpticalPhoton()

Referenced by G4OpticalPhysics::ConstructProcess().

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 147 of file G4OpBoundaryProcess.cc.

149{
150 theStatus = Undefined;
153
154 // Get hyperStep from G4ParallelWorldProcess
155 // NOTE: PostSetpDoIt of this process to be invoked after
156 // G4ParallelWorldProcess!
157 const G4Step* pStep = &aStep;
159 if(hStep)
160 pStep = hStep;
161
162 G4bool isOnBoundary =
164 if(isOnBoundary)
165 {
166 Material1 = pStep->GetPreStepPoint()->GetMaterial();
167 Material2 = pStep->GetPostStepPoint()->GetMaterial();
168 }
169 else
170 {
171 theStatus = NotAtBoundary;
172 if(verboseLevel > 1)
173 BoundaryProcessVerbose();
174 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
175 }
176
179
180 if(verboseLevel > 1)
181 {
182 G4cout << " Photon at Boundary! " << G4endl;
183 if(thePrePV)
184 G4cout << " thePrePV: " << thePrePV->GetName() << G4endl;
185 if(thePostPV)
186 G4cout << " thePostPV: " << thePostPV->GetName() << G4endl;
187 }
188
189 if(aTrack.GetStepLength() <= kCarTolerance)
190 {
191 theStatus = StepTooSmall;
192 if(verboseLevel > 1)
193 BoundaryProcessVerbose();
194 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
195 }
196
197 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
198
199 thePhotonMomentum = aParticle->GetTotalMomentum();
200 OldMomentum = aParticle->GetMomentumDirection();
201 OldPolarization = aParticle->GetPolarization();
202
203 if(verboseLevel > 1)
204 {
205 G4cout << " Old Momentum Direction: " << OldMomentum << G4endl
206 << " Old Polarization: " << OldPolarization << G4endl;
207 }
208
209 G4ThreeVector theGlobalPoint = pStep->GetPostStepPoint()->GetPosition();
210 G4bool valid;
211
212 // ID of Navigator which limits step
216 theGlobalNormal = (iNav[hNavId])->GetGlobalExitNormal(theGlobalPoint, &valid);
217
218 if(valid)
219 {
220 theGlobalNormal = -theGlobalNormal;
221 }
222 else
223 {
225 ed << " G4OpBoundaryProcess/PostStepDoIt(): "
226 << " The Navigator reports that it returned an invalid normal" << G4endl;
228 "G4OpBoundaryProcess::PostStepDoIt", "OpBoun01", EventMustBeAborted, ed,
229 "Invalid Surface Normal - Geometry must return valid surface normal");
230 }
231
232 if(OldMomentum * theGlobalNormal > 0.0)
233 {
234#ifdef G4OPTICAL_DEBUG
236 ed << " G4OpBoundaryProcess/PostStepDoIt(): theGlobalNormal points in a "
237 "wrong direction. "
238 << G4endl
239 << " The momentum of the photon arriving at interface (oldMomentum)"
240 << " must exit the volume cross in the step. " << G4endl
241 << " So it MUST have dot < 0 with the normal that Exits the new "
242 "volume (globalNormal)."
243 << G4endl << " >> The dot product of oldMomentum and global Normal is "
244 << OldMomentum * theGlobalNormal << G4endl
245 << " Old Momentum (during step) = " << OldMomentum << G4endl
246 << " Global Normal (Exiting New Vol) = " << theGlobalNormal << G4endl
247 << G4endl;
248 G4Exception("G4OpBoundaryProcess::PostStepDoIt", "OpBoun02",
249 EventMustBeAborted, // Or JustWarning to see if it happens
250 // repeatedly on one ray
251 ed,
252 "Invalid Surface Normal - Geometry must return valid surface "
253 "normal pointing in the right direction");
254#else
255 theGlobalNormal = -theGlobalNormal;
256#endif
257 }
258
259 G4MaterialPropertyVector* RindexMPV = nullptr;
261 if(MPT)
262 {
263 RindexMPV = MPT->GetProperty(kRINDEX);
264 }
265 else
266 {
267 theStatus = NoRINDEX;
268 if(verboseLevel > 1)
269 BoundaryProcessVerbose();
272 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
273 }
274
275 if(RindexMPV)
276 {
277 Rindex1 = RindexMPV->Value(thePhotonMomentum, idx_rindex1);
278 }
279 else
280 {
281 theStatus = NoRINDEX;
282 if(verboseLevel > 1)
283 BoundaryProcessVerbose();
286 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
287 }
288
289 theReflectivity = 1.;
290 theEfficiency = 0.;
291 theTransmittance = 0.;
292 theSurfaceRoughness = 0.;
293 theModel = glisur;
294 theFinish = polished;
296
297 RindexMPV = nullptr;
298 OpticalSurface = nullptr;
299 G4LogicalSurface* Surface =
300 G4LogicalBorderSurface::GetSurface(thePrePV, thePostPV);
301
302 if(Surface == nullptr)
303 {
304 if(thePostPV->GetMotherLogical() == thePrePV->GetLogicalVolume())
305 {
307 if(Surface == nullptr)
308 {
309 Surface =
311 }
312 }
313 else
314 {
316 if(Surface == nullptr)
317 {
318 Surface =
320 }
321 }
322 }
323
324 if(Surface)
325 {
326 OpticalSurface =
327 dynamic_cast<G4OpticalSurface*>(Surface->GetSurfaceProperty());
328 }
329 if(OpticalSurface)
330 {
331 type = OpticalSurface->GetType();
332 theModel = OpticalSurface->GetModel();
333 theFinish = OpticalSurface->GetFinish();
334
336 OpticalSurface->GetMaterialPropertiesTable();
337 if(sMPT)
338 {
339 if(theFinish == polishedbackpainted || theFinish == groundbackpainted)
340 {
341 RindexMPV = sMPT->GetProperty(kRINDEX);
342 if(RindexMPV)
343 {
344 Rindex2 = RindexMPV->Value(thePhotonMomentum, idx_rindex_surface);
345 }
346 else
347 {
348 theStatus = NoRINDEX;
349 if(verboseLevel > 1)
350 BoundaryProcessVerbose();
353 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
354 }
355 }
356
357 fRealRIndexMPV = sMPT->GetProperty(kREALRINDEX);
358 fImagRIndexMPV = sMPT->GetProperty(kIMAGINARYRINDEX);
359 iTE = iTM = 1;
360
362 if((pp = sMPT->GetProperty(kREFLECTIVITY)))
363 {
364 theReflectivity = pp->Value(thePhotonMomentum, idx_reflect);
365 }
366 else if(fRealRIndexMPV && fImagRIndexMPV)
367 {
368 CalculateReflectivity();
369 }
370
371 if((pp = sMPT->GetProperty(kEFFICIENCY)))
372 {
373 theEfficiency = pp->Value(thePhotonMomentum, idx_eff);
374 }
375 if((pp = sMPT->GetProperty(kTRANSMITTANCE)))
376 {
377 theTransmittance = pp->Value(thePhotonMomentum, idx_trans);
378 }
380 {
381 theSurfaceRoughness = sMPT->GetConstProperty(kSURFACEROUGHNESS);
382 }
383
384 if(theModel == unified)
385 {
386 prob_sl = (pp = sMPT->GetProperty(kSPECULARLOBECONSTANT))
387 ? pp->Value(thePhotonMomentum, idx_lobe)
388 : 0.;
389 prob_ss = (pp = sMPT->GetProperty(kSPECULARSPIKECONSTANT))
390 ? pp->Value(thePhotonMomentum, idx_spike)
391 : 0.;
392 prob_bs = (pp = sMPT->GetProperty(kBACKSCATTERCONSTANT))
393 ? pp->Value(thePhotonMomentum, idx_back)
394 : 0.;
395 }
396 } // end of if(sMPT)
397 else if(theFinish == polishedbackpainted || theFinish == groundbackpainted)
398 {
401 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
402 }
403 } // end of if(OpticalSurface)
404
405 // DIELECTRIC-DIELECTRIC
406 if(type == dielectric_dielectric)
407 {
408 if(theFinish == polished || theFinish == ground)
409 {
410 if(Material1 == Material2)
411 {
412 theStatus = SameMaterial;
413 if(verboseLevel > 1)
414 BoundaryProcessVerbose();
415 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
416 }
417 MPT = Material2->GetMaterialPropertiesTable();
418 if(MPT)
419 {
420 RindexMPV = MPT->GetProperty(kRINDEX);
421 }
422 if(RindexMPV)
423 {
424 Rindex2 = RindexMPV->Value(thePhotonMomentum, idx_rindex2);
425 }
426 else
427 {
428 theStatus = NoRINDEX;
429 if(verboseLevel > 1)
430 BoundaryProcessVerbose();
433 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
434 }
435 }
436 if(theFinish == polishedbackpainted || theFinish == groundbackpainted)
437 {
438 DielectricDielectric();
439 }
440 else
441 {
442 G4double rand = G4UniformRand();
443 if(rand > theReflectivity + theTransmittance)
444 {
445 DoAbsorption();
446 }
447 else if(rand > theReflectivity)
448 {
449 theStatus = Transmission;
450 NewMomentum = OldMomentum;
451 NewPolarization = OldPolarization;
452 }
453 else
454 {
455 if(theFinish == polishedfrontpainted)
456 {
457 DoReflection();
458 }
459 else if(theFinish == groundfrontpainted)
460 {
461 theStatus = LambertianReflection;
462 DoReflection();
463 }
464 else
465 {
466 DielectricDielectric();
467 }
468 }
469 }
470 }
471 else if(type == dielectric_metal)
472 {
473 DielectricMetal();
474 }
475 else if(type == dielectric_LUT)
476 {
477 DielectricLUT();
478 }
479 else if(type == dielectric_LUTDAVIS)
480 {
481 DielectricLUTDAVIS();
482 }
483 else if(type == dielectric_dichroic)
484 {
485 DielectricDichroic();
486 }
487 else
488 {
490 ed << " PostStepDoIt(): Illegal boundary type." << G4endl;
491 G4Exception("G4OpBoundaryProcess", "OpBoun04", JustWarning, ed);
492 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
493 }
494
495 NewMomentum = NewMomentum.unit();
496 NewPolarization = NewPolarization.unit();
497
498 if(verboseLevel > 1)
499 {
500 G4cout << " New Momentum Direction: " << NewMomentum << G4endl
501 << " New Polarization: " << NewPolarization << G4endl;
502 BoundaryProcessVerbose();
503 }
504
506 aParticleChange.ProposePolarization(NewPolarization);
507
508 if(theStatus == FresnelRefraction || theStatus == Transmission)
509 {
510 G4MaterialPropertyVector* groupvel =
512 if(groupvel)
513 {
515 groupvel->Value(thePhotonMomentum, idx_groupvel));
516 }
517 }
518
519 if(theStatus == Detection && fInvokeSD)
520 InvokeSD(pStep);
521 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
522}
@ JustWarning
@ EventMustBeAborted
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
@ kBACKSCATTERCONSTANT
@ kSPECULARLOBECONSTANT
@ kSPECULARSPIKECONSTANT
@ NotAtBoundary
@ Transmission
@ NoRINDEX
@ SameMaterial
@ StepTooSmall
@ LambertianReflection
@ FresnelRefraction
@ Detection
@ unified
@ groundfrontpainted
@ polishedbackpainted
@ ground
@ polishedfrontpainted
@ groundbackpainted
@ fGeomBoundary
Definition: G4StepStatus.hh:43
G4SurfaceType
@ dielectric_metal
@ dielectric_LUT
@ dielectric_dielectric
@ dielectric_LUTDAVIS
@ dielectric_dichroic
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
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, G4bool warning=false)
G4bool ConstPropertyExists(const G4String &key) const
G4double GetConstProperty(const G4String &key) const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:254
G4OpticalSurfaceModel GetModel() const
G4OpticalSurfaceFinish GetFinish() const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
static const G4Step * GetHyperStep()
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, std::size_t &lastidx) const
G4StepStatus GetStepStatus() const
G4Material * GetMaterial() const
const G4ThreeVector & GetPosition() const
G4VPhysicalVolume * GetPhysicalVolume() const
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
const G4SurfaceType & GetType() const
G4double GetVelocity() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetStepLength() const
std::vector< G4Navigator * >::iterator GetActiveNavigatorsIterator()
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:327

◆ PreparePhysicsTable()

void G4OpBoundaryProcess::PreparePhysicsTable ( const G4ParticleDefinition )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 133 of file G4OpBoundaryProcess.cc.

134{
135 Initialise();
136}

◆ SetInvokeSD()

void G4OpBoundaryProcess::SetInvokeSD ( G4bool  flag)
inlinevirtual

Definition at line 265 of file G4OpBoundaryProcess.hh.

265{ fInvokeSD = flag; }

Referenced by Initialise().


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