G4AdjointSimManager.hh

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id$
00027 //
00029 //      Class Name:     G4AdjointSimManager.hh
00030 //      Author:         L. Desorgher
00031 //      Organisation:   SpaceIT GmbH
00032 //      Contract:       ESA contract 21435/08/NL/AT
00033 //      Customer:       ESA/ESTEC
00035 //
00036 // CHANGE HISTORY
00037 // --------------
00038 //      ChangeHistory: 
00039 //              -15-01-2007 creation by L. Desorgher
00040 //              -March 2008 Redesigned as a non RunManager.  L. Desorgher
00041 //              -01-11-2009 Add the possibility to use user defined run, event, tracking, stepping, 
00042 //                              and stacking actions during the adjoint tracking phase. L. Desorgher
00043 //              
00044 //                              
00045 //
00046 //-------------------------------------------------------------
00047 //      Documentation:
00048 //              This class represents the Manager of an adjoint/reverse MC simulation. 
00049 //              An adjoint run is divided in a serie of alternative adjoint and forward tracking 
00050 //              of adjoint and normal particles.
00051 //
00052 //              Reverse tracking phase:
00053 //              -----------------------
00054 //              An adjoint particle of a given type  (adjoint_e-, adjoint_gamma,...) is first generated on the so called adjoint source
00055 //              with a random energy (1/E distribution) and direction. The adjoint source is the 
00056 //              external surface of a user defined volume or of a user defined sphere. The adjoint
00057 //              source should contain one or several sensitive volumes and should be small 
00058 //              compared to the entire geometry.
00059 //              The user can set the min and max energy of the adjoint source. After its
00060 //              generation the adjoint primary  particle is tracked 
00061 //              bacward in the geometry till a user  defined external surface (spherical or boundary of a volume)
00062 //              or is killed before if it reaches a user defined  upper energy limit that represents
00063 //              the maximum energy of the external source. During the reverse tracking, reverse
00064 //              processes take place where  the adjoint particle being tracked can be either scattered
00065 //              or transformed in another type of adjoint paticle. During the reverse tracking the
00066 //              G4SimulationManager replaces the user defined Primary, Run, ... actions, by its own actions.
00067 //              
00068 //              Forward tracking phase
00069 //              -----------------------
00070 //              When an adjoint particle reaches the external surface its weight,type,  position, 
00071 //              and directions are registered and a  normal primary particle with a type equivalent to the last generated primary adjoint is
00072 //              generated with the same energy, position but opposite direction  and is  tracked normally in the sensitive region as in a fwd MC simulation. 
00073 //              During this forward tracking phase the  
00074 //              event, stacking, stepping, tracking actions defined by the user for its general fwd application are used. By this clear separation between
00075 //              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 
00076 //              simulation. Indeed  the computation of the signal is done by the same actions or classes that the one used in the fwd simulation mode.  
00077 // 
00078 //              Modification to brought in a existing G4 application to use the ReverseMC method
00079 //              -------------------------------
00080 //              In order to be able to use the ReverseMC method in his simulation, the user should modify its code as such:
00081 //                      1) Adapt its physics list to use ReverseProcesses for adjoint particles. An example of such physics list is provided in an extended 
00082 //                        example. 
00083 //                      2) Create an instance of    G4AdjointSimManager somewhere in the main code.
00084 //                      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
00085 //                       that reaches the external surface. This is done by using the following method of G4AdjointSimManager.
00086 //                                       
00087 //                                       G4int GetIDOfLastAdjParticleReachingExtSource()                                     
00088 //                                       G4ThreeVector GetPositionAtEndOfLastAdjointTrack(){ return last_pos;}
00089 //                                       G4ThreeVector GetDirectionAtEndOfLastAdjointTrack(){ return last_direction;}
00090 //                                       G4double GetEkinAtEndOfLastAdjointTrack(){ return last_ekin;} 
00091 //                                       G4double GetEkinNucAtEndOfLastAdjointTrack(){ return last_ekin_nuc;}
00092 //                                       G4double GetWeightAtEndOfLastAdjointTrack(){return last_weight;}
00093 //                                       G4double GetCosthAtEndOfLastAdjointTrack(){return last_cos_th;}
00094 //                                       G4String GetFwdParticleNameAtEndOfLastAdjointTrack(){return last_fwd_part_name;}
00095 //                                       G4int GetFwdParticlePDGEncodingAtEndOfLastAdjointTrack(){return last_fwd_part_PDGEncoding;}
00096 //                                       G4int GetFwdParticleIndexAtEndOfLastAdjointTrack().
00097 //
00098 //                        In orther to have a code working for both forward and adjoint simulation mode, the extra code needed in user actions for the adjoint
00099 //                        simulation mode can be seperated to the code needed only for the normal forward  simulation by using the following method
00100 //                                   
00101 //                                      G4bool GetAdjointSimMode() that return true if   an adjoint simulation is running and false if not!
00102 //                                
00103 //                        Example of modification in the analysis part of the code:
00104 //                       -------------------------------------------------------------
00105 //                              Let say that in the forward simulation a G4 application computes the energy deposited in a volume.
00106 //                              The user wants to normalise its results for an external isotropic source of e- with differential spectrum given by f(E).
00107 //                              A possible modification of the code where the deposited energy Edep during an event  is registered would be the following  
00108 //                                      
00109 //                              G4AdjointSimManager* theAdjSimManager =    G4AdjointSimManager::GetInstance();  
00110 //                              if (theAdjSimManager->GetAdjointSimMode()) {
00111 //                                      //code of the user that should be consider only for forwrad simulation
00112 //                                      G4double normalised_edep = 0.;
00113 //                                      if (theAdjSimManager->GetFwdParticleNameAtEndOfLastAdjointTrack() == "e-"){
00114 //                                              G4double ekin_prim = theAdjSimManager->GetEkinAtEndOfLastAdjointTrack();
00115 //                                              G4double weight_prim = theAdjSimManager->GetWeightAtEndOfLastAdjointTrack();
00116 //                                              normalised_edep = weight_prim*f(ekin_prim);
00117 //                                      } 
00118 //                                      //then follow the code where normalised_edep is printed, or registered or whatever ....
00119 //                              }       
00120 //                              
00121 //                              else { //code of the user that should be consider only for forward simulation
00122 //                              }
00123 //                              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
00124 //                              adapted for a normalisatin to several spectra and several type of primary particles  in the same simulation. 
00125 //
00126 
00127 #ifndef G4AdjointSimManager_h
00128 #define G4AdjointSimManager_h 1
00129 #include "globals.hh"
00130 #include "G4ThreeVector.hh"
00131 #include <vector>
00132 
00133 
00134 class G4UserEventAction;
00135 class G4VUserPrimaryGeneratorAction;
00136 class G4UserTrackingAction;
00137 class G4UserSteppingAction;
00138 class G4UserStackingAction;
00139 class G4UserRunAction;
00140 class G4AdjointRunAction;
00141 class G4AdjointPrimaryGeneratorAction;
00142 class G4AdjointSteppingAction;
00143 class G4AdjointEventAction;
00144 class G4AdjointStackingAction;
00145 class G4ParticleDefinition;
00146 class G4AdjointSimMessenger;
00147 class G4PhysicsLogVector;
00148 
00149 class G4AdjointSimManager  
00150 {
00151   public:
00152     
00153     static G4AdjointSimManager* GetInstance();
00154 
00155   public: //publich methods
00156     
00157     void RunAdjointSimulation(G4int nb_evt);
00158     
00159     inline G4int GetNbEvtOfLastRun(){return nb_evt_of_last_run;}
00160     
00161     void SetAdjointTrackingMode(G4bool aBool);
00162     inline G4bool GetAdjointTrackingMode(){return  adjoint_tracking_mode;} //true if an adjoint track is being processed
00163     inline G4bool GetAdjointSimMode(){return  adjoint_sim_mode;} //true if an adjoint simulation is running
00164 
00165     G4bool GetDidAdjParticleReachTheExtSource();
00166     void RegisterAtEndOfAdjointTrack();
00167     void RegisterAdjointPrimaryWeight(G4double aWeight);
00168 
00169     inline G4int GetIDOfLastAdjParticleReachingExtSource(){return ID_of_last_particle_that_reach_the_ext_source;};                                   
00170     inline G4ThreeVector GetPositionAtEndOfLastAdjointTrack(){ return last_pos;}
00171     inline G4ThreeVector GetDirectionAtEndOfLastAdjointTrack(){ return last_direction;}
00172     inline G4double GetEkinAtEndOfLastAdjointTrack(){ return last_ekin;} 
00173     inline G4double GetEkinNucAtEndOfLastAdjointTrack(){ return last_ekin_nuc;}
00174     inline G4double GetWeightAtEndOfLastAdjointTrack(){return last_weight;}
00175     inline G4double GetCosthAtEndOfLastAdjointTrack(){return last_cos_th;}
00176     inline const G4String& GetFwdParticleNameAtEndOfLastAdjointTrack(){return last_fwd_part_name;}
00177     inline G4int GetFwdParticlePDGEncodingAtEndOfLastAdjointTrack(){return last_fwd_part_PDGEncoding;}
00178     inline G4int GetFwdParticleIndexAtEndOfLastAdjointTrack(){return last_fwd_part_index;}
00179     
00180     std::vector<G4ParticleDefinition*> GetListOfPrimaryFwdParticles();
00181      
00182     G4bool DefineSphericalExtSource(G4double radius, G4ThreeVector pos);
00183     G4bool DefineSphericalExtSourceWithCentreAtTheCentreOfAVolume(G4double radius, const G4String& volume_name);
00184     G4bool DefineExtSourceOnTheExtSurfaceOfAVolume(const G4String& volume_name);
00185     void SetExtSourceEmax(G4double Emax);
00186     
00187     //Definition of adjoint source
00188     //---------------------------- 
00189      
00190     G4bool DefineSphericalAdjointSource(G4double radius, G4ThreeVector pos);
00191     G4bool DefineSphericalAdjointSourceWithCentreAtTheCentreOfAVolume(G4double radius, const G4String& volume_name);
00192     G4bool DefineAdjointSourceOnTheExtSurfaceOfAVolume(const G4String& volume_name);
00193     void SetAdjointSourceEmin(G4double Emin);
00194     void SetAdjointSourceEmax(G4double Emax);
00195     inline G4double GetAdjointSourceArea(){return area_of_the_adjoint_source;} 
00196     void ConsiderParticleAsPrimary(const G4String& particle_name);
00197     void NeglectParticleAsPrimary(const G4String& particle_name);
00198     void SetPrimaryIon(G4ParticleDefinition* adjointIon, G4ParticleDefinition* fwdIon);
00199     const G4String& GetPrimaryIonName();
00200     
00201     inline void SetNormalisationMode(G4int n){normalisation_mode=n;};
00202     G4int GetNormalisationMode(){return normalisation_mode;};
00203     G4double GetNumberNucleonsInIon(){return nb_nuc;};
00204 
00205     //Definition of user actions for the adjoint tracking phase
00206     //---------------------------- 
00207     void SetAdjointEventAction(G4UserEventAction* anAction);
00208     void SetAdjointSteppingAction(G4UserSteppingAction* anAction);
00209     void SetAdjointStackingAction(G4UserStackingAction* anAction);
00210     void SetAdjointTrackingAction(G4UserTrackingAction* anAction);
00211     void SetAdjointRunAction(G4UserRunAction* anAction); 
00212     
00213     //Set methods for user run actions
00214     //--------------------------------
00215     inline void UseUserStackingActionInFwdTrackingPhase(G4bool aBool){use_user_StackingAction=aBool;}
00216 
00217     //Convergence test
00218     //-----------------------
00219    /*
00220     void RegisterSignalForConvergenceTest(G4double aSignal);
00221     void DefineExponentialPrimarySpectrumForConvergenceTest(G4ParticleDefinition* aPartDef, G4double E0);
00222     void DefinePowerLawPrimarySpectrumForConvergenceTest(G4ParticleDefinition* aPartDef, G4double alpha);
00223      
00224    */ 
00225 
00226   private: 
00227   
00228     static G4AdjointSimManager* instance;
00229   
00230   private: // methods
00231     
00232     void SetRestOfAdjointActions(); 
00233     void SetAdjointPrimaryRunAndStackingActions();
00234     void ResetRestOfUserActions(); 
00235     void ResetUserPrimaryRunAndStackingActions(); 
00236     void DefineUserActions();
00237 
00238   private: //constructor and destructor
00239   
00240     G4AdjointSimManager();
00241    ~G4AdjointSimManager();
00242 
00243   private ://attributes
00244   
00245   //Messenger
00246   //----------
00247     G4AdjointSimMessenger* theMessenger;
00248   
00249   //user defined actions for the normal fwd simulation. Taken from the G4RunManager
00250   //-------------------------------------------------
00251     bool user_action_already_defined;
00252     G4UserRunAction*            fUserRunAction;
00253     G4UserEventAction*          fUserEventAction;
00254     G4VUserPrimaryGeneratorAction*      fUserPrimaryGeneratorAction;
00255     G4UserTrackingAction*       fUserTrackingAction;
00256     G4UserSteppingAction*       fUserSteppingAction;
00257     G4UserStackingAction*       fUserStackingAction;
00258     bool use_user_StackingAction; //only for fwd part of the adjoint simulation
00259     
00260   //action for adjoint simulation
00261   //-----------------------------
00262     G4UserRunAction*            theAdjointRunAction;
00263     G4UserEventAction*          theAdjointEventAction;
00264     G4AdjointPrimaryGeneratorAction* theAdjointPrimaryGeneratorAction;
00265     G4UserTrackingAction*       theAdjointTrackingAction;
00266     G4AdjointSteppingAction*    theAdjointSteppingAction;
00267     G4AdjointStackingAction*    theAdjointStackingAction;
00268       
00269   //adjoint mode
00270   //-------------
00271     G4bool adjoint_tracking_mode;
00272     G4bool adjoint_sim_mode;   
00273 
00274   //adjoint particle information on the external surface
00275   //----------------------------- 
00276     G4ThreeVector last_pos;
00277     G4ThreeVector last_direction;
00278     G4double last_ekin,last_ekin_nuc; //last_ekin_nuc=last_ekin/nuc, nuc is 1 if not a nucleus
00279     G4double last_cos_th;
00280     G4String last_fwd_part_name;
00281     G4int  last_fwd_part_PDGEncoding;
00282     G4int  last_fwd_part_index;
00283     G4double last_weight;
00284     G4int ID_of_last_particle_that_reach_the_ext_source;
00285       
00286     G4int nb_evt_of_last_run;
00287     G4int normalisation_mode;
00288      
00289   //Adjoint source
00290   //--------------
00291     G4double area_of_the_adjoint_source; 
00292     G4double nb_nuc;
00293     G4double theAdjointPrimaryWeight;
00294 
00295     //Weight Analysis
00296     //----------
00297     G4PhysicsLogVector* electron_last_weight_vector;
00298     G4PhysicsLogVector* proton_last_weight_vector;
00299     G4PhysicsLogVector* gamma_last_weight_vector;
00300     
00301     G4bool welcome_message;
00302     
00303 /* For the future    
00304     //Convergence test
00305     //----------------
00306       
00307      G4double normalised_signal;
00308      G4double error_signal;
00309      G4bool convergence_test_is_used;
00310      G4bool power_law_spectrum_for_convergence_test; // true  PowerLaw, ;
00311      G4ParticleDefinition* the_par_def_for_convergence_test;
00312 */     
00313  
00314 };
00315 
00316 #endif
00317 

Generated on Mon May 27 17:47:38 2013 for Geant4 by  doxygen 1.4.7