Geant4 11.4.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 ()
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 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 95 of file G4OpBoundaryProcess.cc.

97 : G4VDiscreteProcess(processName, ptype)
98{
99 Initialise();
100
101 if(verboseLevel > 0)
102 {
103 G4cout << GetProcessName() << " is created " << G4endl;
104 }
106
107 fStatus = Undefined;
108 fModel = glisur;
109 fFinish = polished;
110 fReflectivity = 1.;
111 fEfficiency = 0.;
112 fTransmittance = 0.;
113 fSurfaceRoughness = 0.;
114 fProb_sl = 0.;
115 fProb_ss = 0.;
116 fProb_bs = 0.;
117
118 fRealRIndexMPV = nullptr;
119 fImagRIndexMPV = nullptr;
120 fMaterial1 = nullptr;
121 fMaterial2 = nullptr;
122 fOpticalSurface = nullptr;
124
125 f_iTE = f_iTM = 0;
126 fPhotonMomentum = 0.;
127 fRindex1 = fRindex2 = 1.;
128 fSint1 = 0.;
129 fDichroicVector = nullptr;
130}
@ fOpBoundary
@ polished
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
G4VDiscreteProcess(const G4String &aName, G4ProcessType aType=fNotDefined)
G4int verboseLevel
void SetProcessSubType(G4int)
const G4String & GetProcessName() const

Referenced by ~G4OpBoundaryProcess().

◆ ~G4OpBoundaryProcess()

G4OpBoundaryProcess::~G4OpBoundaryProcess ( )
virtualdefault

Member Function Documentation

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 1344 of file G4OpBoundaryProcess.cc.

1346{
1347 *condition = Forced;
1348 return DBL_MAX;
1349}
G4double condition(const G4ErrorSymMatrix &m)
@ Forced
#define DBL_MAX
Definition templates.hh:62

◆ GetStatus()

G4OpBoundaryProcessStatus G4OpBoundaryProcess::GetStatus ( ) const
inlinevirtual

Definition at line 287 of file G4OpBoundaryProcess.hh.

288{
289 return fStatus;
290}

◆ Initialise()

void G4OpBoundaryProcess::Initialise ( )
virtual

Definition at line 142 of file G4OpBoundaryProcess.cc.

143{
144 G4OpticalParameters* params = G4OpticalParameters::Instance();
147}
virtual void SetInvokeSD(G4bool)
G4bool GetBoundaryInvokeSD() const
G4int GetBoundaryVerboseLevel() const
static G4OpticalParameters * Instance()

Referenced by G4OpBoundaryProcess(), and PreparePhysicsTable().

◆ IsApplicable()

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

Reimplemented from G4VProcess.

Definition at line 281 of file G4OpBoundaryProcess.hh.

283{
284 return (&aParticleType == G4OpticalPhoton::OpticalPhoton());
285}
static G4OpticalPhoton * OpticalPhoton()

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 150 of file G4OpBoundaryProcess.cc.

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

137{
138 Initialise();
139}

◆ SetInvokeSD()

void G4OpBoundaryProcess::SetInvokeSD ( G4bool flag)
inlinevirtual

Definition at line 1470 of file G4OpBoundaryProcess.cc.

1471{
1472 fInvokeSD = flag;
1474}

Referenced by Initialise().

◆ SetVerboseLevel()

void G4OpBoundaryProcess::SetVerboseLevel ( G4int verbose)

Definition at line 1477 of file G4OpBoundaryProcess.cc.

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


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