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