Geant4 9.6.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// $Id$
27//
28/////////////////////////////////////////////////////////////////////////////
29// Class Name: G4AdjointPrimaryGeneratorAction
30// Author: L. Desorgher
31// Organisation: SpaceIT GmbH
32// Contract: ESA contract 21435/08/NL/AT
33// Customer: ESA/ESTEC
34/////////////////////////////////////////////////////////////////////////////
35
37
39#include "G4Event.hh"
40#include "G4ParticleTable.hh"
42#include "G4ParticleTable.hh"
43#include "G4AdjointSimManager.hh"
45/////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
46//
48 : Emin(0.), Emax(0.), EminIon(0.), EmaxIon(0.), NbOfAdjointPrimaryTypes(0),
49 index_particle(100000), last_generated_part_was_adjoint(false),
50 radius_spherical_source(0.), fwd_ion(0), adj_ion(0),
51 ion_name("not_defined")
52{
53 theAdjointPrimaryGenerator= new G4AdjointPrimaryGenerator();
54
55 PrimariesConsideredInAdjointSim[G4String("e-")]=false;
56 PrimariesConsideredInAdjointSim[G4String("gamma")]=false;
57 PrimariesConsideredInAdjointSim[G4String("proton")]=false;
58 PrimariesConsideredInAdjointSim[G4String("ion")]=false;
59
60 ListOfPrimaryFwdParticles.clear();
61 ListOfPrimaryAdjParticles.clear();
62}
63/////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
64//
66{
67 delete theAdjointPrimaryGenerator;
68}
69/////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
70//
72{
73 if ( !last_generated_part_was_adjoint ) {
74
75 index_particle++;
76 if (index_particle >= ListOfPrimaryAdjParticles.size()) index_particle =0;
77
78
79 G4double E1=Emin;
80 G4double E2=Emax;
81 if (!ListOfPrimaryAdjParticles[index_particle]) UpdateListOfPrimaryParticles();//ion has not been created yet
82
83 if (ListOfPrimaryAdjParticles[index_particle]->GetParticleName() == "adj_proton") {
84 E1=EminIon;
85 E2=EmaxIon;
86 }
87 if (ListOfPrimaryAdjParticles[index_particle]->GetParticleType() == "adjoint_nucleus") {
88 G4int A= ListOfPrimaryAdjParticles[index_particle]->GetAtomicMass();
89 E1=EminIon*A;
90 E2=EmaxIon*A;
91 }
92 theAdjointPrimaryGenerator->GenerateAdjointPrimaryVertex(anEvent,
93 ListOfPrimaryAdjParticles[index_particle],
94 E1,E2);
95 G4PrimaryVertex* aPrimVertex = anEvent->GetPrimaryVertex();
96
97
98 p=aPrimVertex->GetPrimary()->GetMomentum();
99 pos=aPrimVertex->GetPosition();
100 G4double pmag=p.mag();
101
102 G4double m0=ListOfPrimaryAdjParticles[index_particle]->GetPDGMass();
103 G4double ekin=std::sqrt( m0*m0 + pmag*pmag) -m0;
104
105
106 //The factor pi is to normalise the weight to the directional flux
108 G4double adjoint_weight = ComputeEnergyDistWeight(ekin,E1,E2)*adjoint_source_area*pi;
109
110 aPrimVertex->SetWeight(adjoint_weight);
111
112 last_generated_part_was_adjoint =true;
115 }
116 else {
117 //fwd particle equivalent to the last generated adjoint particle ios generated
118 G4PrimaryVertex* aPrimVertex = new G4PrimaryVertex();
119 aPrimVertex->SetPosition(pos.x(),pos.y(),pos.z());
120 aPrimVertex->SetT0(0.);
121 G4PrimaryParticle* aPrimParticle = new G4PrimaryParticle(ListOfPrimaryFwdParticles[index_particle],
122 -p.x(),-p.y(),-p.z());
123
124 aPrimVertex->SetPrimary(aPrimParticle);
125 anEvent->AddPrimaryVertex(aPrimVertex);
126 last_generated_part_was_adjoint =false;
128 }
129}
130/////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
131//
133{
134 Emin=val;
135 EminIon=val;
136}
137/////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
138//
140{
141 Emax=val;
142 EmaxIon=val;
143}
144/////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
145//
147{
148 EminIon=val;
149}
150/////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
151//
153{
154 EmaxIon=val;
155}
156/////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
157//
158G4double G4AdjointPrimaryGeneratorAction::ComputeEnergyDistWeight(G4double E ,G4double E1, G4double E2)
159{
160 // We generate N numbers of primaries with a 1/E energy law distribution.
161 // We have therefore an energy distribution function
162 // f(E)=C/E (1)
163 // with C a constant that is such that
164 // N=Integral(f(E),E1,E2)=C.std::log(E2/E1) (2)
165 // Therefore from (2) we get
166 // C=N/ std::log(E2/E1) (3)
167 // and
168 // f(E)=N/ std::log(E2/E1)/E (4)
169 //For the adjoint simulation we need a energy distribution f'(E)=1..
170 //To get that we need therefore to apply a weight to the primary
171 // W=1/f(E)=E*std::log(E2/E1)/N
172 //
173 return std::log(E2/E1)*E/G4AdjointSimManager::GetInstance()->GetNbEvtOfLastRun();
174
175}
176/////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
177//
179{
180 radius_spherical_source = radius;
181 center_spherical_source = center_pos;
182 type_of_adjoint_source ="Spherical";
183 theAdjointPrimaryGenerator->SetSphericalAdjointPrimarySource(radius,center_pos);
184}
185/////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
186//
188{
189 type_of_adjoint_source ="ExternalSurfaceOfAVolume";
190 theAdjointPrimaryGenerator->SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume(volume_name);
191}
192/////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
193//
195{
196 if (PrimariesConsideredInAdjointSim.find(particle_name) != PrimariesConsideredInAdjointSim.end()){
197 PrimariesConsideredInAdjointSim[particle_name]=true;
198 }
200}
201/////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
202//
204{
205 if (PrimariesConsideredInAdjointSim.find(particle_name) != PrimariesConsideredInAdjointSim.end()){
206 PrimariesConsideredInAdjointSim[particle_name]= false;
207 }
209}
210/////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
211//
213{
215 ListOfPrimaryFwdParticles.clear();
216 ListOfPrimaryAdjParticles.clear();
217 std::map<G4String, G4bool>::iterator iter;
218 for( iter = PrimariesConsideredInAdjointSim.begin(); iter != PrimariesConsideredInAdjointSim.end(); ++iter ) {
219 if(iter->second) {
220 G4String fwd_particle_name = iter->first;
221 if ( fwd_particle_name != "ion") {
222 G4String adj_particle_name = G4String("adj_") + fwd_particle_name;
223 ListOfPrimaryFwdParticles.push_back(theParticleTable->FindParticle(fwd_particle_name));
224 ListOfPrimaryAdjParticles.push_back(theParticleTable->FindParticle(adj_particle_name));
225 }
226 else {
227 if (fwd_ion ){
228 ion_name=fwd_ion->GetParticleName();
229 G4String adj_ion_name=G4String("adj_") +ion_name;
230 ListOfPrimaryFwdParticles.push_back(fwd_ion);
231 ListOfPrimaryAdjParticles.push_back(adj_ion);
232 }
233 else {
234 ListOfPrimaryFwdParticles.push_back(0);
235 ListOfPrimaryAdjParticles.push_back(0);
236
237 }
238 }
239 }
240 }
241}
242/////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
243//
245{
246 fwd_ion = fwdIon;
247 adj_ion = adjointIon;
249}
250
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
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 GenerateAdjointPrimaryVertex(G4Event *anEvt, G4ParticleDefinition *adj_part, G4double E1, G4double E2)
void SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume(const G4String &volume_name)
void SetAdjointTrackingMode(G4bool aBool)
void RegisterAdjointPrimaryWeight(G4double aWeight)
static G4AdjointSimManager * GetInstance()
G4PrimaryVertex * GetPrimaryVertex(G4int i=0) const
Definition: G4Event.hh:155
void AddPrimaryVertex(G4PrimaryVertex *aPrimaryVertex)
Definition: G4Event.hh:142
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
G4int first(char) const