Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Scheduler.hh
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// Author: Mathieu Karamitros
28
29// The code is developed in the framework of the ESA AO7146
30//
31// We would be very happy hearing from you, send us your feedback! :)
32//
33// In order for Geant4-DNA to be maintained and still open-source,
34// article citations are crucial.
35// If you use Geant4-DNA chemistry and you publish papers about your software,
36// in addition to the general paper on Geant4-DNA:
37//
38// Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178
39//
40// we would be very happy if you could please also cite the following
41// reference papers on chemistry:
42//
43// J. Comput. Phys. 274 (2014) 841-882
44// Prog. Nucl. Sci. Tec. 2 (2011) 503-508
45
46#ifndef G4Scheduler_h
47#define G4Scheduler_h
48
49#include "G4ITModelHandler.hh"
50#include "G4ITReaction.hh"
51#include "G4ITStepStatus.hh"
52#include "G4ITTrackHolder.hh"
54#include "G4VStateDependent.hh"
55#include "globals.hh"
56
57#include <G4VScheduler.hh>
58
59#include <map>
60#include <memory>
61#include <vector>
62
66class G4Track;
70class G4ITGun;
71
72#ifndef compTrackPerID__
73# define compTrackPerID__
74struct compTrackPerID
75{
76 G4bool operator()(G4Track* rhs, G4Track* lhs) const
77 {
78 return rhs->GetTrackID() < lhs->GetTrackID();
79 }
80};
81#endif
82
83/**
84 * G4Scheduler synchronizes (in time) track stepping
85 */
87{
88 protected:
89 ~G4Scheduler() override;
90
91 public:
92 G4Scheduler(const G4Scheduler&) = delete;
94
95 static G4Scheduler* Instance();
96 /** DeleteInstance should be used instead
97 * of the destructor
98 */
99 static void DeleteInstance();
100 G4bool Notify(G4ApplicationState requestedState) override;
101
102 void RegisterModel(G4VITStepModel*, G4double) override;
103
104 void Initialize() override;
106 inline G4bool IsInitialized();
107 inline G4bool IsRunning() override { return fRunning; }
108 void Reset() override;
109 void Process() override;
110 void ClearList();
111
112 inline void SetGun(G4ITGun*) override;
113 inline G4ITGun* GetGun();
114
115 inline void Stop();
116 void Clear();
117
118 // To be called only in UserReactionAction::EndProcessing()
119 // after fRunning flag has been turned off.
120 // This is not done automatically before UserReactionAction::EndProcessing()
121 // is called in case one would like to access some track information
122 void EndTracking();
123
124 void SetEndTime(const G4double) override;
125
126 /* Two tracks below the time tolerance are supposed to be
127 * in the same time slice
128 */
129 inline void SetTimeTolerance(G4double) override;
130 inline G4double GetTimeTolerance() const override;
131
132 inline void SetMaxZeroTimeAllowed(G4int) override;
133 inline G4int GetMaxZeroTimeAllowed() const override;
134
135 inline G4ITModelHandler* GetModelHandler() override;
136
137 inline void SetTimeSteps(std::map<G4double, G4double>*) override;
138 inline void AddTimeStep(G4double, G4double) override;
139 inline void SetDefaultTimeStep(G4double) override;
140 G4double GetLimitingTimeStep() const override;
141 inline G4int GetNbSteps() const override;
142 inline void SetMaxNbSteps(G4int) override;
143 inline G4int GetMaxNbSteps() const override;
144 inline G4double GetStartTime() const override;
145 inline G4double GetEndTime() const override;
146 inline G4double GetTimeStep() const override;
147 inline G4double GetPreviousTimeStep() const override;
148 inline G4double GetGlobalTime() const override;
149 inline void SetUserAction(G4UserTimeStepAction*) override;
150 inline G4UserTimeStepAction* GetUserTimeStepAction() const override;
151
152 // To use with transportation only, no reactions
153 inline void UseDefaultTimeSteps(G4bool);
155
156 inline G4ITStepStatus GetStatus() const;
157
158 /* 1 : Reaction information
159 * 2 : (1) + time step information
160 * 3 : (2) + step info for individual tracks
161 * 4 : (2) + trackList processing info + pushed and killed track info
162 */
163 inline void SetVerbose(G4int) override;
164
165 inline G4int GetVerbose() const;
166
167 inline void WhyDoYouStop();
168
171
172 virtual size_t GetNTracks();
173
174 void GetCollisionType(G4String& interactionType);
175
176 void AddWatchedTime(G4double time) { fWatchedTimes.insert(time); }
177
179
180 inline void SetMaxTimeStep(G4double maxTimeStep) { fMaxTimeStep = maxTimeStep; }
181
182 inline G4double GetMaxTimeStep() const { return fMaxTimeStep; }
183
184 inline G4VScavengerMaterial* GetScavengerMaterial() const { return fpUserScavenger.get(); }
185 inline void SetScavengerMaterial(std::unique_ptr<G4VScavengerMaterial> scavengerMaterial)
186 {
187 fpUserScavenger = std::move(scavengerMaterial);
188 }
189
190 protected:
191 void DoProcess();
192 void SynchronizeTracks();
193 void Stepping();
194
196
198
199 void PrintWhyDoYouStop();
200
201 private:
202 G4Scheduler();
203 void Create();
204
205 G4SchedulerMessenger* fpMessenger = nullptr;
206
207 static G4ThreadLocal G4Scheduler* fgScheduler;
208 G4int fVerbose;
209 G4bool fWhyDoYouStop;
210 G4bool fInitialized;
211 G4bool fRunning;
212 G4bool fContinue;
213
214 G4int fNbSteps;
215 G4int fMaxSteps;
216
217 G4ITStepStatus fITStepStatus;
218
219 // Time members
220 G4bool fUseDefaultTimeSteps;
221 G4double fTimeTolerance;
222 G4double fGlobalTime;
223 G4double fStartTime;
224 G4double fStopTime;
225 G4double fEndTime;
226 G4double fPreviousTimeStep;
227 G4int fZeroTimeCount;
228 G4int fMaxNZeroTimeStepsAllowed;
229
230 G4double fTimeStep; // The selected minimum time step
231 G4double fMaxTimeStep;
232
233 // User steps
234 G4bool fUsePreDefinedTimeSteps;
235 G4double fDefaultMinTimeStep;
236 std::map<G4double, G4double>* fpUserTimeSteps = nullptr;
237 // One can give time steps in respect to the global time
238 mutable G4double fUserUpperTimeLimit;
239 G4double fDefinedMinTimeStep;
240 // selected user time step in respect to the global time
241 G4bool fReachedUserTimeLimit; // if fMinTimeStep == the user time step
242
243 std::set<G4double> fWatchedTimes;
244
245 G4UserTimeStepAction* fpUserTimeStepAction;
246
247 std::unique_ptr<G4VScavengerMaterial> fpUserScavenger;
248
249 // ==========================================
250 // TO BE REMOVED
251 G4ITStepProcessor* fpStepProcessor = nullptr;
252 G4ITModelProcessor* fpModelProcessor = nullptr;
253 G4ITTrackingManager* fpTrackingManager = nullptr;
254 G4ITTrackingInteractivity* fpTrackingInteractivity = nullptr;
255 G4ITReactionSet* fReactionSet = nullptr;
256 G4ITTrackHolder& fTrackContainer;
257 G4ITModelHandler* fpModelHandler = nullptr;
258 // ==========================================
259
260 G4double fTSTimeStep;
261 // Time calculated by the time stepper in CalculateMinTimeStep()
262 G4double fILTimeStep;
263 // Time calculated by the interaction length methods
264 // in ComputeInteractionLength()
265
266 G4bool fInteractionStep;
267 // Flag : if the step is driven by the interaction with the matter and
268 // NOT by the reaction between tracks
269
270 G4ITGun* fpGun;
271
272 // ==========================================
273 // Hoang
274 G4bool fResetScavenger;
275
276 public:
277 void ResetScavenger(bool);
278};
279
281{
282 return fInitialized;
283}
284
286{
287 return fpModelHandler;
288}
289
290inline void G4Scheduler::SetEndTime(const G4double __endtime)
291{
292 fEndTime = __endtime;
293}
294
295inline void G4Scheduler::SetTimeSteps(std::map<G4double, G4double>* steps)
296{
297 fUsePreDefinedTimeSteps = true;
298 fpUserTimeSteps = steps;
299}
300
301inline void G4Scheduler::AddTimeStep(G4double startingTime, G4double timeStep)
302{
303 if (fpUserTimeSteps == nullptr) {
304 fpUserTimeSteps = new std::map<G4double, G4double>();
305 fUsePreDefinedTimeSteps = true;
306 }
307
308 (*fpUserTimeSteps)[startingTime] = timeStep;
309}
310
312{
313 return fNbSteps;
314}
315
317{
318 fMaxSteps = maxSteps;
319}
320
322{
323 return fMaxSteps;
324}
325
327{
328 return fStartTime;
329}
330
332{
333 return fEndTime;
334}
335
337{
338 return fTimeStep;
339}
340
342{
343 fDefaultMinTimeStep = timeStep;
344}
345
347{
348 return fGlobalTime;
349}
350
352{
353 fpUserTimeStepAction = userITAction;
354}
355
357{
358 return fpUserTimeStepAction;
359}
360
361inline void G4Scheduler::SetVerbose(G4int verbose)
362{
363 fVerbose = verbose;
364}
365
367{
368 return fVerbose;
369}
370
371inline void G4Scheduler::SetMaxZeroTimeAllowed(G4int maxTimeStepAllowed)
372{
373 fMaxNZeroTimeStepsAllowed = maxTimeStepAllowed;
374}
375
377{
378 return fMaxNZeroTimeStepsAllowed;
379}
380
382{
383 fTimeTolerance = time;
384}
385
387{
388 return fTimeTolerance;
389}
390
392{
393 return fPreviousTimeStep;
394}
395
397{
398 return fITStepStatus;
399}
400
401inline void G4Scheduler::Stop()
402{
403 fContinue = false;
404}
405
407{
408 return fpTrackingInteractivity;
409}
410
412{
413 fpGun = gun;
414}
415
417{
418 return fpGun;
419}
420
422{
423 fWhyDoYouStop = true;
424}
425
427{
428 fUseDefaultTimeSteps = flag;
429}
430
432{
433 return (!fUseDefaultTimeSteps && !fUsePreDefinedTimeSteps);
434}
435
436inline void G4Scheduler::ResetScavenger(bool value)
437{
438 fResetScavenger = value;
439}
440
441#endif
G4ApplicationState
G4ITStepStatus
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4double GetStartTime() const override
void SetEndTime(const G4double) override
void SynchronizeTracks()
G4bool CanICarryOn()
void SetMaxTimeStep(G4double maxTimeStep)
G4double GetEndTime() const override
void SetDefaultTimeStep(G4double) override
void SetUserAction(G4UserTimeStepAction *) override
void EndTracking()
G4double GetTimeStep() const override
void FindUserPreDefinedTimeStep()
G4bool IsRunning() override
void ForceReinitialization()
G4Scheduler & operator=(const G4Scheduler &)=delete
G4bool Notify(G4ApplicationState requestedState) override
G4int GetNbSteps() const override
G4UserTimeStepAction * GetUserTimeStepAction() const override
G4double GetNextWatchedTime() const
G4double GetMaxTimeStep() const
void SetVerbose(G4int) override
G4double GetTimeTolerance() const override
G4bool IsInitialized()
void SetMaxZeroTimeAllowed(G4int) override
G4ITModelHandler * GetModelHandler() override
void Reset() override
void SetGun(G4ITGun *) override
void Stepping()
void SetMaxNbSteps(G4int) override
static G4Scheduler * Instance()
G4ITStepStatus GetStatus() const
void AddTimeStep(G4double, G4double) override
G4int GetMaxNbSteps() const override
G4double GetGlobalTime() const override
void Process() override
G4int GetMaxZeroTimeAllowed() const override
void ClearList()
void Initialize() override
G4bool AreDefaultTimeStepsUsed()
void SetTimeTolerance(G4double) override
G4double GetPreviousTimeStep() const override
void WhyDoYouStop()
virtual size_t GetNTracks()
void AddWatchedTime(G4double time)
G4Scheduler(const G4Scheduler &)=delete
G4ITGun * GetGun()
G4ITTrackingInteractivity * GetInteractivity() override
void GetCollisionType(G4String &interactionType)
void SetInteractivity(G4ITTrackingInteractivity *) override
G4double GetLimitingTimeStep() const override
void ResetScavenger(bool)
void DoProcess()
void PrintWhyDoYouStop()
~G4Scheduler() override
void RegisterModel(G4VITStepModel *, G4double) override
G4int GetVerbose() const
static void DeleteInstance()
void SetScavengerMaterial(std::unique_ptr< G4VScavengerMaterial > scavengerMaterial)
void SetTimeSteps(std::map< G4double, G4double > *) override
void UseDefaultTimeSteps(G4bool)
G4VScavengerMaterial * GetScavengerMaterial() const
G4int GetTrackID() const
G4VStateDependent(G4bool bottom=false)
G4bool operator()(G4Track *rhs, G4Track *lhs) const
#define G4ThreadLocal
Definition tls.hh:77