Geant4 9.6.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 ()
 
virtual ~G4ITModelProcessor ()
 
void SetModelHandler (G4ITModelHandler *)
 
void Initialize ()
 
void CleanProcessor ()
 
void InitializeStepper (const G4double &currentGlobalTime, const G4double &userMinTime)
 
void CalculateTimeStep (const G4Track *, const G4double)
 
void DoCalculateStep ()
 
void FindReaction (std::map< G4Track *, G4TrackVectorHandle > *, const double currentStepTime, const double previousStepTime, const bool reachedUserStepTimeLimit)
 
const std::vector< std::vector< G4VITModel * > > * GetCurrentModel ()
 
std::vector< G4ITReactionChange * > * GetReactionInfo ()
 
const G4TrackGetTrack () const
 

Protected Member Functions

void SetTrack (const G4Track *)
 
 G4ITModelProcessor (const G4ITModelProcessor &other)
 
G4ITModelProcessoroperator= (const G4ITModelProcessor &other)
 

Protected Attributes

G4bool fInitialized
 
G4ITModelHandlerfpModelHandler
 
const G4TrackfpTrack
 
G4double fUserMinTimeStep
 
std::vector< std::vector< G4VITModel * > > fCurrentModel
 
G4VITModelfpModel
 
G4ITModelManagerfpModelManager
 
G4ITType fCurrentType1
 
G4ITType fCurrentType2
 
std::vector< G4ITReactionChange * > fReactionInfo
 

Static Protected Attributes

static std::map< const G4Track *, G4boolfHasReacted
 

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.

Definition at line 62 of file G4ITModelProcessor.hh.

Constructor & Destructor Documentation

◆ G4ITModelProcessor() [1/2]

G4ITModelProcessor::G4ITModelProcessor ( )

Default constructor

Definition at line 42 of file G4ITModelProcessor.cc.

43{
44 //ctor
45 fpTrack = 0;
47 fpModel = 0;
48 fInitialized = false;
50 fCurrentModel.assign(G4ITType::size(), std::vector<G4VITModel*>());
51
52 for(int i = 0 ; i < (int) G4ITType::size() ; i++)
53 {
54 fCurrentModel[i].assign(G4ITType::size(),0);
55 }
56 fUserMinTimeStep = -1.;
57}
const G4Track * fpTrack
G4ITModelManager * fpModelManager
std::vector< std::vector< G4VITModel * > > fCurrentModel
G4ITModelHandler * fpModelHandler
static size_t size()
Definition: G4ITType.cc:41

◆ ~G4ITModelProcessor()

G4ITModelProcessor::~G4ITModelProcessor ( )
virtual

Default destructor

Definition at line 59 of file G4ITModelProcessor.cc.

60{
61 //dtor
62// if(fpModelHandler) delete fpModelHandler; deleted by G4ITStepManager
63 fCurrentModel.clear();
64 fReactionInfo.clear();
65}
std::vector< G4ITReactionChange * > fReactionInfo

◆ G4ITModelProcessor() [2/2]

G4ITModelProcessor::G4ITModelProcessor ( const G4ITModelProcessor other)
protected

Copy constructor

Parameters
otherObject to copy from

Definition at line 68 of file G4ITModelProcessor.cc.

69{
70 //copy ctorr
71 fpTrack = 0;
73 fpModel = 0;
74 fInitialized = false;
76 fUserMinTimeStep = -1.;
77}

Member Function Documentation

◆ CalculateTimeStep()

void G4ITModelProcessor::CalculateTimeStep ( const G4Track track,
const G4double  userMinTimeStep 
)

Definition at line 156 of file G4ITModelProcessor.cc.

157{
158 // G4cout << "G4ITModelProcessor::CalculateStep" << G4endl;
160 if(track == 0)
161 {
162 G4ExceptionDescription exceptionDescription ;
163 exceptionDescription << "No track was passed to the method (track == 0).";
164 G4Exception("G4ITModelProcessor::CalculateStep","ITModelProcessor004",
165 FatalErrorInArgument,exceptionDescription);
166 }
167 SetTrack(track);
168 fUserMinTimeStep = userMinTimeStep ;
169
171}
@ FatalErrorInArgument
void SetTrack(const G4Track *)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

