Geant4 11.4.0
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
51 PrimariesConsideredInAdjointSim[G4String("e-")] = false;
52 PrimariesConsideredInAdjointSim[G4String("gamma")] = false;
53 PrimariesConsideredInAdjointSim[G4String("proton")] = false;
54 PrimariesConsideredInAdjointSim[G4String("ion")] = false;
55
56 ListOfPrimaryFwdParticles.clear();
57 ListOfPrimaryAdjParticles.clear();
58}
59
60// --------------------------------------------------------------------
61//
65
66// --------------------------------------------------------------------
67//
69{
70 G4int evt_id = anEvent->GetEventID();
71 std::size_t n = ListOfPrimaryAdjParticles.size();
72 index_particle = std::size_t(evt_id) - n * (std::size_t(evt_id) / n);
73
74 G4double E1 = Emin;
75 G4double E2 = Emax;
76 if (ListOfPrimaryAdjParticles[index_particle] == nullptr)
77 UpdateListOfPrimaryParticles(); // ion has not been created yet
78
79 if (ListOfPrimaryAdjParticles[index_particle]->GetParticleName() == "adj_proton") {
80 E1 = EminIon;
81 E2 = EmaxIon;
82 }
83 if (ListOfPrimaryAdjParticles[index_particle]->GetParticleType() == "adjoint_nucleus") {
84 G4int A = ListOfPrimaryAdjParticles[index_particle]->GetAtomicMass();
85 E1 = EminIon * A;
86 E2 = EmaxIon * A;
87 }
88 // Generate first the forwrad primaries
89 G4AdjointPrimaryGenerator::GetInstance()->GenerateFwdPrimaryVertex(anEvent, ListOfPrimaryFwdParticles[index_particle], E1, E2);
90
91 G4PrimaryVertex* fwdPrimVertex = anEvent->GetPrimaryVertex();
92
93 p = fwdPrimVertex->GetPrimary()->GetMomentum();
94 pos = fwdPrimVertex->GetPosition();
95 G4double pmag = p.mag();
96 G4double m0 = ListOfPrimaryFwdParticles[index_particle]->GetPDGMass();
97 G4double ekin = std::sqrt(m0 * m0 + pmag * pmag) - m0;
98
99 G4double weight_correction = 1.;
100 // For gamma generate the particle along the backward ray
101 G4ThreeVector dir = -p / p.mag();
102
103 weight_correction = 1.;
104
105 if (ListOfPrimaryFwdParticles[index_particle] == G4Gamma::Gamma() && nb_fwd_gammas_per_event > 1)
106 {
107 G4double weight = (1. / nb_fwd_gammas_per_event);
108 fwdPrimVertex->SetWeight(weight);
109 for (G4int i = 0; i < nb_fwd_gammas_per_event - 1; ++i) {
110 auto newFwdPrimVertex = new G4PrimaryVertex();
111 newFwdPrimVertex->SetPosition(pos.x(), pos.y(), pos.z());
112 newFwdPrimVertex->SetT0(0.);
113 auto aPrimParticle =
114 new G4PrimaryParticle(ListOfPrimaryFwdParticles[index_particle], p.x(), p.y(), p.z());
115 newFwdPrimVertex->SetPrimary(aPrimParticle);
116 newFwdPrimVertex->SetWeight(weight);
117 anEvent->AddPrimaryVertex(newFwdPrimVertex);
118 }
119 }
120
121 // Now generate the adjoint primaries
122 auto adjPrimVertex = new G4PrimaryVertex();
123 adjPrimVertex->SetPosition(pos.x(), pos.y(), pos.z());
124 adjPrimVertex->SetT0(0.);
125 auto aPrimParticle =
126 new G4PrimaryParticle(ListOfPrimaryAdjParticles[index_particle], -p.x(), -p.y(), -p.z());
127
128 adjPrimVertex->SetPrimary(aPrimParticle);
129 anEvent->AddPrimaryVertex(adjPrimVertex);
130
131 // The factor pi is to normalise the weight to the directional flux
133 G4double adjoint_weight =
134 weight_correction * ComputeEnergyDistWeight(ekin, E1, E2) * adjoint_source_area * pi;
135 if (ListOfPrimaryAdjParticles[index_particle]->GetParticleName() == "adj_gamma") {
136 // The weight will be corrected at the end of the track if splitted tracks
137 // are used
138 adjoint_weight = adjoint_weight / nb_adj_primary_gammas_per_event;
139 for (G4int i = 0; i < nb_adj_primary_gammas_per_event - 1; ++i) {
140 auto newAdjPrimVertex = new G4PrimaryVertex();
141 newAdjPrimVertex->SetPosition(pos.x(), pos.y(), pos.z());
142 newAdjPrimVertex->SetT0(0.);
143 aPrimParticle =
144 new G4PrimaryParticle(ListOfPrimaryAdjParticles[index_particle], -p.x(), -p.y(), -p.z());
145 newAdjPrimVertex->SetPrimary(aPrimParticle);
146 newAdjPrimVertex->SetWeight(adjoint_weight);
147 anEvent->AddPrimaryVertex(newAdjPrimVertex);
148 }
149 }
150 else if (ListOfPrimaryAdjParticles[index_particle]->GetParticleName() == "adj_electron") {
151 // The weight will be corrected at the end of the track if splitted tracks
152 // are used
153 adjoint_weight = adjoint_weight / nb_adj_primary_electrons_per_event;
154 for (G4int i = 0; i < nb_adj_primary_electrons_per_event - 1; ++i) {
155 auto newAdjPrimVertex = new G4PrimaryVertex();
156 newAdjPrimVertex->SetPosition(pos.x(), pos.y(), pos.z());
157 newAdjPrimVertex->SetT0(0.);
158 aPrimParticle =
159 new G4PrimaryParticle(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 adjPrimVertex->SetWeight(adjoint_weight);
166
167 // Call some methods of G4AdjointSimManager
171}
172
173// --------------------------------------------------------------------
174//
176{
177 Emin = val;
178 EminIon = val;
179}
180
181// --------------------------------------------------------------------
182//
184{
185 Emax = val;
186 EmaxIon = val;
187}
188
189// --------------------------------------------------------------------
190//
192{
193 EminIon = val;
194}
195
196// --------------------------------------------------------------------
197//
199{
200 EmaxIon = val;
201}
202
203// --------------------------------------------------------------------
204//
205G4double G4AdjointPrimaryGeneratorAction::ComputeEnergyDistWeight(G4double E, G4double E1,
206 G4double E2)
207{
208 // We generate N numbers of primaries with a 1/E energy law distribution.
209 // We have therefore an energy distribution function
210 // f(E)=C/E (1)
211 // with C a constant that is such that
212 // N=Integral(f(E),E1,E2)=C.std::log(E2/E1) (2)
213 // Therefore from (2) we get
214 // C=N/ std::log(E2/E1) (3)
215 // and
216 // f(E)=N/ std::log(E2/E1)/E (4)
217 // For the adjoint simulation we need a energy distribution f'(E)=1..
218 // To get that we need therefore to apply a weight to the primary
219 // W=1/f(E)=E*std::log(E2/E1)/N
220 //
221 return std::log(E2 / E1) * E / G4AdjointSimManager::GetInstance()->GetNbEvtOfLastRun();
222}
223
224// --------------------------------------------------------------------
225//
227 G4ThreeVector center_pos)
228{
229 radius_spherical_source = radius;
230 center_spherical_source = center_pos;
231 type_of_adjoint_source = "Spherical";
233}
234
235// --------------------------------------------------------------------
236//
238 const G4String& volume_name)
239{
240 type_of_adjoint_source = "ExternalSurfaceOfAVolume";
242}
243
244// --------------------------------------------------------------------
245//
247{
248 if (PrimariesConsideredInAdjointSim.find(particle_name) != PrimariesConsideredInAdjointSim.end())
249 {
250 PrimariesConsideredInAdjointSim[particle_name] = true;
251 }
253}
254
255// --------------------------------------------------------------------
256//
258{
259 if (PrimariesConsideredInAdjointSim.find(particle_name) != PrimariesConsideredInAdjointSim.end())
260 {
261 PrimariesConsideredInAdjointSim[particle_name] = false;
262 }
264}
265
266// --------------------------------------------------------------------
267//
269{
271 ListOfPrimaryFwdParticles.clear();
272 ListOfPrimaryAdjParticles.clear();
273 for (const auto& iter : PrimariesConsideredInAdjointSim) {
274 if (iter.second) {
275 G4String fwd_particle_name = iter.first;
276 if (fwd_particle_name != "ion") {
277 G4String adj_particle_name = G4String("adj_") + fwd_particle_name;
278 ListOfPrimaryFwdParticles.push_back(theParticleTable->FindParticle(fwd_particle_name));
279 ListOfPrimaryAdjParticles.push_back(theParticleTable->FindParticle(adj_particle_name));
280 }
281 else {
282 if (fwd_ion != nullptr) {
283 ion_name = fwd_ion->GetParticleName();
284 G4String adj_ion_name = G4String("adj_") + ion_name;
285 ListOfPrimaryFwdParticles.push_back(fwd_ion);
286 ListOfPrimaryAdjParticles.push_back(adj_ion);
287 }
288 else {
289 ListOfPrimaryFwdParticles.push_back(nullptr);
290 ListOfPrimaryAdjParticles.push_back(nullptr);
291 }
292 }
293 }
294 }
295}
296
297// --------------------------------------------------------------------
298//
300 G4ParticleDefinition* fwdIon)
301{
302 fwd_ion = fwdIon;
303 adj_ion = adjointIon;
305}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
const G4double A[17]
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)
static G4AdjointPrimaryGenerator * GetInstance()
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:145
G4int GetEventID() const
Definition G4Event.hh:126
void AddPrimaryVertex(G4PrimaryVertex *aPrimaryVertex)
Definition G4Event.hh:129
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
G4ThreeVector GetMomentum() const
G4ThreeVector GetPosition() const
void SetWeight(G4double w)
G4PrimaryParticle * GetPrimary(G4int i=0) const