Geant4 11.2.2
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 ()
 
void SetVerboseLevel (G4int)
 
- 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 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 &)
 
virtual G4double GetCrossSection (const G4double, const G4MaterialCutsCouple *)
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
- 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
 
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 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 const G4VProcessGetCreatorProcess () const
 
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)
 
- Protected Member Functions inherited from G4VDiscreteProcess
- 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 119 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, ptype)
94{
95 Initialise();
96
97 if(verboseLevel > 0)
98 {
99 G4cout << GetProcessName() << " is created " << G4endl;
100 }
102
103 fStatus = Undefined;
104 fModel = glisur;
105 fFinish = polished;
106 fReflectivity = 1.;
107 fEfficiency = 0.;
108 fTransmittance = 0.;
109 fSurfaceRoughness = 0.;
110 fProb_sl = 0.;
111 fProb_ss = 0.;
112 fProb_bs = 0.;
113
114 fRealRIndexMPV = nullptr;
115 fImagRIndexMPV = nullptr;
116 fMaterial1 = nullptr;
117 fMaterial2 = nullptr;
118 fOpticalSurface = nullptr;
120
121 f_iTE = f_iTM = 0;
122 fPhotonMomentum = 0.;
123 fRindex1 = fRindex2 = 1.;
124 fSint1 = 0.;
125 fDichroicVector = nullptr;
126}
@ fOpBoundary
@ polished
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
G4int verboseLevel
void SetProcessSubType(G4int)
const G4String & GetProcessName() const

◆ ~G4OpBoundaryProcess()

G4OpBoundaryProcess::~G4OpBoundaryProcess ( )
virtualdefault

Member Function Documentation

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 1352 of file G4OpBoundaryProcess.cc.

1354{
1355 *condition = Forced;
1356 return DBL_MAX;
1357}
G4double condition(const G4ErrorSymMatrix &m)
@ Forced
#define DBL_MAX
Definition templates.hh:62

◆ GetStatus()

G4OpBoundaryProcessStatus G4OpBoundaryProcess::GetStatus ( ) const
inlinevirtual

Definition at line 276 of file G4OpBoundaryProcess.hh.

277{
278 return fStatus;
279}

◆ Initialise()

void G4OpBoundaryProcess::Initialise ( )
virtual

◆ IsApplicable()

G4bool G4OpBoundaryProcess::IsApplicable ( const G4ParticleDefinition & aParticleType)
inlineoverridevirtual

Reimplemented from G4VProcess.

Definition at line 270 of file G4OpBoundaryProcess.hh.

272{
273 return (&aParticleType == G4OpticalPhoton::OpticalPhoton());
274}
static G4OpticalPhoton * OpticalPhoton()

Referenced by G4OpticalPhysics::ConstructProcess().

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 146 of file G4OpBoundaryProcess.cc.

