Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes
G4ITModelProcessor Class Reference

#include <G4ITModelProcessor.hh>

Public Member Functions

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

Protected Member Functions

void SetTrack (const G4Track *)
 

Protected Attributes

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

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 55 of file G4ITModelProcessor.cc.

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

References DBL_MAX, fComputeReaction, fComputeTimeStep, fInitialized, fpActiveModelWithMinTimeStep, fpModelHandler, fpTrack, fpTrackContainer, fpTrackingManager, fReactionSet, fTSTimeStep, and fUserMinTimeStep.

◆ 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 95 of file G4ITModelProcessor.cc.

97{
98
99#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
100 MemStat mem_first, mem_second, mem_diff;
101 mem_first = MemoryUsage();
102#endif
103
106
107 InitializeStepper(currentGlobalTime, definedMinTimeStep);
108
109#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
110 mem_second = MemoryUsage();
111 mem_diff = mem_second-mem_first;
112 G4cout << "\t || MEM || G4Scheduler::CalculateMinTimeStep || After "
113 "computing fpModelProcessor -> InitializeStepper, diff is : "
114 << mem_diff
115 << G4endl;
116#endif
117
118 for (auto& pStepModel : fActiveModels)
119 {
121 pStepModel->GetTimeStepper()->CalculateMinTimeStep(
122 currentGlobalTime,
123 definedMinTimeStep);
124
125 fpActiveModelWithMinTimeStep = pStepModel;
126
127 if(fTSTimeStep == -1){
129 if(fReactionSet->Empty()) return DBL_MAX;
130 auto fReactionSetInTime = fReactionSet->GetReactionsPerTime();
131 fTSTimeStep = fReactionSetInTime.begin()->get()->GetTime() - currentGlobalTime;
132 }
133 }
134
135#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
136 mem_second = MemoryUsage();
137 mem_diff = mem_second-mem_first;
138 G4cout << "\t || MEM || G4Scheduler::CalculateMinTimeStep || "
139 "After looping on tracks, diff is : " << mem_diff << G4endl;
140#endif
141 return fTSTimeStep;
142}
#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

References DBL_MAX, G4ITReactionSet::Empty(), fActiveModels, fpActiveModelWithMinTimeStep, fReactionSet, fTSTimeStep, G4cout, G4endl, G4VITStepModel::GetReactionProcess(), G4ITReactionSet::GetReactionsPerTime(), G4VITReactionProcess::Initialize(), InitializeStepper(), and G4MemStat::MemoryUsage().

Referenced by G4Scheduler::Stepping().

◆ CleanProcessor()

void G4ITModelProcessor::CleanProcessor ( )

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

Definition at line 375 of file G4ITModelProcessor.cc.

376{
377 fpTrack = nullptr;
378}

References fpTrack.

◆ ComputeTrackReaction()

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

Definition at line 173 of file G4ITModelProcessor.cc.

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

References eCollisionBetweenTracks, G4ITReactionSet::Empty(), G4ITTrackingManager::EndTracking(), FatalErrorInArgument, G4VITReactionProcess::FindReaction(), fpActiveModelWithMinTimeStep, fpTrackContainer, fpTrackingManager, fReactionInfo, fReactionSet, fStopAndKill, G4BestUnit, G4cout, G4endl, G4Exception(), G4Track::GetGlobalTime(), GetIT(), G4IT::GetName(), G4Track::GetParticleDefinition(), G4ParticleDefinition::GetParticleName(), G4VITStepModel::GetReactionProcess(), G4Track::GetStep(), G4Track::GetTrackID(), G4ITTrackHolder::KillTracks(), G4ITTrackHolder::MergeSecondariesWithMainList(), G4IT::SetParentID(), and G4UserTimeStepAction::UserReactionAction().

Referenced by G4Scheduler::Stepping().

◆ GetComputeTimeStep()

bool G4ITModelProcessor::GetComputeTimeStep ( ) const

Definition at line 380 of file G4ITModelProcessor.cc.

381{
382 return fComputeTimeStep;
383}

References fComputeTimeStep.

Referenced by G4Scheduler::Stepping().

◆ GetTrack()

const G4Track * G4ITModelProcessor::GetTrack ( ) const

Definition at line 385 of file G4ITModelProcessor.cc.

386{
387 return fpTrack;
388}

References fpTrack.

◆ Initialize()

void G4ITModelProcessor::Initialize ( )

◆ InitializeStepper()

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

Definition at line 146 of file G4ITModelProcessor.cc.

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

References fActiveModels, fpModelHandler, G4cout, G4endl, G4ITModelHandler::GetActiveModels(), G4MemStat::MemoryUsage(), and G4VITTimeStepComputer::SetTimes().

Referenced by CalculateMinTimeStep().

◆ operator=()

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

◆ RegisterModel()

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

Definition at line 72 of file G4ITModelProcessor.cc.

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

References fpModelHandler, and G4ITModelHandler::RegisterModel().

◆ SetModelHandler()

void G4ITModelProcessor::SetModelHandler ( G4ITModelHandler pModelHandler)

Definition at line 362 of file G4ITModelProcessor.cc.

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

References FatalErrorInArgument, fInitialized, fpModelHandler, and G4Exception().

Referenced by G4Scheduler::Initialize().

◆ SetTrack()

void G4ITModelProcessor::SetTrack ( const G4Track track)
protected

Definition at line 357 of file G4ITModelProcessor.cc.

358{
359 fpTrack = track;
360}

References fpTrack.

◆ SetTrackingManager()

void G4ITModelProcessor::SetTrackingManager ( G4ITTrackingManager trackingManager)

Definition at line 390 of file G4ITModelProcessor.cc.

391{
392 fpTrackingManager = pTrackingManager;
393}

References fpTrackingManager.

Referenced by G4Scheduler::Initialize().

Field 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: