Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4AdjointForcedInteractionForGamma.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: G4AdjointForcedInteractionForGamma.cc 87443 2014-12-04 12:26:31Z gunter $
27 //
29 #include "G4SystemOfUnits.hh"
30 #include "G4AdjointCSManager.hh"
31 #include "G4AdjointCSMatrix.hh"
32 #include "G4VEmAdjointModel.hh"
33 #include "G4MaterialCutsCouple.hh"
34 #include "G4ParticleChange.hh"
35 #include "G4AdjointGamma.hh"
36 
37 
40  G4VContinuousDiscreteProcess(process_name),theAdjointComptonModel(0),theAdjointBremModel(0)
41 { theAdjointCSManager = G4AdjointCSManager::GetAdjointCSManager();
42  fParticleChange=new G4ParticleChange();
43  lastAdjCS=0.;
44  trackid = nstep = 0;
45  is_free_flight_gamma = false;
46  copy_gamma_for_forced_interaction = false;
47  last_free_flight_trackid=1000;
48 
49  theAdjointComptonModel =0;
50  theAdjointBremModel=0;
51 
52  acc_track_length=0.;
53  acc_nb_adj_interaction_length=0.;
54  acc_nb_fwd_interaction_length=0.;
55  total_acc_nb_adj_interaction_length=0.;
56  total_acc_nb_fwd_interaction_length=0.;
57  continue_gamma_as_new_free_flight =false;
58 
59 
60 
61 }
63 //
66 { if (fParticleChange) delete fParticleChange;
67 }
69 //
71 {;
72 }
74 //
76 {
77  theAdjointCSManager->BuildCrossSectionMatrices(); //do not worry it will be done just once
78  theAdjointCSManager->BuildTotalSigmaTables();
79 }
80 //Note on weight correction for forced interaction
81 //For the forced interaction applied here we do use a truncated exponential law for the probability of survival
82 //over a fixed total length. This is done by using a linear transformation of the non biased probability survival
83 //In mathematic this writes
84 //P'(x)=C1P(x)+C2
85 //With P(x)=exp(-sum(sigma_ixi)) x and L can cross different volumes with different cross section sigma.
86 //For forced interaction we get the following limit conditions
87 //P'(L)=0 P'(0)=1 (L can be used over different volumes)
88 //From simple solving of linear equation we get
89 //C1=1/(1-P(L)) et C2=-P(L)/(1-P(L))
90 //P'(x)=(P(x)-P(L))/(1-P(L))
91 //For the probability over a step x1 to x2
92 //P'(x1->x2)=P'(x2)/P'(x1)
93 //The effective cross section is defined -d(P'(x))/dx/P'(x)
94 //We get therefore
95 //sigma_eff=C1sigmaP(x)/(C1P(x)+C2)=sigmaP(x)/(P(x)+C2/C1)=sigmaP(x)/(P(x)-P(L))=sigma/(1-P(L)/P(x))
97 //
99 { fParticleChange->Initialize(track);
100  //For the free flight gamma no interaction occur but a gamma with same property is
101  //produces for further forced interaction
102  //It is done at the very beginning of the track such that the weight can be the same
103  if (copy_gamma_for_forced_interaction) {
104  G4ThreeVector theGammaMomentum = track.GetMomentum();
105  fParticleChange->AddSecondary(new G4DynamicParticle(G4AdjointGamma::AdjointGamma(),theGammaMomentum));
106  fParticleChange->SetParentWeightByProcess(false);
107  fParticleChange->SetSecondaryWeightByProcess(false);
108  }
109  else { //Occurrence of forced interaction
110 
111  //Selection of the model to be called
112  G4VEmAdjointModel* theSelectedModel =0;
113  G4bool is_scat_proj_to_proj_case=false;
114  if (!theAdjointComptonModel && !theAdjointBremModel) return fParticleChange;
115  if (!theAdjointComptonModel) {
116  theSelectedModel = theAdjointBremModel;
117  is_scat_proj_to_proj_case=false;
118  //This is needed because the results of it will be used in the post step do it weight correction inside the model
119  theAdjointBremModel->AdjointCrossSection(
120  track.GetMaterialCutsCouple(),track.GetKineticEnergy(), false);
121 
122  }
123  else if (!theAdjointBremModel) {
124  theSelectedModel = theAdjointComptonModel;
125  is_scat_proj_to_proj_case=true;
126  }
127  else { //Choose the model according to cross sections
128 
129  G4double bremAdjCS = theAdjointBremModel->AdjointCrossSection(
130  track.GetMaterialCutsCouple(),track.GetKineticEnergy(), false);
131  if (G4UniformRand()*lastAdjCS<bremAdjCS) {
132  theSelectedModel = theAdjointBremModel;
133  is_scat_proj_to_proj_case=false;
134  }
135  else {
136  theSelectedModel = theAdjointComptonModel;
137  is_scat_proj_to_proj_case=true;
138  }
139  }
140 
141  //Compute the weight correction factor
142  G4double one_over_effectiveAdjointCS= (1.-std::exp(acc_nb_adj_interaction_length-total_acc_nb_adj_interaction_length))/lastAdjCS;
143  G4double weight_correction_factor = lastAdjCS*one_over_effectiveAdjointCS;
144  //G4cout<<"Weight correction factor start "<<weight_correction_factor<<std::endl;
145  //Call the selected model without correction of the weight in the model
146  theSelectedModel->SetCorrectWeightForPostStepInModel(false);
147  theSelectedModel->SetAdditionalWeightCorrectionFactorForPostStepOutsideModel(weight_correction_factor);
148  theSelectedModel->SampleSecondaries(track,is_scat_proj_to_proj_case,fParticleChange);
149  theSelectedModel->SetCorrectWeightForPostStepInModel(true);
150 
151  continue_gamma_as_new_free_flight =true;
152  }
153  return fParticleChange;
154 }
156 //
158 { fParticleChange->Initialize(track);
159  //Compute nb of interactions length over step length
161  G4double stepLength = track.GetStep()->GetStepLength();
162  G4double ekin = track.GetKineticEnergy();
163  G4double nb_fwd_interaction_length_over_step=0.;
164  G4double nb_adj_interaction_length_over_step=0.;
167  ekin,track.GetMaterialCutsCouple());
168  nb_fwd_interaction_length_over_step = stepLength*lastFwdCS;
169  nb_adj_interaction_length_over_step = stepLength*lastAdjCS;
170  G4double fwd_survival_probability=std::exp(-nb_fwd_interaction_length_over_step);
171  G4double mc_induced_survival_probability=1.;
172 
173  if (is_free_flight_gamma) { //for free_flight survival probability stays 1
174  //Accumulate the number of interaction lengths during free flight of gamma
175  total_acc_nb_fwd_interaction_length+=nb_fwd_interaction_length_over_step;
176  total_acc_nb_adj_interaction_length+=nb_adj_interaction_length_over_step;
177  acc_track_length+=stepLength;
178  }
179  else {
180  G4double previous_acc_nb_adj_interaction_length =acc_nb_adj_interaction_length;
181  acc_nb_fwd_interaction_length+=nb_fwd_interaction_length_over_step;
182  acc_nb_adj_interaction_length+=nb_adj_interaction_length_over_step;
183  theNumberOfInteractionLengthLeft-=nb_adj_interaction_length_over_step;
184 
185  //Following condition to remove very rare FPE issue
186  if (total_acc_nb_adj_interaction_length <= 1.e-50 && theNumberOfInteractionLengthLeft<=1.e-50) { //condition added to avoid FPE issue
187  mc_induced_survival_probability = 1.e50;
188  }
189  else {
190  mc_induced_survival_probability= std::exp(-acc_nb_adj_interaction_length)-std::exp(-total_acc_nb_adj_interaction_length);
191  mc_induced_survival_probability=mc_induced_survival_probability/(std::exp(-previous_acc_nb_adj_interaction_length)-std::exp(-total_acc_nb_adj_interaction_length));
192  }
193  }
194  G4double weight_correction = fwd_survival_probability/mc_induced_survival_probability;
195  //weight_correction = 1.;
196  //Caution!!!
197  // It is important to select the weight of the post_step_point
198  // as the current weight and not the weight of the track, as t
199  // the weight of the track is changed after having applied all
200  // the along_step_do_it.
201  G4double new_weight=weight_correction*track.GetStep()->GetPostStepPoint()->GetWeight();
202 /*
203  G4cout<<"New weight "<<new_weight<<std::endl;
204  G4cout<<"Weight correction "<<weight_correction<<std::endl;
205  */
206 
207  fParticleChange->SetParentWeightByProcess(false);
208  fParticleChange->SetSecondaryWeightByProcess(false);
209  fParticleChange->ProposeParentWeight(new_weight);
210 
211  return fParticleChange;
212 }
214 //
216  const G4Track& track,
217  G4double ,
219 { G4int step_id = track.GetCurrentStepNumber();
220  *condition = NotForced;
221  copy_gamma_for_forced_interaction = false;
222  G4int track_id=track.GetTrackID();
223  is_free_flight_gamma = (track_id != last_free_flight_trackid+1 || continue_gamma_as_new_free_flight);
224  if (is_free_flight_gamma) {
225  if (step_id == 1 || continue_gamma_as_new_free_flight) {
226  *condition=Forced;
227  //A gamma with same conditions will be generate at next post_step do it for the forced interaction
228  copy_gamma_for_forced_interaction = true;
229  last_free_flight_trackid = track_id;
230  acc_track_length=0.;
231  total_acc_nb_adj_interaction_length=0.;
232  total_acc_nb_fwd_interaction_length=0.;
233  continue_gamma_as_new_free_flight=false;
234  return 1.e-90;
235  }
236  else {
237  //Computation of accumulated length for
238  return DBL_MAX;
239  }
240  }
241  else { //compute the interaction length for forced interaction
242  if (step_id ==1) {
243  G4double min_val= std::exp(-total_acc_nb_adj_interaction_length);
244  theNumberOfInteractionLengthLeft = -std::log( min_val+G4UniformRand()*(1.-min_val));
246  acc_nb_adj_interaction_length=0.;
247  acc_nb_fwd_interaction_length=0.;
248  }
249  G4VPhysicalVolume* thePostPhysVolume = track.GetStep()->GetPreStepPoint()->GetPhysicalVolume();
250  G4double ekin =track.GetKineticEnergy();
251  G4double postCS=0.;
252  if (thePostPhysVolume){
254  ekin,thePostPhysVolume->GetLogicalVolume()->GetMaterialCutsCouple());
255  }
256  if (postCS>0.) return theNumberOfInteractionLengthLeft/postCS;
257  else return DBL_MAX;
258  }
259 }
261 //
263  G4double ,
264  G4double ,
265  G4double& )
266 {return DBL_MAX;
267 }
269 //Not used in this process but should be implemented as virtual method
271  G4double ,
273 { return 0.;
274 }
275 
276 
G4double condition(const G4ErrorSymMatrix &m)
static G4AdjointGamma * AdjointGamma()
G4double GetTotalAdjointCS(G4ParticleDefinition *aPartDef, G4double Ekin, const G4MaterialCutsCouple *aCouple)
G4ParticleDefinition * GetDefinition() const
G4double GetWeight() const
G4double GetStepLength() const
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
void ProposeParentWeight(G4double finalWeight)
const G4ThreeVector & GetPosition() const
void SetCorrectWeightForPostStepInModel(G4bool aBool)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
virtual G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
void SetAdditionalWeightCorrectionFactorForPostStepOutsideModel(G4double factor)
const G4Step * GetStep() const
int G4int
Definition: G4Types.hh:78
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
virtual void SampleSecondaries(const G4Track &aTrack, G4bool IsScatProjToProjCase, G4ParticleChange *fParticleChange)=0
G4StepPoint * GetPreStepPoint() const
void SetSecondaryWeightByProcess(G4bool)
void SetParentWeightByProcess(G4bool)
G4double GetKineticEnergy() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4int GetCurrentStepNumber() const
G4VPhysicalVolume * GetPhysicalVolume() const
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &step)
bool G4bool
Definition: G4Types.hh:79
Definition: G4Step.hh:76
G4int GetTrackID() const
G4double GetTotalForwardCS(G4ParticleDefinition *aPartDef, G4double Ekin, const G4MaterialCutsCouple *aCouple)
G4ThreeVector GetMomentum() const
virtual void Initialize(const G4Track &)
virtual G4double AdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)
G4LogicalVolume * GetLogicalVolume() const
G4StepPoint * GetPostStepPoint() const
virtual G4double GetContinuousStepLimit(const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
void AddSecondary(G4Track *aSecondary)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4double theInitialNumberOfInteractionLength
Definition: G4VProcess.hh:300
double G4double
Definition: G4Types.hh:76
G4ForceCondition
#define DBL_MAX
Definition: templates.hh:83
static G4AdjointCSManager * GetAdjointCSManager()
void BuildPhysicsTable(const G4ParticleDefinition &)
void PreparePhysicsTable(const G4ParticleDefinition &)