Geant4  10.02.p03
G4AdjointPrimaryGenerator Class Reference

#include <G4AdjointPrimaryGenerator.hh>

Collaboration diagram for G4AdjointPrimaryGenerator:

Public Member Functions

 G4AdjointPrimaryGenerator ()
 
 ~G4AdjointPrimaryGenerator ()
 
void GenerateAdjointPrimaryVertex (G4Event *anEvt, G4ParticleDefinition *adj_part, G4double E1, G4double E2)
 
void GenerateFwdPrimaryVertex (G4Event *anEvt, G4ParticleDefinition *adj_part, G4double E1, G4double E2)
 
void SetSphericalAdjointPrimarySource (G4double radius, G4ThreeVector pos)
 
void SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume (const G4String &volume_name)
 
void ComputeAccumulatedDepthVectorAlongBackRay (G4ThreeVector glob_pos, G4ThreeVector direction, G4double ekin, G4ParticleDefinition *aPartDef)
 
G4double SampleDistanceAlongBackRayAndComputeWeightCorrection (G4double &weight_corr)
 

Private Member Functions

 G4AdjointPrimaryGenerator (const G4AdjointPrimaryGenerator &)
 
G4AdjointPrimaryGeneratoroperator= (const G4AdjointPrimaryGenerator &)
 

Private Attributes

G4AdjointPosOnPhysVolGeneratortheG4AdjointPosOnPhysVolGenerator
 
G4SingleParticleSourcetheSingleParticleSource
 
G4String type_of_adjoint_source
 
G4double radius_spherical_source
 
G4ThreeVector center_spherical_source
 
G4NavigatorfLinearNavigator
 
G4PhysicsOrderedFreeVectortheAccumulatedDepthVector
 

Detailed Description

Definition at line 67 of file G4AdjointPrimaryGenerator.hh.

Constructor & Destructor Documentation

◆ G4AdjointPrimaryGenerator() [1/2]

G4AdjointPrimaryGenerator::G4AdjointPrimaryGenerator ( )

Definition at line 53 of file G4AdjointPrimaryGenerator.cc.

55 {
57  type_of_adjoint_source="Spherical";
59 
64 
66 
67 }
G4SPSPosDistribution * GetPosDist() const
CLHEP::Hep3Vector G4ThreeVector
G4AdjointPosOnPhysVolGenerator * theG4AdjointPosOnPhysVolGenerator
G4SPSEneDistribution * GetEneDist() const
G4SPSAngDistribution * GetAngDist() const
G4PhysicsOrderedFreeVector * theAccumulatedDepthVector
static G4AdjointPosOnPhysVolGenerator * GetInstance()
G4SingleParticleSource * theSingleParticleSource
Here is the call graph for this function:

◆ ~G4AdjointPrimaryGenerator()

G4AdjointPrimaryGenerator::~G4AdjointPrimaryGenerator ( )

Definition at line 70 of file G4AdjointPrimaryGenerator.cc.

71 {
73 }
G4SingleParticleSource * theSingleParticleSource

◆ G4AdjointPrimaryGenerator() [2/2]

G4AdjointPrimaryGenerator::G4AdjointPrimaryGenerator ( const G4AdjointPrimaryGenerator )
private

Member Function Documentation

◆ ComputeAccumulatedDepthVectorAlongBackRay()

void G4AdjointPrimaryGenerator::ComputeAccumulatedDepthVectorAlongBackRay ( G4ThreeVector  glob_pos,
G4ThreeVector  direction,
G4double  ekin,
G4ParticleDefinition aPartDef 
)

Definition at line 154 of file G4AdjointPrimaryGenerator.cc.

