Geant4  10.02.p03
G4AdjointPrimaryGeneratorAction Class Reference

#include <G4AdjointPrimaryGeneratorAction.hh>

Inheritance diagram for G4AdjointPrimaryGeneratorAction:
Collaboration diagram for G4AdjointPrimaryGeneratorAction:

Public Member Functions

 G4AdjointPrimaryGeneratorAction ()
 
 ~G4AdjointPrimaryGeneratorAction ()
 
void GeneratePrimaries (G4Event *)
 
void SetRndmFlag (const G4String &val)
 
void SetEmin (G4double val)
 
void SetEmax (G4double val)
 
void SetEminIon (G4double val)
 
void SetEmaxIon (G4double val)
 
void SetSphericalAdjointPrimarySource (G4double radius, G4ThreeVector pos)
 
void SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume (const G4String &volume_name)
 
void ConsiderParticleAsPrimary (const G4String &particle_name)
 
void NeglectParticleAsPrimary (const G4String &particle_name)
 
void SetPrimaryIon (G4ParticleDefinition *adjointIon, G4ParticleDefinition *fwdIon)
 
void UpdateListOfPrimaryParticles ()
 
size_t GetNbOfAdjointPrimaryTypes ()
 
std::vector< G4ParticleDefinition * > * GetListOfPrimaryFwdParticles ()
 
const G4StringGetPrimaryIonName ()
 
void SetNbPrimaryFwdGammasPerEvent (G4int nb)
 
- Public Member Functions inherited from G4VUserPrimaryGeneratorAction
 G4VUserPrimaryGeneratorAction ()
 
virtual ~G4VUserPrimaryGeneratorAction ()
 

Private Member Functions

G4double ComputeEnergyDistWeight (G4double energy, G4double E1, G4double E2)
 
 G4AdjointPrimaryGeneratorAction (const G4AdjointPrimaryGeneratorAction &)
 
G4AdjointPrimaryGeneratorActionoperator= (const G4AdjointPrimaryGeneratorAction &)
 

Private Attributes

G4String rndmFlag
 
G4AdjointPrimaryGeneratortheAdjointPrimaryGenerator
 
G4double Emin
 
G4double Emax
 
G4double EminIon
 
G4double EmaxIon
 
std::vector< G4ParticleDefinition * > ListOfPrimaryFwdParticles
 
std::vector< G4ParticleDefinition * > ListOfPrimaryAdjParticles
 
std::map< G4String, G4boolPrimariesConsideredInAdjointSim
 
size_t index_particle
 
G4ThreeVector pos
 
G4ThreeVector direction
 
G4ThreeVector p
 
G4String type_of_adjoint_source
 
G4double radius_spherical_source
 
G4ThreeVector center_spherical_source
 
G4int nb_fwd_gammas_per_event
 
G4ParticleDefinitionfwd_ion
 
G4ParticleDefinitionadj_ion
 
G4String ion_name
 

Detailed Description

Definition at line 78 of file G4AdjointPrimaryGeneratorAction.hh.

Constructor & Destructor Documentation

◆ G4AdjointPrimaryGeneratorAction() [1/2]

G4AdjointPrimaryGeneratorAction::G4AdjointPrimaryGeneratorAction ( )

Definition at line 49 of file G4AdjointPrimaryGeneratorAction.cc.

50  : Emin(0.), Emax(0.), EminIon(0.), EmaxIon(0.),
51  index_particle(100000),
53  ion_name("not_defined")
54 {
56 
61 
65 }
std::vector< G4ParticleDefinition * > ListOfPrimaryAdjParticles
std::vector< G4ParticleDefinition * > ListOfPrimaryFwdParticles
std::map< G4String, G4bool > PrimariesConsideredInAdjointSim
G4AdjointPrimaryGenerator * theAdjointPrimaryGenerator

◆ ~G4AdjointPrimaryGeneratorAction()

G4AdjointPrimaryGeneratorAction::~G4AdjointPrimaryGeneratorAction ( )

Definition at line 68 of file G4AdjointPrimaryGeneratorAction.cc.