◆ CleanProcessor()

void G4ITModelProcessor::CleanProcessor ( )
inline

Restaure original state of the modelProcessor. This method should be call only by the ITStepManager

Definition at line 178 of file G4ITModelProcessor.hh.

179{
180 fpTrack = 0;
181}

Referenced by CalculateTimeStep().

◆ DoCalculateStep()

void G4ITModelProcessor::DoCalculateStep ( )

Definition at line 174 of file G4ITModelProcessor.cc.

175{
176 if(fpModel) // ie only one model has been declared and will be used
177 {
178 fpModel -> GetTimeStepper()->CalculateStep(*fpTrack, fUserMinTimeStep);
179 }
180 else // ie many models have been declared and will be used
181 {
182 std::vector<G4VITModel*>& model = fCurrentModel[GetIT(fpTrack)->GetITType()];
183
184 for(int i =0 ; i < (int) model.size() ; i++)
185 {
186 if(model[i] == 0) continue;
187 model[i]->GetTimeStepper()->CalculateStep(*fpTrack, fUserMinTimeStep);
188 }
189 }
190}
G4IT * GetIT(const G4Track *track)
Definition: G4IT.cc:48
virtual const G4ITType GetITType() const =0

Referenced by CalculateTimeStep().

◆ FindReaction()

void G4ITModelProcessor::FindReaction ( std::map< G4Track *, G4TrackVectorHandle > *  tracks,
const double  currentStepTime,
const double  previousStepTime,
const bool  reachedUserStepTimeLimit 
)

Get track A

Definition at line 193 of file G4ITModelProcessor.cc.

197{
198 // DEBUG
199 // G4cout << "G4ITReactionManager::FindReaction" << G4endl;
200 if(tracks == 0) return ;
201
202 if(fpModelHandler->GetAllModelManager()->empty()) return ;
203
204 std::map<G4Track*, G4TrackVectorHandle>::iterator tracks_i = tracks->begin();;
205
206 for(tracks_i = tracks->begin() ; tracks_i != tracks-> end() ; tracks_i ++)
207 {
208 /// Get track A
209 G4Track* trackA = tracks_i->first;
210
211 if(trackA == 0) continue;
212
213 std::map<const G4Track*, G4bool>::iterator it_hasReacted = fHasReacted.find(trackA);
214 if(it_hasReacted != fHasReacted.end()) continue;
215 if(trackA->GetTrackStatus() == fStopAndKill) continue;
216
217 G4IT* ITA = GetIT(trackA);
218 G4ITType ITypeA = ITA -> GetITType();
219
220 const std::vector<G4VITModel*> model = fCurrentModel[ITypeA];
221
222 G4TrackVectorHandle& trackB_vector = tracks_i->second ;
223 std::vector<G4Track*>::iterator trackB_i = trackB_vector->begin();
224
225 G4Track* trackB = 0 ;
226 G4ITType ITypeB(-1);
227 G4VITReactionProcess* process = 0;
228 G4ITReactionChange* changes = 0;
229
230 for(; trackB_i != trackB_vector->end() ; trackB_i++)
231 {
232 trackB = *trackB_i;
233
234 if(trackB == 0) continue;
235 it_hasReacted = fHasReacted.find(trackB);
236 if(it_hasReacted != fHasReacted.end()) continue;
237 if(trackB->GetTrackStatus() == fStopAndKill) continue;
238
239 // DEBUG
240 // G4cout << "Couple : " << trackA->GetParticleDefinition->GetParticleName() << " ("
241 // << trackA->GetTrackID() << ") "
242 // << trackB->GetParticleDefinition->GetParticleName() << " ("
243 // << trackB->GetTrackID() << ")"
244 // << G4endl;
245
246 if(trackB == trackA)
247 {
248 G4ExceptionDescription exceptionDescription ;
249 exceptionDescription << "The IT reaction process sent back a reaction between trackA and trackB. ";
250 exceptionDescription << "The problem is trackA == trackB";
251 G4Exception("G4ITModelProcessor::FindReaction","ITModelProcessor005",
252 FatalErrorInArgument,exceptionDescription);
253 }
254
255 G4IT* ITB = GetIT(trackB);
256 G4ITType ITypeBtmp = ITB -> GetITType();
257
258 if(ITypeB != ITypeBtmp)
259 {
260 ITypeB = ITypeBtmp ;
261
262 if(model[ITypeB])
263 process = model[ITypeB]->GetReactionProcess();
264 }
265
266 if(process && process -> TestReactibility(*trackA, *trackB,
267 currentStepTime, previousStepTime,
268 reachedUserStepTimeLimit))
269 {
270 changes = process->MakeReaction(*trackA, *trackB);
271 }
272
273 if(changes)
274 {
275 fHasReacted[trackA] = true;
276 fHasReacted[trackB] = true;
277 changes -> GetTrackA();
278 changes -> GetTrackB();
279
280 fReactionInfo.push_back(changes);
281
282 process->ResetChanges();
283 changes = 0;
284
285 break;
286 }
287 }
288 }
289
290 fHasReacted.clear();
291}
@ fStopAndKill
const std::vector< std::vector< G4ITModelManager * > > * GetAllModelManager()
static std::map< const G4Track *, G4bool > fHasReacted
Definition: G4IT.hh:83
G4TrackStatus GetTrackStatus() const
virtual G4ITReactionChange * MakeReaction(const G4Track &, const G4Track &)=0

