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

#include <G4AdjointForcedInteractionForGamma.hh>

+ Inheritance diagram for G4AdjointForcedInteractionForGamma:

Public Member Functions

 G4AdjointForcedInteractionForGamma (const G4String &process_name)
 
 ~G4AdjointForcedInteractionForGamma () override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &step) override
 
void RegisterAdjointComptonModel (G4VEmAdjointModel *adjModel)
 
void RegisterAdjointBremModel (G4VEmAdjointModel *adjModel)
 
void ProcessDescription (std::ostream &) const override
 
void DumpInfo () const override
 
 G4AdjointForcedInteractionForGamma (G4AdjointForcedInteractionForGamma &)=delete
 
G4AdjointForcedInteractionForGammaoperator= (const G4AdjointForcedInteractionForGamma &right)=delete
 
- Public Member Functions inherited from G4VContinuousDiscreteProcess
 G4VContinuousDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VContinuousDiscreteProcess (G4VContinuousDiscreteProcess &)
 
virtual ~G4VContinuousDiscreteProcess ()
 
G4VContinuousDiscreteProcessoperator= (const G4VContinuousDiscreteProcess &)=delete
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual const G4VProcessGetCreatorProcess () const
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4double GetContinuousStepLimit (const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety) override
 
- Protected Member Functions inherited from G4VContinuousDiscreteProcess
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 47 of file G4AdjointForcedInteractionForGamma.hh.

Constructor & Destructor Documentation

◆ G4AdjointForcedInteractionForGamma() [1/2]

G4AdjointForcedInteractionForGamma::G4AdjointForcedInteractionForGamma ( const G4String & process_name)
explicit

Definition at line 36 of file G4AdjointForcedInteractionForGamma.cc.

38 : G4VContinuousDiscreteProcess(process_name)
39 , fAdjointComptonModel(nullptr)
40 , fAdjointBremModel(nullptr)
41{
43 fParticleChange = new G4ParticleChange();
44}
static G4AdjointCSManager * GetAdjointCSManager()
G4VContinuousDiscreteProcess(const G4String &, G4ProcessType aType=fNotDefined)

Referenced by G4AdjointForcedInteractionForGamma(), and operator=().

◆ ~G4AdjointForcedInteractionForGamma()

G4AdjointForcedInteractionForGamma::~G4AdjointForcedInteractionForGamma ( )
override

Definition at line 47 of file G4AdjointForcedInteractionForGamma.cc.

48{
49 if(fParticleChange)
50 delete fParticleChange;
51}

◆ G4AdjointForcedInteractionForGamma() [2/2]

G4AdjointForcedInteractionForGamma::G4AdjointForcedInteractionForGamma ( G4AdjointForcedInteractionForGamma & )
delete

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4AdjointForcedInteractionForGamma::AlongStepDoIt ( const G4Track & track,
const G4Step & step )
overridevirtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 160 of file G4AdjointForcedInteractionForGamma.cc.

162{
163 fParticleChange->Initialize(track);
164 // Compute nb of interactions length over step length
165 G4ThreeVector position = track.GetPosition();
166 G4double stepLength = track.GetStep()->GetStepLength();
167 G4double ekin = track.GetKineticEnergy();
168 fLastAdjCS = fCSManager->GetTotalAdjointCS(track.GetDefinition(), ekin,
169 track.GetMaterialCutsCouple());
170 G4double nb_fwd_interaction_length_over_step =
171 stepLength * fCSManager->GetTotalForwardCS(G4AdjointGamma::AdjointGamma(),
172 ekin,
173 track.GetMaterialCutsCouple());
174
175 G4double nb_adj_interaction_length_over_step = stepLength * fLastAdjCS;
176 G4double fwd_survival_probability =
177 std::exp(-nb_fwd_interaction_length_over_step);
178 G4double mc_induced_survival_probability = 1.;
179
180 if(fFreeFlightGamma)
181 { // for free_flight survival probability stays 1
182 // Accumulate the number of interaction lengths during free flight of gamma
183 fTotNbAdjIntLength += nb_adj_interaction_length_over_step;
184 fAccTrackLength += stepLength;
185 }
186 else
187 {
188 G4double previous_acc_nb_adj_interaction_length = fNbAdjIntLength;
189 fNbAdjIntLength += fCSBias*nb_adj_interaction_length_over_step;
190 theNumberOfInteractionLengthLeft -= fCSBias*nb_adj_interaction_length_over_step;
191
192 // protection against rare race condition
193 if(std::abs(fTotNbAdjIntLength - previous_acc_nb_adj_interaction_length) <=
194 1.e-15)
195 {
196 mc_induced_survival_probability = 1.e50;
197 }
198 else
199 {
200 mc_induced_survival_probability =
201 std::exp(-fNbAdjIntLength) - std::exp(-fTotNbAdjIntLength);
202 mc_induced_survival_probability /=
203 (std::exp(-previous_acc_nb_adj_interaction_length) -
204 std::exp(-fTotNbAdjIntLength));
205 }
206 }
207 G4double weight_correction =
208 fwd_survival_probability / mc_induced_survival_probability;
209
210 // Caution!!!
211 // It is important to select the weight of the post_step_point as the
212 // current weight and not the weight of the track, as the weight of the track
213 // is changed after having applied all the along_step_do_it.
214 G4double new_weight =
215 weight_correction * track.GetStep()->GetPostStepPoint()->GetWeight();
216
217 fParticleChange->SetParentWeightByProcess(false);
218 fParticleChange->SetSecondaryWeightByProcess(false);
219 fParticleChange->ProposeParentWeight(new_weight);
220
221 return fParticleChange;
222}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
static G4AdjointGamma * AdjointGamma()
G4double GetWeight() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
const G4ThreeVector & GetPosition() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
const G4Step * GetStep() const
G4double theNumberOfInteractionLengthLeft

