Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4AdjointPrimaryGenerator Class Reference

#include <G4AdjointPrimaryGenerator.hh>

Public Member Functions

 G4AdjointPrimaryGenerator ()
 
 ~G4AdjointPrimaryGenerator ()
 
 G4AdjointPrimaryGenerator (const G4AdjointPrimaryGenerator &)=delete
 
G4AdjointPrimaryGeneratoroperator= (const G4AdjointPrimaryGenerator &)=delete
 
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 &v_name)
 
void ComputeAccumulatedDepthVectorAlongBackRay (G4ThreeVector glob_pos, G4ThreeVector direction, G4double ekin, G4ParticleDefinition *aPDef)
 
G4double SampleDistanceAlongBackRayAndComputeWeightCorrection (G4double &weight_corr)
 

Detailed Description

Definition at line 58 of file G4AdjointPrimaryGenerator.hh.

Constructor & Destructor Documentation

◆ G4AdjointPrimaryGenerator() [1/2]

G4AdjointPrimaryGenerator::G4AdjointPrimaryGenerator ( )

Definition at line 47 of file G4AdjointPrimaryGenerator.cc.

48{
49 center_spherical_source = G4ThreeVector(0.,0.,0.);
50 type_of_adjoint_source="Spherical";
51 theSingleParticleSource = new G4SingleParticleSource();
52
53 theSingleParticleSource->GetEneDist()->SetEnergyDisType("Pow");
54 theSingleParticleSource->GetEneDist()->SetAlpha(-1.);
55 theSingleParticleSource->GetPosDist()->SetPosDisType("Point");
56 theSingleParticleSource->GetAngDist()->SetAngDistType("planar");
57
58 theG4AdjointPosOnPhysVolGenerator = G4AdjointPosOnPhysVolGenerator::GetInstance();
59}
CLHEP::Hep3Vector G4ThreeVector
static G4AdjointPosOnPhysVolGenerator * GetInstance()
void SetAngDistType(const G4String &)
void SetEnergyDisType(const G4String &)
void SetPosDisType(const G4String &)
G4SPSAngDistribution * GetAngDist() const
G4SPSEneDistribution * GetEneDist() const
G4SPSPosDistribution * GetPosDist() const

◆ ~G4AdjointPrimaryGenerator()

G4AdjointPrimaryGenerator::~G4AdjointPrimaryGenerator ( )

Definition at line 63 of file G4AdjointPrimaryGenerator.cc.

64{
65 delete theSingleParticleSource;
66}

◆ G4AdjointPrimaryGenerator() [2/2]

G4AdjointPrimaryGenerator::G4AdjointPrimaryGenerator ( const G4AdjointPrimaryGenerator )
delete

Member Function Documentation

◆ ComputeAccumulatedDepthVectorAlongBackRay()

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

Definition at line 161 of file G4AdjointPrimaryGenerator.cc.

165{
166 if (fLinearNavigator == nullptr)
167 {
170 }
171 G4ThreeVector position = glob_pos;
172 G4double safety=1.;
173 G4VPhysicalVolume* thePhysVolume =
174 fLinearNavigator->LocateGlobalPointAndSetup(position);
175 G4double newStep = fLinearNavigator->ComputeStep(position,direction,1.e50,
176 safety);
177 if (theAccumulatedDepthVector != nullptr) {delete theAccumulatedDepthVector;}
178 theAccumulatedDepthVector = new G4PhysicsOrderedFreeVector();
179
180 G4double acc_depth=0.;
181 G4double acc_length=0.;
182 theAccumulatedDepthVector->InsertValues(acc_length,acc_depth);
183
184 while (newStep > 0. && thePhysVolume != nullptr)
185 {
186 acc_length+=newStep;
187 acc_depth+=newStep*thePhysVolume->GetLogicalVolume()
188 ->GetMaterial()->GetDensity();
189 theAccumulatedDepthVector->InsertValues(acc_length,acc_depth);
190 position=position+newStep*direction;
191 thePhysVolume = fLinearNavigator
192 ->LocateGlobalPointAndSetup(position,nullptr,false);
193 newStep = fLinearNavigator->ComputeStep(position,direction,1.e50,safety);
194 }
195}
double G4double
Definition: G4Types.hh:83
G4Material * GetMaterial() const
G4double GetDensity() const
Definition: G4Material.hh:178
virtual G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
Definition: G4Navigator.cc:764
virtual G4VPhysicalVolume * LocateGlobalPointAndSetup(const G4ThreeVector &point, const G4ThreeVector *direction=nullptr, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
Definition: G4Navigator.cc:126
void InsertValues(G4double energy, G4double value)
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
G4LogicalVolume * GetLogicalVolume() const

◆ GenerateAdjointPrimaryVertex()

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

Definition at line 70 of file G4AdjointPrimaryGenerator.cc.

73{
74 if (type_of_adjoint_source == "ExternalSurfaceOfAVolume")
75 {
76 // Generate position and direction relative to the external surface
77 // of sensitive volume
78
79 G4double costh_to_normal=1.;
81 G4ThreeVector direction = G4ThreeVector(0.,0.,1.);
82 theG4AdjointPosOnPhysVolGenerator
84 costh_to_normal);
85 if (costh_to_normal <1.e-4) { costh_to_normal = 1.e-4; }
86
87 // compute now the position along the ray backward direction
88 //
89 theSingleParticleSource->GetAngDist()
90 ->SetParticleMomentumDirection(-direction);
91 theSingleParticleSource->GetPosDist()->SetCentreCoords(pos);
92 }
93
94 theSingleParticleSource->GetEneDist()->SetEmin(E1);
95 theSingleParticleSource->GetEneDist()->SetEmax(E2);
96
97 theSingleParticleSource->SetParticleDefinition(adj_part);
98 theSingleParticleSource->GeneratePrimaryVertex(anEvent);
99}
void GenerateAPositionOnTheExtSurfaceOfThePhysicalVolume(G4ThreeVector &p, G4ThreeVector &direction)
void SetParticleMomentumDirection(const G4ParticleMomentum &aMomDirection)
void SetCentreCoords(const G4ThreeVector &)
void SetParticleDefinition(G4ParticleDefinition *aParticleDefinition)
void GeneratePrimaryVertex(G4Event *evt)

