Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4AdjointSimManager.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: G4AdjointSimManager.hh 102566 2017-02-09 08:35:24Z gcosmo $
27 //
29 // Class Name: G4AdjointSimManager.hh
30 // Author: L. Desorgher
31 // Organisation: SpaceIT GmbH
32 // Contract: ESA contract 21435/08/NL/AT
33 // Customer: ESA/ESTEC
35 //
36 // CHANGE HISTORY
37 // --------------
38 // ChangeHistory:
39 // -15-01-2007 creation by L. Desorgher
40 // -March 2008 Redesigned as a non RunManager. L. Desorgher
41 // -01-11-2009 Add the possibility to use user defined run, event, tracking, stepping,
42 // and stacking actions during the adjoint tracking phase. L. Desorgher
43 //
44 //
45 //
46 //-------------------------------------------------------------
47 // Documentation:
48 // This class represents the Manager of an adjoint/reverse MC simulation.
49 // An adjoint run is divided in a serie of alternative adjoint and forward tracking
50 // of adjoint and normal particles.
51 //
52 // Reverse tracking phase:
53 // -----------------------
54 // An adjoint particle of a given type (adjoint_e-, adjoint_gamma,...) is first generated on the so called adjoint source
55 // with a random energy (1/E distribution) and direction. The adjoint source is the
56 // external surface of a user defined volume or of a user defined sphere. The adjoint
57 // source should contain one or several sensitive volumes and should be small
58 // compared to the entire geometry.
59 // The user can set the min and max energy of the adjoint source. After its
60 // generation the adjoint primary particle is tracked
61 // bacward in the geometry till a user defined external surface (spherical or boundary of a volume)
62 // or is killed before if it reaches a user defined upper energy limit that represents
63 // the maximum energy of the external source. During the reverse tracking, reverse
64 // processes take place where the adjoint particle being tracked can be either scattered
65 // or transformed in another type of adjoint paticle. During the reverse tracking the
66 // G4SimulationManager replaces the user defined Primary, Run, ... actions, by its own actions.
67 //
68 // Forward tracking phase
69 // -----------------------
70 // When an adjoint particle reaches the external surface its weight,type, position,
71 // and directions are registered and a normal primary particle with a type equivalent to the last generated primary adjoint is
72 // generated with the same energy, position but opposite direction and is tracked normally in the sensitive region as in a fwd MC simulation.
73 // During this forward tracking phase the
74 // event, stacking, stepping, tracking actions defined by the user for its general fwd application are used. By this clear separation between
75 // adjoint and fwd tracking phases , the code of the user developed for a fwd simulation should be only slightly modified to adapt it for an adjoint
76 // simulation. Indeed the computation of the signal is done by the same actions or classes that the one used in the fwd simulation mode.
77 //
78 // Modification to brought in a existing G4 application to use the ReverseMC method
79 // -------------------------------
80 // In order to be able to use the ReverseMC method in his simulation, the user should modify its code as such:
81 // 1) Adapt its physics list to use ReverseProcesses for adjoint particles. An example of such physics list is provided in an extended
82 // example.
83 // 2) Create an instance of G4AdjointSimManager somewhere in the main code.
84 // 3) Modify the analysis part of the code to normalise the signal computed during the fwd phase to the weight of the last adjoint particle
85 // that reaches the external surface. This is done by using the following method of G4AdjointSimManager.
86 //
87 // G4int GetIDOfLastAdjParticleReachingExtSource()
88 // G4ThreeVector GetPositionAtEndOfLastAdjointTrack(){ return last_pos;}
89 // G4ThreeVector GetDirectionAtEndOfLastAdjointTrack(){ return last_direction;}
90 // G4double GetEkinAtEndOfLastAdjointTrack(){ return last_ekin;}
91 // G4double GetEkinNucAtEndOfLastAdjointTrack(){ return last_ekin_nuc;}
92 // G4double GetWeightAtEndOfLastAdjointTrack(){return last_weight;}
93 // G4double GetCosthAtEndOfLastAdjointTrack(){return last_cos_th;}
94 // G4String GetFwdParticleNameAtEndOfLastAdjointTrack(){return last_fwd_part_name;}
95 // G4int GetFwdParticlePDGEncodingAtEndOfLastAdjointTrack(){return last_fwd_part_PDGEncoding;}
96 // G4int GetFwdParticleIndexAtEndOfLastAdjointTrack().
97 //
98 // In orther to have a code working for both forward and adjoint simulation mode, the extra code needed in user actions for the adjoint
99 // simulation mode can be seperated to the code needed only for the normal forward simulation by using the following method
100 //
101 // G4bool GetAdjointSimMode() that return true if an adjoint simulation is running and false if not!
102 //
103 // Example of modification in the analysis part of the code:
104 // -------------------------------------------------------------
105 // Let say that in the forward simulation a G4 application computes the energy deposited in a volume.
106 // The user wants to normalise its results for an external isotropic source of e- with differential spectrum given by f(E).
107 // A possible modification of the code where the deposited energy Edep during an event is registered would be the following
108 //
109 // G4AdjointSimManager* theAdjSimManager = G4AdjointSimManager::GetInstance();
110 // if (theAdjSimManager->GetAdjointSimMode()) {
111 // //code of the user that should be consider only for forwrad simulation
112 // G4double normalised_edep = 0.;
113 // if (theAdjSimManager->GetFwdParticleNameAtEndOfLastAdjointTrack() == "e-"){
114 // G4double ekin_prim = theAdjSimManager->GetEkinAtEndOfLastAdjointTrack();
115 // G4double weight_prim = theAdjSimManager->GetWeightAtEndOfLastAdjointTrack();
116 // normalised_edep = weight_prim*f(ekin_prim);
117 // }
118 // //then follow the code where normalised_edep is printed, or registered or whatever ....
119 // }
120 //
121 // else { //code of the user that should be consider only for forward simulation
122 // }
123 // Note that in this example a normalisation to only primary e- with only one spectrum f(E) is considered. The example code could be easily
124 // adapted for a normalisatin to several spectra and several type of primary particles in the same simulation.
125 //
126 
127 #ifndef G4AdjointSimManager_h
128 #define G4AdjointSimManager_h 1
129 #include "globals.hh"
130 #include "G4ThreeVector.hh"
131 #include <vector>
132 #include "G4UserRunAction.hh"
133 
134 class G4UserEventAction;
139 class G4AdjointRunAction;
142 class G4AdjointEventAction;
147 class G4PhysicsLogVector;
148 class G4Run;
149 
151 {
152  public:
153 
155 
156  public: //public methods
157 
158  virtual void BeginOfRunAction(const G4Run* aRun);
159  virtual void EndOfRunAction(const G4Run* aRun);
160  void RunAdjointSimulation(G4int nb_evt);
161 
162  inline G4int GetNbEvtOfLastRun(){return nb_evt_of_last_run;}
163 
164  void SetAdjointTrackingMode(G4bool aBool);
165  G4bool GetAdjointTrackingMode(); //true if an adjoint track is being processed
166  inline G4bool GetAdjointSimMode(){return adjoint_sim_mode;} //true if an adjoint simulation is running
167 
172  //to continue here
173  inline G4int GetIDOfLastAdjParticleReachingExtSource(){return ID_of_last_particle_that_reach_the_ext_source;};
186 
187 
188 
189 
190  std::vector<G4ParticleDefinition*>* GetListOfPrimaryFwdParticles();
192 
197 
198  //Definition of adjoint source
199  //----------------------------
200 
206  inline G4double GetAdjointSourceArea(){return area_of_the_adjoint_source;}
207  void ConsiderParticleAsPrimary(const G4String& particle_name);
208  void NeglectParticleAsPrimary(const G4String& particle_name);
209  void SetPrimaryIon(G4ParticleDefinition* adjointIon, G4ParticleDefinition* fwdIon);
210  const G4String& GetPrimaryIonName();
211 
212  inline void SetNormalisationMode(G4int n){normalisation_mode=n;};
213  G4int GetNormalisationMode(){return normalisation_mode;};
214  G4double GetNumberNucleonsInIon(){return nb_nuc;};
215 
216  //Definition of user actions for the adjoint tracking phase
217  //----------------------------
221  void SetAdjointRunAction(G4UserRunAction* anAction);
222 
223  //Set methods for user run actions
224  //--------------------------------
225  inline void UseUserStackingActionInFwdTrackingPhase(G4bool aBool){use_user_StackingAction=aBool;}
226  inline void UseUserTrackingActionInFwdTrackingPhase(G4bool aBool){use_user_TrackingAction=aBool;}
227 
228 
229  //Set nb of primary fwd gamma
230  //---------------------------
232 
233 
234  //Set nb of adjoint primaries for reverse splitting
235  //-------------------------------------------------
238 
239  //Convergence test
240  //-----------------------
241  /*
242  void RegisterSignalForConvergenceTest(G4double aSignal);
243  void DefineExponentialPrimarySpectrumForConvergenceTest(G4ParticleDefinition* aPartDef, G4double E0);
244  void DefinePowerLawPrimarySpectrumForConvergenceTest(G4ParticleDefinition* aPartDef, G4double alpha);
245 
246  */
247 
248  private:
249 
250  static G4ThreadLocal G4AdjointSimManager* instance;
251 
252 
253  private: // methods
254 
255  void SetRestOfAdjointActions();
256  void SetAdjointPrimaryRunAndStackingActions();
257  void SetAdjointActions();
258  void ResetRestOfUserActions();
259  void ResetUserPrimaryRunAndStackingActions();
260  void ResetUserActions();
261  void DefineUserActions();
262  public:
265 
266 
267  private: //constructor and destructor
268 
271 
272  private ://attributes
273 
274  //Messenger
275  //----------
276  G4AdjointSimMessenger* theMessenger;
277 
278  //user defined actions for the normal fwd simulation. Taken from the G4RunManager
279  //-------------------------------------------------
280  bool user_action_already_defined;
281  G4UserRunAction* fUserRunAction;
282  G4UserEventAction* fUserEventAction;
283  G4VUserPrimaryGeneratorAction* fUserPrimaryGeneratorAction;
284  G4UserTrackingAction* fUserTrackingAction;
285  G4UserSteppingAction* fUserSteppingAction;
286  G4UserStackingAction* fUserStackingAction;
287  bool use_user_StackingAction; //only for fwd part of the adjoint simulation
288  bool use_user_TrackingAction;
289 
290  //action for adjoint simulation
291  //-----------------------------
292  G4UserRunAction* theAdjointRunAction;
293  G4UserEventAction* theAdjointEventAction;
294  G4AdjointPrimaryGeneratorAction* theAdjointPrimaryGeneratorAction;
295  G4AdjointTrackingAction* theAdjointTrackingAction;
296  G4AdjointSteppingAction* theAdjointSteppingAction;
297  G4AdjointStackingAction* theAdjointStackingAction;
298 
299  //adjoint mode
300  //-------------
301  G4bool adjoint_tracking_mode;
302  G4bool adjoint_sim_mode;
303 
304  //adjoint particle information on the external surface
305  //-----------------------------
306  std::vector<G4ThreeVector> last_pos_vec;
307  std::vector<G4ThreeVector> last_direction_vec;
308  std::vector<G4double> last_ekin_vec;
309  std::vector<G4double> last_ekin_nuc_vec;
310  std::vector<G4double> last_cos_th_vec;
311  std::vector<G4double> last_weight_vec;
312  std::vector<G4int> last_fwd_part_PDGEncoding_vec;
313  std::vector<G4int> last_fwd_part_index_vec;
314  std::vector<G4int> ID_of_last_particle_that_reach_the_ext_source_vec;
315 
316 
317 
318 
319  G4ThreeVector last_pos;
320  G4ThreeVector last_direction;
321  G4double last_ekin,last_ekin_nuc; //last_ekin_nuc=last_ekin/nuc, nuc is 1 if not a nucleus
322  G4double last_cos_th;
323  G4String last_fwd_part_name;
324  G4int last_fwd_part_PDGEncoding;
325  G4int last_fwd_part_index;
326  G4double last_weight;
327  G4int ID_of_last_particle_that_reach_the_ext_source;
328 
329  G4int nb_evt_of_last_run;
330  G4int normalisation_mode;
331 
332  //Adjoint source
333  //--------------
334  G4double area_of_the_adjoint_source;
335  G4double nb_nuc;
336  G4double theAdjointPrimaryWeight;
337 
338  //Weight Analysis
339  //----------
340  /*G4PhysicsLogVector* electron_last_weight_vector;
341  G4PhysicsLogVector* proton_last_weight_vector;
342  G4PhysicsLogVector* gamma_last_weight_vector;*/
343 
344  G4bool welcome_message;
345 
346 /* For the future
347  //Convergence test
348  //----------------
349 
350  G4double normalised_signal;
351  G4double error_signal;
352  G4bool convergence_test_is_used;
353  G4bool power_law_spectrum_for_convergence_test; // true PowerLaw, ;
354  G4ParticleDefinition* the_par_def_for_convergence_test;
355 */
356 
357 };
358 
359 #endif
360 
static G4AdjointSimManager * GetInstance()
G4double GetCosthAtEndOfLastAdjointTrack(size_t i=0)
void UseUserTrackingActionInFwdTrackingPhase(G4bool aBool)
G4int GetFwdParticleIndexAtEndOfLastAdjointTrack(size_t i=0)
std::vector< G4ParticleDefinition * > * GetListOfPrimaryFwdParticles()
G4bool GetDidAdjParticleReachTheExtSource()
G4int GetIDOfLastAdjParticleReachingExtSource()
void SetAdjointSourceEmax(G4double Emax)
void SetAdjointStackingAction(G4UserStackingAction *anAction)
G4ThreeVector GetDirectionAtEndOfLastAdjointTrack(size_t i=0)
G4bool DefineSphericalAdjointSourceWithCentreAtTheCentreOfAVolume(G4double radius, const G4String &volume_name)
G4bool DefineSphericalAdjointSource(G4double radius, G4ThreeVector pos)
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
virtual void EndOfRunAction(const G4Run *aRun)
void SetAdjointSourceEmin(G4double Emin)
void SetAdjointTrackingMode(G4bool aBool)
virtual void BeginOfRunAction(const G4Run *aRun)
G4bool DefineAdjointSourceOnTheExtSurfaceOfAVolume(const G4String &volume_name)
bool G4bool
Definition: G4Types.hh:79
void SetNbOfPrimaryFwdGammasPerEvent(G4int)
G4double GetEkinNucAtEndOfLastAdjointTrack(size_t i=0)
Definition: G4Run.hh:46
void RunAdjointSimulation(G4int nb_evt)
const G4int n
G4double GetWeightAtEndOfLastAdjointTrack(size_t i=0)
G4bool DefineSphericalExtSourceWithCentreAtTheCentreOfAVolume(G4double radius, const G4String &volume_name)
const G4String & GetFwdParticleNameAtEndOfLastAdjointTrack()
void SetPrimaryIon(G4ParticleDefinition *adjointIon, G4ParticleDefinition *fwdIon)
const G4String & GetPrimaryIonName()
G4bool DefineSphericalExtSource(G4double radius, G4ThreeVector pos)
void SetAdjointRunAction(G4UserRunAction *anAction)
void SetAdjointEventAction(G4UserEventAction *anAction)
G4ThreeVector GetPositionAtEndOfLastAdjointTrack(size_t i=0)
void ConsiderParticleAsPrimary(const G4String &particle_name)
void SetNormalisationMode(G4int n)
static const G4double Emin
G4int GetFwdParticlePDGEncodingAtEndOfLastAdjointTrack(size_t i=0)
void RegisterAdjointPrimaryWeight(G4double aWeight)
G4bool DefineExtSourceOnTheExtSurfaceOfAVolume(const G4String &volume_name)
static const G4double Emax
void SetNbAdjointPrimaryElectronsPerEvent(G4int)
void SetNbAdjointPrimaryGammasPerEvent(G4int)
G4ParticleDefinition * GetLastGeneratedFwdPrimaryParticle()
void UseUserStackingActionInFwdTrackingPhase(G4bool aBool)
double G4double
Definition: G4Types.hh:76
void SetExtSourceEmax(G4double Emax)
void SetAdjointSteppingAction(G4UserSteppingAction *anAction)
G4double GetEkinAtEndOfLastAdjointTrack(size_t i=0)
static const G4double pos
size_t GetNbOfAdointTracksReachingTheExternalSurface()
void NeglectParticleAsPrimary(const G4String &particle_name)
void ResetDidOneAdjPartReachExtSourceDuringEvent()