◆ GetCurrentModel()

const std::vector< std::vector< G4VITModel * > > * G4ITModelProcessor::GetCurrentModel ( )
inline

Definition at line 161 of file G4ITModelProcessor.hh.

162{
163 return &fCurrentModel ;
164}

◆ GetReactionInfo()

std::vector< G4ITReactionChange * > * G4ITModelProcessor::GetReactionInfo ( )
inline

Definition at line 104 of file G4ITModelProcessor.hh.

105 {
106 return &fReactionInfo;
107 }

◆ GetTrack()

const G4Track * G4ITModelProcessor::GetTrack ( ) const
inline

Definition at line 109 of file G4ITModelProcessor.hh.

110 {
111 return fpTrack;
112 }

◆ Initialize()

void G4ITModelProcessor::Initialize ( )

Definition at line 87 of file G4ITModelProcessor.cc.

88{
90 fInitialized = true;
91}

◆ InitializeStepper()

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

Definition at line 94 of file G4ITModelProcessor.cc.

96{
97 // G4cout << "G4ITModelProcessor::InitializeStepper" << G4endl;
98 if(fpModelHandler==0)
99 {
100 G4ExceptionDescription exceptionDescription ;
101 exceptionDescription << "No G4ITModelHandler was passed to the modelProcessor.";
102 G4Exception("G4ITModelProcessor::InitializeStepper","ITModelProcessor002",
103 FatalErrorInArgument,exceptionDescription);
104 }
105 const std::vector<std::vector<G4ITModelManager*> >* modelManager = fpModelHandler
107
108 if(modelManager==0)
109 {
110 G4ExceptionDescription exceptionDescription ;
111 exceptionDescription << "No G4ITModelManager was register to G4ITModelHandler.";
112 G4Exception("G4ITModelProcessor::InitializeStepper","ITModelProcessor003",
113 FatalErrorInArgument,exceptionDescription);
114 }
115
116 int nbModels1 = modelManager->size() ;
117
118 G4VITTimeStepper::SetTimes(currentGlobalTime, userMinTime) ;
119
120 // TODO !!!
121 // if( nbModels1 != 1 || (nbModels1 == 1 && !fpModelManager) )
122 {
123 int nbModels2 = -1;
124 G4VITModel* model = 0;
125 G4ITModelManager* modman = 0;
126
127 for(int i = 0 ; i < nbModels1 ; i++)
128 {
129 nbModels2 = (*modelManager)[i].size();
130
131 for(int j = 0 ; j <= i ; j++)
132 {
133 modman = (*modelManager)[i][j];
134
135 if(modman == 0) continue ;
136
137 model = modman -> GetModel(currentGlobalTime);
138 G4VITTimeStepper* stepper = model->GetTimeStepper() ;
139
140// stepper -> PrepareForAllProcessors() ;
141 stepper -> Prepare() ;
142 fCurrentModel[i][j] = model;
143 }
144 }
145
146 if(nbModels1 == 1 && nbModels2 ==1)
147 {
148 fpModelManager = modman;
149 fpModel = model;
150 }
151 else fpModel = 0;
152 }
153}
G4VITTimeStepper * GetTimeStepper()
Definition: G4VITModel.hh:123
static void SetTimes(const G4double &, const G4double &)

◆ operator=()

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

Assignment operator

Parameters
otherObject to assign from
Returns
A reference to this

Definition at line 80 of file G4ITModelProcessor.cc.

81{
82 if (this == &rhs) return *this; // handle self assignment
83 //assignment operator
84 return *this;
85}

◆ SetModelHandler()

void G4ITModelProcessor::SetModelHandler ( G4ITModelHandler modelHandler)
inline

Definition at line 166 of file G4ITModelProcessor.hh.

167{
168 if(fInitialized == 1)
169 {
170 G4ExceptionDescription exceptionDescription ;
171 exceptionDescription << "You are trying to set a new model while the model processor has alreaday be initialized";
172 G4Exception("G4ITModelProcessor::SetModelHandler","ITModelProcessor001",
173 FatalErrorInArgument,exceptionDescription);
174 }
175 fpModelHandler = modelHandler;
176}

◆ SetTrack()

void G4ITModelProcessor::SetTrack ( const G4Track track)
inlineprotected

Definition at line 156 of file G4ITModelProcessor.hh.

157{
158 fpTrack = track;
159}

Referenced by CalculateTimeStep().

Member Data Documentation

◆ fCurrentModel

std::vector<std::vector<G4VITModel*> > G4ITModelProcessor::fCurrentModel
protected

◆ fCurrentType1

G4ITType G4ITModelProcessor::fCurrentType1
protected

Definition at line 144 of file G4ITModelProcessor.hh.

◆ fCurrentType2

G4ITType G4ITModelProcessor::fCurrentType2
protected

Definition at line 145 of file G4ITModelProcessor.hh.

◆ fHasReacted

std::map< const G4Track *, G4bool > G4ITModelProcessor::fHasReacted
staticprotected

Definition at line 149 of file G4ITModelProcessor.hh.

Referenced by FindReaction().

◆ fInitialized

G4bool G4ITModelProcessor::fInitialized
protected

Definition at line 128 of file G4ITModelProcessor.hh.

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

◆ fpModel

G4VITModel* G4ITModelProcessor::fpModel
protected

Definition at line 140 of file G4ITModelProcessor.hh.

Referenced by DoCalculateStep(), G4ITModelProcessor(), and InitializeStepper().

◆ fpModelHandler

G4ITModelHandler* G4ITModelProcessor::fpModelHandler
protected

◆ fpModelManager

G4ITModelManager* G4ITModelProcessor::fpModelManager
protected

Definition at line 141 of file G4ITModelProcessor.hh.

Referenced by G4ITModelProcessor(), and InitializeStepper().

◆ fpTrack

const G4Track* G4ITModelProcessor::fpTrack
protected

◆ fReactionInfo

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

Definition at line 148 of file G4ITModelProcessor.hh.

Referenced by FindReaction(), GetReactionInfo(), and ~G4ITModelProcessor().

◆ fUserMinTimeStep

G4double G4ITModelProcessor::fUserMinTimeStep
protected

Definition at line 132 of file G4ITModelProcessor.hh.

Referenced by CalculateTimeStep(), DoCalculateStep(), and G4ITModelProcessor().


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