source: trunk/source/run/include/G4AdjointSimManager.hh @ 1202

Last change on this file since 1202 was 1197, checked in by garnier, 15 years ago

nvx fichiers dans CVS

File size: 15.1 KB
Line 
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,v 1.2 2009/11/18 18:02:06 gcosmo Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
28//
29/////////////////////////////////////////////////////////////////////////////////
30//      Class Name:     G4AdjointSimManager.hh
31//      Author:         L. Desorgher
32//      Organisation:   SpaceIT GmbH
33//      Contract:       ESA contract 21435/08/NL/AT
34//      Customer:       ESA/ESTEC
35/////////////////////////////////////////////////////////////////////////////////
36//
37// CHANGE HISTORY
38// --------------
39//      ChangeHistory:
40//              -15-01-2007 creation by L. Desorgher
41//              -March 2008 Redesigned as a non RunManager.  L. Desorgher
42//              -01-11-2009 Add the possibility to use user defined run, event, tracking, stepping,
43//                              and stacking actions during the adjoint tracking phase. L. Desorgher
44//             
45//                             
46//
47//-------------------------------------------------------------
48//      Documentation:
49//              This class represents the Manager of an adjoint/reverse MC simulation.
50//              An adjoint run is divided in a serie of alternative adjoint and forward tracking
51//              of adjoint and normal particles.
52//
53//              Reverse tracking phase:
54//              -----------------------
55//              An adjoint particle of a given type  (adjoint_e-, adjoint_gamma,...) is first generated on the so called adjoint source
56//              with a random energy (1/E distribution) and direction. The adjoint source is the
57//              external surface of a user defined volume or of a user defined sphere. The adjoint
58//              source should contain one or several sensitive volumes and should be small
59//              compared to the entire geometry.
60//              The user can set the min and max energy of the adjoint source. After its
61//              generation the adjoint primary  particle is tracked
62//              bacward in the geometry till a user  defined external surface (spherical or boundary of a volume)
63//              or is killed before if it reaches a user defined  upper energy limit that represents
64//              the maximum energy of the external source. During the reverse tracking, reverse
65//              processes take place where  the adjoint particle being tracked can be either scattered
66//              or transformed in another type of adjoint paticle. During the reverse tracking the
67//              G4SimulationManager replaces the user defined Primary, Run, ... actions, by its own actions.
68//             
69//              Forward tracking phase
70//              -----------------------
71//              When an adjoint particle reaches the external surface its weight,type,  position,
72//              and directions are registered and a  normal primary particle with a type equivalent to the last generated primary adjoint is
73//              generated with the same energy, position but opposite direction  and is  tracked normally in the sensitive region as in a fwd MC simulation.
74//              During this forward tracking phase the 
75//              event, stacking, stepping, tracking actions defined by the user for its general fwd application are used. By this clear separation between
76//              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
77//              simulation. Indeed  the computation of the signal is done by the same actions or classes that the one used in the fwd simulation mode. 
78//
79//              Modification to brought in a existing G4 application to use the ReverseMC method
80//              -------------------------------
81//              In order to be able to use the ReverseMC method in his simulation, the user should modify its code as such:
82//                      1) Adapt its physics list to use ReverseProcesses for adjoint particles. An example of such physics list is provided in an extended
83//                        example.
84//                      2) Create an instance of    G4AdjointSimManager somewhere in the main code.
85//                      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
86//                       that reaches the external surface. This is done by using the following method of G4AdjointSimManager.
87//                                       
88//                                       G4int GetIDOfLastAdjParticleReachingExtSource()                                     
89//                                       G4ThreeVector GetPositionAtEndOfLastAdjointTrack(){ return last_pos;}
90//                                       G4ThreeVector GetDirectionAtEndOfLastAdjointTrack(){ return last_direction;}
91//                                       G4double GetEkinAtEndOfLastAdjointTrack(){ return last_ekin;}
92//                                       G4double GetEkinNucAtEndOfLastAdjointTrack(){ return last_ekin_nuc;}
93//                                       G4double GetWeightAtEndOfLastAdjointTrack(){return last_weight;}
94//                                       G4double GetCosthAtEndOfLastAdjointTrack(){return last_cos_th;}
95//                                       G4String GetFwdParticleNameAtEndOfLastAdjointTrack(){return last_fwd_part_name;}
96//                                       G4int GetFwdParticlePDGEncodingAtEndOfLastAdjointTrack(){return last_fwd_part_PDGEncoding;}
97//                                       G4int GetFwdParticleIndexAtEndOfLastAdjointTrack().
98//
99//                        In orther to have a code working for both forward and adjoint simulation mode, the extra code needed in user actions for the adjoint
100//                        simulation mode can be seperated to the code needed only for the normal forward  simulation by using the following method
101//                                   
102//                                      G4bool GetAdjointSimMode() that return true if   an adjoint simulation is running and false if not!
103//                               
104//                        Example of modification in the analysis part of the code:
105//                       -------------------------------------------------------------
106//                              Let say that in the forward simulation a G4 application computes the energy deposited in a volume.
107//                              The user wants to normalise its results for an external isotropic source of e- with differential spectrum given by f(E).
108//                              A possible modification of the code where the deposited energy Edep during an event  is registered would be the following 
109//                                     
110//                              G4AdjointSimManager* theAdjSimManager =    G4AdjointSimManager::GetInstance(); 
111//                              if (theAdjSimManager->GetAdjointSimMode()) {
112//                                      //code of the user that should be consider only for forwrad simulation
113//                                      G4double normalised_edep = 0.;
114//                                      if (theAdjSimManager->GetFwdParticleNameAtEndOfLastAdjointTrack() == "e-"){
115//                                              G4double ekin_prim = theAdjSimManager->GetEkinAtEndOfLastAdjointTrack();
116//                                              G4double weight_prim = theAdjSimManager->GetWeightAtEndOfLastAdjointTrack();
117//                                              normalised_edep = weight_prim*f(ekin_prim);
118//                                      }
119//                                      //then follow the code where normalised_edep is printed, or registered or whatever ....
120//                              }       
121//                             
122//                              else { //code of the user that should be consider only for forward simulation
123//                              }
124//                              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
125//                              adapted for a normalisatin to several spectra and several type of primary particles  in the same simulation.
126//
127
128#ifndef G4AdjointSimManager_h
129#define G4AdjointSimManager_h 1
130#include "globals.hh"
131#include "G4ThreeVector.hh"
132#include <vector>
133
134
135class G4UserEventAction;
136class G4VUserPrimaryGeneratorAction;
137class G4UserTrackingAction;
138class G4UserSteppingAction;
139class G4UserStackingAction;
140class G4UserRunAction;
141class G4AdjointRunAction;
142class G4AdjointPrimaryGeneratorAction;
143class G4AdjointSteppingAction;
144class G4AdjointEventAction;
145class G4AdjointStackingAction;
146class G4ParticleDefinition;
147class G4AdjointSimMessenger;
148class G4PhysicsLogVector;
149
150class G4AdjointSimManager 
151{
152  public:
153   
154    static G4AdjointSimManager* GetInstance();
155
156  public: //publich methods
157   
158    void RunAdjointSimulation(G4int nb_evt);
159   
160    inline G4int GetNbEvtOfLastRun(){return nb_evt_of_last_run;}
161   
162    void SetAdjointTrackingMode(G4bool aBool);
163    inline G4bool GetAdjointTrackingMode(){return  adjoint_tracking_mode;} //true if an adjoint track is being processed
164    inline G4bool GetAdjointSimMode(){return  adjoint_sim_mode;} //true if an adjoint simulation is running
165
166    G4bool GetDidAdjParticleReachTheExtSource();
167    void RegisterAtEndOfAdjointTrack();
168    void RegisterAdjointPrimaryWeight(G4double aWeight);
169
170    inline G4int GetIDOfLastAdjParticleReachingExtSource(){return ID_of_last_particle_that_reach_the_ext_source;};                                   
171    inline G4ThreeVector GetPositionAtEndOfLastAdjointTrack(){ return last_pos;}
172    inline G4ThreeVector GetDirectionAtEndOfLastAdjointTrack(){ return last_direction;}
173    inline G4double GetEkinAtEndOfLastAdjointTrack(){ return last_ekin;} 
174    inline G4double GetEkinNucAtEndOfLastAdjointTrack(){ return last_ekin_nuc;}
175    inline G4double GetWeightAtEndOfLastAdjointTrack(){return last_weight;}
176    inline G4double GetCosthAtEndOfLastAdjointTrack(){return last_cos_th;}
177    inline const G4String& GetFwdParticleNameAtEndOfLastAdjointTrack(){return last_fwd_part_name;}
178    inline G4int GetFwdParticlePDGEncodingAtEndOfLastAdjointTrack(){return last_fwd_part_PDGEncoding;}
179    inline G4int GetFwdParticleIndexAtEndOfLastAdjointTrack(){return last_fwd_part_index;}
180   
181    std::vector<G4ParticleDefinition*> GetListOfPrimaryFwdParticles();
182     
183    G4bool DefineSphericalExtSource(G4double radius, G4ThreeVector pos);
184    G4bool DefineSphericalExtSourceWithCentreAtTheCentreOfAVolume(G4double radius, const G4String& volume_name);
185    G4bool DefineExtSourceOnTheExtSurfaceOfAVolume(const G4String& volume_name);
186    void SetExtSourceEmax(G4double Emax);
187   
188    //Definition of adjoint source
189    //----------------------------
190     
191    G4bool DefineSphericalAdjointSource(G4double radius, G4ThreeVector pos);
192    G4bool DefineSphericalAdjointSourceWithCentreAtTheCentreOfAVolume(G4double radius, const G4String& volume_name);
193    G4bool DefineAdjointSourceOnTheExtSurfaceOfAVolume(const G4String& volume_name);
194    void SetAdjointSourceEmin(G4double Emin);
195    void SetAdjointSourceEmax(G4double Emax);
196    inline G4double GetAdjointSourceArea(){return area_of_the_adjoint_source;} 
197    void ConsiderParticleAsPrimary(const G4String& particle_name);
198    void NeglectParticleAsPrimary(const G4String& particle_name);
199    void SetPrimaryIon(G4ParticleDefinition* adjointIon, G4ParticleDefinition* fwdIon);
200    const G4String& GetPrimaryIonName();
201   
202    inline void SetNormalisationMode(G4int n){normalisation_mode=n;};
203    G4int GetNormalisationMode(){return normalisation_mode;};
204    G4double GetNumberNucleonsInIon(){return nb_nuc;};
205
206    //Definition of user actions for the adjoint tracking phase
207    //----------------------------
208    void SetAdjointEventAction(G4UserEventAction* anAction);
209    void SetAdjointSteppingAction(G4UserSteppingAction* anAction);
210    void SetAdjointStackingAction(G4UserStackingAction* anAction);
211    void SetAdjointTrackingAction(G4UserTrackingAction* anAction);
212    void SetAdjointRunAction(G4UserRunAction* anAction); 
213   
214    //Set methods for user run actions
215    //--------------------------------
216    inline void UseUserStackingActionInFwdTrackingPhase(G4bool aBool){use_user_StackingAction=aBool;}
217
218    //Convergence test
219    //-----------------------
220   /*
221    void RegisterSignalForConvergenceTest(G4double aSignal);
222    void DefineExponentialPrimarySpectrumForConvergenceTest(G4ParticleDefinition* aPartDef, G4double E0);
223    void DefinePowerLawPrimarySpectrumForConvergenceTest(G4ParticleDefinition* aPartDef, G4double alpha);
224     
225   */ 
226
227  private: 
228 
229    static G4AdjointSimManager* instance;
230 
231  private: // methods
232   
233    void SetRestOfAdjointActions(); 
234    void SetAdjointPrimaryRunAndStackingActions();
235    void ResetRestOfUserActions(); 
236    void ResetUserPrimaryRunAndStackingActions(); 
237    void DefineUserActions();
238
239  private: //constructor and destructor
240 
241    G4AdjointSimManager();
242   ~G4AdjointSimManager();
243
244  private ://attributes
245 
246  //Messenger
247  //----------
248    G4AdjointSimMessenger* theMessenger;
249 
250  //user defined actions for the normal fwd simulation. Taken from the G4RunManager
251  //-------------------------------------------------
252    bool user_action_already_defined;
253    G4UserRunAction*            fUserRunAction;
254    G4UserEventAction*          fUserEventAction;
255    G4VUserPrimaryGeneratorAction*      fUserPrimaryGeneratorAction;
256    G4UserTrackingAction*       fUserTrackingAction;
257    G4UserSteppingAction*       fUserSteppingAction;
258    G4UserStackingAction*       fUserStackingAction;
259    bool use_user_StackingAction; //only for fwd part of the adjoint simulation
260   
261  //action for adjoint simulation
262  //-----------------------------
263    G4UserRunAction*            theAdjointRunAction;
264    G4UserEventAction*          theAdjointEventAction;
265    G4AdjointPrimaryGeneratorAction* theAdjointPrimaryGeneratorAction;
266    G4UserTrackingAction*       theAdjointTrackingAction;
267    G4AdjointSteppingAction*    theAdjointSteppingAction;
268    G4AdjointStackingAction*    theAdjointStackingAction;
269     
270  //adjoint mode
271  //-------------
272    G4bool adjoint_tracking_mode;
273    G4bool adjoint_sim_mode;   
274
275  //adjoint particle information on the external surface
276  //-----------------------------
277    G4ThreeVector last_pos;
278    G4ThreeVector last_direction;
279    G4double last_ekin,last_ekin_nuc; //last_ekin_nuc=last_ekin/nuc, nuc is 1 if not a nucleus
280    G4double last_cos_th;
281    G4String last_fwd_part_name;
282    G4int  last_fwd_part_PDGEncoding;
283    G4int  last_fwd_part_index;
284    G4double last_weight;
285    G4int ID_of_last_particle_that_reach_the_ext_source;
286     
287    G4int nb_evt_of_last_run;
288    G4int normalisation_mode;
289     
290  //Adjoint source
291  //--------------
292    G4double area_of_the_adjoint_source; 
293    G4double nb_nuc;
294    G4double theAdjointPrimaryWeight;
295
296    //Weight Analysis
297    //----------
298    G4PhysicsLogVector* electron_last_weight_vector;
299    G4PhysicsLogVector* proton_last_weight_vector;
300    G4PhysicsLogVector* gamma_last_weight_vector;
301   
302    G4bool welcome_message;
303   
304/* For the future   
305    //Convergence test
306    //----------------
307     
308     G4double normalised_signal;
309     G4double error_signal;
310     G4bool convergence_test_is_used;
311     G4bool power_law_spectrum_for_convergence_test; // true  PowerLaw, ;
312     G4ParticleDefinition* the_par_def_for_convergence_test;
313*/     
314 
315};
316
317#endif
318
Note: See TracBrowser for help on using the repository browser.