Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ITModelProcessor Class Reference

#include <G4ITModelProcessor.hh>

Public Member Functions

 G4ITModelProcessor ()
 
 G4ITModelProcessor (const G4ITModelProcessor &other)=delete
 
G4ITModelProcessoroperator= (const G4ITModelProcessor &other)=delete
 
virtual ~G4ITModelProcessor ()
 
void SetModelHandler (G4ITModelHandler *)
 
void SetTrackingManager (G4ITTrackingManager *trackingManager)
 
void Initialize ()
 
void RegisterModel (double time, G4VITStepModel *)
 
void CleanProcessor ()
 
G4double CalculateMinTimeStep (G4double currentGlobalTime, G4double definedMinTimeStep)
 
void ComputeTrackReaction (G4ITStepStatus fITStepStatus, G4double fGlobalTime, G4double currentTimeStep, G4double previousTimeStep, G4bool reachedUserTimeLimit, G4double fTimeTolerance, G4UserTimeStepAction *fpUserTimeStepAction, G4int fVerbose)
 
void InitializeStepper (G4double currentGlobalTime, G4double userMinTime)
 
bool GetComputeTimeStep () const
 
const G4TrackGetTrack () const
 

Protected Member Functions

void SetTrack (const G4Track *)
 

Protected Attributes

G4double fTSTimeStep
 
G4ITReactionSetfReactionSet
 
G4ITTrackingManagerfpTrackingManager
 
G4ITTrackHolderfpTrackContainer
 
G4bool fInitialized
 
G4ITModelHandlerfpModelHandler
 
const G4TrackfpTrack
 
G4double fUserMinTimeStep
 
std::vector< G4VITStepModel * > fActiveModels
 
G4VITStepModelfpActiveModelWithMinTimeStep
 
std::vector< std::unique_ptr< G4ITReactionChange > > fReactionInfo
 
bool fComputeTimeStep
 
bool fComputeReaction
 

Detailed Description

The G4ITModelProcessor will call the two processes defined in G4VITModel. This processes act at the beginning and end of each step. The first one, the TimeStepper will calculate a time step to propagate all the track and eventually it can return some tracks that can likely react at the end of the step. The second one, the ReactionProcess will make the tracks reacting.

Deprecated:
This class will be removed

Definition at line 71 of file G4ITModelProcessor.hh.

Constructor & Destructor Documentation

◆ G4ITModelProcessor() [1/2]

G4ITModelProcessor::G4ITModelProcessor ( )

Definition at line 53 of file G4ITModelProcessor.cc.

54{
55 fpTrack = nullptr;
56 fInitialized = false;
57 fUserMinTimeStep = -1.;
59 fpTrackingManager = nullptr;
60 fReactionSet = nullptr;
61 fpTrackContainer = nullptr;
62 fpModelHandler = nullptr;
64 fComputeTimeStep = false;
65 fComputeReaction = false;
66}
G4ITTrackHolder * fpTrackContainer
G4ITReactionSet * fReactionSet
G4VITStepModel * fpActiveModelWithMinTimeStep
const G4Track * fpTrack
G4ITTrackingManager * fpTrackingManager
G4ITModelHandler * fpModelHandler
#define DBL_MAX
Definition: templates.hh:62

◆ G4ITModelProcessor() [2/2]

G4ITModelProcessor::G4ITModelProcessor ( const G4ITModelProcessor other)
delete

◆ ~G4ITModelProcessor()

G4ITModelProcessor::~G4ITModelProcessor ( )
virtualdefault

Member Function Documentation

◆ CalculateMinTimeStep()

G4double G4ITModelProcessor::CalculateMinTimeStep ( G4double  currentGlobalTime,
G4double  definedMinTimeStep 
)

Definition at line 93 of file G4ITModelProcessor.cc.

95{
96
97#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
98 MemStat mem_first, mem_second, mem_diff;
99 mem_first = MemoryUsage();
100#endif
101
104
105 InitializeStepper(currentGlobalTime, definedMinTimeStep);
106
107#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
108 mem_second = MemoryUsage();
109 mem_diff = mem_second-mem_first;
110 G4cout << "\t || MEM || G4Scheduler::CalculateMinTimeStep || After "
111 "computing fpModelProcessor -> InitializeStepper, diff is : "
112 << mem_diff
113 << G4endl;
114#endif
115
116 for (auto& pStepModel : fActiveModels)
117 {
119 pStepModel->GetTimeStepper()->CalculateMinTimeStep(
120 currentGlobalTime,
121 definedMinTimeStep);
122
123 fpActiveModelWithMinTimeStep = pStepModel;
124
125 if(fTSTimeStep == -1){
127 if(fReactionSet->Empty()) return DBL_MAX;
128 auto fReactionSetInTime = fReactionSet->GetReactionsPerTime();
129 fTSTimeStep = fReactionSetInTime.begin()->get()->GetTime() - currentGlobalTime;
130 }
131 }
132
133#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
134 mem_second = MemoryUsage();
135 mem_diff = mem_second-mem_first;
136 G4cout << "\t || MEM || G4Scheduler::CalculateMinTimeStep || "
137 "After looping on tracks, diff is : " << mem_diff << G4endl;
138#endif
139 return fTSTimeStep;
140}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void InitializeStepper(G4double currentGlobalTime, G4double userMinTime)
std::vector< G4VITStepModel * > fActiveModels
G4ITReactionPerTime & GetReactionsPerTime()
G4VITReactionProcess * GetReactionProcess()
MemStat MemoryUsage()
Definition: G4MemStat.cc:55

