Geant4 11.4.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 }
190 return fInteractionStep;
191 }
192 inline void SetInteractionStep(G4bool InteractionStep){
193 fInteractionStep = InteractionStep;
194 }
195
196 protected:
197 void DoProcess();
198 void SynchronizeTracks();
199 void Stepping();
200
202
204
205 void PrintWhyDoYouStop();
206
207 private:
208 G4Scheduler();
209 void Create();
210
211 G4SchedulerMessenger* fpMessenger = nullptr;
212
213 static G4ThreadLocal G4Scheduler* fgScheduler;
214 G4int fVerbose;
215 G4bool fWhyDoYouStop;
216 G4bool fInitialized;
217 G4bool fRunning;
218 G4bool fContinue;
219
220 G4int fNbSteps;
221 G4int fMaxSteps;
222
223 G4ITStepStatus fITStepStatus;
224
225 // Time members
226 G4bool fUseDefaultTimeSteps;
227 G4double fTimeTolerance;
228 G4double fGlobalTime;
229 G4double fTmpGlobalTime;
230 G4double fStartTime;
231 G4double fStopTime;
232 G4double fEndTime;
233 G4double fPreviousTimeStep;
234 G4int fZeroTimeCount;
235 G4int fMaxNZeroTimeStepsAllowed;
236
237 G4double fTimeStep; // The selected minimum time step
238 G4double fMaxTimeStep;
239
240 // User steps
241 G4bool fUsePreDefinedTimeSteps;
242 G4double fDefaultMinTimeStep;
243 std::map<G4double, G4double>* fpUserTimeSteps = nullptr;
244 // One can give time steps in respect to the global time
245 mutable G4double fUserUpperTimeLimit;
246 G4double fDefinedMinTimeStep;
247 // selected user time step in respect to the global time
248 G4bool fReachedUserTimeLimit; // if fMinTimeStep == the user time step
249
250 std::set<G4double> fWatchedTimes;
251
252 G4UserTimeStepAction* fpUserTimeStepAction;
253
254 std::unique_ptr<G4VScavengerMaterial> fpUserScavenger;
255
256 // ==========================================
257 // TO BE REMOVED
258 G4ITStepProcessor* fpStepProcessor = nullptr;
259 G4ITModelProcessor* fpModelProcessor = nullptr;
260 G4ITTrackingManager* fpTrackingManager = nullptr;
261 G4ITTrackingInteractivity* fpTrackingInteractivity = nullptr;
262 G4ITReactionSet* fReactionSet = nullptr;
263 G4ITTrackHolder& fTrackContainer;
264 G4ITModelHandler* fpModelHandler = nullptr;
265 // ==========================================
266
267 G4double fTSTimeStep;
268 // Time calculated by the time stepper in CalculateMinTimeStep()
269 G4double fILTimeStep;
270 // Time calculated by the interaction length methods
271 // in ComputeInteractionLength()
272
273 G4bool fInteractionStep;
274 // Flag : if the step is driven by the interaction with the matter and
275 // NOT by the reaction between tracks
276
277 G4ITGun* fpGun;
278
279 // ==========================================
280 // Hoang
281 G4bool fResetScavenger;
282
283 public:
284 void ResetScavenger(bool);
285};
286
288{
289 return fInitialized;
290}
291
293{
294 return fpModelHandler;
295}
296
297inline void G4Scheduler::SetEndTime(const G4double __endtime)
298{
299 fEndTime = __endtime;
300}
301
302inline void G4Scheduler::SetTimeSteps(std::map<G4double, G4double>* steps)
303{
304 fUsePreDefinedTimeSteps = true;
305 fpUserTimeSteps = steps;
306}
307
308inline void G4Scheduler::AddTimeStep(G4double startingTime, G4double timeStep)
309{
310 if (fpUserTimeSteps == nullptr) {
311 fpUserTimeSteps = new std::map<G4double, G4double>();
312 fUsePreDefinedTimeSteps = true;
313 }
314
315 (*fpUserTimeSteps)[startingTime] = timeStep;
316}
317
319{
320 return fNbSteps;
321}
322
324{
325 fMaxSteps = maxSteps;
326}
327
329{
330 return fMaxSteps;
331}
332
334{
335 return fStartTime;
336}
337
339{
340 return fEndTime;
341}
342
344{
345 return fTimeStep;
346}
347
349{
350 fDefaultMinTimeStep = timeStep;
351}
352
354{
355 return fGlobalTime;
356}
357
359{
360 fpUserTimeStepAction = userITAction;
361}
362
364{
365 return fpUserTimeStepAction;
366}
367
368inline void G4Scheduler::SetVerbose(G4int verbose)
369{
370 fVerbose = verbose;
371}
372
374{
375 return fVerbose;
376}
377
378inline void G4Scheduler::SetMaxZeroTimeAllowed(G4int maxTimeStepAllowed)
379{
380 fMaxNZeroTimeStepsAllowed = maxTimeStepAllowed;
381}
382
384{
385 return fMaxNZeroTimeStepsAllowed;
386}
387
389{
390 fTimeTolerance = time;
391}
392
394{
395 return fTimeTolerance;
396}
397
399{
400 return fPreviousTimeStep;
401}
402
404{
405 return fITStepStatus;
406}
407
408inline void G4Scheduler::Stop()
409{
410 fContinue = false;
411}
412
414{
415 return fpTrackingInteractivity;
416}
417
419{
420 fpGun = gun;
421}
422
424{
425 return fpGun;
426}
427
429{
430 fWhyDoYouStop = true;
431}
432
434{
435 fUseDefaultTimeSteps = flag;
436}
437
439{
440 return (!fUseDefaultTimeSteps && !fUsePreDefinedTimeSteps);
441}
442
443inline void G4Scheduler::ResetScavenger(bool value)
444{
445 fResetScavenger = value;
446}
447
448#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
void SetInteractionStep(G4bool InteractionStep)
G4bool IsInitialized()
G4bool IsInteractionStep()
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