162  G4ThreeVector position = glob_pos;
163  G4double safety=1.;
164  G4VPhysicalVolume* thePhysVolume =
166  G4double newStep =fLinearNavigator->ComputeStep(position,direction,1.e50,
167  safety);
170  //if (theAccumulatedCSDepthVector) delete theAccumulatedCSDepthVector;
171  //theAccumulatedCSDepthVector = new G4PhysicsOrderedFreeVector();
172 
173  G4double acc_depth=0.;
174  G4double acc_length=0.;
175  //G4double acc_cs_depth=0.;
176  //theAccumulatedCSDepthVector->InsertValues(acc_cs_depth, acc_length);
177  theAccumulatedDepthVector->InsertValues(acc_length,acc_depth);
178 
179  while (newStep > 0. && thePhysVolume) {
180  acc_length+=newStep;
181  /*
182  const G4MaterialCutsCouple* theMatCutsCouple=
183  thePhysVolume->GetLogicalVolume()->GetMaterialCutsCouple();
184 
185 
186  acc_cs_depth+=newStep*G4AdjointCSManager::GetAdjointCSManager()->GetTotalAdjointCS(aPartDef,
187  ekin,
188  theMatCutsCouple);
189  theAccumulatedCSDepthVector->InsertValues(acc_cs_depth, acc_length);*/
190 
191  acc_depth+=newStep*thePhysVolume->GetLogicalVolume()->GetMaterial()->GetDensity();
192  theAccumulatedDepthVector->InsertValues(acc_length,acc_depth);
193  position=position+newStep*direction;
194  thePhysVolume =
195  fLinearNavigator->LocateGlobalPointAndSetup(position,0,false);
196  newStep =fLinearNavigator->ComputeStep(position,direction,1.e50,
197  safety);
198  }
199 
200 
201 }
G4Material * GetMaterial() const
virtual G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
Definition: G4Navigator.cc:747
G4double GetDensity() const
Definition: G4Material.hh:180
void InsertValues(G4double energy, G4double value)
G4Navigator * GetNavigatorForTracking() const
G4PhysicsOrderedFreeVector * theAccumulatedDepthVector
static G4TransportationManager * GetTransportationManager()
virtual G4VPhysicalVolume * LocateGlobalPointAndSetup(const G4ThreeVector &point, const G4ThreeVector *direction=0, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
Definition: G4Navigator.cc:125
double G4double
Definition: G4Types.hh:76
G4LogicalVolume * GetLogicalVolume() const
Here is the call graph for this function:

◆ GenerateAdjointPrimaryVertex()

void G4AdjointPrimaryGenerator::GenerateAdjointPrimaryVertex ( G4Event anEvt,
G4ParticleDefinition adj_part,
G4double  E1,
G4double  E2 
)

Definition at line 76 of file G4AdjointPrimaryGenerator.cc.

77 {
78  if (type_of_adjoint_source == "ExternalSurfaceOfAVolume") {
79 
80  //Generate position and direction relative to the external surface of sensitive volume
81  //-------------------------------------------------------------
82 
83  G4double costh_to_normal=1.;
84  G4ThreeVector pos =G4ThreeVector(0.,0.,0.);
85  G4ThreeVector direction = G4ThreeVector(0.,0.,1.);
87  if (costh_to_normal <1.e-4) costh_to_normal =1.e-4;
88  //compute now the position along the ray backward direction
89 
92  }
93 
96 
99 
100 
101 
102 }
G4SPSPosDistribution * GetPosDist() const
CLHEP::Hep3Vector G4ThreeVector
G4AdjointPosOnPhysVolGenerator * theG4AdjointPosOnPhysVolGenerator
G4SPSEneDistribution * GetEneDist() const
void SetCentreCoords(G4ThreeVector)
G4SPSAngDistribution * GetAngDist() const
void SetParticleDefinition(G4ParticleDefinition *aParticleDefinition)
void GenerateAPositionOnTheExtSurfaceOfThePhysicalVolume(G4ThreeVector &p, G4ThreeVector &direction)
void SetParticleMomentumDirection(G4ParticleMomentum aMomentumDirection)
void GeneratePrimaryVertex(G4Event *evt)
double G4double
Definition: G4Types.hh:76
G4SingleParticleSource * theSingleParticleSource
static const G4double pos
Here is the call graph for this function:

◆ GenerateFwdPrimaryVertex()

void G4AdjointPrimaryGenerator::GenerateFwdPrimaryVertex ( G4Event anEvt,
G4ParticleDefinition adj_part,
G4double  E1,
G4double  E2 
)

Definition at line 105 of file G4AdjointPrimaryGenerator.cc.

106 {
107  if (type_of_adjoint_source == "ExternalSurfaceOfAVolume") {
108 
109  //Generate position and direction relative to the external surface of sensitive volume
110  //-------------------------------------------------------------
111 
112  G4double costh_to_normal=1.;
113  G4ThreeVector pos =G4ThreeVector(0.,0.,0.);
114  G4ThreeVector direction = G4ThreeVector(0.,0.,1.);
116  if (costh_to_normal <1.e-4) costh_to_normal =1.e-4;
119  }
120 
123 
126 }
G4SPSPosDistribution * GetPosDist() const
CLHEP::Hep3Vector G4ThreeVector
G4AdjointPosOnPhysVolGenerator * theG4AdjointPosOnPhysVolGenerator
G4SPSEneDistribution * GetEneDist() const
void SetCentreCoords(G4ThreeVector)
G4SPSAngDistribution * GetAngDist() const
void SetParticleDefinition(G4ParticleDefinition *aParticleDefinition)
void GenerateAPositionOnTheExtSurfaceOfThePhysicalVolume(G4ThreeVector &p, G4ThreeVector &direction)
void SetParticleMomentumDirection(G4ParticleMomentum aMomentumDirection)
void GeneratePrimaryVertex(G4Event *evt)
double G4double
Definition: G4Types.hh:76
G4SingleParticleSource * theSingleParticleSource
static const G4double pos
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator=()

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

◆ SampleDistanceAlongBackRayAndComputeWeightCorrection()

G4double G4AdjointPrimaryGenerator::SampleDistanceAlongBackRayAndComputeWeightCorrection ( G4double weight_corr)

Definition at line 204 of file G4AdjointPrimaryGenerator.cc.

205 {G4double rand = G4UniformRand();
207  /*
208  G4double acc_cs_depth=theAccumulatedCSDepthVector->GetEnergy(distance);
209  weight_corr=std::exp(-acc_cs_depth);*/
210  weight_corr=1.;
211  return distance;
212 }
#define G4UniformRand()
Definition: Randomize.hh:97
G4PhysicsOrderedFreeVector * theAccumulatedDepthVector
double G4double
Definition: G4Types.hh:76
G4double FindLinearEnergy(G4double rand) const
Here is the call graph for this function:

◆ SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume()

void G4AdjointPrimaryGenerator::SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume ( const G4String volume_name)

Definition at line 144 of file G4AdjointPrimaryGenerator.cc.

145 {
147  type_of_adjoint_source ="ExternalSurfaceOfAVolume";
150 }
G4SPSPosDistribution * GetPosDist() const
G4AdjointPosOnPhysVolGenerator * theG4AdjointPosOnPhysVolGenerator
G4SPSAngDistribution * GetAngDist() const
void DefinePhysicalVolume1(const G4String &aName)
G4SingleParticleSource * theSingleParticleSource
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetSphericalAdjointPrimarySource()

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

Definition at line 129 of file G4AdjointPrimaryGenerator.cc.

130 {
132  center_spherical_source = center_pos;
133  type_of_adjoint_source ="Spherical";
141 }
G4SPSPosDistribution * GetPosDist() const
static const double halfpi
Definition: G4SIunits.hh:76
void SetCentreCoords(G4ThreeVector)
G4SPSAngDistribution * GetAngDist() const
static const double pi
Definition: G4SIunits.hh:74
G4SingleParticleSource * theSingleParticleSource
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ center_spherical_source

G4ThreeVector G4AdjointPrimaryGenerator::center_spherical_source
private

Definition at line 97 of file G4AdjointPrimaryGenerator.hh.

◆ fLinearNavigator

G4Navigator* G4AdjointPrimaryGenerator::fLinearNavigator
private

Definition at line 98 of file G4AdjointPrimaryGenerator.hh.

◆ radius_spherical_source

G4double G4AdjointPrimaryGenerator::radius_spherical_source
private

Definition at line 96 of file G4AdjointPrimaryGenerator.hh.

◆ theAccumulatedDepthVector

G4PhysicsOrderedFreeVector* G4AdjointPrimaryGenerator::theAccumulatedDepthVector
private

Definition at line 99 of file G4AdjointPrimaryGenerator.hh.

◆ theG4AdjointPosOnPhysVolGenerator

G4AdjointPosOnPhysVolGenerator* G4AdjointPrimaryGenerator::theG4AdjointPosOnPhysVolGenerator
private

Definition at line 89 of file G4AdjointPrimaryGenerator.hh.

◆ theSingleParticleSource

G4SingleParticleSource* G4AdjointPrimaryGenerator::theSingleParticleSource
private

Definition at line 91 of file G4AdjointPrimaryGenerator.hh.

◆ type_of_adjoint_source

G4String G4AdjointPrimaryGenerator::type_of_adjoint_source
private

Definition at line 95 of file G4AdjointPrimaryGenerator.hh.


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