69 {
71 }
G4AdjointPrimaryGenerator * theAdjointPrimaryGenerator

◆ G4AdjointPrimaryGeneratorAction() [2/2]

G4AdjointPrimaryGeneratorAction::G4AdjointPrimaryGeneratorAction ( const G4AdjointPrimaryGeneratorAction )
private

Member Function Documentation

◆ ComputeEnergyDistWeight()

G4double G4AdjointPrimaryGeneratorAction::ComputeEnergyDistWeight ( G4double  energy,
G4double  E1,
G4double  E2 
)
private

Definition at line 245 of file G4AdjointPrimaryGeneratorAction.cc.

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 }
static G4AdjointSimManager * GetInstance()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConsiderParticleAsPrimary()

void G4AdjointPrimaryGeneratorAction::ConsiderParticleAsPrimary ( const G4String particle_name)

Definition at line 281 of file G4AdjointPrimaryGeneratorAction.cc.

282 {
283  if (PrimariesConsideredInAdjointSim.find(particle_name) != PrimariesConsideredInAdjointSim.end()){
284  PrimariesConsideredInAdjointSim[particle_name]=true;
285  }
287 }
std::map< G4String, G4bool > PrimariesConsideredInAdjointSim
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GeneratePrimaries()

void G4AdjointPrimaryGeneratorAction::GeneratePrimaries ( G4Event anEvent)
virtual

Implements G4VUserPrimaryGeneratorAction.

Definition at line 74 of file G4AdjointPrimaryGeneratorAction.cc.

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 
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 }
static G4AdjointSimManager * GetInstance()
G4ThreeVector GetPosition() const
void SetWeight(G4double w)
void AddPrimaryVertex(G4PrimaryVertex *aPrimaryVertex)
Definition: G4Event.hh:154
G4PrimaryParticle * GetPrimary(G4int i=0) const
TDirectory * dir
double weight
Definition: plottest35.C:25
int G4int
Definition: G4Types.hh:78
std::vector< G4ParticleDefinition * > ListOfPrimaryAdjParticles
Char_t n[5]
G4double ComputeEnergyDistWeight(G4double energy, G4double E1, G4double E2)
double A(double temperature)
void SetAdjointTrackingMode(G4bool aBool)
std::vector< G4ParticleDefinition * > ListOfPrimaryFwdParticles
double mag() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void GenerateFwdPrimaryVertex(G4Event *anEvt, G4ParticleDefinition *adj_part, G4double E1, G4double E2)
double x() const
G4PrimaryVertex * GetPrimaryVertex(G4int i=0) const
Definition: G4Event.hh:167
static const double pi
Definition: G4SIunits.hh:74
double y() const
G4int GetEventID() const
Definition: G4Event.hh:151
void SetT0(G4double t0)
double z() const
void SetPrimary(G4PrimaryParticle *pp)
void SetPosition(G4double x0, G4double y0, G4double z0)
G4ThreeVector GetMomentum() const
G4AdjointPrimaryGenerator * theAdjointPrimaryGenerator
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ GetListOfPrimaryFwdParticles()

std::vector<G4ParticleDefinition*>* G4AdjointPrimaryGeneratorAction::GetListOfPrimaryFwdParticles ( )
inline

Definition at line 100 of file G4AdjointPrimaryGeneratorAction.hh.

100  {
101  return &ListOfPrimaryFwdParticles;}
std::vector< G4ParticleDefinition * > ListOfPrimaryFwdParticles
Here is the caller graph for this function:

◆ GetNbOfAdjointPrimaryTypes()

size_t G4AdjointPrimaryGeneratorAction::GetNbOfAdjointPrimaryTypes ( )
inline

Definition at line 99 of file G4AdjointPrimaryGeneratorAction.hh.

99 {return ListOfPrimaryAdjParticles.size();}
std::vector< G4ParticleDefinition * > ListOfPrimaryAdjParticles
Here is the caller graph for this function:

◆ GetPrimaryIonName()

