Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FastSimulationManagerProcess Class Reference

#include <G4FastSimulationManagerProcess.hh>

+ Inheritance diagram for G4FastSimulationManagerProcess:

Public Member Functions

 G4FastSimulationManagerProcess (const G4String &processName="G4FastSimulationManagerProcess", G4ProcessType theType=fParameterisation)
 
 G4FastSimulationManagerProcess (const G4String &processName, const G4String &worldVolumeName, G4ProcessType theType=fParameterisation)
 
 G4FastSimulationManagerProcess (const G4String &processName, G4VPhysicalVolume *worldVolume, G4ProcessType theType=fParameterisation)
 
virtual ~G4FastSimulationManagerProcess ()
 
G4VPhysicalVolumeGetWorldVolume () const
 
void SetWorldVolume (G4String)
 
void SetWorldVolume (G4VPhysicalVolume *)
 
void StartTracking (G4Track *)
 
void EndTracking ()
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &step)
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
void Verbose () const
 
- 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)
 
- 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 77 of file G4FastSimulationManagerProcess.hh.

Constructor & Destructor Documentation

◆ G4FastSimulationManagerProcess() [1/3]

G4FastSimulationManagerProcess::G4FastSimulationManagerProcess ( const G4String processName = "G4FastSimulationManagerProcess",
G4ProcessType  theType = fParameterisation 
)

Definition at line 52 of file G4FastSimulationManagerProcess.cc.

54 :
55 G4VProcess(processName,theType),
56 fWorldVolume ( nullptr ),
57 fIsTrackingTime ( false ),
58 fIsFirstStep ( false ),
59 fGhostNavigator ( nullptr ),
60 fGhostNavigatorIndex ( -1 ),
61 fIsGhostGeometry ( false ),
62 fGhostSafety ( -1.0 ),
63 fFieldTrack ( '0' ),
64 fFastSimulationManager( nullptr ),
65 fFastSimulationTrigger( false )
66{
67 // -- set Process Sub Type
69
70
71 fPathFinder = G4PathFinder::GetInstance();
73
74 SetWorldVolume(fTransportationManager->GetNavigatorForTracking()->GetWorldVolume()->GetName());
75 if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
76 << "' is created, and will message geometry with world volume `"
77 << fWorldVolume->GetName() << "'." << G4endl;
79}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void AddFSMP(G4FastSimulationManagerProcess *)
static G4GlobalFastSimulationManager * GetGlobalFastSimulationManager()
G4VPhysicalVolume * GetWorldVolume() const
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:52
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
const G4String & GetName() const
G4int verboseLevel
Definition: G4VProcess.hh:356
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

◆ G4FastSimulationManagerProcess() [2/3]

G4FastSimulationManagerProcess::G4FastSimulationManagerProcess ( const G4String processName,
const G4String worldVolumeName,
G4ProcessType  theType = fParameterisation 
)

Definition at line 82 of file G4FastSimulationManagerProcess.cc.

85 :
86 G4VProcess(processName,theType),
87 fWorldVolume ( nullptr ),
88 fIsTrackingTime ( false ),
89 fIsFirstStep ( false ),
90 fGhostNavigator ( nullptr ),
91 fGhostNavigatorIndex ( -1 ),
92 fIsGhostGeometry ( false ),
93 fGhostSafety ( -1.0 ),
94 fFieldTrack ( '0' ),
95 fFastSimulationManager( nullptr ),
96 fFastSimulationTrigger( false )
97{
98 // -- set Process Sub Type
100
101
102 fPathFinder = G4PathFinder::GetInstance();
103 fTransportationManager = G4TransportationManager::GetTransportationManager();
104
105 SetWorldVolume(worldVolumeName);
106 if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
107 << "' is created, and will message geometry with world volume `"
108 << fWorldVolume->GetName() << "'." << G4endl;
110}

◆ G4FastSimulationManagerProcess() [3/3]

G4FastSimulationManagerProcess::G4FastSimulationManagerProcess ( const G4String processName,
G4VPhysicalVolume worldVolume,
G4ProcessType  theType = fParameterisation 
)

Definition at line 113 of file G4FastSimulationManagerProcess.cc.