Referenced by G4Scheduler::Stepping().

◆ CleanProcessor()

void G4ITModelProcessor::CleanProcessor ( )

Restore the original state. This method should be called only by G4Scheduler

Definition at line 373 of file G4ITModelProcessor.cc.

374{
375 fpTrack = nullptr;
376}

◆ ComputeTrackReaction()

void G4ITModelProcessor::ComputeTrackReaction ( G4ITStepStatus  fITStepStatus,
G4double  fGlobalTime,
G4double  currentTimeStep,
G4double  previousTimeStep,
G4bool  reachedUserTimeLimit,
G4double  fTimeTolerance,
G4UserTimeStepAction fpUserTimeStepAction,
G4int  fVerbose 
)

Definition at line 171 of file G4ITModelProcessor.cc.

183{
184 if (fReactionSet->Empty())
185 {
186 return;
187 }
188
189 if (fITStepStatus == eCollisionBetweenTracks)
190 {
192 fReactionInfo = pReactionProcess->FindReaction(fReactionSet,
193 currentTimeStep,
194 fGlobalTime,
195 reachedUserTimeLimit);
196
197 // TODO
198 // A ne faire uniquement si le temps choisis est celui calculé par le time stepper
199 // Sinon utiliser quelque chose comme : fModelProcessor->FindReaction(&fMainList);
200
201 for (auto& pChanges : fReactionInfo)
202 {
203 auto pTrackA = const_cast<G4Track*>(pChanges->GetTrackA());
204 auto pTrackB = const_cast<G4Track*>(pChanges->GetTrackB());
205
206 if (pTrackA == nullptr
207 || pTrackB == nullptr
208 || pTrackA->GetTrackStatus() == fStopAndKill
209 || pTrackB->GetTrackStatus() == fStopAndKill)
210 {
211 continue;
212 }
213
214 G4int nbSecondaries = pChanges->GetNumberOfSecondaries();
215 const std::vector<G4Track*>* productsVector = pChanges->GetfSecondary();
216
217 if (fpUserTimeStepAction)
218 {
219 fpUserTimeStepAction->UserReactionAction(*pTrackA,
220 *pTrackB,
221 productsVector);
222 }
223
224#ifdef G4VERBOSE
225 if (fVerbose)
226 {
227 G4cout << "At time : " << std::setw(7) << G4BestUnit(fGlobalTime, "Time")
228 << " Reaction : " << GetIT(pTrackA)->GetName() << " ("
229 << pTrackA->GetTrackID() << ") + " << GetIT(pTrackB)->GetName() << " ("
230 << pTrackB->GetTrackID() << ") -> ";
231 }
232#endif
233
234 if (nbSecondaries > 0)
235 {
236 for (int i = 0; i < nbSecondaries; ++i)
237 {
238#ifdef G4VERBOSE
239 if (fVerbose && i != 0)
240 {
241 G4cout << " + ";
242 }
243#endif
244
245 G4Track* secondary = (*productsVector)[i]; //changes->GetSecondary(i);
246// fpTrackContainer->_PushTrack(secondary);
247 GetIT(secondary)->SetParentID(pTrackA->GetTrackID(),
248 pTrackB->GetTrackID());
249
250 if (secondary->GetGlobalTime() - fGlobalTime > fTimeTolerance)
251 {
252 G4ExceptionDescription exceptionDescription;
253 exceptionDescription << "The time of the secondary should not be bigger than the"
254 " current global time."
255 << " This may cause synchronization problem. If the process you"
256 " are using required "
257 << "such feature please contact the developers." << G4endl
258 << "The global time in the step manager : "
259 << G4BestUnit(fGlobalTime, "Time")
260 << G4endl
261 << "The global time of the track : "
262 << G4BestUnit(secondary->GetGlobalTime(), "Time")
263 << G4endl;
264
265 G4Exception("G4Scheduler::ComputeInteractionBetweenTracks",
266 "ITScheduler010",
268 exceptionDescription);
269 }
270
271#ifdef G4VERBOSE
272 if (fVerbose)
273 {
274 G4cout << GetIT(secondary)->GetName() << " ("
275 << secondary->GetTrackID() << ")";
276 }
277#endif
278 }
279 }
280 else
281 {
282#ifdef G4VERBOSE
283 if (fVerbose)
284 {
285 G4cout << "No product";
286 }
287#endif
288 }
289#ifdef G4VERBOSE
290 if (fVerbose)
291 {
292 G4cout << G4endl;
293 }
294#endif
295 if (pTrackA->GetTrackID() == 0 || pTrackB->GetTrackID() == 0)
296 {
297 G4Track* pTrack = nullptr;
298 if (pTrackA->GetTrackID() == 0)
299 {
300 pTrack = pTrackA;
301 }
302 else
303 {
304 pTrack = pTrackB;
305 }
306
307 G4ExceptionDescription exceptionDescription;
308 exceptionDescription
309 << "The problem was found for the reaction between tracks :"
310 << pTrackA->GetParticleDefinition()->GetParticleName() << " ("
311 << pTrackA->GetTrackID() << ") & "
312 << pTrackB->GetParticleDefinition()->GetParticleName() << " ("
313 << pTrackB->GetTrackID() << "). \n";
314
315 if (pTrack->GetStep() == nullptr)
316 {
317 exceptionDescription << "Also no step was found"
318 << " ie track->GetStep() == 0 \n";
319 }
320
321 exceptionDescription << "Parent ID of trackA : "
322 << pTrackA->GetParentID() << "\n";
323 exceptionDescription << "Parent ID of trackB : "
324 << pTrackB->GetParentID() << "\n";
325
326 exceptionDescription
327 << "The ID of one of the reaction track was not setup.";
328 G4Exception("G4Scheduler::ComputeInteractionBetweenTracks",
329 "ITScheduler011",
331 exceptionDescription);
332 }
333
334 if (pChanges->WereParentsKilled())
335 {
336 pTrackA->SetTrackStatus(fStopAndKill);
337 pTrackB->SetTrackStatus(fStopAndKill);
338
341 }
342
343 pChanges.reset(nullptr);
344 }
345
346 fReactionInfo.clear();
347 }
348
349// fReactionSet->CleanAllReaction();
350
353}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
@ eCollisionBetweenTracks
G4IT * GetIT(const G4Track *track)
Definition: G4IT.cc:48
#define G4BestUnit(a, b)
@ fStopAndKill
int G4int
Definition: G4Types.hh:85
std::vector< std::unique_ptr< G4ITReactionChange > > fReactionInfo
void MergeSecondariesWithMainList()
void EndTracking(G4Track *)
void SetParentID(int, int)
Definition: G4IT.hh:228
virtual const G4String & GetName() const =0
const G4String & GetParticleName() const
G4int GetTrackID() const
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetGlobalTime() const
const G4Step * GetStep() const
virtual void UserReactionAction(const G4Track &, const G4Track &, const std::vector< G4Track * > *)
virtual std::vector< std::unique_ptr< G4ITReactionChange > > FindReaction(G4ITReactionSet *, const double, const double, const bool)=0