148{
149 fStatus = Undefined;
152
153 // Get hyperStep from G4ParallelWorldProcess
154 // NOTE: PostSetpDoIt of this process to be invoked after
155 // G4ParallelWorldProcess!
156 const G4Step* pStep = &aStep;
158 if(hStep != nullptr)
159 pStep = hStep;
160
162 {
163 fMaterial1 = pStep->GetPreStepPoint()->GetMaterial();
164 fMaterial2 = pStep->GetPostStepPoint()->GetMaterial();
165 }
166 else
167 {
168 fStatus = NotAtBoundary;
169 if(verboseLevel > 1)
170 BoundaryProcessVerbose();
171 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
172 }
173
176
177 if(verboseLevel > 1)
178 {
179 G4cout << " Photon at Boundary! " << G4endl;
180 if(thePrePV != nullptr)
181 G4cout << " thePrePV: " << thePrePV->GetName() << G4endl;
182 if(thePostPV != nullptr)
183 G4cout << " thePostPV: " << thePostPV->GetName() << G4endl;
184 }
185
186 G4double stepLength = aTrack.GetStepLength();
187 if(stepLength <= fCarTolerance)
188 {
189 fStatus = StepTooSmall;
190 if(verboseLevel > 1)
191 BoundaryProcessVerbose();
192
193 G4MaterialPropertyVector* groupvel = nullptr;
195 if(aMPT != nullptr)
196 {
197 groupvel = aMPT->GetProperty(kGROUPVEL);
198 }
199
200 if(groupvel != nullptr)
201 {
203 groupvel->Value(fPhotonMomentum, idx_groupvel));
204 }
205 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
206 }
207 else if (stepLength <= 10.*fCarTolerance && fNumSmallStepWarnings < 10)
208 { // see bug 2510
209 ++fNumSmallStepWarnings;
210 if(verboseLevel > 0)
211 {
213 ed << "G4OpBoundaryProcess: "
214 << "Opticalphoton step length: " << stepLength/mm << " mm." << G4endl
215 << "This is larger than the threshold " << fCarTolerance/mm << " mm "
216 "to set status StepTooSmall." << G4endl
217 << "Boundary scattering may be incorrect. ";
218 if(fNumSmallStepWarnings == 10)
219 {
220 ed << G4endl << "*** Step size warnings stopped.";
221 }
222 G4Exception("G4OpBoundaryProcess", "OpBoun06", JustWarning, ed, "");
223 }
224 }
225
226 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
227
228 fPhotonMomentum = aParticle->GetTotalMomentum();
229 fOldMomentum = aParticle->GetMomentumDirection();
230 fOldPolarization = aParticle->GetPolarization();
231
232 if(verboseLevel > 1)
233 {
234 G4cout << " Old Momentum Direction: " << fOldMomentum << G4endl
235 << " Old Polarization: " << fOldPolarization << G4endl;
236 }
237
238 G4ThreeVector theGlobalPoint = pStep->GetPostStepPoint()->GetPosition();
239 G4bool valid;
240
241 // ID of Navigator which limits step
245 fGlobalNormal = (iNav[hNavId])->GetGlobalExitNormal(theGlobalPoint, &valid);
246
247 if(valid)
248 {
249 fGlobalNormal = -fGlobalNormal;
250 }
251 else
252 {
254 ed << " G4OpBoundaryProcess/PostStepDoIt(): "
255 << " The Navigator reports that it returned an invalid normal" << G4endl;
257 "G4OpBoundaryProcess::PostStepDoIt", "OpBoun01", EventMustBeAborted, ed,
258 "Invalid Surface Normal - Geometry must return valid surface normal");
259 }
260
261 if(fOldMomentum * fGlobalNormal > 0.0)
262 {
263#ifdef G4OPTICAL_DEBUG
265 ed << " G4OpBoundaryProcess/PostStepDoIt(): fGlobalNormal points in a "
266 "wrong direction. "
267 << G4endl
268 << " The momentum of the photon arriving at interface (oldMomentum)"
269 << " must exit the volume cross in the step. " << G4endl
270 << " So it MUST have dot < 0 with the normal that Exits the new "
271 "volume (globalNormal)."
272 << G4endl << " >> The dot product of oldMomentum and global Normal is "
273 << fOldMomentum * fGlobalNormal << G4endl
274 << " Old Momentum (during step) = " << fOldMomentum << G4endl
275 << " Global Normal (Exiting New Vol) = " << fGlobalNormal << G4endl
276 << G4endl;
277 G4Exception("G4OpBoundaryProcess::PostStepDoIt", "OpBoun02",
278 EventMustBeAborted, // Or JustWarning to see if it happens
279 // repeatedly on one ray
280 ed,
281 "Invalid Surface Normal - Geometry must return valid surface "
282 "normal pointing in the right direction");
283#else
284 fGlobalNormal = -fGlobalNormal;
285#endif
286 }
287
288 G4MaterialPropertyVector* rIndexMPV = nullptr;
290 if(MPT != nullptr)
291 {
292 rIndexMPV = MPT->GetProperty(kRINDEX);
293 }
294 if(rIndexMPV != nullptr)
295 {
296 fRindex1 = rIndexMPV->Value(fPhotonMomentum, idx_rindex1);
297 }
298 else
299 {
300 fStatus = NoRINDEX;
301 if(verboseLevel > 1)
302 BoundaryProcessVerbose();
305 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
306 }
307
308 fReflectivity = 1.;
309 fEfficiency = 0.;
310 fTransmittance = 0.;
311 fSurfaceRoughness = 0.;
312 fModel = glisur;
313 fFinish = polished;
315
316 rIndexMPV = nullptr;
317 fOpticalSurface = nullptr;
318
319 G4LogicalSurface* surface =
320 G4LogicalBorderSurface::GetSurface(thePrePV, thePostPV);
321 if(surface == nullptr)
322 {
323 if(thePostPV->GetMotherLogical() == thePrePV->GetLogicalVolume())
324 {
326 if(surface == nullptr)
327 {
328 surface =
330 }
331 }
332 else
333 {
335 if(surface == nullptr)
336 {
337 surface =
339 }
340 }
341 }
342
343 if(surface != nullptr)
344 {
345 fOpticalSurface =
346 dynamic_cast<G4OpticalSurface*>(surface->GetSurfaceProperty());
347 }
348 if(fOpticalSurface != nullptr)
349 {
350 type = fOpticalSurface->GetType();
351 fModel = fOpticalSurface->GetModel();
352 fFinish = fOpticalSurface->GetFinish();
353
355 fOpticalSurface->GetMaterialPropertiesTable();
356 if(sMPT != nullptr)
357 {
358 if(fFinish == polishedbackpainted || fFinish == groundbackpainted)
359 {
360 rIndexMPV = sMPT->GetProperty(kRINDEX);
361 if(rIndexMPV != nullptr)
362 {
363 fRindex2 = rIndexMPV->Value(fPhotonMomentum, idx_rindex_surface);
364 }
365 else
366 {
367 fStatus = NoRINDEX;
368 if(verboseLevel > 1)
369 BoundaryProcessVerbose();
372 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
373 }
374 }
375
376 fRealRIndexMPV = sMPT->GetProperty(kREALRINDEX);
377 fImagRIndexMPV = sMPT->GetProperty(kIMAGINARYRINDEX);
378 f_iTE = f_iTM = 1;
379
381 if((pp = sMPT->GetProperty(kREFLECTIVITY)))
382 {
383 fReflectivity = pp->Value(fPhotonMomentum, idx_reflect);
384 }
385 else if(fRealRIndexMPV && fImagRIndexMPV)
386 {
387 CalculateReflectivity();
388 }
389
390 if((pp = sMPT->GetProperty(kEFFICIENCY)))
391 {
392 fEfficiency = pp->Value(fPhotonMomentum, idx_eff);
393 }
394 if((pp = sMPT->GetProperty(kTRANSMITTANCE)))
395 {
396 fTransmittance = pp->Value(fPhotonMomentum, idx_trans);
397 }
399 {
400 fSurfaceRoughness = sMPT->GetConstProperty(kSURFACEROUGHNESS);
401 }
402
403 if(fModel == unified)
404 {
405 fProb_sl = (pp = sMPT->GetProperty(kSPECULARLOBECONSTANT))
406 ? pp->Value(fPhotonMomentum, idx_lobe)
407 : 0.;
408 fProb_ss = (pp = sMPT->GetProperty(kSPECULARSPIKECONSTANT))
409 ? pp->Value(fPhotonMomentum, idx_spike)
410 : 0.;
411 fProb_bs = (pp = sMPT->GetProperty(kBACKSCATTERCONSTANT))
412 ? pp->Value(fPhotonMomentum, idx_back)
413 : 0.;
414 }
415 } // end of if(sMPT)
416 else if(fFinish == polishedbackpainted || fFinish == groundbackpainted)
417 {
420 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
421 }
422 } // end of if(fOpticalSurface)
423
424 // DIELECTRIC-DIELECTRIC
425 if(type == dielectric_dielectric)
426 {
427 if(fFinish == polished || fFinish == ground)
428 {
429 if(fMaterial1 == fMaterial2)
430 {
431 fStatus = SameMaterial;
432 if(verboseLevel > 1)
433 BoundaryProcessVerbose();
434 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
435 }
436 MPT = fMaterial2->GetMaterialPropertiesTable();
437 rIndexMPV = nullptr;
438 if(MPT != nullptr)
439 {
440 rIndexMPV = MPT->GetProperty(kRINDEX);
441 }
442 if(rIndexMPV != nullptr)
443 {
444 fRindex2 = rIndexMPV->Value(fPhotonMomentum, idx_rindex2);
445 }
446 else
447 {
448 fStatus = NoRINDEX;
449 if(verboseLevel > 1)
450 BoundaryProcessVerbose();
453 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
454 }
455 }
456 if(fFinish == polishedbackpainted || fFinish == groundbackpainted)
457 {
458 DielectricDielectric();
459 }
460 else
461 {
462 G4double rand = G4UniformRand();
463 if(rand > fReflectivity + fTransmittance)
464 {
465 DoAbsorption();
466 }
467 else if(rand > fReflectivity)
468 {
469 fStatus = Transmission;
470 fNewMomentum = fOldMomentum;
471 fNewPolarization = fOldPolarization;
472 }
473 else
474 {
475 if(fFinish == polishedfrontpainted)
476 {
477 DoReflection();
478 }
479 else if(fFinish == groundfrontpainted)
480 {
481 fStatus = LambertianReflection;
482 DoReflection();
483 }
484 else
485 {
486 DielectricDielectric();
487 }
488 }
489 }
490 }
491 else if(type == dielectric_metal)
492 {
493 DielectricMetal();
494 }
495 else if(type == dielectric_LUT)
496 {
497 DielectricLUT();
498 }
499 else if(type == dielectric_LUTDAVIS)
500 {
501 DielectricLUTDAVIS();
502 }
503 else if(type == dielectric_dichroic)
504 {
505 DielectricDichroic();
506 }
507 else if(type == coated)
508 {
509 CoatedDielectricDielectric();
510 }
511 else
512 {
513 if(fNumBdryTypeWarnings <= 10)
514 {
515 ++fNumBdryTypeWarnings;
516 if(verboseLevel > 0)
517 {
519 ed << " PostStepDoIt(): Illegal boundary type." << G4endl;
520 if(fNumBdryTypeWarnings == 10)
521 {
522 ed << "** Boundary type warnings stopped." << G4endl;
523 }
524 G4Exception("G4OpBoundaryProcess", "OpBoun04", JustWarning, ed);
525 }
526 }
527 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
528 }
529
530 fNewMomentum = fNewMomentum.unit();
531 fNewPolarization = fNewPolarization.unit();
532
533 if(verboseLevel > 1)
534 {
535 G4cout << " New Momentum Direction: " << fNewMomentum << G4endl
536 << " New Polarization: " << fNewPolarization << G4endl;
537 BoundaryProcessVerbose();
538 }
539
541 aParticleChange.ProposePolarization(fNewPolarization);
542
543 if(fStatus == FresnelRefraction || fStatus == Transmission)
544 {
545 // not all surface types check that fMaterial2 has an MPT
547 G4MaterialPropertyVector* groupvel = nullptr;
548 if(aMPT != nullptr)
549 {
550 groupvel = aMPT->GetProperty(kGROUPVEL);
551 }
552 if(groupvel != nullptr)
553 {
555 groupvel->Value(fPhotonMomentum, idx_groupvel));
556 }
557 }
558
559 if(fStatus == Detection && fInvokeSD)
560 InvokeSD(pStep);
561 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
562}
@ JustWarning
@ EventMustBeAborted
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
@ NotAtBoundary
@ Transmission
@ SameMaterial
@ StepTooSmall
@ LambertianReflection
@ FresnelRefraction
@ unified
@ groundfrontpainted
@ polishedbackpainted
@ ground
@ polishedfrontpainted
@ groundbackpainted
@ fGeomBoundary
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
G4bool ConstPropertyExists(const G4String &key) const
G4double GetConstProperty(const G4String &key) const
G4MaterialPropertyVector * GetProperty(const char *key) const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
G4OpticalSurfaceModel GetModel() const
G4OpticalSurfaceFinish GetFinish() const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
static const G4Step * GetHyperStep()
void ProposePolarization(G4double Px, G4double Py, G4double Pz)
void Initialize(const G4Track &) override
void ProposeVelocity(G4double finalVelocity)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double Value(const G4double energy, std::size_t &lastidx) const
G4StepStatus GetStepStatus() const
G4Material * GetMaterial() 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
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

◆ PreparePhysicsTable()

void G4OpBoundaryProcess::PreparePhysicsTable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 132 of file G4OpBoundaryProcess.cc.

133{
134 Initialise();
135}

◆ SetInvokeSD()

void G4OpBoundaryProcess::SetInvokeSD ( G4bool flag)
inlinevirtual

Definition at line 1511 of file G4OpBoundaryProcess.cc.

1512{
1513 fInvokeSD = flag;
1515}

Referenced by Initialise().

◆ SetVerboseLevel()

void G4OpBoundaryProcess::SetVerboseLevel ( G4int verbose)

Definition at line 1518 of file G4OpBoundaryProcess.cc.

Referenced by LBE::ConstructOp(), and Initialise().


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