Geant4  10.02.p01
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: G4AdjointPrimaryGeneratorAction.cc 81773 2014-06-05 08:35:38Z gcosmo $
27 //
29 // Class Name: G4AdjointPrimaryGeneratorAction
30 // Author: L. Desorgher
31 // Organisation: SpaceIT GmbH
32 // Contract: ESA contract 21435/08/NL/AT
33 // Customer: ESA/ESTEC
35 
37 
38 #include "G4PhysicalConstants.hh"
39 #include "G4Event.hh"
40 #include "G4ParticleTable.hh"
41 #include "G4ParticleDefinition.hh"
42 #include "G4ParticleTable.hh"
43 #include "G4AdjointSimManager.hh"
45 #include "G4Gamma.hh"
46 
48 //
50  : Emin(0.), Emax(0.), EminIon(0.), EmaxIon(0.),
51  index_particle(100000),
52  radius_spherical_source(0.), fwd_ion(0), adj_ion(0),
53  ion_name("not_defined")
54 {
56 
61 
65 }
67 //
69 {
71 }
73 //
75 { G4int evt_id=anEvent->GetEventID();
76  size_t n=ListOfPrimaryAdjParticles.size();
77  index_particle=size_t(evt_id)-n*(size_t(evt_id)/n);
78 
79 
80  G4double E1=Emin;
81  G4double E2=Emax;
82  if (!ListOfPrimaryAdjParticles[index_particle]) UpdateListOfPrimaryParticles();//ion has not been created yet
83 
84  if (ListOfPrimaryAdjParticles[index_particle]->GetParticleName() == "adj_proton") {
85  E1=EminIon;
86  E2=EmaxIon;
87  }
88  if (ListOfPrimaryAdjParticles[index_particle]->GetParticleType() == "adjoint_nucleus") {
90  E1=EminIon*A;
91  E2=EmaxIon*A;
92  }
94  G4PrimaryVertex* fwdPrimVertex = anEvent->GetPrimaryVertex();
95 
96  p=fwdPrimVertex->GetPrimary()->GetMomentum();
97  pos=fwdPrimVertex->GetPosition();
98  G4double pmag=p.mag();
100  G4double ekin=std::sqrt( m0*m0 + pmag*pmag) -m0;
101 
102  G4double weight_correction=1.;
103  //For gamma generate the particle along the backward ray
104  G4ThreeVector dir=-p/p.mag();
105 
106  /*if (ListOfPrimaryAdjParticles[index_particle]->GetParticleName() == "adj_gamma"){
107 
108  theAdjointPrimaryGenerator
109  ->ComputeAccumulatedDepthVectorAlongBackRay(pos,dir,ekin,ListOfPrimaryAdjParticles[index_particle]);
110 
111 
112  G4double distance = theAdjointPrimaryGenerator
113  ->SampleDistanceAlongBackRayAndComputeWeightCorrection(weight_correction);
114 
115  //pos=pos+dir*distance;
116  //fwdPrimVertex->SetPosition(pos[0],pos[1],pos[2]);
117  }
118  */
119  weight_correction=1.;
120 
122  G4double weight = (1./nb_fwd_gammas_per_event);
123  fwdPrimVertex->SetWeight(weight);
124  for (int i=0;i<nb_fwd_gammas_per_event-1;i++){
125  G4PrimaryVertex* newFwdPrimVertex = new G4PrimaryVertex();
126  newFwdPrimVertex->SetPosition(pos.x(),pos.y(),pos.z());
127  newFwdPrimVertex->SetT0(0.);
128  G4PrimaryParticle* aPrimParticle = new G4PrimaryParticle(ListOfPrimaryFwdParticles[index_particle],
129  p.x(),p.y(),p.z());
130  newFwdPrimVertex->SetPrimary(aPrimParticle);
131  newFwdPrimVertex->SetWeight(weight);
132  anEvent->AddPrimaryVertex(newFwdPrimVertex);
133  }
134  }
135 
136 
137  G4PrimaryVertex* adjPrimVertex = new G4PrimaryVertex();
138  adjPrimVertex->SetPosition(pos.x(),pos.y(),pos.z());
139  adjPrimVertex->SetT0(0.);
140  G4PrimaryParticle* aPrimParticle = new G4PrimaryParticle(ListOfPrimaryAdjParticles[index_particle],
141  -p.x(),-p.y(),-p.z());
142 
143  adjPrimVertex->SetPrimary(aPrimParticle);
144  anEvent->AddPrimaryVertex(adjPrimVertex);
145 
146  //The factor pi is to normalise the weight to the directional flux
148  G4double adjoint_weight = weight_correction*ComputeEnergyDistWeight(ekin,E1,E2)*adjoint_source_area*pi;
149  if (ListOfPrimaryFwdParticles[index_particle] ==G4Gamma::Gamma()) adjoint_weight = adjoint_weight/3.;
150  adjPrimVertex->SetWeight(adjoint_weight);
151 
152 
153 
154 
156 
157 
158 
159  /* if ( !last_generated_part_was_adjoint ) {
160 
161  index_particle++;
162  if (index_particle >= ListOfPrimaryAdjParticles.size()) index_particle =0;
163 
164 
165  G4double E1=Emin;
166  G4double E2=Emax;
167  if (!ListOfPrimaryAdjParticles[index_particle]) UpdateListOfPrimaryParticles();//ion has not been created yet
168 
169  if (ListOfPrimaryAdjParticles[index_particle]->GetParticleName() == "adj_proton") {
170  E1=EminIon;
171  E2=EmaxIon;
172  }
173  if (ListOfPrimaryAdjParticles[index_particle]->GetParticleType() == "adjoint_nucleus") {
174  G4int A= ListOfPrimaryAdjParticles[index_particle]->GetAtomicMass();
175  E1=EminIon*A;
176  E2=EmaxIon*A;
177  }
178  theAdjointPrimaryGenerator->GenerateAdjointPrimaryVertex(anEvent,
179  ListOfPrimaryAdjParticles[index_particle],
180  E1,E2);
181  G4PrimaryVertex* aPrimVertex = anEvent->GetPrimaryVertex();
182 
183 
184  p=aPrimVertex->GetPrimary()->GetMomentum();
185  pos=aPrimVertex->GetPosition();
186  G4double pmag=p.mag();
187 
188  G4double m0=ListOfPrimaryAdjParticles[index_particle]->GetPDGMass();
189  G4double ekin=std::sqrt( m0*m0 + pmag*pmag) -m0;
190 
191 
192  //The factor pi is to normalise the weight to the directional flux
193  G4double adjoint_source_area = G4AdjointSimManager::GetInstance()->GetAdjointSourceArea();
194  G4double adjoint_weight = ComputeEnergyDistWeight(ekin,E1,E2)*adjoint_source_area*pi;
195 
196  aPrimVertex->SetWeight(adjoint_weight);
197 
198  last_generated_part_was_adjoint =true;
199  G4AdjointSimManager::GetInstance()->SetAdjointTrackingMode(true);
200  G4AdjointSimManager::GetInstance()->RegisterAdjointPrimaryWeight(adjoint_weight);
201  }
202  else {
203  //fwd particle equivalent to the last generated adjoint particle ios generated
204  G4PrimaryVertex* aPrimVertex = new G4PrimaryVertex();
205  aPrimVertex->SetPosition(pos.x(),pos.y(),pos.z());
206  aPrimVertex->SetT0(0.);
207  G4PrimaryParticle* aPrimParticle = new G4PrimaryParticle(ListOfPrimaryFwdParticles[index_particle],
208  -p.x(),-p.y(),-p.z());
209 
210  aPrimVertex->SetPrimary(aPrimParticle);
211  anEvent->AddPrimaryVertex(aPrimVertex);
212  last_generated_part_was_adjoint =false;
213  G4AdjointSimManager::GetInstance()->SetAdjointTrackingMode(false);
214  */
215 
216 }
218 //
220 {
221  Emin=val;
222  EminIon=val;
223 }
225 //
227 {
228  Emax=val;
229  EmaxIon=val;
230 }
232 //
234 {
235  EminIon=val;
236 }
238 //
240 {
241  EmaxIon=val;
242 }
244 //
246 {
247  // We generate N numbers of primaries with a 1/E energy law distribution.
248  // We have therefore an energy distribution function
249  // f(E)=C/E (1)
250  // with C a constant that is such that
251  // N=Integral(f(E),E1,E2)=C.std::log(E2/E1) (2)
252  // Therefore from (2) we get
253  // C=N/ std::log(E2/E1) (3)
254  // and
255  // f(E)=N/ std::log(E2/E1)/E (4)
256  //For the adjoint simulation we need a energy distribution f'(E)=1..
257  //To get that we need therefore to apply a weight to the primary
258  // W=1/f(E)=E*std::log(E2/E1)/N
259  //
260  return std::log(E2/E1)*E/G4AdjointSimManager::GetInstance()->GetNbEvtOfLastRun();
261 
262 }
264 //
266 {
267  radius_spherical_source = radius;
268  center_spherical_source = center_pos;
269  type_of_adjoint_source ="Spherical";
271 }
273 //
275 {
276  type_of_adjoint_source ="ExternalSurfaceOfAVolume";
278 }
280 //
282 {
283  if (PrimariesConsideredInAdjointSim.find(particle_name) != PrimariesConsideredInAdjointSim.end()){
284  PrimariesConsideredInAdjointSim[particle_name]=true;
285  }
287 }
289 //
291 {
292  if (PrimariesConsideredInAdjointSim.find(particle_name) != PrimariesConsideredInAdjointSim.end()){
293  PrimariesConsideredInAdjointSim[particle_name]= false;
294  }
296 }
298 //
300 {
304  std::map<G4String, G4bool>::iterator iter;
305  for( iter = PrimariesConsideredInAdjointSim.begin(); iter != PrimariesConsideredInAdjointSim.end(); ++iter ) {
306  if(iter->second) {
307  G4String fwd_particle_name = iter->first;
308  if ( fwd_particle_name != "ion") {
309  G4String adj_particle_name = G4String("adj_") + fwd_particle_name;
310  ListOfPrimaryFwdParticles.push_back(theParticleTable->FindParticle(fwd_particle_name));
311  ListOfPrimaryAdjParticles.push_back(theParticleTable->FindParticle(adj_particle_name));
312  if ( fwd_particle_name == "gamma") {
313  for (G4int i=0;i<2;i++){
314  ListOfPrimaryFwdParticles.push_back(theParticleTable->FindParticle(fwd_particle_name));
315  ListOfPrimaryAdjParticles.push_back(theParticleTable->FindParticle(adj_particle_name));
316  }
317  }
318  }
319  else {
320  if (fwd_ion ){
322  G4String adj_ion_name=G4String("adj_") +ion_name;
325  }
326  else {
327  ListOfPrimaryFwdParticles.push_back(0);
328  ListOfPrimaryAdjParticles.push_back(0);
329 
330  }
331  }
332  }
333  }
334 }
336 //
338 {
339  fwd_ion = fwdIon;
340  adj_ion = adjointIon;
342 }
343 
static G4AdjointSimManager * GetInstance()
void SetSphericalAdjointPrimarySource(G4double radius, G4ThreeVector pos)
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
G4ThreeVector GetMomentum() const
G4ThreeVector GetPosition() const
CLHEP::Hep3Vector G4ThreeVector
G4int first(char) const
void SetWeight(G4double w)
void AddPrimaryVertex(G4PrimaryVertex *aPrimaryVertex)
Definition: G4Event.hh:154
void SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume(const G4String &volume_name)
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4int GetEventID() const
Definition: G4Event.hh:151
void SetPrimaryIon(G4ParticleDefinition *adjointIon, G4ParticleDefinition *fwdIon)
G4PrimaryParticle * GetPrimary(G4int i=0) const
std::vector< G4ParticleDefinition * > ListOfPrimaryAdjParticles
G4double ComputeEnergyDistWeight(G4double energy, G4double E1, G4double E2)
void ConsiderParticleAsPrimary(const G4String &particle_name)
double A(double temperature)
void SetAdjointTrackingMode(G4bool aBool)
std::vector< G4ParticleDefinition * > ListOfPrimaryFwdParticles
const G4int n
G4PrimaryVertex * GetPrimaryVertex(G4int i=0) const
Definition: G4Event.hh:167
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void GenerateFwdPrimaryVertex(G4Event *anEvt, G4ParticleDefinition *adj_part, G4double E1, G4double E2)
std::map< G4String, G4bool > PrimariesConsideredInAdjointSim
static const double pi
Definition: G4SIunits.hh:74
static G4ParticleTable * GetParticleTable()
void NeglectParticleAsPrimary(const G4String &particle_name)
void SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume(const G4String &volume_name)
void SetT0(G4double t0)
static const G4double Emin
void SetPrimary(G4PrimaryParticle *pp)
void SetPosition(G4double x0, G4double y0, G4double z0)
static const G4double Emax
G4AdjointPrimaryGenerator * theAdjointPrimaryGenerator
double G4double
Definition: G4Types.hh:76
void SetSphericalAdjointPrimarySource(G4double radius, G4ThreeVector pos)