Referenced by G4Scheduler::Stepping().

◆ GetComputeTimeStep()

bool G4ITModelProcessor::GetComputeTimeStep ( ) const

Definition at line 378 of file G4ITModelProcessor.cc.

379{
380 return fComputeTimeStep;
381}

Referenced by G4Scheduler::Stepping().

◆ GetTrack()

const G4Track * G4ITModelProcessor::GetTrack ( ) const

Definition at line 383 of file G4ITModelProcessor.cc.

384{
385 return fpTrack;
386}

◆ Initialize()

void G4ITModelProcessor::Initialize ( )

Definition at line 75 of file G4ITModelProcessor.cc.

76{
80 fInitialized = true;
81 fComputeTimeStep = false;
82 fComputeReaction = false;
84 {
85 fComputeTimeStep = true;
86 }
88 {
89 fComputeReaction = true;
90 }
91}
static G4ITReactionSet * Instance()
static G4ITTrackHolder * Instance()

Referenced by G4Scheduler::Process().

◆ InitializeStepper()

void G4ITModelProcessor::InitializeStepper ( G4double  currentGlobalTime,
G4double  userMinTime 
)

Definition at line 144 of file G4ITModelProcessor.cc.

146{
147 G4VITTimeStepComputer::SetTimes(currentGlobalTime, userMinTime);
148
149#if defined (DEBUG_MEM)
150 MemStat mem_first, mem_second, mem_diff;
151 mem_first = MemoryUsage();
152#endif
153
154 fActiveModels = fpModelHandler->GetActiveModels(currentGlobalTime);
155
156 for (auto& pModel : fActiveModels)
157 {
158 pModel->PrepareNewTimeStep();
159 }
160
161#if defined (DEBUG_MEM)
162 mem_second = MemoryUsage();
163 mem_diff = mem_second-mem_first;
164 G4cout << "\t || MEM || G4ITModelProcessor::InitializeStepper || After computing stepper -> Prepare(), diff is : " << mem_diff << G4endl;
165#endif
166
167}
std::vector< G4VITStepModel * > GetActiveModels(G4double globalTime) const
static void SetTimes(const G4double &, const G4double &)