const G4String& G4AdjointPrimaryGeneratorAction::GetPrimaryIonName ( )
inline

Definition at line 102 of file G4AdjointPrimaryGeneratorAction.hh.

Here is the caller graph for this function:

◆ NeglectParticleAsPrimary()

void G4AdjointPrimaryGeneratorAction::NeglectParticleAsPrimary ( const G4String particle_name)

Definition at line 290 of file G4AdjointPrimaryGeneratorAction.cc.

291 {
292  if (PrimariesConsideredInAdjointSim.find(particle_name) != PrimariesConsideredInAdjointSim.end()){
293  PrimariesConsideredInAdjointSim[particle_name]= false;
294  }
296 }
std::map< G4String, G4bool > PrimariesConsideredInAdjointSim
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator=()

G4AdjointPrimaryGeneratorAction& G4AdjointPrimaryGeneratorAction::operator= ( const G4AdjointPrimaryGeneratorAction )
private

◆ SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume()

void G4AdjointPrimaryGeneratorAction::SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume ( const G4String volume_name)

Definition at line 274 of file G4AdjointPrimaryGeneratorAction.cc.

275 {
276  type_of_adjoint_source ="ExternalSurfaceOfAVolume";
278 }
void SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume(const G4String &volume_name)
G4AdjointPrimaryGenerator * theAdjointPrimaryGenerator
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetEmax()

void G4AdjointPrimaryGeneratorAction::SetEmax ( G4double  val)

Definition at line 226 of file G4AdjointPrimaryGeneratorAction.cc.

Here is the caller graph for this function:

◆ SetEmaxIon()

void G4AdjointPrimaryGeneratorAction::SetEmaxIon ( G4double  val)

Definition at line 239 of file G4AdjointPrimaryGeneratorAction.cc.

Here is the caller graph for this function:

◆ SetEmin()

void G4AdjointPrimaryGeneratorAction::SetEmin ( G4double  val)

Definition at line 219 of file G4AdjointPrimaryGeneratorAction.cc.

Here is the caller graph for this function:

◆ SetEminIon()

void G4AdjointPrimaryGeneratorAction::SetEminIon ( G4double  val)

Definition at line 233 of file G4AdjointPrimaryGeneratorAction.cc.

Here is the caller graph for this function:

◆ SetNbPrimaryFwdGammasPerEvent()

void G4AdjointPrimaryGeneratorAction::SetNbPrimaryFwdGammasPerEvent ( G4int  nb)
inline

Definition at line 103 of file G4AdjointPrimaryGeneratorAction.hh.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetPrimaryIon()

void G4AdjointPrimaryGeneratorAction::SetPrimaryIon ( G4ParticleDefinition adjointIon,
G4ParticleDefinition fwdIon 
)

Definition at line 337 of file G4AdjointPrimaryGeneratorAction.cc.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetRndmFlag()

void G4AdjointPrimaryGeneratorAction::SetRndmFlag ( const G4String val)
inline

Definition at line 88 of file G4AdjointPrimaryGeneratorAction.hh.

Here is the call graph for this function:

◆ SetSphericalAdjointPrimarySource()

void G4AdjointPrimaryGeneratorAction::SetSphericalAdjointPrimarySource ( G4double  radius,
G4ThreeVector  pos 
)

Definition at line 265 of file G4AdjointPrimaryGeneratorAction.cc.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ UpdateListOfPrimaryParticles()

void G4AdjointPrimaryGeneratorAction::UpdateListOfPrimaryParticles ( )

Definition at line 299 of file G4AdjointPrimaryGeneratorAction.cc.

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 }
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
G4int first(char) const
int G4int
Definition: G4Types.hh:78
std::vector< G4ParticleDefinition * > ListOfPrimaryAdjParticles
const G4String & GetParticleName() const
std::vector< G4ParticleDefinition * > ListOfPrimaryFwdParticles
std::map< G4String, G4bool > PrimariesConsideredInAdjointSim
static G4ParticleTable * GetParticleTable()
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ adj_ion

G4ParticleDefinition* G4AdjointPrimaryGeneratorAction::adj_ion
private

