Geant4  10.01.p02
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 90769 2015-06-09 10:33:41Z gcosmo $
27 //
28 // Author: Mathieu Karamitros, kara@cenbg.in2p3.fr
29 
30 // The code is developed in the framework of the ESA AO7146
31 //
32 // We would be very happy hearing from you, send us your feedback! :)
33 //
34 // In order for Geant4-DNA to be maintained and still open-source,
35 // article citations are crucial.
36 // If you use Geant4-DNA chemistry and you publish papers about your software,
37 // in addition to the general paper on Geant4-DNA:
38 //
39 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178
40 //
41 // we would be very happy if you could please also cite the following
42 // reference papers on chemistry:
43 //
44 // J. Comput. Phys. 274 (2014) 841-882
45 // Prog. Nucl. Sci. Tec. 2 (2011) 503-508
46 
47 #ifndef G4ITMODELPROCESSOR_H
48 #define G4ITMODELPROCESSOR_H
49 
50 #include <vector>
51 #include "G4ITReactionChange.hh"
52 #include "G4ITType.hh"
53 #include "G4ITModelHandler.hh"
54 
57 class G4ITModelHandler;
58 class G4ITReactionSet;
59 
60 //#ifndef compTrackPerID__
61 //#define compTrackPerID__
62 // struct compTrackPerID
63 // {
64 // bool operator()(G4Track* rhs, G4Track* lhs) const
65 // {
66 // return rhs->GetTrackID() < lhs->GetTrackID();
67 // }
68 // };
69 //#endif
70 
80 {
81 public:
83  virtual ~G4ITModelProcessor();
84 
85  inline void SetModelHandler(G4ITModelHandler*);
86  void Initialize();
87 
92  inline void CleanProcessor();
93 
94  //____________________________________________________________
95  // Time stepper part
96  void InitializeStepper(const G4double& currentGlobalTime,
97  const G4double& userMinTime);
98 protected:
99 
100  inline void SetTrack(const G4Track*);
101 
102 public:
103  void CalculateTimeStep(const G4Track*, const G4double);
104  void DoCalculateStep();
105 
106  //____________________________________________________________
107  // Reaction process part
108  void FindReaction(G4ITReactionSet* reactionSet,
109  const double currentStepTime,
110  const double previousStepTime,
111  const bool reachedUserStepTimeLimit);
112 
113  //____________________________________________________________
114  // Get results
115  inline const std::vector<std::vector<G4VITStepModel*> >* GetCurrentModel();
116 
117  inline std::vector<G4ITReactionChange*>* GetReactionInfo()
118  {
119  return &fReactionInfo;
120  }
121 
122  const G4Track* GetTrack() const
123  {
124  return fpTrack;
125  }
126 
127 protected:
137 
138  //_____________________________
139  // Members
142 
143  const G4Track* fpTrack;
145 
146  // Attributes for interaction between many IT types
147  // eg : electron/proton
148  std::vector<std::vector<G4VITStepModel*> > fCurrentModel;
149 
150  // Attributes for interaction between one type of IT
151  // eg : molecule/molecule or electron/electron
154 // G4double fNextTimeChangeModel ;
155 
158 
159  // Atribute for reactions
160  std::vector<G4ITReactionChange*> fReactionInfo;
161  static G4ThreadLocal std::map<const G4Track*, G4bool> *fHasReacted;
162 };
163 
165 // Inline methods
167 
168 inline void G4ITModelProcessor::SetTrack(const G4Track* track)
169 {
170  fpTrack = track;
171 }
172 
173 inline const std::vector<std::vector<G4VITStepModel*> >* G4ITModelProcessor::GetCurrentModel()
174 {
175  return &fCurrentModel;
176 }
177 
179 {
180  if (fInitialized == 1)
181  {
182  G4ExceptionDescription exceptionDescription;
183  exceptionDescription
184  << "You are trying to set a new model while the model processor has alreaday be initialized";
185  G4Exception("G4ITModelProcessor::SetModelHandler", "ITModelProcessor001",
186  FatalErrorInArgument, exceptionDescription);
187  }
188  fpModelHandler = modelHandler;
189 }
190 
192 {
193  fpTrack = 0;
194 }
195 #endif // G4ITMODELPROCESSOR_H
G4ITModelHandler holds for two IT types the corresponding model manager.
static G4ThreadLocal std::map< const G4Track *, G4bool > * fHasReacted
The G4ITModelProcessor will call the two processes defined in G4VITModel.
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void SetTrack(const G4Track *)
G4VITStepModel * fpModel
std::vector< G4ITReactionChange * > fReactionInfo
Define what to do before stepping and after stepping.
std::vector< std::vector< G4VITStepModel * > > fCurrentModel
void CalculateTimeStep(const G4Track *, const G4double)
G4VITReactionProcess defines the reaction between two G4IT.
Tag the G4IT Should be automatically setup by G4IT using : ITDef(MyIT) and ITImp(MyIT) ...
Definition: G4ITType.hh:60
#define G4ThreadLocal
Definition: tls.hh:89
void FindReaction(G4ITReactionSet *reactionSet, const double currentStepTime, const double previousStepTime, const bool reachedUserStepTimeLimit)
bool G4bool
Definition: G4Types.hh:79
const std::vector< std::vector< G4VITStepModel * > > * GetCurrentModel()
const G4Track * GetTrack() const
std::vector< G4ITReactionChange * > * GetReactionInfo()
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4ITModelProcessor & operator=(const G4ITModelProcessor &other)
Assignment operator.
G4ITModelManager * fpModelManager
void InitializeStepper(const G4double &currentGlobalTime, const G4double &userMinTime)
const G4Track * fpTrack
void SetModelHandler(G4ITModelHandler *)
double G4double
Definition: G4Types.hh:76
Before stepping all tracks G4Scheduler calls all the G4VITModel which may contain a G4VITTimeStepper ...
G4ITModelHandler * fpModelHandler
G4ITModelManager chooses which model to use according to the global simulation time.
void CleanProcessor()
Restaure original state of the modelProcessor.