Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes
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 G4ThreadLocal std::map
< const G4Track *, G4bool > * 
fHasReacted = 0
 

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::G4ITModelProcessor ( )

Definition at line 42 of file G4ITModelProcessor.cc.

References fCurrentModel, fHasReacted, fInitialized, fpModel, fpModelHandler, fpModelManager, fpTrack, fUserMinTimeStep, int(), and G4ITType::size().

43 {
44  //ctor
45  if (!fHasReacted) fHasReacted = new std::map<const G4Track*, G4bool>;
46  fpTrack = 0;
47  fpModelHandler = 0;
48  fpModel = 0;
49  fInitialized = false;
50  fpModelManager = 0;
51  fCurrentModel.assign(G4ITType::size(), std::vector<G4VITModel*>());
52 
53  for(int i = 0 ; i < (int) G4ITType::size() ; i++)
54  {
55  fCurrentModel[i].assign(G4ITType::size(),0);
56  }
57  fUserMinTimeStep = -1.;
58 }
static G4ThreadLocal std::map< const G4Track *, G4bool > * fHasReacted
typedef int(XMLCALL *XML_NotStandaloneHandler)(void *userData)
static size_t size()
Definition: G4ITType.cc:46
G4ITModelManager * fpModelManager
const G4Track * fpTrack
std::vector< std::vector< G4VITModel * > > fCurrentModel
G4ITModelHandler * fpModelHandler
G4ITModelProcessor::~G4ITModelProcessor ( )
virtual

Definition at line 60 of file G4ITModelProcessor.cc.

References fCurrentModel, and fReactionInfo.

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

Copy constructor

Parameters
otherObject to copy from

Definition at line 69 of file G4ITModelProcessor.cc.

References fInitialized, fpModel, fpModelHandler, fpModelManager, fpTrack, and fUserMinTimeStep.

70 {
71  //copy ctorr
72  fpTrack = 0;
73  fpModelHandler = 0;
74  fpModel = 0;
75  fInitialized = false;
76  fpModelManager = 0;
77  fUserMinTimeStep = -1.;
78 }
G4ITModelManager * fpModelManager
const G4Track * fpTrack
G4ITModelHandler * fpModelHandler

Member Function Documentation

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

Definition at line 157 of file G4ITModelProcessor.cc.

References CleanProcessor(), DoCalculateStep(), FatalErrorInArgument, fUserMinTimeStep, G4Exception(), and SetTrack().

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

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

Definition at line 174 of file G4ITModelProcessor.hh.

References fpTrack.

Referenced by CalculateTimeStep().

175 {
176  fpTrack = 0;
177 }
const G4Track * fpTrack
void G4ITModelProcessor::DoCalculateStep ( )

Definition at line 175 of file G4ITModelProcessor.cc.

References fCurrentModel, fpModel, fpTrack, fUserMinTimeStep, GetIT(), G4IT::GetITType(), and int().

Referenced by CalculateTimeStep().

176 {
177  if(fpModel) // ie only one model has been declared and will be used
178  {
179  fpModel -> GetTimeStepper()->CalculateStep(*fpTrack, fUserMinTimeStep);
180  }
181  else // ie many models have been declared and will be used
182  {
183  std::vector<G4VITModel*>& model = fCurrentModel[GetIT(fpTrack)->GetITType()];
184 
185  for(int i =0 ; i < (int) model.size() ; i++)
186  {
187  if(model[i] == 0) continue;
188  model[i]->GetTimeStepper()->CalculateStep(*fpTrack, fUserMinTimeStep);
189  }
190  }
191 }
typedef int(XMLCALL *XML_NotStandaloneHandler)(void *userData)
G4IT * GetIT(const G4Track *track)
Definition: G4IT.cc:48
const XML_Char XML_Content * model
virtual const G4ITType GetITType() const =0
const G4Track * fpTrack
std::vector< std::vector< G4VITModel * > > fCurrentModel
void G4ITModelProcessor::FindReaction ( std::map< G4Track *, G4TrackVectorHandle > *  tracks,
const double  currentStepTime,
const double  previousStepTime,
const bool  reachedUserStepTimeLimit 
)

Get track A

Definition at line 194 of file G4ITModelProcessor.cc.