◆ BuildPhysicsTable()

void G4AdjointForcedInteractionForGamma::BuildPhysicsTable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 61 of file G4AdjointForcedInteractionForGamma.cc.

63{
64 fCSManager->BuildCrossSectionMatrices(); // it will be done just once
65 fCSManager->BuildTotalSigmaTables();
66}

◆ DumpInfo()

void G4AdjointForcedInteractionForGamma::DumpInfo ( ) const
inlineoverridevirtual

Reimplemented from G4VProcess.

Definition at line 76 of file G4AdjointForcedInteractionForGamma.hh.

G4GLOB_DLL std::ostream G4cout
void ProcessDescription(std::ostream &) const override

◆ GetContinuousStepLimit()

G4double G4AdjointForcedInteractionForGamma::GetContinuousStepLimit ( const G4Track & aTrack,
G4double previousStepSize,
G4double currentMinimumStep,
G4double & currentSafety )
overrideprotectedvirtual

Implements G4VContinuousDiscreteProcess.

Definition at line 284 of file G4AdjointForcedInteractionForGamma.cc.

286{
287 return DBL_MAX;
288}
#define DBL_MAX
Definition templates.hh:62

◆ GetMeanFreePath()

G4double G4AdjointForcedInteractionForGamma::GetMeanFreePath ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overrideprotectedvirtual

Implements G4VContinuousDiscreteProcess.

Definition at line 292 of file G4AdjointForcedInteractionForGamma.cc.

295{
296 return 0.;
297}

◆ operator=()

G4AdjointForcedInteractionForGamma & G4AdjointForcedInteractionForGamma::operator= ( const G4AdjointForcedInteractionForGamma & right)
delete

◆ PostStepDoIt()

G4VParticleChange * G4AdjointForcedInteractionForGamma::PostStepDoIt ( const G4Track & track,
const G4Step &  )
overridevirtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 86 of file G4AdjointForcedInteractionForGamma.cc.

88{
89 fParticleChange->Initialize(track);
90 // For the free flight gamma no interaction occurs but a gamma with same
91 // properties is produced for further forced interaction. It is done at the
92 // very beginning of the track so that the weight can be the same
93 if(fCopyGammaForForced)
94 {
95 G4ThreeVector theGammaMomentum = track.GetMomentum();
96 fParticleChange->AddSecondary(
97 new G4DynamicParticle(G4AdjointGamma::AdjointGamma(), theGammaMomentum));
98 fParticleChange->SetParentWeightByProcess(false);
99 fParticleChange->SetSecondaryWeightByProcess(false);
100 }
101 else
102 { // Occurrence of forced interaction
103 // Selection of the model to be called
104 G4VEmAdjointModel* theSelectedModel = nullptr;
105 G4bool is_scat_proj_to_proj_case = false;
106 G4double factor=1.;
107 if(!fAdjointComptonModel && !fAdjointBremModel)
108 return fParticleChange;
109 if(!fAdjointComptonModel)
110 {
111 theSelectedModel = fAdjointBremModel;
112 is_scat_proj_to_proj_case = false;
113 // This is needed because the results of it will be used in the post step
114 // do it weight correction inside the model
115 fAdjointBremModel->AdjointCrossSection(track.GetMaterialCutsCouple(),
116 track.GetKineticEnergy(), false);
117 }
118 else if(!fAdjointBremModel)
119 {
120 theSelectedModel = fAdjointComptonModel;
121 is_scat_proj_to_proj_case = true;
122 }
123 else
124 { // Choose the model according to a 50-50 % probability
125 G4double bremAdjCS = fAdjointBremModel->AdjointCrossSection(
126 track.GetMaterialCutsCouple(), track.GetKineticEnergy(), false);
127 if(G4UniformRand() < 0.5)
128 {
129 theSelectedModel = fAdjointBremModel;
130 is_scat_proj_to_proj_case = false;
131 factor=bremAdjCS/fLastAdjCS/0.5;
132 }
133 else
134 {
135 theSelectedModel = fAdjointComptonModel;
136 is_scat_proj_to_proj_case = true;
137 factor=(fLastAdjCS-bremAdjCS)/fLastAdjCS/0.5;
138 }
139 }
140
141 // Compute the weight correction factor
142 G4double invEffectiveAdjointCS =
143 (1. - std::exp(fNbAdjIntLength - fTotNbAdjIntLength)) / fLastAdjCS/fCSBias;
144
145 // Call the selected model without correction of the weight in the model
146 theSelectedModel->SetCorrectWeightForPostStepInModel(false);
147 theSelectedModel
149 factor*fLastAdjCS * invEffectiveAdjointCS);
150 theSelectedModel->SampleSecondaries(track, is_scat_proj_to_proj_case,
151 fParticleChange);
152 theSelectedModel->SetCorrectWeightForPostStepInModel(true);
153
154 fContinueGammaAsNewFreeFlight = true;
155 }
156 return fParticleChange;
157}
bool G4bool
Definition G4Types.hh:86
#define G4UniformRand()
Definition Randomize.hh:52
G4ThreeVector GetMomentum() const
virtual void SampleSecondaries(const G4Track &aTrack, G4bool isScatProjToProj, G4ParticleChange *fParticleChange)=0
void SetAdditionalWeightCorrectionFactorForPostStepOutsideModel(G4double factor)
virtual G4double AdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool isScatProjToProj)
void SetCorrectWeightForPostStepInModel(G4bool aBool)

