Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ITModelProcessor.hh
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.hh 64057 2012-10-30 15:04:49Z gcosmo $
27//
28// Author: Mathieu Karamitros (kara (AT) cenbg . in2p3 . fr)
29//
30// WARNING : This class is released as a prototype.
31// It might strongly evolve or even disapear in the next releases.
32//
33// History:
34// -----------
35// 10 Oct 2011 M.Karamitros created
36//
37// -------------------------------------------------------------------
38
39#ifndef G4ITMODELPROCESSOR_H
40#define G4ITMODELPROCESSOR_H
41
42#include <vector>
43#include "G4ITReactionChange.hh"
44#include "G4ITType.hh"
45#include "G4ITModelHandler.hh"
46
50
52
53 /**
54 * The G4ITModelProcessor will call the two processes defined in G4VITModel.
55 * This processes act at the beginning and end of each step.
56 * The first one, the TimeStepper will calculate a time step to propagate all
57 * the track and eventually it can return some tracks that can likely react
58 * at the end of the step.
59 * The second one, the ReactionProcess will make the tracks reacting.
60 */
61
63{
64public:
65 /** Default constructor */
67 /** Default destructor */
68 virtual ~G4ITModelProcessor();
69
70
72
73 void Initialize();
74
75 /**
76 * Restaure original state of the modelProcessor.
77 * This method should be call only by the ITStepManager
78 */
79 inline void CleanProcessor();
80
81 //____________________________________________________________
82 // Time stepper part
83 void InitializeStepper(const G4double& currentGlobalTime,
84 const G4double& userMinTime);
85protected :
86
87 inline void SetTrack(const G4Track*);
88
89public :
90 void CalculateTimeStep(const G4Track*, const G4double);
91 void DoCalculateStep();
92
93 //____________________________________________________________
94 // Reaction process part
95 void FindReaction(std::map<G4Track*, G4TrackVectorHandle>*,
96 const double currentStepTime,
97 const double previousStepTime,
98 const bool reachedUserStepTimeLimit) ;
99
100 //____________________________________________________________
101 // Get results
102 inline const std::vector<std::vector<G4VITModel*> >* GetCurrentModel();
103
104 inline std::vector<G4ITReactionChange*>* GetReactionInfo()
105 {
106 return &fReactionInfo;
107 }
108
109 const G4Track* GetTrack() const
110 {
111 return fpTrack;
112 }
113
114
115protected:
116 /** Copy constructor
117 * \param other Object to copy from
118 */
120 /** Assignment operator
121 * \param other Object to assign from
122 * \return A reference to this
123 */
125
126 //_____________________________
127 // Members
130
133
134 // Attributes for interaction between many IT types
135 // eg : electron/proton
136 std::vector<std::vector<G4VITModel*> > fCurrentModel;
137
138 // Attributes for interaction between one type of IT
139 // eg : molecule/molecule or electron/electron
142// G4double fNextTimeChangeModel ;
143
146
147 // Atribute for reactions
148 std::vector<G4ITReactionChange*> fReactionInfo ;
149 static std::map<const G4Track*, G4bool> fHasReacted;
150};
151
152///
153// Inline methods
154///
155
156inline void G4ITModelProcessor::SetTrack(const G4Track* track)
157{
158 fpTrack = track;
159}
160
161inline const std::vector<std::vector<G4VITModel*> >* G4ITModelProcessor::GetCurrentModel()
162{
163 return &fCurrentModel ;
164}
165
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}
177
179{
180 fpTrack = 0;
181}
182#endif // G4ITMODELPROCESSOR_H
@ FatalErrorInArgument
G4ReferenceCountedHandle< std::vector< G4Track * > > G4TrackVectorHandle
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
static std::map< const G4Track *, G4bool > fHasReacted
void SetModelHandler(G4ITModelHandler *)
void CalculateTimeStep(const G4Track *, const G4double)
const std::vector< std::vector< G4VITModel * > > * GetCurrentModel()
G4ITModelProcessor & operator=(const G4ITModelProcessor &other)
std::vector< G4ITReactionChange * > * GetReactionInfo()
const G4Track * GetTrack() const
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
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76