Geant4 11.2.2
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
47
48#ifndef G4Scheduler_h
49#define G4Scheduler_h
50
51#include <G4VScheduler.hh>
52#include <vector>
53#include <map>
54#include <memory>
55
56#include "globals.hh"
57
58#include "G4ITModelHandler.hh"
59#include "G4ITStepStatus.hh"
60#include "G4ITTrackHolder.hh"
61#include "G4VStateDependent.hh"
62#include "G4ITReaction.hh"
64
68class G4Track;
72class G4ITGun;
73
74#ifndef compTrackPerID__
75#define compTrackPerID__
76 struct compTrackPerID
77 {
78 G4bool operator()(G4Track* rhs, G4Track* lhs) const
79 {
80 return rhs->GetTrackID() < lhs->GetTrackID();
81 }
82 };
83#endif
84
85/**
86 * G4Scheduler synchronizes (in time) track stepping
87 */
89 public G4VScheduler,
91{
92protected:
93 ~G4Scheduler() override;
94
95public:
96 G4Scheduler(const G4Scheduler&) = delete;
98
99 static G4Scheduler* Instance();
100 /** DeleteInstance should be used instead
101 * of the destructor
102 */
103 static void DeleteInstance();
104 G4bool Notify(G4ApplicationState requestedState) override;
105
106 void RegisterModel(G4VITStepModel*, G4double) override;
107
108 void Initialize() override;
110 inline G4bool IsInitialized();
111 inline G4bool IsRunning() override{return fRunning;}
112 void Reset() override;
113 void Process() override;
114 void ClearList();
115
116 inline void SetGun(G4ITGun*) override;
117 inline G4ITGun* GetGun();
118
119 inline void Stop();
120 void Clear();
121
122 // To be called only in UserReactionAction::EndProcessing()
123 // after fRunning flag has been turned off.
124 // This is not done automatically before UserReactionAction::EndProcessing()
125 // is called in case one would like to access some track information
126 void EndTracking();
127
128 void SetEndTime(const G4double) override;
129
130 /* Two tracks below the time tolerance are supposed to be
131 * in the same time slice
132 */
133 inline void SetTimeTolerance(G4double) override;
134 inline G4double GetTimeTolerance() const override;
135
136 inline void SetMaxZeroTimeAllowed(G4int) override;
137 inline G4int GetMaxZeroTimeAllowed() const override;
138
139 inline G4ITModelHandler* GetModelHandler() override;
140
141 inline void SetTimeSteps(std::map<G4double, G4double>*) override;
142 inline void AddTimeStep(G4double, G4double) override;
143 inline void SetDefaultTimeStep(G4double) override;
144 G4double GetLimitingTimeStep() const override;
145 inline G4int GetNbSteps() const override;
146 inline void SetMaxNbSteps(G4int) override;
147 inline G4int GetMaxNbSteps() const override;
148 inline G4double GetStartTime() const override;
149 inline G4double GetEndTime() const override;
150 inline G4double GetTimeStep() const override;
151 inline G4double GetPreviousTimeStep() const override;
152 inline G4double GetGlobalTime() const override;
153 inline void SetUserAction(G4UserTimeStepAction*) override;
154 inline G4UserTimeStepAction* GetUserTimeStepAction() const override;
155
156 // To use with transportation only, no reactions
157 inline void UseDefaultTimeSteps(G4bool);
159
160 inline G4ITStepStatus GetStatus() const;
161
162 /* 1 : Reaction information
163 * 2 : (1) + time step information
164 * 3 : (2) + step info for individual tracks
165 * 4 : (2) + trackList processing info + pushed and killed track info
166 */
167 inline void SetVerbose(G4int) override;
168
169 inline G4int GetVerbose() const;
170
171 inline void WhyDoYouStop();
172
175
176 virtual size_t GetNTracks();
177
178 void GetCollisionType(G4String& interactionType);
179
181 {
182 fWatchedTimes.insert(time);
183 }
184
186
187 inline void SetMaxTimeStep(G4double maxTimeStep)
188 {
189 fMaxTimeStep = maxTimeStep;
190 }
191
193 {
194 return fMaxTimeStep;
195 }
196
198 {
199 return fpUserScavenger.get();
200 }
201 inline void SetScavengerMaterial(std::unique_ptr<G4VScavengerMaterial> scavengerMaterial)
202 {
203 fpUserScavenger = std::move(scavengerMaterial);
204 }
205
206protected:
207
208 void DoProcess();
209 void SynchronizeTracks();
210 void Stepping();
211
213
215
216 void PrintWhyDoYouStop();
217
218private:
219 G4Scheduler();
220 void Create();
221
222 G4SchedulerMessenger* fpMessenger;
223
224 static G4ThreadLocal G4Scheduler* fgScheduler;
225 G4int fVerbose;
226 G4bool fWhyDoYouStop;
227 G4bool fInitialized;
228 G4bool fRunning;
229 G4bool fContinue;
230
231 G4int fNbSteps;
232 G4int fMaxSteps;
233
234 G4ITStepStatus fITStepStatus;
235
236 // Time members
237 G4bool fUseDefaultTimeSteps;
238 G4double fTimeTolerance;
239 G4double fGlobalTime;
240 G4double fTmpGlobalTime;
241 G4double fStartTime;
242 G4double fStopTime;
243 G4double fEndTime;
244 G4double fPreviousTimeStep;
245 G4int fZeroTimeCount;
246 G4int fMaxNZeroTimeStepsAllowed;
247
248 G4double fTimeStep; // The selected minimum time step
249 G4double fMaxTimeStep;
250
251 // User steps
252 G4bool fUsePreDefinedTimeSteps;
253 G4double fDefaultMinTimeStep;
254 std::map<G4double, G4double>* fpUserTimeSteps;
255 // One can give time steps in respect to the global time
256 mutable G4double fUserUpperTimeLimit;
257 G4double fDefinedMinTimeStep;
258 // selected user time step in respect to the global time
259 G4bool fReachedUserTimeLimit; // if fMinTimeStep == the user time step
260
261 std::set<G4double> fWatchedTimes;
262
263 G4UserTimeStepAction* fpUserTimeStepAction;
264
265 std::unique_ptr<G4VScavengerMaterial> fpUserScavenger;
266
267 // ==========================================
268 // TO BE REMOVED
269 G4ITStepProcessor* fpStepProcessor;
270 G4ITModelProcessor* fpModelProcessor;
271 G4ITTrackingManager* fpTrackingManager;
272 G4ITTrackingInteractivity* fpTrackingInteractivity;
273 G4ITReactionSet* fReactionSet;
274 G4ITTrackHolder& fTrackContainer;
275 G4ITModelHandler* fpModelHandler;
276 // ==========================================
277
278 G4double fTSTimeStep;
279 // Time calculated by the time stepper in CalculateMinTimeStep()
280 G4double fILTimeStep;
281 // Time calculated by the interaction length methods
282 // in ComputeInteractionLength()
283
284 G4bool fInteractionStep;
285 // Flag : if the step is driven by the interaction with the matter and
286 // NOT by the reaction between tracks
287
288 G4ITGun* fpGun;
289
290 // ==========================================
291 //Hoang
292 bool fResetScavenger;
293public:
294 void ResetScavenger(bool);
295};
296
298{
299 return fInitialized;
300}
301
303{
304 return fpModelHandler;
305}
306
307inline void G4Scheduler::SetEndTime(const G4double __endtime)
308{
309 fEndTime = __endtime;
310}
311
312inline
313void G4Scheduler::SetTimeSteps(std::map<G4double, G4double>* steps)
314{
315 fUsePreDefinedTimeSteps = true;
316 fpUserTimeSteps = steps;
317}
318
319inline void G4Scheduler::AddTimeStep(G4double startingTime, G4double timeStep)
320{
321 if (fpUserTimeSteps == nullptr)
322 {
323 fpUserTimeSteps = new std::map<G4double, G4double>();
324 fUsePreDefinedTimeSteps = true;
325 }
326
327 (*fpUserTimeSteps)[startingTime] = timeStep;
328}
329
331{
332 return fNbSteps;
333}
334
336{
337 fMaxSteps = maxSteps;
338}
339
341{
342 return fMaxSteps;
343}
344
346{
347 return fStartTime;
348}
349
351{
352 return fEndTime;
353}
354
356{
357 return fTimeStep;
358}
359
361{
362 fDefaultMinTimeStep = timeStep;
363}
364
366{
367 return fGlobalTime;
368}
369
370inline
372{
373 fpUserTimeStepAction = userITAction;
374}
375
377{
378 return fpUserTimeStepAction;
379}
380
381inline void G4Scheduler::SetVerbose(G4int verbose)
382{
383 fVerbose = verbose;
384}
385
387{
388 return fVerbose;
389}
390
391inline
393{
394 fMaxNZeroTimeStepsAllowed = maxTimeStepAllowed;
395}
396
398{
399 return fMaxNZeroTimeStepsAllowed;
400}
401
403{
404 fTimeTolerance = time;
405}
406
408{
409 return fTimeTolerance;
410}
411
413{
414 return fPreviousTimeStep;
415}
416
418{
419 return fITStepStatus;
420}
421
422inline void G4Scheduler::Stop()
423{
424 fContinue = false;
425}
426
428{
429 return fpTrackingInteractivity;
430}
431
433{
434 fpGun = gun;
435}
436
438{
439 return fpGun;
440}
441
443{
444 fWhyDoYouStop = true;
445}
446
448{
449 fUseDefaultTimeSteps = flag;
450}
451
453{
454 return (!fUseDefaultTimeSteps && !fUsePreDefinedTimeSteps);
455}
456
457inline void G4Scheduler::ResetScavenger(bool value)
458{
459 fResetScavenger = value;
460}
461
462#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
G4bool operator()(G4Track *rhs, G4Track *lhs) const
#define G4ThreadLocal
Definition tls.hh:77