Referenced by CalculateMinTimeStep().

◆ operator=()

G4ITModelProcessor & G4ITModelProcessor::operator= ( const G4ITModelProcessor other)
delete

◆ RegisterModel()

void G4ITModelProcessor::RegisterModel ( double  time,
G4VITStepModel model 
)

Definition at line 70 of file G4ITModelProcessor.cc.

71{
72 fpModelHandler->RegisterModel(model, time);
73}
void RegisterModel(G4VITStepModel *pModel, G4double globalTime)

◆ SetModelHandler()

void G4ITModelProcessor::SetModelHandler ( G4ITModelHandler pModelHandler)

Definition at line 360 of file G4ITModelProcessor.cc.

361{
362 if (fInitialized)
363 {
364 G4ExceptionDescription exceptionDescription;
365 exceptionDescription
366 << "You are trying to set a new model while the model processor has alreaday be initialized";
367 G4Exception("G4ITModelProcessor::SetModelHandler", "ITModelProcessor001",
368 FatalErrorInArgument, exceptionDescription);
369 }
370 fpModelHandler = pModelHandler;
371}

Referenced by G4Scheduler::Initialize().

◆ SetTrack()

void G4ITModelProcessor::SetTrack ( const G4Track track)
protected

Definition at line 355 of file G4ITModelProcessor.cc.

356{
357 fpTrack = track;
358}

◆ SetTrackingManager()

void G4ITModelProcessor::SetTrackingManager ( G4ITTrackingManager trackingManager)

Definition at line 388 of file G4ITModelProcessor.cc.

389{
390 fpTrackingManager = pTrackingManager;
391}

Referenced by G4Scheduler::Initialize().

Member Data Documentation

◆ fActiveModels

std::vector<G4VITStepModel*> G4ITModelProcessor::fActiveModels
protected

Definition at line 123 of file G4ITModelProcessor.hh.

Referenced by CalculateMinTimeStep(), and InitializeStepper().

◆ fComputeReaction

bool G4ITModelProcessor::fComputeReaction
protected

Definition at line 129 of file G4ITModelProcessor.hh.

Referenced by G4ITModelProcessor(), and Initialize().

◆ fComputeTimeStep

bool G4ITModelProcessor::fComputeTimeStep
protected

Definition at line 128 of file G4ITModelProcessor.hh.

Referenced by G4ITModelProcessor(), GetComputeTimeStep(), and Initialize().

◆ fInitialized

G4bool G4ITModelProcessor::fInitialized
protected

Definition at line 117 of file G4ITModelProcessor.hh.

Referenced by G4ITModelProcessor(), Initialize(), and SetModelHandler().

◆ fpActiveModelWithMinTimeStep

G4VITStepModel* G4ITModelProcessor::fpActiveModelWithMinTimeStep
protected

◆ fpModelHandler

G4ITModelHandler* G4ITModelProcessor::fpModelHandler
protected

◆ fpTrack

const G4Track* G4ITModelProcessor::fpTrack
protected

Definition at line 120 of file G4ITModelProcessor.hh.

Referenced by CleanProcessor(), G4ITModelProcessor(), GetTrack(), and SetTrack().

◆ fpTrackContainer

G4ITTrackHolder* G4ITModelProcessor::fpTrackContainer
protected

Definition at line 115 of file G4ITModelProcessor.hh.

Referenced by ComputeTrackReaction(), G4ITModelProcessor(), and Initialize().

◆ fpTrackingManager

G4ITTrackingManager* G4ITModelProcessor::fpTrackingManager
protected

◆ fReactionInfo

std::vector<std::unique_ptr<G4ITReactionChange> > G4ITModelProcessor::fReactionInfo
protected

Definition at line 126 of file G4ITModelProcessor.hh.

Referenced by ComputeTrackReaction().

◆ fReactionSet

G4ITReactionSet* G4ITModelProcessor::fReactionSet
protected

◆ fTSTimeStep

G4double G4ITModelProcessor::fTSTimeStep
protected

Definition at line 112 of file G4ITModelProcessor.hh.

Referenced by CalculateMinTimeStep(), and G4ITModelProcessor().

◆ fUserMinTimeStep

G4double G4ITModelProcessor::fUserMinTimeStep
protected

Definition at line 121 of file G4ITModelProcessor.hh.

Referenced by G4ITModelProcessor().


The documentation for this class was generated from the following files: