Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4ITModelProcessor.cc
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// $Id: G4ITModelProcessor.cc 64057 2012-10-30 15:04:49Z gcosmo $
27//
28// Author: Mathieu Karamitros (kara (AT) cenbg . in2p3 . fr)
29//
30// History:
31// -----------
32// 10 Oct 2011 M.Karamitros created
33//
34// -------------------------------------------------------------------
35
36#include "G4ITModelProcessor.hh"
37#include "G4VITTimeStepper.hh"
39
40std::map<const G4Track*, G4bool> G4ITModelProcessor::fHasReacted ;
41
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}
58
60{
61 //dtor
62// if(fpModelHandler) delete fpModelHandler; deleted by G4ITStepManager
63 fCurrentModel.clear();
64 fReactionInfo.clear();
65}
66
67// Should not be used
69{
70 //copy ctorr
71 fpTrack = 0;
73 fpModel = 0;
74 fInitialized = false;
76 fUserMinTimeStep = -1.;
77}
78
79// Should not be used
81{
82 if (this == &rhs) return *this; // handle self assignment
83 //assignment operator
84 return *this;
85}
86//______________________________________________________________________________
88{
90 fInitialized = true;
91}
92
93//______________________________________________________________________________
94void G4ITModelProcessor::InitializeStepper(const G4double& currentGlobalTime,
95 const G4double& userMinTime)
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}
154
155//______________________________________________________________________________
156void G4ITModelProcessor::CalculateTimeStep(const G4Track* track, const G4double userMinTimeStep)
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}
172
173//______________________________________________________________________________
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}
191
192//______________________________________________________________________________
193void G4ITModelProcessor::FindReaction(std::map<G4Track*, G4TrackVectorHandle>* tracks,
194 const double currentStepTime,
195 const double previousStepTime,
196 const bool reachedUserStepTimeLimit)
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}
@ FatalErrorInArgument
G4IT * GetIT(const G4Track *track)
Definition: G4IT.cc:48
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
const std::vector< std::vector< G4ITModelManager * > > * GetAllModelManager()
static std::map< const G4Track *, G4bool > fHasReacted
void CalculateTimeStep(const G4Track *, const G4double)
G4ITModelProcessor & operator=(const G4ITModelProcessor &other)
void FindReaction(std::map< G4Track *, G4TrackVectorHandle > *, const double currentStepTime, const double previousStepTime, const bool reachedUserStepTimeLimit)
void InitializeStepper(const G4double &currentGlobalTime, const G4double &userMinTime)
const G4Track * fpTrack
G4ITModelManager * fpModelManager
void SetTrack(const G4Track *)
std::vector< G4ITReactionChange * > fReactionInfo
std::vector< std::vector< G4VITModel * > > fCurrentModel
G4ITModelHandler * fpModelHandler
Definition: G4IT.hh:83
virtual const G4ITType GetITType() const =0
G4TrackStatus GetTrackStatus() const
G4VITTimeStepper * GetTimeStepper()
Definition: G4VITModel.hh:123
virtual G4ITReactionChange * MakeReaction(const G4Track &, const G4Track &)=0
static void SetTimes(const G4double &, const G4double &)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
static size_t size()
Definition: G4ITType.cc:41