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
G4AdjointPrimaryGeneratorAction.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// G4AdjointPrimaryGeneratorAction
27//
28// Class description:
29//
30// This class represents the PrimaryGeneratorAction that is used during
31// the entire adjoint simulation.
32// It uses the class G4AdjointPrimaryGenerator to generate randomly
33// adjoint primary particles on a user selected adjoint source
34// (External surface of a volume or Sphere).
35// The spectrum of the primary adjoint particles is set as 1/E with
36// user defined max and min energy.
37// The weight of the primary is set according to ReverseMC theory as
38// w=log(Emax/Emin)*E*adjoint_source_area*pi/n, with E the energy of
39// the particle, n the number of adjoint primary particles of same type
40// that will be generated during the simulation.
41// Different types of adjoint particles are generated event after event
42// in order to cover all the type of primaries and secondaries needed
43// for the simulation. For example if reverse e- ionisation, brem,
44// photo electric effect, and compton are considered both adjoint gamma
45// and adjoint e- will be considered alternatively as adjoint primary.
46// The user can decide to consider/neglect some type of particle by
47// using the macro commands "/adjoint/ConsiderAsPrimary" and
48// "/adjoint/NeglectAsPrimary". If an adjoint primary or its secondary
49// has reached the external surface, in the next event a fwd primary
50// particle equivalent to the last generated adjoint primary is
51// generated with the same position, energy but opposite direction
52// and the forward tracking phase starts.
53
54// --------------------------------------------------------------------
55// Class Name: G4AdjointPrimaryGeneratorAction
56// Author: L. Desorgher, 2007-2009
57// Organisation: SpaceIT GmbH
58// Contract: ESA contract 21435/08/NL/AT
59// Customer: ESA/ESTEC
60// --------------------------------------------------------------------
61#ifndef G4AdjointPrimaryGeneratorAction_hh
62#define G4AdjointPrimaryGeneratorAction_hh 1
63
64#include <iterator>
65#include <map>
66#include <vector>
67
68#include "globals.hh"
69#include "G4ThreeVector.hh"
71
73class G4ParticleGun;
74class G4Event;
77
78// --------------------------------------------------------------------
79
81{
82 public:
83
86
88 const G4AdjointPrimaryGeneratorAction&) = delete;
90 const G4AdjointPrimaryGeneratorAction&) = delete;
91
93 void SetEmin(G4double val);
94 void SetEmax(G4double val);
95 void SetEminIon(G4double val);
96 void SetEmaxIon(G4double val);
98 G4ThreeVector pos);
100 const G4String& volume_name);
101 void ConsiderParticleAsPrimary(const G4String& particle_name);
102 void NeglectParticleAsPrimary(const G4String& particle_name);
103 void SetPrimaryIon(G4ParticleDefinition* adjointIon,
104 G4ParticleDefinition* fwdIon);
106
107 inline void SetRndmFlag(const G4String& val)
108 {
109 rndmFlag = val;
110 }
112 {
113 return ListOfPrimaryAdjParticles.size();
114 }
115 inline std::vector<G4ParticleDefinition*>* GetListOfPrimaryFwdParticles()
116 {
117 return &ListOfPrimaryFwdParticles;
118 }
120 {
121 return ion_name;
122 }
124 {
125 nb_fwd_gammas_per_event = nb;
126 }
128 {
129 nb_adj_primary_gammas_per_event = nb;
130 }
132 {
133 nb_adj_primary_electrons_per_event = nb;
134 }
136 {
137 return ListOfPrimaryFwdParticles[index_particle];
138 }
139
140 private: // methods
141
142 G4double ComputeEnergyDistWeight(G4double energy, G4double E1, G4double E2);
143
144 private: // attributes
145
146 G4String rndmFlag; // flag for a rndm impact point
147
148 // The generator of primary vertex except for weight
149 G4AdjointPrimaryGenerator* theAdjointPrimaryGenerator = nullptr;
150
151 // Emin and Emax energies of the adjoint source
152 //---------------------------------------------
153 G4double Emin = 0.0;
154 G4double Emax = 0.0;
155 G4double EminIon = 0.0;
156 G4double EmaxIon = 0.0;
157
158 // List of type of primary adjoint and forward particle used in the
159 // simulation
160 //------------------------------------------------------------------
161 std::vector<G4ParticleDefinition*> ListOfPrimaryFwdParticles;
162 std::vector<G4ParticleDefinition*> ListOfPrimaryAdjParticles;
163 std::map<G4String, G4bool> PrimariesConsideredInAdjointSim;
164 // if true considered if false not considered
165
166 std::size_t index_particle = 100000;
167
168 G4ThreeVector pos, direction, p;
169
170 G4String type_of_adjoint_source; // Spherical ExtSurfaceOfAVolume
171 G4double radius_spherical_source = 0.0;
172 G4ThreeVector center_spherical_source;
173 G4int nb_fwd_gammas_per_event = 1;
174 G4int nb_adj_primary_gammas_per_event = 1;
175 G4int nb_adj_primary_electrons_per_event = 1;
176
177 // For simulation with ions
178 //--------------------------
179 G4ParticleDefinition* fwd_ion = nullptr;
180 G4ParticleDefinition* adj_ion = nullptr;
181 G4String ion_name = "not_defined";
182};
183
184#endif
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4AdjointPrimaryGeneratorAction & operator=(const G4AdjointPrimaryGeneratorAction &)=delete
std::vector< G4ParticleDefinition * > * GetListOfPrimaryFwdParticles()
void ConsiderParticleAsPrimary(const G4String &particle_name)
G4AdjointPrimaryGeneratorAction(const G4AdjointPrimaryGeneratorAction &)=delete
void SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume(const G4String &volume_name)
void SetSphericalAdjointPrimarySource(G4double radius, G4ThreeVector pos)
G4ParticleDefinition * GetLastGeneratedFwdPrimaryParticle()
void NeglectParticleAsPrimary(const G4String &particle_name)
void SetPrimaryIon(G4ParticleDefinition *adjointIon, G4ParticleDefinition *fwdIon)