116 :
117 G4VProcess(processName,theType),
118 fWorldVolume ( nullptr ),
119 fIsTrackingTime ( false ),
120 fIsFirstStep ( false ),
121 fGhostNavigator ( nullptr ),
122 fGhostNavigatorIndex ( -1 ),
123 fIsGhostGeometry ( false ),
124 fGhostSafety ( -1.0 ),
125 fFieldTrack ( '0' ),
126 fFastSimulationManager( nullptr ),
127 fFastSimulationTrigger( false )
128{
129 // -- set Process Sub Type
130 SetProcessSubType(static_cast<int>(FASTSIM_ManagerProcess));
131
132
133 fPathFinder = G4PathFinder::GetInstance();
134 fTransportationManager = G4TransportationManager::GetTransportationManager();
135
136 SetWorldVolume(worldVolume);
137 if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
138 << "' is created, and will message geometry with world volume `"
139 << fWorldVolume->GetName() << "'." << G4endl;
141}

◆ ~G4FastSimulationManagerProcess()

G4FastSimulationManagerProcess::~G4FastSimulationManagerProcess ( )
virtual

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4FastSimulationManagerProcess::AlongStepDoIt ( const G4Track track,
const G4Step step 
)
virtual

Implements G4VProcess.

Definition at line 361 of file G4FastSimulationManagerProcess.cc.

364{
365 fDummyParticleChange.Initialize(track);
366 return &fDummyParticleChange;
367}
virtual void Initialize(const G4Track &)

◆ AlongStepGetPhysicalInteractionLength()

G4double G4FastSimulationManagerProcess::AlongStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double proposedSafety,
G4GPILSelection selection 
)
virtual

Implements G4VProcess.

Definition at line 292 of file G4FastSimulationManagerProcess.cc.

298{
299
300 *selection = NotCandidateForSelection;
301 G4double returnedStep = DBL_MAX;
302
303 // ---------------------------------------------------
304 // -- Below code valid for ghost geometry, otherwise
305 // -- useless for fast simulation attached to mass
306 // -- geometry. Warn user in case along used for
307 // -- mass geometry ?
308 // --------------------------------------------------
309 if ( fIsGhostGeometry )
310 {
311 static G4ThreadLocal G4FieldTrack *endTrack_G4MT_TLS_ = 0 ;
312 if (!endTrack_G4MT_TLS_) endTrack_G4MT_TLS_ = new G4FieldTrack ('0') ;
313 G4FieldTrack &endTrack = *endTrack_G4MT_TLS_;
314
315 static G4ThreadLocal ELimited *eLimited_G4MT_TLS_ = 0 ;
316 if (!eLimited_G4MT_TLS_) eLimited_G4MT_TLS_ = new ELimited ;
317 ELimited &eLimited = *eLimited_G4MT_TLS_;
318
319 if (previousStepSize > 0.) fGhostSafety -= previousStepSize;
320 if (fGhostSafety < 0.) fGhostSafety = 0.0;
321
322 // ------------------------------------------
323 // Determination of the proposed step length:
324 // ------------------------------------------
325 if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
326 {
327 // -- No chance to limit the step, as proposed move inside safety
328 returnedStep = currentMinimumStep;
329 proposedSafety = fGhostSafety - currentMinimumStep;
330 }
331 else
332 {
333 // -- Proposed move exceeds safety, need to state
334 G4FieldTrackUpdator::Update(&fFieldTrack, &track);
335 returnedStep = fPathFinder->ComputeStep(fFieldTrack,
336 currentMinimumStep,
337 fGhostNavigatorIndex,
338 track.GetCurrentStepNumber(),
339 fGhostSafety,
340 eLimited,
341 endTrack,
342 track.GetVolume());
343
344 if(eLimited == kDoNot) fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition()); // -- step no limited by ghost
345 proposedSafety = fGhostSafety;
346 if (eLimited == kUnique || eLimited == kSharedOther) *selection = CandidateForSelection;
347 else if (eLimited == kSharedTransport) returnedStep *= (1.0 + 1.0e-9); // -- Expand to disable its selection in Step Manager comparison
348 }
349 }
350
351
352 // ----------------------------------------------
353 // Returns the fGhostSafety as the proposedSafety
354 // The SteppingManager will take care of keeping
355 // the smallest one.
356 // ----------------------------------------------
357 return returnedStep;
358}
@ CandidateForSelection
@ NotCandidateForSelection
ELimited
@ kDoNot
@ kUnique
@ kSharedOther
@ kSharedTransport
double G4double
Definition: G4Types.hh:83
static void Update(G4FieldTrack *, const G4Track *)
G4ThreeVector GetPosition() const
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
G4VPhysicalVolume * GetVolume() const
G4int GetCurrentStepNumber() const
#define DBL_MAX
Definition: templates.hh:62
#define G4ThreadLocal
Definition: tls.hh:77