Definition at line 142 of file G4AdjointPrimaryGeneratorAction.hh.

◆ center_spherical_source

G4ThreeVector G4AdjointPrimaryGeneratorAction::center_spherical_source
private

Definition at line 136 of file G4AdjointPrimaryGeneratorAction.hh.

◆ direction

G4ThreeVector G4AdjointPrimaryGeneratorAction::direction
private

Definition at line 132 of file G4AdjointPrimaryGeneratorAction.hh.

◆ Emax

G4double G4AdjointPrimaryGeneratorAction::Emax
private

Definition at line 119 of file G4AdjointPrimaryGeneratorAction.hh.

◆ EmaxIon

G4double G4AdjointPrimaryGeneratorAction::EmaxIon
private

Definition at line 121 of file G4AdjointPrimaryGeneratorAction.hh.

◆ Emin

G4double G4AdjointPrimaryGeneratorAction::Emin
private

Definition at line 118 of file G4AdjointPrimaryGeneratorAction.hh.

◆ EminIon

G4double G4AdjointPrimaryGeneratorAction::EminIon
private

Definition at line 120 of file G4AdjointPrimaryGeneratorAction.hh.

◆ fwd_ion

G4ParticleDefinition* G4AdjointPrimaryGeneratorAction::fwd_ion
private

Definition at line 141 of file G4AdjointPrimaryGeneratorAction.hh.

◆ index_particle

size_t G4AdjointPrimaryGeneratorAction::index_particle
private

Definition at line 130 of file G4AdjointPrimaryGeneratorAction.hh.

◆ ion_name

G4String G4AdjointPrimaryGeneratorAction::ion_name
private

Definition at line 143 of file G4AdjointPrimaryGeneratorAction.hh.

◆ ListOfPrimaryAdjParticles

std::vector<G4ParticleDefinition*> G4AdjointPrimaryGeneratorAction::ListOfPrimaryAdjParticles
private

Definition at line 126 of file G4AdjointPrimaryGeneratorAction.hh.

◆ ListOfPrimaryFwdParticles

std::vector<G4ParticleDefinition*> G4AdjointPrimaryGeneratorAction::ListOfPrimaryFwdParticles
private

Definition at line 125 of file G4AdjointPrimaryGeneratorAction.hh.

◆ nb_fwd_gammas_per_event

G4int G4AdjointPrimaryGeneratorAction::nb_fwd_gammas_per_event
private

Definition at line 137 of file G4AdjointPrimaryGeneratorAction.hh.

◆ p

G4ThreeVector G4AdjointPrimaryGeneratorAction::p
private

Definition at line 132 of file G4AdjointPrimaryGeneratorAction.hh.

◆ pos

G4ThreeVector G4AdjointPrimaryGeneratorAction::pos
private

Definition at line 132 of file G4AdjointPrimaryGeneratorAction.hh.

◆ PrimariesConsideredInAdjointSim

std::map<G4String, G4bool> G4AdjointPrimaryGeneratorAction::PrimariesConsideredInAdjointSim
private

Definition at line 127 of file G4AdjointPrimaryGeneratorAction.hh.

◆ radius_spherical_source

G4double G4AdjointPrimaryGeneratorAction::radius_spherical_source
private

Definition at line 135 of file G4AdjointPrimaryGeneratorAction.hh.

◆ rndmFlag

G4String G4AdjointPrimaryGeneratorAction::rndmFlag
private

Definition at line 111 of file G4AdjointPrimaryGeneratorAction.hh.

◆ theAdjointPrimaryGenerator

G4AdjointPrimaryGenerator* G4AdjointPrimaryGeneratorAction::theAdjointPrimaryGenerator
private

Definition at line 114 of file G4AdjointPrimaryGeneratorAction.hh.

◆ type_of_adjoint_source

G4String G4AdjointPrimaryGeneratorAction::type_of_adjoint_source
private

Definition at line 134 of file G4AdjointPrimaryGeneratorAction.hh.


The documentation for this class was generated from the following files: