Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ParallelGeometriesLimiterProcess.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//
28//
29
30#include "G4ios.hh"
33#include "G4ProcessManager.hh"
35#include "G4PathFinder.hh"
37
38#include "G4SystemOfUnits.hh"
39
41 G4VProcess(processName, fParallel),
42 fParallelWorldSafety( 0.0 ),
43 fIsTrackingTime ( false ),
44 fFieldTrack ( '0' )
45{
46 // -- Process Sub Type ? §§
47
48 fPathFinder = G4PathFinder::GetInstance();
50}
51
52
53// ----------------------------
54// -- Add/Remove world volumes:
55// ----------------------------
57{
58
59 // -- Refuse adding parallel geometry during tracking time:
60 if (fIsTrackingTime)
61 {
63 ed << "G4ParallelGeometriesLimiterProcess `" << GetProcessName()
64 << "': adding a parallel world volume at tracking time is not allowed." << G4endl;
65 G4Exception("G4ParallelGeometriesLimiterProcess::AddParallelWorld(const G4String& parallelWorldName)",
66 "BIAS.GEN.21",
67 JustWarning, ed,
68 "Call ignored.");
69 return;
70 }
71
72 else
73
74 {
75 G4VPhysicalVolume* newWorld = fTransportationManager->IsWorldExisting( parallelWorldName );
76
77 // -- Fatal exception if requested world does not exist:
78 if (newWorld == 0)
79 {
80 G4ExceptionDescription tellWhatIsWrong;
81 tellWhatIsWrong << "Volume `" << parallelWorldName
82 << "' is not a parallel world nor the mass world volume."
83 << G4endl;
84 G4Exception("G4ParallelGeometriesLimiterProcess::SetWorldVolume(const G4String)",
85 "BIAS.GEN.22",
87 tellWhatIsWrong);
88 }
89
90 // -- Protection against adding the mass geometry world as parallel world:
91 if ( newWorld == fTransportationManager->GetNavigatorForTracking()->GetWorldVolume() )
92 {
94 ed << "G4ParallelGeometriesLimiterProcess `" << GetProcessName()
95 << "': trying to add the world volume for tracking as a parallel world." << G4endl;
96 G4Exception("G4ParallelGeometriesLimiterProcess::AddParallelWorld(const G4String& parallelWorldName)",
97 "BIAS.GEN.23",
98 JustWarning, ed,
99 "Call ignored.");
100 return;
101 }
102
103 // -- Add parallel world, taking care it is not in the list yet:
104 G4bool isNew = true;
105 for ( auto knownWorld : fParallelWorlds )
106 {
107 if ( knownWorld == newWorld ) isNew = false;
108 }
109 if ( isNew ) fParallelWorlds.push_back( newWorld );
110 else
111 {
113 ed << "G4ParallelGeometriesLimiterProcess `" << GetProcessName()
114 << "': trying to re-add the parallel world volume `" << parallelWorldName << "'." << G4endl;
115 G4Exception("G4ParallelGeometriesLimiterProcess::AddParallelWorld(const G4String& parallelWorldName)",
116 "BIAS.GEN.24",
117 JustWarning, ed,
118 "Call ignored.");
119 return;
120 }
121 }
122
123}
124
125
127{
128
129 // -- Refuse refuse removing parallel geometry during tracking time:
130 if (fIsTrackingTime)
131 {
133 ed << "G4ParallelGeometriesLimiterProcess `" << GetProcessName()
134 << "': removing a parallel world volume at tracking time is not allowed." << G4endl;
135 G4Exception("G4ParallelGeometriesLimiterProcess::RemoveParallelWorld(const G4String& parallelWorldName)",
136 "BIAS.GEN.25",
137 JustWarning, ed,
138 "Call ignored.");
139 return;
140 }
141
142 else
143
144 {
145 G4VPhysicalVolume* newWorld = fTransportationManager->IsWorldExisting( parallelWorldName );
146
147 if (newWorld == 0)
148 {
149
151 ed << "G4ParallelGeometriesLimiterProcess `" << GetProcessName()
152 << "': trying to remove an inexisting parallel world '" << parallelWorldName << "'." << G4endl;
153 G4Exception("G4ParallelGeometriesLimiterProcess::RemoveParallelWorld(const G4String& parallelWorldName)",
154 "BIAS.GEN.26",
155 JustWarning, ed,
156 "Call ignored.");
157 return;
158 }
159
160 // -- get position of world volume in list:
161 size_t iWorld = 0;
162 for ( auto knownWorld : fParallelWorlds )
163 {
164 if ( knownWorld == newWorld ) break;
165 iWorld++;
166 }
167
168 if ( iWorld == fParallelWorlds.size() )
169 {
171 ed << "G4ParallelGeometriesLimiterProcess `" << GetProcessName()
172 << "': trying to remove an non-registerered parallel world '" << parallelWorldName << "'." << G4endl;
173 G4Exception("G4ParallelGeometriesLimiterProcess::RemoveParallelWorld(const G4String& parallelWorldName)",
174 "BIAS.GEN.27",
175 JustWarning, ed,
176 "Call ignored.");
177 return;
178 }
179
180 // -- remove from vector:
181 fParallelWorlds.erase( fParallelWorlds.begin() + iWorld );
182
183 }
184
185
186
187}
188
189
190// --------------------
191// Start/End tracking:
192// --------------------
194{
195 fIsTrackingTime = true;
196
197 // -- fetch the navigators, their indeces, and activate:
198 fParallelWorldNavigators .clear();
199 fParallelWorldNavigatorIndeces.clear();
200 fParallelWorldSafeties .clear();
201 fParallelWorldIsLimiting .clear();
202 fParallelWorldWasLimiting .clear();
203 fCurrentVolumes .clear();
204 fPreviousVolumes .clear();
205 for ( auto parallelWorld : fParallelWorlds )
206 {
207 fParallelWorldNavigators .push_back( fTransportationManager-> GetNavigator( parallelWorld ) );
208 fParallelWorldNavigatorIndeces.push_back( fTransportationManager->ActivateNavigator( fParallelWorldNavigators.back() ) );
209 fParallelWorldSafeties .push_back( 0.0 );
210 fParallelWorldIsLimiting .push_back( false );
211 fParallelWorldWasLimiting .push_back( false );
212 }
213
214 fPathFinder->PrepareNewTrack( track->GetPosition(), track->GetMomentumDirection() );
215 // -- §§ does it work at this level, after "PrepareNewTrack" above ?
216 for ( auto navigatorIndex : fParallelWorldNavigatorIndeces )
217 {
218 fPreviousVolumes.push_back( nullptr );
219 fCurrentVolumes .push_back( fPathFinder->GetLocatedVolume( navigatorIndex ) );
220 }
221
222 // -- will force updating safety:
223 fParallelWorldSafety = 0.0;
224 for ( size_t i = 0 ; i < fParallelWorldNavigatorIndeces.size() ; i++ ) fParallelWorldSafeties[i] = 0.0;
225}
226
227
229{
230 fIsTrackingTime = false;
231 for ( auto parallelWorldNavigator : fParallelWorldNavigators )
232 fTransportationManager->DeActivateNavigator( parallelWorldNavigator );
233}
234
235
237{
238
239 // -- push previous step limitation flags and volumes:
240 // -- §§ consider switching pointers insteads of making copies of std::vector's:
241 fParallelWorldWasLimiting = fParallelWorldIsLimiting;
242 fPreviousVolumes = fCurrentVolumes;
243
244 // -- update volumes:
245 size_t i = 0;
246 for ( auto navigatorIndex : fParallelWorldNavigatorIndeces ) fCurrentVolumes[i++] = fPathFinder->GetLocatedVolume( navigatorIndex );
247
249 return DBL_MAX;
250}
251
252
254 G4double previousStepSize,
255 G4double currentMinimumStep,
256 G4double& proposedSafety,
257 G4GPILSelection* selection)
258{
259
260 // -- Init:
261 // -- Note that the returnedStep must be physically meaningful, even if we return NotCandidateForSelection as condition;
262 // -- the reason is that the stepping manager always takes the smallest alongstep among the returned ones (point related
263 // -- to geometry step length wrt to true path length).
264 *selection = NotCandidateForSelection;
265 G4double returnedStep = DBL_MAX;
266
267 // -- G4FieldTrack and ELimited:
268 static G4ThreadLocal G4FieldTrack *endTrack_G4MT_TLS_ = 0 ;
269 if (!endTrack_G4MT_TLS_) endTrack_G4MT_TLS_ = new G4FieldTrack ('0') ;
270 G4FieldTrack &endTrack = *endTrack_G4MT_TLS_;
271
272 static G4ThreadLocal ELimited *eLimited_G4MT_TLS_ = 0 ;
273 if (!eLimited_G4MT_TLS_) eLimited_G4MT_TLS_ = new ELimited ;
274 ELimited &eLimited = *eLimited_G4MT_TLS_;
275
276
277 // -------------------
278 // -- Update safeties:
279 // -------------------
280 if ( previousStepSize > 0.0 )
281 {
282 for ( auto& parallelWorldSafety : fParallelWorldSafeties )
283 {
284 parallelWorldSafety -= previousStepSize;
285 if ( parallelWorldSafety < 0. ) parallelWorldSafety = 0.0;
286 fParallelWorldSafety = parallelWorldSafety < fParallelWorldSafety ? parallelWorldSafety : fParallelWorldSafety ;
287 }
288 }
289
290
291 // ------------------------------------------
292 // Determination of the proposed step length:
293 // ------------------------------------------
294 if ( ( currentMinimumStep <= fParallelWorldSafety ) && ( currentMinimumStep > 0. ) )
295 {
296 // -- No chance to limit the step, as proposed move inside safety
297
298 returnedStep = currentMinimumStep;
299 proposedSafety = fParallelWorldSafety - currentMinimumStep;
300 }
301 else
302 {
303 // -- Proposed move exceeds common safety, need to state
304 G4double smallestReturnedStep = -1.0;
305 ELimited eLimitedForSmallestStep = kDoNot;
306 for ( size_t i = 0 ; i < fParallelWorldNavigatorIndeces.size() ; i++ )
307 {
308 // -- Update safety of geometries having safety smaller than current minimum step
309 if ( currentMinimumStep >= fParallelWorldSafeties[i] )
310 {
311 G4FieldTrackUpdator::Update(&fFieldTrack, &track);
312 G4double tmpReturnedStep = fPathFinder->ComputeStep(fFieldTrack,
313 currentMinimumStep,
314 fParallelWorldNavigatorIndeces[i],
315 track.GetCurrentStepNumber(),
316 fParallelWorldSafeties[i],
317 eLimited,
318 endTrack,
319 track.GetVolume());
320
321 if ( ( smallestReturnedStep < 0.0 ) || ( tmpReturnedStep <= smallestReturnedStep ) )
322 {
323 smallestReturnedStep = tmpReturnedStep;
324 eLimitedForSmallestStep = eLimited;
325 }
326
327 if (eLimited == kDoNot)
328 {
329 // -- Step not limited by this geometry
330 fParallelWorldSafeties[i] = fParallelWorldNavigators[i]->ComputeSafety(endTrack.GetPosition());
331 fParallelWorldIsLimiting[i] = false;
332 }
333 else
334 {
335 fParallelWorldIsLimiting[i] = true;
336 }
337 }
338
339 // -- update with smallest safety:
340 fParallelWorldSafety = fParallelWorldSafeties[i] < fParallelWorldSafety ? fParallelWorldSafeties[i] : fParallelWorldSafety ;
341 }
342
343 // -- no geometry limitation among all geometries, can return currentMinimumStep (or DBL_MAX):
344 // -- Beware : the returnedStep must be physically meaningful, even if we say "NotCandidateForSelection" !
345 if ( eLimitedForSmallestStep == kDoNot )
346 {
347 returnedStep = currentMinimumStep;
348 }
349 // -- proposed step length of limiting geometry:
350 if ( eLimitedForSmallestStep == kUnique ||
351 eLimitedForSmallestStep == kSharedOther )
352 {
353 *selection = CandidateForSelection;
354 returnedStep = smallestReturnedStep;
355 }
356 else if ( eLimitedForSmallestStep == kSharedTransport)
357 {
358 returnedStep = smallestReturnedStep* (1.0 + 1.0e-9); // -- Expand to disable its selection in Step Manager comparison
359 }
360
361 // -- and smallest safety among geometries:
362 proposedSafety = fParallelWorldSafety ;
363 }
364
365 // -- returns step length, and proposedSafety
366 return returnedStep;
367}
368
369
371 const G4Step& )
372{
373
374 fDummyParticleChange.Initialize(track);
375 return &fDummyParticleChange;
376}
377
378
380{
381 G4BiasingProcessSharedData *sharedData(nullptr);
382
383 // -- initialize sharedData pointer:
384 if ( G4BiasingProcessSharedData::fSharedDataMap.Find(mgr) == G4BiasingProcessSharedData::fSharedDataMap.End() )
385 {
386 sharedData = new G4BiasingProcessSharedData( mgr );
387 G4BiasingProcessSharedData::fSharedDataMap[mgr] = sharedData;
388 }
389 else sharedData = G4BiasingProcessSharedData::fSharedDataMap[mgr] ;
390
391 // -- add itself to the shared data:
392 if ( sharedData->fParallelGeometriesLimiterProcess == nullptr ) sharedData->fParallelGeometriesLimiterProcess = this;
393 else
394 {
396 ed << " Trying to add more than one G4ParallelGeometriesLimiterProcess process to the process manager " << mgr
397 << " (process manager for `" << mgr->GetParticleType()->GetParticleName() << "'). Only one is needed. Call ignored." << G4endl;
398 G4Exception(" G4ParallelGeometriesLimiterProcess::SetProcessManager(...)",
399 "BIAS.GEN.29",
401 ed);
402 }
403}
404
405
407{
408 G4int toReturn = -1;
409 G4int iWorld = 0;
410 for ( auto world : fParallelWorlds )
411 {
412 if ( world == parallelWorld )
413 {
414 toReturn = iWorld;
415 break;
416 }
417 iWorld++;
418 }
419 return toReturn;
420}
421
422
424{
425 G4VPhysicalVolume* aWorld = fTransportationManager->IsWorldExisting( parallelWorldName ); // note aWorld might be nullptr
426 return GetParallelWorldIndex( aWorld );
427}
428
G4double condition(const G4ErrorSymMatrix &m)
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ForceCondition
@ NotForced
G4GPILSelection
@ CandidateForSelection
@ NotCandidateForSelection
ELimited
@ kDoNot
@ kUnique
@ kSharedOther
@ kSharedTransport
@ fParallel
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
static void Update(G4FieldTrack *, const G4Track *)
G4ThreeVector GetPosition() const
G4VPhysicalVolume * GetWorldVolume() const
G4ParallelGeometriesLimiterProcess(const G4String &processName="biasLimiter")
void AddParallelWorld(const G4String &parallelWorldName)
G4int GetParallelWorldIndex(const G4VPhysicalVolume *parallelWorld) const
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double, G4ForceCondition *)
const G4Navigator * GetNavigator(G4int worldIndex) const
void RemoveParallelWorld(const G4String &parallelWorldName)
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
virtual void SetProcessManager(const G4ProcessManager *)
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &step)
virtual void Initialize(const G4Track &track)
const G4String & GetParticleName() 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=nullptr)
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:52
G4ParticleDefinition * GetParticleType() const
Definition: G4Step.hh:62
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)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
#define DBL_MAX
Definition: templates.hh:62
#define G4ThreadLocal
Definition: tls.hh:77