◆ GenerateFwdPrimaryVertex()

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

Definition at line 103 of file G4AdjointPrimaryGenerator.cc.

106{
107 if (type_of_adjoint_source == "ExternalSurfaceOfAVolume")
108 {
109 // Generate position and direction relative to the external surface
110 // of sensitive volume
111
112 G4double costh_to_normal=1.;
114 G4ThreeVector direction = G4ThreeVector(0.,0.,1.);
115 theG4AdjointPosOnPhysVolGenerator
117 costh_to_normal);
118 if (costh_to_normal <1.e-4) { costh_to_normal =1.e-4; }
119 theSingleParticleSource->GetAngDist()
120 ->SetParticleMomentumDirection(direction);
121 theSingleParticleSource->GetPosDist()->SetCentreCoords(pos);
122 }
123
124 theSingleParticleSource->GetEneDist()->SetEmin(E1);
125 theSingleParticleSource->GetEneDist()->SetEmax(E2);
126
127 theSingleParticleSource->SetParticleDefinition(fwd_part);
128 theSingleParticleSource->GeneratePrimaryVertex(anEvent);
129}

Referenced by G4AdjointPrimaryGeneratorAction::GeneratePrimaries().

◆ operator=()

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

◆ SampleDistanceAlongBackRayAndComputeWeightCorrection()

G4double G4AdjointPrimaryGenerator::SampleDistanceAlongBackRayAndComputeWeightCorrection ( G4double weight_corr)

Definition at line 199 of file G4AdjointPrimaryGenerator.cc.

201{
202 G4double rand = G4UniformRand();
203 G4double distance = theAccumulatedDepthVector->FindLinearEnergy(rand);
204 weight_corr=1.;
205 return distance;
206}
#define G4UniformRand()
Definition: Randomize.hh:52
G4double FindLinearEnergy(G4double rand) const

◆ SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume()

void G4AdjointPrimaryGenerator::SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume ( const G4String v_name)

Definition at line 150 of file G4AdjointPrimaryGenerator.cc.

152{
153 theG4AdjointPosOnPhysVolGenerator->DefinePhysicalVolume1(volume_name);
154 type_of_adjoint_source ="ExternalSurfaceOfAVolume";
155 theSingleParticleSource->GetPosDist()->SetPosDisType("Point");
156 theSingleParticleSource->GetAngDist()->SetAngDistType("planar");
157}
void DefinePhysicalVolume1(const G4String &aName)

Referenced by G4AdjointPrimaryGeneratorAction::SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume().

◆ SetSphericalAdjointPrimarySource()

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

Definition at line 133 of file G4AdjointPrimaryGenerator.cc.

135{
136 radius_spherical_source = radius;
137 center_spherical_source = center_pos;
138 type_of_adjoint_source = "Spherical";
139 theSingleParticleSource->GetPosDist()->SetPosDisType("Surface");
140 theSingleParticleSource->GetPosDist()->SetPosDisShape("Sphere");
141 theSingleParticleSource->GetPosDist()->SetCentreCoords(center_pos);
142 theSingleParticleSource->GetPosDist()->SetRadius(radius);
143 theSingleParticleSource->GetAngDist()->SetAngDistType("cos");
144 theSingleParticleSource->GetAngDist()->SetMaxTheta(pi);
145 theSingleParticleSource->GetAngDist()->SetMinTheta(halfpi);
146}
void SetPosDisShape(const G4String &)

Referenced by G4AdjointPrimaryGeneratorAction::SetSphericalAdjointPrimarySource().


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