◆ PostStepGetPhysicalInteractionLength()

G4double G4AdjointForcedInteractionForGamma::PostStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overridevirtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 226 of file G4AdjointForcedInteractionForGamma.cc.

228{
229 G4int step_id = track.GetCurrentStepNumber();
231 fCopyGammaForForced = false;
232 G4int track_id = track.GetTrackID();
233 fFreeFlightGamma =
234 (track_id != fLastFreeFlightTrackId + 1 || fContinueGammaAsNewFreeFlight);
235 if(fFreeFlightGamma)
236 {
237 if(step_id == 1 || fContinueGammaAsNewFreeFlight)
238 {
239 *condition = Forced;
240 // A gamma with same conditions will be generate at next post_step do it
241 // for the forced interaction
242 fCopyGammaForForced = true;
243 fLastFreeFlightTrackId = track_id;
244 fAccTrackLength = 0.;
245 fTotNbAdjIntLength = 0.;
246 fContinueGammaAsNewFreeFlight = false;
247 return 1.e-90;
248 }
249 else
250 {
251 return DBL_MAX;
252 }
253 }
254 else
255 { // compute the interaction length for forced interaction
256 if(step_id == 1)
257 {
258 fCSBias=0.000001/fTotNbAdjIntLength;
259 fTotNbAdjIntLength*=fCSBias;
260 G4double min_val = std::exp(-fTotNbAdjIntLength);
262 -std::log(min_val + G4UniformRand() * (1. - min_val));
264 fNbAdjIntLength = 0.;
265 }
266 G4VPhysicalVolume* thePostPhysVolume =
268 G4double ekin = track.GetKineticEnergy();
269 G4double postCS = 0.;
270 if(thePostPhysVolume)
271 {
272 postCS = fCSManager->GetTotalAdjointCS(
274 thePostPhysVolume->GetLogicalVolume()->GetMaterialCutsCouple());
275 }
276 if(postCS > 0.)
277 return theNumberOfInteractionLengthLeft / postCS /fCSBias;
278 else
279 return DBL_MAX;
280 }
281}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
@ Forced
int G4int
Definition G4Types.hh:85
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4VPhysicalVolume * GetPhysicalVolume() const
G4StepPoint * GetPreStepPoint() const
G4int GetTrackID() const
G4int GetCurrentStepNumber() const
G4LogicalVolume * GetLogicalVolume() const
G4double theInitialNumberOfInteractionLength

◆ ProcessDescription()

void G4AdjointForcedInteractionForGamma::ProcessDescription ( std::ostream & out) const
overridevirtual

Reimplemented from G4VProcess.

Definition at line 54 of file G4AdjointForcedInteractionForGamma.cc.

56{
57 out << "Forced interaction for gamma.\n";
58}

Referenced by DumpInfo().

◆ RegisterAdjointBremModel()

void G4AdjointForcedInteractionForGamma::RegisterAdjointBremModel ( G4VEmAdjointModel * adjModel)
inline

Definition at line 70 of file G4AdjointForcedInteractionForGamma.hh.

71 {
72 fAdjointBremModel = adjModel;
73 }

◆ RegisterAdjointComptonModel()

void G4AdjointForcedInteractionForGamma::RegisterAdjointComptonModel ( G4VEmAdjointModel * adjModel)
inline

Definition at line 65 of file G4AdjointForcedInteractionForGamma.hh.

66 {
67 fAdjointComptonModel = adjModel;
68 }

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