◆ AtRestDoIt()

G4VParticleChange * G4FastSimulationManagerProcess::AtRestDoIt ( const G4Track ,
const G4Step  
)
virtual

Implements G4VProcess.

Definition at line 404 of file G4FastSimulationManagerProcess.cc.

405{
406 return fFastSimulationManager->InvokeAtRestDoIt();
407}
G4VParticleChange * InvokeAtRestDoIt()

◆ AtRestGetPhysicalInteractionLength()

G4double G4FastSimulationManagerProcess::AtRestGetPhysicalInteractionLength ( const G4Track track,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 376 of file G4FastSimulationManagerProcess.cc.

379{
380 const G4VPhysicalVolume* currentVolume(0);
381 if ( fIsGhostGeometry ) currentVolume = fPathFinder->GetLocatedVolume(fGhostNavigatorIndex);
382 else currentVolume = track.GetVolume();
383 fFastSimulationManager = currentVolume->GetLogicalVolume()->GetFastSimulationManager();
384 if( fFastSimulationManager )
385 {
386 // Ask for trigger:
387 fFastSimulationTrigger = fFastSimulationManager->AtRestGetFastSimulationManagerTrigger(track, fGhostNavigator);
388 if( fFastSimulationTrigger )
389 {
390 // Dirty trick to take control over stepping. Does anyone will ever use that ?
392 return -1.0;
393 }
394 }
395
396 // -- no fast simulation occuring there:
398 return DBL_MAX;
399}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
G4bool AtRestGetFastSimulationManagerTrigger(const G4Track &, const G4Navigator *a=0)
G4VPhysicalVolume * GetLocatedVolume(G4int navId) const

◆ EndTracking()

void G4FastSimulationManagerProcess::EndTracking ( )
virtual

Reimplemented from G4VProcess.

Definition at line 227 of file G4FastSimulationManagerProcess.cc.

229{
230 fIsTrackingTime = false;
231 if ( fIsGhostGeometry ) fTransportationManager->DeActivateNavigator(fGhostNavigator);
232}
void DeActivateNavigator(G4Navigator *aNavigator)

◆ GetWorldVolume()

G4VPhysicalVolume * G4FastSimulationManagerProcess::GetWorldVolume ( ) const
inline

Definition at line 103 of file G4FastSimulationManagerProcess.hh.

103{return fWorldVolume;}

Referenced by G4FastSimulationPhysics::ConstructProcess().

◆ PostStepDoIt()

G4VParticleChange * G4FastSimulationManagerProcess::PostStepDoIt ( const G4Track ,
const G4Step  
)
virtual

Implements G4VProcess.

Definition at line 278 of file G4FastSimulationManagerProcess.cc.

281{
282 G4VParticleChange* finalState = fFastSimulationManager->InvokePostStepDoIt();
283
284 // If the particle is still alive, suspend it to force physics re-initialisation:
285 if (finalState->GetTrackStatus() != fStopAndKill) finalState->ProposeTrackStatus(fSuspend);
286
287 return finalState;
288}
@ fSuspend
@ fStopAndKill
G4VParticleChange * InvokePostStepDoIt()
void ProposeTrackStatus(G4TrackStatus status)
G4TrackStatus GetTrackStatus() const

◆ PostStepGetPhysicalInteractionLength()

G4double G4FastSimulationManagerProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 239 of file G4FastSimulationManagerProcess.cc.

243{
244 // -- Get current volume, and check for presence of fast simulation manager.
245 // -- For the case of the navigator for tracking (fGhostNavigatorIndex == 0)
246 // -- we use the track volume. This allows the code to be valid for both
247 // -- cases where the PathFinder is used (G4CoupledTranportation) or not
248 // -- (G4Transportation).
249 const G4VPhysicalVolume* currentVolume(0);
250 if ( fIsGhostGeometry ) currentVolume = fPathFinder->GetLocatedVolume(fGhostNavigatorIndex);
251 else currentVolume = track.GetVolume();
252
253 if ( currentVolume )
254 {
255 fFastSimulationManager = currentVolume->GetLogicalVolume()->GetFastSimulationManager();
256 if( fFastSimulationManager )
257 {
258 // Ask for trigger:
259 fFastSimulationTrigger = fFastSimulationManager->PostStepGetFastSimulationManagerTrigger(track, fGhostNavigator);
260 if( fFastSimulationTrigger )
261 {
262 // Take control over stepping:
264 return 0.0;
265 }
266 }
267 }
268
269 // -- no fast simulation occuring there:
271 return DBL_MAX;
272}
@ ExclusivelyForced
G4bool PostStepGetFastSimulationManagerTrigger(const G4Track &, const G4Navigator *a=0)

◆ SetWorldVolume() [1/2]

void G4FastSimulationManagerProcess::SetWorldVolume ( G4String  newWorldName)

Definition at line 153 of file G4FastSimulationManagerProcess.cc.

154{
155 if (fIsTrackingTime)
156 {
158 ed << "G4FastSimulationManagerProcess `" << GetProcessName()
159 << "': changing of world volume at tracking time is not allowed." << G4endl;
160 G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4String)",
161 "FastSim002",
162 JustWarning, ed,
163 "Call ignored.");
164 }
165 else
166 {
167 G4VPhysicalVolume* newWorld = fTransportationManager->IsWorldExisting(newWorldName);
168 if (newWorld == 0)
169 {
170 G4ExceptionDescription tellWhatIsWrong;
171 tellWhatIsWrong << "Volume newWorldName = `" << newWorldName
172 << "' is not a parallel world nor the mass world volume."
173 << G4endl;
174 G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4String)",
175 "FastSim003",
177 tellWhatIsWrong);
178 }
179 if (verboseLevel>0)
180 {
181 if (fWorldVolume) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
182 << "': changing world volume from '" << fWorldVolume->GetName()
183 << "' to `" << newWorld << "'." << G4endl;
184 else G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
185 << "': setting world volume from to `"<< newWorld->GetName() << "'." << G4endl;
186 }
187 fWorldVolume = newWorld;
188 }
189}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4VPhysicalVolume * IsWorldExisting(const G4String &worldName)

