Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4AdjointPrimaryGeneratorAction.cc
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 implementation
27//
28// --------------------------------------------------------------------
29// Class Name: G4AdjointPrimaryGeneratorAction
30// Author: L. Desorgher, 2007-2009
31// Organisation: SpaceIT GmbH
32// Contract: ESA contract 21435/08/NL/AT
33// Customer: ESA/ESTEC
34// --------------------------------------------------------------------
35
37
40#include "G4Event.hh"
41#include "G4Gamma.hh"
43#include "G4ParticleTable.hh"
45
46// --------------------------------------------------------------------
47//
49{
50 theAdjointPrimaryGenerator = new G4AdjointPrimaryGenerator();
51
52 PrimariesConsideredInAdjointSim[G4String("e-")] = false;
53 PrimariesConsideredInAdjointSim[G4String("gamma")] = false;
54 PrimariesConsideredInAdjointSim[G4String("proton")] = false;
55 PrimariesConsideredInAdjointSim[G4String("ion")] = false;
56
57 ListOfPrimaryFwdParticles.clear();
58 ListOfPrimaryAdjParticles.clear();
59}
60
61// --------------------------------------------------------------------
62//
64{
65 delete theAdjointPrimaryGenerator;
66}
67
68// --------------------------------------------------------------------
69//
71{
72 G4int evt_id = anEvent->GetEventID();
73 std::size_t n = ListOfPrimaryAdjParticles.size();
74 index_particle = std::size_t(evt_id) - n * (std::size_t(evt_id) / n);
75
76 G4double E1 = Emin;
77 G4double E2 = Emax;
78 if(ListOfPrimaryAdjParticles[index_particle] == nullptr)
79 UpdateListOfPrimaryParticles(); // ion has not been created yet
80
81 if(ListOfPrimaryAdjParticles[index_particle]->GetParticleName() ==
82 "adj_proton")
83 {
84 E1 = EminIon;
85 E2 = EmaxIon;
86 }
87 if(ListOfPrimaryAdjParticles[index_particle]->GetParticleType() ==
88 "adjoint_nucleus")
89 {
90 G4int A = ListOfPrimaryAdjParticles[index_particle]->GetAtomicMass();
91 E1 = EminIon * A;
92 E2 = EmaxIon * A;
93 }
94 // Generate first the forwrad primaries
95 theAdjointPrimaryGenerator->GenerateFwdPrimaryVertex(
96 anEvent, ListOfPrimaryFwdParticles[index_particle], E1, E2);
97 G4PrimaryVertex* fwdPrimVertex = anEvent->GetPrimaryVertex();
98
99 p = fwdPrimVertex->GetPrimary()->GetMomentum();
100 pos = fwdPrimVertex->GetPosition();
101 G4double pmag = p.mag();
102 G4double m0 = ListOfPrimaryFwdParticles[index_particle]->GetPDGMass();
103 G4double ekin = std::sqrt(m0 * m0 + pmag * pmag) - m0;
104
105 G4double weight_correction = 1.;
106 // For gamma generate the particle along the backward ray
107 G4ThreeVector dir = -p / p.mag();
108
109 weight_correction = 1.;
110
111 if(ListOfPrimaryFwdParticles[index_particle] == G4Gamma::Gamma() &&
112 nb_fwd_gammas_per_event > 1)
113 {
114 G4double weight = (1. / nb_fwd_gammas_per_event);
115 fwdPrimVertex->SetWeight(weight);
116 for(G4int i = 0; i < nb_fwd_gammas_per_event - 1; ++i)
117 {
118 G4PrimaryVertex* newFwdPrimVertex = new G4PrimaryVertex();
119 newFwdPrimVertex->SetPosition(pos.x(), pos.y(), pos.z());
120 newFwdPrimVertex->SetT0(0.);
121 G4PrimaryParticle* aPrimParticle = new G4PrimaryParticle(
122 ListOfPrimaryFwdParticles[index_particle], p.x(), p.y(), p.z());
123 newFwdPrimVertex->SetPrimary(aPrimParticle);
124 newFwdPrimVertex->SetWeight(weight);
125 anEvent->AddPrimaryVertex(newFwdPrimVertex);
126 }
127 }
128
129 // Now generate the adjoint primaries
130 G4PrimaryVertex* adjPrimVertex = new G4PrimaryVertex();
131 adjPrimVertex->SetPosition(pos.x(), pos.y(), pos.z());
132 adjPrimVertex->SetT0(0.);
133 G4PrimaryParticle* aPrimParticle = new G4PrimaryParticle(
134 ListOfPrimaryAdjParticles[index_particle], -p.x(), -p.y(), -p.z());
135
136 adjPrimVertex->SetPrimary(aPrimParticle);
137 anEvent->AddPrimaryVertex(adjPrimVertex);
138
139 // The factor pi is to normalise the weight to the directional flux
140 G4double adjoint_source_area =
142 G4double adjoint_weight = weight_correction *
143 ComputeEnergyDistWeight(ekin, E1, E2) *
144 adjoint_source_area * pi;
145 // if (ListOfPrimaryFwdParticles[index_particle] ==G4Gamma::Gamma())
146 // adjoint_weight = adjoint_weight/3.;
147 if(ListOfPrimaryAdjParticles[index_particle]->GetParticleName() ==
148 "adj_gamma")
149 {
150 // The weight will be corrected at the end of the track if splitted tracks
151 // are used
152 adjoint_weight = adjoint_weight / nb_adj_primary_gammas_per_event;
153 for(G4int i = 0; i < nb_adj_primary_gammas_per_event - 1; ++i)
154 {
155 G4PrimaryVertex* newAdjPrimVertex = new G4PrimaryVertex();
156 newAdjPrimVertex->SetPosition(pos.x(), pos.y(), pos.z());
157 newAdjPrimVertex->SetT0(0.);
158 aPrimParticle = new G4PrimaryParticle(
159 ListOfPrimaryAdjParticles[index_particle], -p.x(), -p.y(), -p.z());
160 newAdjPrimVertex->SetPrimary(aPrimParticle);
161 newAdjPrimVertex->SetWeight(adjoint_weight);
162 anEvent->AddPrimaryVertex(newAdjPrimVertex);
163 }
164 }
165 else if(ListOfPrimaryAdjParticles[index_particle]->GetParticleName() ==
166 "adj_electron")
167 {
168 // The weight will be corrected at the end of the track if splitted tracks
169 // are used
170 adjoint_weight = adjoint_weight / nb_adj_primary_electrons_per_event;
171 for(G4int i = 0; i < nb_adj_primary_electrons_per_event - 1; ++i)
172 {
173 G4PrimaryVertex* newAdjPrimVertex = new G4PrimaryVertex();
174 newAdjPrimVertex->SetPosition(pos.x(), pos.y(), pos.z());
175 newAdjPrimVertex->SetT0(0.);
176 aPrimParticle = new G4PrimaryParticle(
177 ListOfPrimaryAdjParticles[index_particle], -p.x(), -p.y(), -p.z());
178 newAdjPrimVertex->SetPrimary(aPrimParticle);
179 newAdjPrimVertex->SetWeight(adjoint_weight);
180 anEvent->AddPrimaryVertex(newAdjPrimVertex);
181 }
182 }
183 adjPrimVertex->SetWeight(adjoint_weight);
184
185 // Call some methods of G4AdjointSimManager
190
191}
192
193// --------------------------------------------------------------------
194//
196{
197 Emin = val;
198 EminIon = val;
199}
200
201// --------------------------------------------------------------------
202//
204{
205 Emax = val;
206 EmaxIon = val;
207}
208
209// --------------------------------------------------------------------
210//
212{
213 EminIon = val;
214}
215
216// --------------------------------------------------------------------
217//
219{
220 EmaxIon = val;
221}
222
223// --------------------------------------------------------------------
224//
225G4double G4AdjointPrimaryGeneratorAction::ComputeEnergyDistWeight(G4double E,
226 G4double E1,
227 G4double E2)
228{
229 // We generate N numbers of primaries with a 1/E energy law distribution.
230 // We have therefore an energy distribution function
231 // f(E)=C/E (1)
232 // with C a constant that is such that
233 // N=Integral(f(E),E1,E2)=C.std::log(E2/E1) (2)
234 // Therefore from (2) we get
235 // C=N/ std::log(E2/E1) (3)
236 // and
237 // f(E)=N/ std::log(E2/E1)/E (4)
238 // For the adjoint simulation we need a energy distribution f'(E)=1..
239 // To get that we need therefore to apply a weight to the primary
240 // W=1/f(E)=E*std::log(E2/E1)/N
241 //
242 return std::log(E2 / E1) * E /
244}
245
246// --------------------------------------------------------------------
247//
249 G4double radius, G4ThreeVector center_pos)
250{
251 radius_spherical_source = radius;
252 center_spherical_source = center_pos;
253 type_of_adjoint_source = "Spherical";
254 theAdjointPrimaryGenerator->SetSphericalAdjointPrimarySource(radius,
255 center_pos);
256}
257
258// --------------------------------------------------------------------
259//
262{
263 type_of_adjoint_source = "ExternalSurfaceOfAVolume";
264 theAdjointPrimaryGenerator->SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume(
265 volume_name);
266}
267
268// --------------------------------------------------------------------
269//
271 const G4String& particle_name)
272{
273 if(PrimariesConsideredInAdjointSim.find(particle_name) !=
274 PrimariesConsideredInAdjointSim.end())
275 {
276 PrimariesConsideredInAdjointSim[particle_name] = true;
277 }
279}
280
281// --------------------------------------------------------------------
282//
284 const G4String& particle_name)
285{
286 if(PrimariesConsideredInAdjointSim.find(particle_name) !=
287 PrimariesConsideredInAdjointSim.end())
288 {
289 PrimariesConsideredInAdjointSim[particle_name] = false;
290 }
292}
293
294// --------------------------------------------------------------------
295//
297{
299 ListOfPrimaryFwdParticles.clear();
300 ListOfPrimaryAdjParticles.clear();
301 for(auto iter = PrimariesConsideredInAdjointSim.cbegin();
302 iter != PrimariesConsideredInAdjointSim.cend(); ++iter)
303 {
304 if(iter->second)
305 {
306 G4String fwd_particle_name = iter->first;
307 if(fwd_particle_name != "ion")
308 {
309 G4String adj_particle_name = G4String("adj_") + fwd_particle_name;
310 ListOfPrimaryFwdParticles.push_back(
311 theParticleTable->FindParticle(fwd_particle_name));
312 ListOfPrimaryAdjParticles.push_back(
313 theParticleTable->FindParticle(adj_particle_name));
314 }
315 else
316 {
317 if(fwd_ion)
318 {
319 ion_name = fwd_ion->GetParticleName();
320 G4String adj_ion_name = G4String("adj_") + ion_name;
321 ListOfPrimaryFwdParticles.push_back(fwd_ion);
322 ListOfPrimaryAdjParticles.push_back(adj_ion);
323 }
324 else
325 {
326 ListOfPrimaryFwdParticles.push_back(nullptr);
327 ListOfPrimaryAdjParticles.push_back(nullptr);
328 }
329 }
330 }
331 }
332}
333
334// --------------------------------------------------------------------
335//
337 G4ParticleDefinition* adjointIon, G4ParticleDefinition* fwdIon)
338{
339 fwd_ion = fwdIon;
340 adj_ion = adjointIon;
342}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4double A[17]
double z() const
double x() const
double y() const
double mag() const
void ConsiderParticleAsPrimary(const G4String &particle_name)
void SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume(const G4String &volume_name)
void SetSphericalAdjointPrimarySource(G4double radius, G4ThreeVector pos)
void NeglectParticleAsPrimary(const G4String &particle_name)
void SetPrimaryIon(G4ParticleDefinition *adjointIon, G4ParticleDefinition *fwdIon)
void SetSphericalAdjointPrimarySource(G4double radius, G4ThreeVector pos)
void SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume(const G4String &v_name)
void GenerateFwdPrimaryVertex(G4Event *anEvt, G4ParticleDefinition *adj_part, G4double E1, G4double E2)
void SetAdjointTrackingMode(G4bool aBool)
void ResetDidOneAdjPartReachExtSourceDuringEvent()
static G4AdjointSimManager * GetInstance()
G4PrimaryVertex * GetPrimaryVertex(G4int i=0) const
Definition: G4Event.hh:137
G4int GetEventID() const
Definition: G4Event.hh:118
void AddPrimaryVertex(G4PrimaryVertex *aPrimaryVertex)
Definition: G4Event.hh:121
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
const G4String & GetParticleName() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
G4ThreeVector GetMomentum() const
void SetPosition(G4double x0, G4double y0, G4double z0)
G4ThreeVector GetPosition() const
void SetPrimary(G4PrimaryParticle *pp)
void SetWeight(G4double w)
void SetT0(G4double t0)
G4PrimaryParticle * GetPrimary(G4int i=0) const