References FatalErrorInArgument, fCurrentModel, fHasReacted, fpModelHandler, fReactionInfo, fStopAndKill, G4Exception(), G4ITModelHandler::GetAllModelManager(), GetIT(), G4Track::GetTrackStatus(), G4VITReactionProcess::MakeReaction(), and G4VITReactionProcess::ResetChanges().

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

Definition at line 157 of file G4ITModelProcessor.hh.

References fCurrentModel.

158 {
159  return &fCurrentModel ;
160 }
std::vector< std::vector< G4VITModel * > > fCurrentModel
std::vector<G4ITReactionChange*>* G4ITModelProcessor::GetReactionInfo ( )
inline

Definition at line 100 of file G4ITModelProcessor.hh.

References fReactionInfo.

101  {
102  return &fReactionInfo;
103  }
std::vector< G4ITReactionChange * > fReactionInfo
const G4Track* G4ITModelProcessor::GetTrack ( ) const
inline

Definition at line 105 of file G4ITModelProcessor.hh.

References fpTrack.

106  {
107  return fpTrack;
108  }
const G4Track * fpTrack
void G4ITModelProcessor::Initialize ( )

Definition at line 88 of file G4ITModelProcessor.cc.

References fInitialized, fpModelHandler, and G4ITModelHandler::Initialize().

89 {
91  fInitialized = true;
92 }
G4ITModelHandler * fpModelHandler
void G4ITModelProcessor::InitializeStepper ( const G4double currentGlobalTime,
const G4double userMinTime 
)

Definition at line 95 of file G4ITModelProcessor.cc.

References FatalErrorInArgument, fCurrentModel, fpModel, fpModelHandler, fpModelManager, G4Exception(), G4ITModelHandler::GetAllModelManager(), G4VITModel::GetTimeStepper(), and G4VITTimeStepper::SetTimes().

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

Assignment operator

Parameters
otherObject to assign from
Returns
A reference to this

Definition at line 81 of file G4ITModelProcessor.cc.

82 {
83  if (this == &rhs) return *this; // handle self assignment
84  //assignment operator
85  return *this;
86 }
void G4ITModelProcessor::SetModelHandler ( G4ITModelHandler modelHandler)
inline

Definition at line 162 of file G4ITModelProcessor.hh.

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

163 {
164  if(fInitialized == 1)
165  {
166  G4ExceptionDescription exceptionDescription ;
167  exceptionDescription << "You are trying to set a new model while the model processor has alreaday be initialized";
168  G4Exception("G4ITModelProcessor::SetModelHandler","ITModelProcessor001",
169  FatalErrorInArgument,exceptionDescription);
170  }
171  fpModelHandler = modelHandler;
172 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4ITModelHandler * fpModelHandler
void G4ITModelProcessor::SetTrack ( const G4Track track)
inlineprotected

Definition at line 152 of file G4ITModelProcessor.hh.

References fpTrack.

Referenced by CalculateTimeStep().

153 {
154  fpTrack = track;
155 }
const G4Track * fpTrack

Field Documentation

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

Definition at line 140 of file G4ITModelProcessor.hh.

G4ITType G4ITModelProcessor::fCurrentType2
protected

Definition at line 141 of file G4ITModelProcessor.hh.

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

Definition at line 145 of file G4ITModelProcessor.hh.

Referenced by FindReaction(), and G4ITModelProcessor().

G4bool G4ITModelProcessor::fInitialized
protected

Definition at line 124 of file G4ITModelProcessor.hh.

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

G4VITModel* G4ITModelProcessor::fpModel
protected

Definition at line 136 of file G4ITModelProcessor.hh.

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

G4ITModelHandler* G4ITModelProcessor::fpModelHandler
protected
G4ITModelManager* G4ITModelProcessor::fpModelManager
protected

Definition at line 137 of file G4ITModelProcessor.hh.

Referenced by G4ITModelProcessor(), and InitializeStepper().

const G4Track* G4ITModelProcessor::fpTrack
protected
std::vector<G4ITReactionChange*> G4ITModelProcessor::fReactionInfo
protected

Definition at line 144 of file G4ITModelProcessor.hh.

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

G4double G4ITModelProcessor::fUserMinTimeStep
protected

Definition at line 128 of file G4ITModelProcessor.hh.

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


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