Referenced by G4FastSimulationManagerProcess(), and SetWorldVolume().

◆ SetWorldVolume() [2/2]

void G4FastSimulationManagerProcess::SetWorldVolume ( G4VPhysicalVolume newWorld)

Definition at line 192 of file G4FastSimulationManagerProcess.cc.

193{
194 if (newWorld) SetWorldVolume(newWorld->GetName());
195 else
196 {
197 G4ExceptionDescription tellWhatIsWrong;
198 tellWhatIsWrong << "Null pointer passed for world volume." << G4endl;
199 G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4VPhysicalVolume* newWorld)",
200 "FastSim004",
202 tellWhatIsWrong);
203 }
204}

◆ StartTracking()

void G4FastSimulationManagerProcess::StartTracking ( G4Track track)
virtual

Reimplemented from G4VProcess.

Definition at line 210 of file G4FastSimulationManagerProcess.cc.

211{
212 fIsTrackingTime = true;
213 fIsFirstStep = true;
214
215 // -- fetch the navigator (and its index) and activate it:
217 fGhostNavigator = transportationManager->GetNavigator(fWorldVolume);
218 fIsGhostGeometry = (fGhostNavigator != transportationManager->GetNavigatorForTracking());
219 if (fIsGhostGeometry) fGhostNavigatorIndex = transportationManager->ActivateNavigator(fGhostNavigator);
220 else fGhostNavigatorIndex = -1;
221
222 fPathFinder->PrepareNewTrack(track->GetPosition(), track->GetMomentumDirection());
223}
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=nullptr)
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
G4int ActivateNavigator(G4Navigator *aNavigator)
G4Navigator * GetNavigator(const G4String &worldName)

◆ Verbose()

void G4FastSimulationManagerProcess::Verbose ( ) const

Definition at line 410 of file G4FastSimulationManagerProcess.cc.

411{
412 /* G4cout << " >>>>> Trigger Status : ";
413 switch(fFastSimulationManager->GetTriggerStatus())
414 {
415 case NoModel:
416 G4cout << "NoModel" << G4endl;
417 break;
418 case OnBoundaryButLeaving:
419 G4cout << "OnBoundaryButLeaving" << G4endl;
420 break;
421 case OneModelTrigger:
422 G4cout << "OneModelTrigger" << G4endl;
423 break;
424 case NoModelTrigger:
425 G4cout << "NoModelTrigger" << G4endl;
426 break;
427 case Undefined:
428 G4cout << "Undefined" << G4endl;
429 break;
430 default:
431 G4cout << " Bizarre..." << G4endl;
432 break;
433 }*/
434}

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