Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FastSimulationManagerProcess.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// $Id$
28//
29//
30//---------------------------------------------------------------
31//
32// G4FastSimulationProcess.cc
33//
34// Description:
35// The process that triggers the parameterised simulations,
36// if any.
37//
38// History:
39// August 97: First implementation. Verderi && MoraDeFreitas.
40// October 06: move to parallel geometry scheme, M. Verderi
41//---------------------------------------------------------------
42
43#include "G4ios.hh"
47#include "G4PathFinder.hh"
48#include "G4ParticleChange.hh"
50
51#define PARANOIA
52
55 G4ProcessType theType) :
56 G4VProcess(processName,theType),
57 fWorldVolume(0),
58 fIsTrackingTime(false),
59 fGhostNavigator(0),
60 fGhostNavigatorIndex(-1),
61 fIsGhostGeometry(false),
62 fGhostSafety(-1.0),
63 fFieldTrack('0'),
64 fFastSimulationManager(0),
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}
80
81
84 const G4String& worldVolumeName,
85 G4ProcessType theType) :
86 G4VProcess(processName,theType),
87 fWorldVolume(0),
88 fIsTrackingTime(false),
89 fGhostNavigator(0),
90 fGhostNavigatorIndex(-1),
91 fIsGhostGeometry(false),
92 fGhostSafety(-1.0),
93 fFieldTrack('0'),
94 fFastSimulationManager(0),
95 fFastSimulationTrigger(false)
96{
97 // -- set Process Sub Type
99
100
101 fPathFinder = G4PathFinder::GetInstance();
102 fTransportationManager = G4TransportationManager::GetTransportationManager();
103
104 SetWorldVolume(worldVolumeName);
105 if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
106 << "' is created, and will message geometry with world volume `"
107 << fWorldVolume->GetName() << "'." << G4endl;
109}
110
111
114 G4VPhysicalVolume* worldVolume,
115 G4ProcessType theType) :
116 G4VProcess(processName,theType),
117 fWorldVolume(0),
118 fIsTrackingTime(false),
119 fGhostNavigator(0),
120 fGhostNavigatorIndex(-1),
121 fIsGhostGeometry(false),
122 fGhostSafety(-1.0),
123 fFieldTrack('0'),
124 fFastSimulationManager(0),
125 fFastSimulationTrigger(false)
126{
127 // -- set Process Sub Type
128 SetProcessSubType(static_cast<int>(FASTSIM_ManagerProcess));
129
130
131 fPathFinder = G4PathFinder::GetInstance();
132 fTransportationManager = G4TransportationManager::GetTransportationManager();
133
134 SetWorldVolume(worldVolume);
135 if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
136 << "' is created, and will message geometry with world volume `"
137 << fWorldVolume->GetName() << "'." << G4endl;
139}
140
141
143{
145}
146
147
148// -----------------------
149// User access methods:
150// -----------------------
152{
153 if (fIsTrackingTime)
154 {
156 ed << "G4FastSimulationManagerProcess `" << GetProcessName()
157 << "': changing of world volume at tracking time is not allowed." << G4endl;
158 G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4String)",
159 "FastSim002",
160 JustWarning, ed,
161 "Call ignored.");
162 }
163 else
164 {
165 G4VPhysicalVolume* newWorld = fTransportationManager->IsWorldExisting(newWorldName);
166 if (newWorld == 0)
167 {
168 G4ExceptionDescription tellWhatIsWrong;
169 tellWhatIsWrong << "Volume newWorldName = `" << newWorldName
170 << "' is not a parallel world nor the mass world volume."
171 << G4endl;
172 G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4String)",
173 "FastSim003",
175 tellWhatIsWrong);
176 }
177 if (verboseLevel>0)
178 {
179 if (fWorldVolume) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
180 << "': changing world volume from '" << fWorldVolume->GetName()
181 << "' to `" << newWorld << "'." << G4endl;
182 else G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
183 << "': setting world volume from to `"<< newWorld->GetName() << "'." << G4endl;
184 }
185 fWorldVolume = newWorld;
186 }
187}
188
189
191{
192 if (newWorld) SetWorldVolume(newWorld->GetName());
193 else
194 {
195 G4ExceptionDescription tellWhatIsWrong;
196 tellWhatIsWrong << "Null pointer passed for world volume." << G4endl;
197 G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4VPhysicalVolume* newWorld)",
198 "FastSim004",
200 tellWhatIsWrong);
201 }
202}
203
204
205// --------------------
206// Start/End tracking:
207// --------------------
208void
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}
224
225
226void
229{
230 fIsTrackingTime = false;
231 if ( fIsGhostGeometry ) fTransportationManager->DeActivateNavigator(fGhostNavigator);
232}
233
234
235// ------------------------------------------
236// PostStepGetPhysicalInteractionLength():
237// ------------------------------------------
241 G4double,
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}
273
274//------------------------------------
275// PostStepDoIt()
276//------------------------------------
279PostStepDoIt(const G4Track&,
280 const G4Step&)
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}
289
290
294 G4double previousStepSize,
295 G4double currentMinimumStep,
296 G4double& proposedSafety,
297 G4GPILSelection* selection)
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 G4FieldTrack endTrack('0');
312 static ELimited eLimited;
313
314 if (previousStepSize > 0.) fGhostSafety -= previousStepSize;
315 if (fGhostSafety < 0.) fGhostSafety = 0.0;
316
317 // ------------------------------------------
318 // Determination of the proposed step length:
319 // ------------------------------------------
320 if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
321 {
322 // -- No chance to limit the step, as proposed move inside safety
323 returnedStep = currentMinimumStep;
324 proposedSafety = fGhostSafety - currentMinimumStep;
325 }
326 else
327 {
328 // -- Proposed move exceeds safety, need to state
329 G4FieldTrackUpdator::Update(&fFieldTrack, &track);
330 returnedStep = fPathFinder->ComputeStep(fFieldTrack,
331 currentMinimumStep,
332 fGhostNavigatorIndex,
333 track.GetCurrentStepNumber(),
334 fGhostSafety,
335 eLimited,
336 endTrack,
337 track.GetVolume());
338
339 if(eLimited == kDoNot) fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition()); // -- step no limited by ghost
340 proposedSafety = fGhostSafety;
341 if (eLimited == kUnique || eLimited == kSharedOther) *selection = CandidateForSelection;
342 else if (eLimited == kSharedTransport) returnedStep *= (1.0 + 1.0e-9); // -- Expand to disable its selection in Step Manager comparison
343 }
344 }
345
346
347 // ----------------------------------------------
348 // Returns the fGhostSafety as the proposedSafety
349 // The SteppingManager will take care of keeping
350 // the smallest one.
351 // ----------------------------------------------
352 return returnedStep;
353}
354
357AlongStepDoIt(const G4Track& track,
358 const G4Step&)
359{
360 fDummyParticleChange.Initialize(track);
361 return &fDummyParticleChange;
362}
363
364
365//--------------------------------------------
366// At Rest parameterisation:
367//--------------------------------------------
368// AtRestGetPhysiscalInteractionLength:
369//--------------------------------------------
374{
375 const G4VPhysicalVolume* currentVolume(0);
376 if ( fIsGhostGeometry ) currentVolume = fPathFinder->GetLocatedVolume(fGhostNavigatorIndex);
377 else currentVolume = track.GetVolume();
378 fFastSimulationManager = currentVolume->GetLogicalVolume()->GetFastSimulationManager();
379 if( fFastSimulationManager )
380 {
381 // Ask for trigger:
382 fFastSimulationTrigger = fFastSimulationManager->AtRestGetFastSimulationManagerTrigger(track, fGhostNavigator);
383 if( fFastSimulationTrigger )
384 {
385 // Dirty trick to take control over stepping. Does anyone will ever use that ?
387 return -1.0;
388 }
389 }
390
391 // -- no fast simulation occuring there:
393 return DBL_MAX;
394}
395
396//-----------------------------------------------
397// AtRestDoIt:
398//-----------------------------------------------
400{
401 return fFastSimulationManager->InvokeAtRestDoIt();
402}
403
404
406{
407 /* G4cout << " >>>>> Trigger Status : ";
408 switch(fFastSimulationManager->GetTriggerStatus())
409 {
410 case NoModel:
411 G4cout << "NoModel" << G4endl;
412 break;
413 case OnBoundaryButLeaving:
414 G4cout << "OnBoundaryButLeaving" << G4endl;
415 break;
416 case OneModelTrigger:
417 G4cout << "OneModelTrigger" << G4endl;
418 break;
419 case NoModelTrigger:
420 G4cout << "NoModelTrigger" << G4endl;
421 break;
422 case Undefined:
423 G4cout << "Undefined" << G4endl;
424 break;
425 default:
426 G4cout << " Bizarre..." << G4endl;
427 break;
428 }*/
429}
430
431
G4double condition(const G4ErrorSymMatrix &m)
@ JustWarning
@ FatalException
G4ForceCondition
@ NotForced
@ ExclusivelyForced
G4GPILSelection
@ CandidateForSelection
@ NotCandidateForSelection
ELimited
@ kDoNot
@ kUnique
@ kSharedOther
@ kSharedTransport
G4ProcessType
@ fSuspend
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
G4FastSimulationManagerProcess(const G4String &processName="G4FastSimulationManagerProcess", G4ProcessType theType=fParameterisation)
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &step)
G4VParticleChange * InvokePostStepDoIt()
G4VParticleChange * InvokeAtRestDoIt()
G4bool PostStepGetFastSimulationManagerTrigger(const G4Track &, const G4Navigator *a=0)
G4bool AtRestGetFastSimulationManagerTrigger(const G4Track &, const G4Navigator *a=0)
static void Update(G4FieldTrack *, const G4Track *)
G4ThreeVector GetPosition() const
void RemoveFSMP(G4FastSimulationManagerProcess *)
void AddFSMP(G4FastSimulationManagerProcess *)
static G4GlobalFastSimulationManager * GetGlobalFastSimulationManager()
G4FastSimulationManager * GetFastSimulationManager() const
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
G4VPhysicalVolume * GetWorldVolume() const
G4VPhysicalVolume * GetLocatedVolume(G4int navId) const
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=0)
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:57
Definition: G4Step.hh:78
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
G4int GetCurrentStepNumber() const
const G4ThreeVector & GetMomentumDirection() const
static G4TransportationManager * GetTransportationManager()
G4VPhysicalVolume * IsWorldExisting(const G4String &worldName)
G4Navigator * GetNavigatorForTracking() const
G4int ActivateNavigator(G4Navigator *aNavigator)
void DeActivateNavigator(G4Navigator *aNavigator)
G4Navigator * GetNavigator(const G4String &worldName)
void ProposeTrackStatus(G4TrackStatus status)
virtual void Initialize(const G4Track &)
G4TrackStatus GetTrackStatus() const
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
G4int verboseLevel
Definition: G4VProcess.hh:368
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
#define DBL_MAX
Definition: templates.hh:83