Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes
G4VEmAdjointModel Class Referenceabstract

#include <G4VEmAdjointModel.hh>

Inheritance diagram for G4VEmAdjointModel:
G4AdjointBremsstrahlungModel G4AdjointComptonModel G4AdjointeIonisationModel G4AdjointhIonisationModel G4AdjointIonIonisationModel G4AdjointPhotoElectricModel

Public Member Functions

 G4VEmAdjointModel (const G4String &nam)
 
virtual ~G4VEmAdjointModel ()
 
virtual void SampleSecondaries (const G4Track &aTrack, G4bool IsScatProjToProjCase, G4ParticleChange *fParticleChange)=0
 
virtual G4double AdjointCrossSection (const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)
 
virtual G4double GetAdjointCrossSection (const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)
 
virtual G4double DiffCrossSectionPerAtomPrimToSecond (G4double kinEnergyProj, G4double kinEnergyProd, G4double Z, G4double A=0.)
 
virtual G4double DiffCrossSectionPerAtomPrimToScatPrim (G4double kinEnergyProj, G4double kinEnergyScatProj, G4double Z, G4double A=0.)
 
virtual G4double DiffCrossSectionPerVolumePrimToSecond (const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyProd)
 
virtual G4double DiffCrossSectionPerVolumePrimToScatPrim (const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyScatProj)
 
virtual G4double GetSecondAdjEnergyMaxForScatProjToProjCase (G4double PrimAdjEnergy)
 
virtual G4double GetSecondAdjEnergyMinForScatProjToProjCase (G4double PrimAdjEnergy, G4double Tcut=0)
 
virtual G4double GetSecondAdjEnergyMaxForProdToProjCase (G4double PrimAdjEnergy)
 
virtual G4double GetSecondAdjEnergyMinForProdToProjCase (G4double PrimAdjEnergy)
 
void DefineCurrentMaterial (const G4MaterialCutsCouple *couple)
 
std::vector< std::vector
< double > * > 
ComputeAdjointCrossSectionVectorPerAtomForSecond (G4double kinEnergyProd, G4double Z, G4double A=0., G4int nbin_pro_decade=10)
 
std::vector< std::vector
< double > * > 
ComputeAdjointCrossSectionVectorPerAtomForScatProj (G4double kinEnergyProd, G4double Z, G4double A=0., G4int nbin_pro_decade=10)
 
std::vector< std::vector
< double > * > 
ComputeAdjointCrossSectionVectorPerVolumeForSecond (G4Material *aMaterial, G4double kinEnergyProd, G4int nbin_pro_decade=10)
 
std::vector< std::vector
< double > * > 
ComputeAdjointCrossSectionVectorPerVolumeForScatProj (G4Material *aMaterial, G4double kinEnergyProd, G4int nbin_pro_decade=10)
 
void SetCSMatrices (std::vector< G4AdjointCSMatrix * > *Vec1CSMatrix, std::vector< G4AdjointCSMatrix * > *Vec2CSMatrix)
 
G4ParticleDefinitionGetAdjointEquivalentOfDirectPrimaryParticleDefinition ()
 
G4ParticleDefinitionGetAdjointEquivalentOfDirectSecondaryParticleDefinition ()
 
G4double GetHighEnergyLimit ()
 
G4double GetLowEnergyLimit ()
 
void SetHighEnergyLimit (G4double aVal)
 
void SetLowEnergyLimit (G4double aVal)
 
void DefineDirectEMModel (G4VEmModel *aModel)
 
void SetAdjointEquivalentOfDirectPrimaryParticleDefinition (G4ParticleDefinition *aPart)
 
void SetAdjointEquivalentOfDirectSecondaryParticleDefinition (G4ParticleDefinition *aPart)
 
void SetSecondPartOfSameType (G4bool aBool)
 
G4bool GetSecondPartOfSameType ()
 
void SetUseMatrix (G4bool aBool)
 
void SetUseMatrixPerElement (G4bool aBool)
 
void SetUseOnlyOneMatrixForAllElements (G4bool aBool)
 
void SetApplyCutInRange (G4bool aBool)
 
G4bool GetUseMatrix ()
 
G4bool GetUseMatrixPerElement ()
 
G4bool GetUseOnlyOneMatrixForAllElements ()
 
G4bool GetApplyCutInRange ()
 
G4String GetName ()
 
virtual void SetCSBiasingFactor (G4double aVal)
 

Protected Member Functions

G4double DiffCrossSectionFunction1 (G4double kinEnergyProj)
 
G4double DiffCrossSectionFunction2 (G4double kinEnergyProj)
 
G4double DiffCrossSectionPerVolumeFunctionForIntegrationOverEkinProj (G4double EkinProd)
 
G4double SampleAdjSecEnergyFromCSMatrix (size_t MatrixIndex, G4double prim_energy, G4bool IsScatProjToProjCase)
 
G4double SampleAdjSecEnergyFromCSMatrix (G4double prim_energy, G4bool IsScatProjToProjCase)
 
void SelectCSMatrix (G4bool IsScatProjToProjCase)
 
virtual G4double SampleAdjSecEnergyFromDiffCrossSectionPerAtom (G4double prim_energy, G4bool IsScatProjToProjCase)
 
virtual void CorrectPostStepWeight (G4ParticleChange *fParticleChange, G4double old_weight, G4double adjointPrimKinEnergy, G4double projectileKinEnergy, G4bool IsScatProjToProjCase)
 

Protected Attributes

G4VEmModeltheDirectEMModel
 
G4VParticleChangepParticleChange
 
const G4String name
 
G4int ASelectedNucleus
 
G4int ZSelectedNucleus
 
G4MaterialSelectedMaterial
 
G4double kinEnergyProdForIntegration
 
G4double kinEnergyScatProjForIntegration
 
G4double kinEnergyProjForIntegration
 
std::vector< G4AdjointCSMatrix * > * pOnCSMatrixForProdToProjBackwardScattering
 
std::vector< G4AdjointCSMatrix * > * pOnCSMatrixForScatProjToProjBackwardScattering
 
std::vector< G4doubleCS_Vs_ElementForScatProjToProjCase
 
std::vector< G4doubleCS_Vs_ElementForProdToProjCase
 
G4double lastCS
 
G4double lastAdjointCSForScatProjToProjCase
 
G4double lastAdjointCSForProdToProjCase
 
G4ParticleDefinitiontheAdjEquivOfDirectPrimPartDef
 
G4ParticleDefinitiontheAdjEquivOfDirectSecondPartDef
 
G4ParticleDefinitiontheDirectPrimaryPartDef
 
G4bool second_part_of_same_type
 
G4double preStepEnergy
 
G4MaterialcurrentMaterial
 
G4MaterialCutsCouplecurrentCouple
 
size_t currentMaterialIndex
 
size_t currentCoupleIndex
 
G4double currentTcutForDirectPrim
 
G4double currentTcutForDirectSecond
 
G4bool ApplyCutInRange
 
G4double mass_ratio_product
 
G4double mass_ratio_projectile
 
G4double HighEnergyLimit
 
G4double LowEnergyLimit
 
G4double CS_biasing_factor
 
G4bool UseMatrix
 
G4bool UseMatrixPerElement
 
G4bool UseOnlyOneMatrixForAllElements
 
size_t indexOfUsedCrossSectionMatrix
 
size_t model_index
 

Detailed Description

Definition at line 72 of file G4VEmAdjointModel.hh.

Constructor & Destructor Documentation

G4VEmAdjointModel::G4VEmAdjointModel ( const G4String nam)

Definition at line 41 of file G4VEmAdjointModel.cc.

References currentCouple, G4AdjointCSManager::GetAdjointCSManager(), mass_ratio_product, mass_ratio_projectile, model_index, G4AdjointCSManager::RegisterEmAdjointModel(), second_part_of_same_type, and theDirectEMModel.

41  :
42 name(nam)
43 // lowLimit(0.1*keV), highLimit(100.0*TeV), fluc(0), name(nam), pParticleChange(0)
44 {
50  currentCouple=0;
51 }
size_t RegisterEmAdjointModel(G4VEmAdjointModel *)
G4VEmModel * theDirectEMModel
G4MaterialCutsCouple * currentCouple
static G4AdjointCSManager * GetAdjointCSManager()
G4VEmAdjointModel::~G4VEmAdjointModel ( )
virtual

Definition at line 54 of file G4VEmAdjointModel.cc.

55 {;}

Member Function Documentation

G4double G4VEmAdjointModel::AdjointCrossSection ( const G4MaterialCutsCouple aCouple,
G4double  primEnergy,
G4bool  IsScatProjToProjCase 
)
virtual

Reimplemented in G4AdjointhIonisationModel, G4AdjointBremsstrahlungModel, G4AdjointComptonModel, and G4AdjointPhotoElectricModel.

Definition at line 58 of file G4VEmAdjointModel.cc.

References G4AdjointCSManager::ComputeAdjointCS(), CS_Vs_ElementForProdToProjCase, CS_Vs_ElementForScatProjToProjCase, currentMaterial, currentTcutForDirectSecond, DefineCurrentMaterial(), G4AdjointCSManager::GetAdjointCSManager(), lastAdjointCSForProdToProjCase, lastAdjointCSForScatProjToProjCase, lastCS, and preStepEnergy.

Referenced by G4AdjointComptonModel::AdjointCrossSection(), G4AdjointBremsstrahlungModel::AdjointCrossSection(), G4AdjointhIonisationModel::AdjointCrossSection(), G4AdjointCSManager::ComputeAdjointCS(), CorrectPostStepWeight(), and GetAdjointCrossSection().

61 {
62  DefineCurrentMaterial(aCouple);
63  preStepEnergy=primEnergy;
64 
65  std::vector<G4double>* CS_Vs_Element = &CS_Vs_ElementForProdToProjCase;
66  if (IsScatProjToProjCase) CS_Vs_Element = &CS_Vs_ElementForScatProjToProjCase;
68  this,
69  primEnergy,
71  IsScatProjToProjCase,
72  *CS_Vs_Element);
73  if (IsScatProjToProjCase) lastAdjointCSForScatProjToProjCase = lastCS;
75 
76 
77 
78  return lastCS;
79 
80 }
G4double lastAdjointCSForScatProjToProjCase
G4Material * currentMaterial
std::vector< G4double > CS_Vs_ElementForProdToProjCase
std::vector< G4double > CS_Vs_ElementForScatProjToProjCase
void DefineCurrentMaterial(const G4MaterialCutsCouple *couple)
G4double currentTcutForDirectSecond
G4double lastAdjointCSForProdToProjCase
static G4AdjointCSManager * GetAdjointCSManager()
G4double ComputeAdjointCS(G4Material *aMaterial, G4VEmAdjointModel *aModel, G4double PrimEnergy, G4double Tcut, G4bool IsScatProjToProjCase, std::vector< G4double > &AdjointCS_for_each_element)
std::vector< std::vector< G4double > * > G4VEmAdjointModel::ComputeAdjointCrossSectionVectorPerAtomForScatProj ( G4double  kinEnergyProd,
G4double  Z,
G4double  A = 0.,
G4int  nbin_pro_decade = 10 
)

Definition at line 291 of file G4VEmAdjointModel.cc.

References ASelectedNucleus, DiffCrossSectionFunction2(), GetLowEnergyLimit(), GetSecondAdjEnergyMaxForScatProjToProjCase(), GetSecondAdjEnergyMinForScatProjToProjCase(), int(), kinEnergyScatProjForIntegration, G4INCL::Math::max(), G4INCL::Math::min(), G4Integrator< T, F >::Simpson(), and ZSelectedNucleus.

299  kinEnergyScatProjForIntegration = kinEnergyScatProj;
300 
301  //compute the vector of integrated cross sections
302  //-------------------
303 
304  G4double minEProj= GetSecondAdjEnergyMinForScatProjToProjCase(kinEnergyScatProj);
305  G4double maxEProj= GetSecondAdjEnergyMaxForScatProjToProjCase(kinEnergyScatProj);
306  G4double dEmax=maxEProj-kinEnergyScatProj;
307  G4double dEmin=GetLowEnergyLimit();
308  G4double dE1=dEmin;
309  G4double dE2=dEmin;
310 
311 
312  std::vector< double>* log_ESec_vector = new std::vector< double>();
313  std::vector< double>* log_Prob_vector = new std::vector< double>();
314  log_ESec_vector->push_back(std::log(dEmin));
315  log_Prob_vector->push_back(-50.);
316  G4int nbins=std::max( int(std::log10(dEmax/dEmin))*nbin_pro_decade,5);
317  G4double fE=std::pow(dEmax/dEmin,1./nbins);
318 
319 
320 
321 
322 
323  G4double int_cross_section=0.;
324 
325  while (dE1 <dEmax*0.9999999999999){
326  dE2=dE1*fE;
327  int_cross_section +=integral.Simpson(this,
328  &G4VEmAdjointModel::DiffCrossSectionFunction2,minEProj+dE1,std::min(minEProj+dE2,maxEProj), 5);
329  //G4cout<<"int_cross_section "<<minEProj+dE1<<'\t'<<int_cross_section<<G4endl;
330  log_ESec_vector->push_back(std::log(std::min(dE2,maxEProj-minEProj)));
331  log_Prob_vector->push_back(std::log(int_cross_section));
332  dE1=dE2;
333 
334  }
335 
336 
337  std::vector< std::vector<G4double> *> res_mat;
338  res_mat.clear();
339  if (int_cross_section >0.) {
340  res_mat.push_back(log_ESec_vector);
341  res_mat.push_back(log_Prob_vector);
342  }
343 
344  return res_mat;
345 }
G4double kinEnergyScatProjForIntegration
typedef int(XMLCALL *XML_NotStandaloneHandler)(void *userData)
int G4int
Definition: G4Types.hh:78
virtual G4double GetSecondAdjEnergyMaxForScatProjToProjCase(G4double PrimAdjEnergy)
G4double GetLowEnergyLimit()
G4double DiffCrossSectionFunction2(G4double kinEnergyProj)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
virtual G4double GetSecondAdjEnergyMinForScatProjToProjCase(G4double PrimAdjEnergy, G4double Tcut=0)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double Simpson(T &typeT, F f, G4double a, G4double b, G4int n)
double G4double
Definition: G4Types.hh:76
std::vector< std::vector< G4double > * > G4VEmAdjointModel::ComputeAdjointCrossSectionVectorPerAtomForSecond ( G4double  kinEnergyProd,
G4double  Z,
G4double  A = 0.,
G4int  nbin_pro_decade = 10 
)

Definition at line 238 of file G4VEmAdjointModel.cc.

References ASelectedNucleus, DiffCrossSectionFunction1(), GetSecondAdjEnergyMaxForProdToProjCase(), GetSecondAdjEnergyMinForProdToProjCase(), int(), kinEnergyProdForIntegration, G4INCL::Math::min(), G4Integrator< T, F >::Simpson(), and ZSelectedNucleus.

243 {
245  ASelectedNucleus= int(A);
247  kinEnergyProdForIntegration = kinEnergyProd;
248 
249  //compute the vector of integrated cross sections
250  //-------------------
251 
252  G4double minEProj= GetSecondAdjEnergyMinForProdToProjCase(kinEnergyProd);
253  G4double maxEProj= GetSecondAdjEnergyMaxForProdToProjCase(kinEnergyProd);
254  G4double E1=minEProj;
255  std::vector< double>* log_ESec_vector = new std::vector< double>();
256  std::vector< double>* log_Prob_vector = new std::vector< double>();
257  log_ESec_vector->clear();
258  log_Prob_vector->clear();
259  log_ESec_vector->push_back(std::log(E1));
260  log_Prob_vector->push_back(-50.);
261 
262  G4double E2=std::pow(10.,double( int(std::log10(minEProj)*nbin_pro_decade)+1)/nbin_pro_decade);
263  G4double fE=std::pow(10.,1./nbin_pro_decade);
264  G4double int_cross_section=0.;
265 
266  if (std::pow(fE,5.)>(maxEProj/minEProj)) fE = std::pow(maxEProj/minEProj,0.2);
267 
268  while (E1 <maxEProj*0.9999999){
269  //G4cout<<E1<<'\t'<<E2<<G4endl;
270 
271  int_cross_section +=integral.Simpson(this,
272  &G4VEmAdjointModel::DiffCrossSectionFunction1,E1,std::min(E2,maxEProj*0.99999999), 5);
273  log_ESec_vector->push_back(std::log(std::min(E2,maxEProj)));
274  log_Prob_vector->push_back(std::log(int_cross_section));
275  E1=E2;
276  E2*=fE;
277 
278  }
279  std::vector< std::vector<G4double>* > res_mat;
280  res_mat.clear();
281  if (int_cross_section >0.) {
282  res_mat.push_back(log_ESec_vector);
283  res_mat.push_back(log_Prob_vector);
284  }
285 
286  return res_mat;
287 }
virtual G4double GetSecondAdjEnergyMaxForProdToProjCase(G4double PrimAdjEnergy)
typedef int(XMLCALL *XML_NotStandaloneHandler)(void *userData)
G4double DiffCrossSectionFunction1(G4double kinEnergyProj)
virtual G4double GetSecondAdjEnergyMinForProdToProjCase(G4double PrimAdjEnergy)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double Simpson(T &typeT, F f, G4double a, G4double b, G4int n)
double G4double
Definition: G4Types.hh:76
G4double kinEnergyProdForIntegration
std::vector< std::vector< G4double > * > G4VEmAdjointModel::ComputeAdjointCrossSectionVectorPerVolumeForScatProj ( G4Material aMaterial,
G4double  kinEnergyProd,
G4int  nbin_pro_decade = 10 
)

Definition at line 399 of file G4VEmAdjointModel.cc.

References DiffCrossSectionFunction2(), GetLowEnergyLimit(), GetSecondAdjEnergyMaxForScatProjToProjCase(), GetSecondAdjEnergyMinForScatProjToProjCase(), kinEnergyScatProjForIntegration, G4INCL::Math::max(), G4INCL::Math::min(), SelectedMaterial, and G4Integrator< T, F >::Simpson().

404  SelectedMaterial= aMaterial;
405  kinEnergyScatProjForIntegration = kinEnergyScatProj;
406 
407  //compute the vector of integrated cross sections
408  //-------------------
409 
410  G4double minEProj= GetSecondAdjEnergyMinForScatProjToProjCase(kinEnergyScatProj);
411  G4double maxEProj= GetSecondAdjEnergyMaxForScatProjToProjCase(kinEnergyScatProj);
412 
413 
414  G4double dEmax=maxEProj-kinEnergyScatProj;
415  G4double dEmin=GetLowEnergyLimit();
416  G4double dE1=dEmin;
417  G4double dE2=dEmin;
418 
419 
420  std::vector< double>* log_ESec_vector = new std::vector< double>();
421  std::vector< double>* log_Prob_vector = new std::vector< double>();
422  log_ESec_vector->push_back(std::log(dEmin));
423  log_Prob_vector->push_back(-50.);
424  G4int nbins=std::max( int(std::log10(dEmax/dEmin))*nbin_pro_decade,5);
425  G4double fE=std::pow(dEmax/dEmin,1./nbins);
426 
427  G4double int_cross_section=0.;
428 
429  while (dE1 <dEmax*0.9999999999999){
430  dE2=dE1*fE;
431  int_cross_section +=integral.Simpson(this,
432  &G4VEmAdjointModel::DiffCrossSectionFunction2,minEProj+dE1,std::min(minEProj+dE2,maxEProj), 5);
433  log_ESec_vector->push_back(std::log(std::min(dE2,maxEProj-minEProj)));
434  log_Prob_vector->push_back(std::log(int_cross_section));
435  dE1=dE2;
436 
437  }
438 
439 
440 
441 
442 
443  std::vector< std::vector<G4double> *> res_mat;
444  res_mat.clear();
445  if (int_cross_section >0.) {
446  res_mat.push_back(log_ESec_vector);
447  res_mat.push_back(log_Prob_vector);
448  }
449 
450  return res_mat;
451 }
G4double kinEnergyScatProjForIntegration
G4Material * SelectedMaterial
int G4int
Definition: G4Types.hh:78
virtual G4double GetSecondAdjEnergyMaxForScatProjToProjCase(G4double PrimAdjEnergy)
G4double GetLowEnergyLimit()
G4double DiffCrossSectionFunction2(G4double kinEnergyProj)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
virtual G4double GetSecondAdjEnergyMinForScatProjToProjCase(G4double PrimAdjEnergy, G4double Tcut=0)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double Simpson(T &typeT, F f, G4double a, G4double b, G4int n)
double G4double
Definition: G4Types.hh:76
std::vector< std::vector< G4double > * > G4VEmAdjointModel::ComputeAdjointCrossSectionVectorPerVolumeForSecond ( G4Material aMaterial,
G4double  kinEnergyProd,
G4int  nbin_pro_decade = 10 
)

Definition at line 348 of file G4VEmAdjointModel.cc.

References DiffCrossSectionFunction1(), GetSecondAdjEnergyMaxForProdToProjCase(), GetSecondAdjEnergyMinForProdToProjCase(), kinEnergyProdForIntegration, G4INCL::Math::min(), SelectedMaterial, and G4Integrator< T, F >::Simpson().

353  SelectedMaterial= aMaterial;
354  kinEnergyProdForIntegration = kinEnergyProd;
355  //compute the vector of integrated cross sections
356  //-------------------
357 
358  G4double minEProj= GetSecondAdjEnergyMinForProdToProjCase(kinEnergyProd);
359  G4double maxEProj= GetSecondAdjEnergyMaxForProdToProjCase(kinEnergyProd);
360  G4double E1=minEProj;
361  std::vector< double>* log_ESec_vector = new std::vector< double>();
362  std::vector< double>* log_Prob_vector = new std::vector< double>();
363  log_ESec_vector->clear();
364  log_Prob_vector->clear();
365  log_ESec_vector->push_back(std::log(E1));
366  log_Prob_vector->push_back(-50.);
367 
368  G4double E2=std::pow(10.,double( int(std::log10(minEProj)*nbin_pro_decade)+1)/nbin_pro_decade);
369  G4double fE=std::pow(10.,1./nbin_pro_decade);
370  G4double int_cross_section=0.;
371 
372  if (std::pow(fE,5.)>(maxEProj/minEProj)) fE = std::pow(maxEProj/minEProj,0.2);
373 
374  while (E1 <maxEProj*0.9999999){
375 
376  int_cross_section +=integral.Simpson(this,
377  &G4VEmAdjointModel::DiffCrossSectionFunction1,E1,std::min(E2,maxEProj*0.99999999), 5);
378  log_ESec_vector->push_back(std::log(std::min(E2,maxEProj)));
379  log_Prob_vector->push_back(std::log(int_cross_section));
380  E1=E2;
381  E2*=fE;
382 
383  }
384  std::vector< std::vector<G4double>* > res_mat;
385  res_mat.clear();
386 
387  if (int_cross_section >0.) {
388  res_mat.push_back(log_ESec_vector);
389  res_mat.push_back(log_Prob_vector);
390  }
391 
392 
393 
394  return res_mat;
395 }
virtual G4double GetSecondAdjEnergyMaxForProdToProjCase(G4double PrimAdjEnergy)
G4Material * SelectedMaterial
G4double DiffCrossSectionFunction1(G4double kinEnergyProj)
virtual G4double GetSecondAdjEnergyMinForProdToProjCase(G4double PrimAdjEnergy)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double Simpson(T &typeT, F f, G4double a, G4double b, G4int n)
double G4double
Definition: G4Types.hh:76
G4double kinEnergyProdForIntegration
void G4VEmAdjointModel::CorrectPostStepWeight ( G4ParticleChange fParticleChange,
G4double  old_weight,
G4double  adjointPrimKinEnergy,
G4double  projectileKinEnergy,
G4bool  IsScatProjToProjCase 
)
protectedvirtual

Reimplemented in G4AdjointIonIonisationModel, and G4AdjointPhotoElectricModel.

Definition at line 624 of file G4VEmAdjointModel.cc.

References AdjointCrossSection(), CS_biasing_factor, currentCouple, G4AdjointCSManager::GetAdjointCSManager(), G4AdjointCSManager::GetPostStepWeightCorrection(), lastAdjointCSForProdToProjCase, lastAdjointCSForScatProjToProjCase, lastCS, preStepEnergy, G4VParticleChange::ProposeParentWeight(), G4VParticleChange::SetParentWeightByProcess(), and G4VParticleChange::SetSecondaryWeightByProcess().

Referenced by G4AdjointComptonModel::SampleSecondaries(), G4AdjointeIonisationModel::SampleSecondaries(), G4AdjointBremsstrahlungModel::SampleSecondaries(), and G4AdjointhIonisationModel::SampleSecondaries().

629 {
630  G4double new_weight=old_weight;
631  G4double w_corr =1./CS_biasing_factor;
633 
634 
636  if ( !IsScatProjToProjCase) lastCS=lastAdjointCSForProdToProjCase;
637  if ((adjointPrimKinEnergy-preStepEnergy)/preStepEnergy>0.001){ //Is that in all cases needed???
638  G4double post_stepCS=AdjointCrossSection(currentCouple, adjointPrimKinEnergy
639  ,IsScatProjToProjCase );
640  if (post_stepCS>0 && lastCS>0) w_corr*=post_stepCS/lastCS;
641  }
642 
643  new_weight*=w_corr;
644 
645  //G4cout<<"Post step "<<new_weight<<'\t'<<w_corr<<'\t'<<old_weight<<G4endl;
646  new_weight*=projectileKinEnergy/adjointPrimKinEnergy;//This is needed due to the biasing of diff CS
647  //by the factor adjointPrimKinEnergy/projectileKinEnergy
648 
649 
650 
651  fParticleChange->SetParentWeightByProcess(false);
652  fParticleChange->SetSecondaryWeightByProcess(false);
653  fParticleChange->ProposeParentWeight(new_weight);
654 }
G4double lastAdjointCSForScatProjToProjCase
G4double GetPostStepWeightCorrection()
void ProposeParentWeight(G4double finalWeight)
void SetSecondaryWeightByProcess(G4bool)
void SetParentWeightByProcess(G4bool)
G4MaterialCutsCouple * currentCouple
virtual G4double AdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)
G4double lastAdjointCSForProdToProjCase
double G4double
Definition: G4Types.hh:76
static G4AdjointCSManager * GetAdjointCSManager()
void G4VEmAdjointModel::DefineCurrentMaterial ( const G4MaterialCutsCouple couple)

Definition at line 683 of file G4VEmAdjointModel.cc.

References G4AdjointElectron::AdjointElectron(), G4AdjointGamma::AdjointGamma(), G4AdjointPositron::AdjointPositron(), currentCouple, currentCoupleIndex, currentMaterial, currentMaterialIndex, currentTcutForDirectSecond, G4ProductionCutsTable::GetEnergyCutsVector(), G4MaterialCutsCouple::GetIndex(), G4Material::GetIndex(), G4MaterialCutsCouple::GetMaterial(), G4ProductionCutsTable::GetProductionCutsTable(), and theAdjEquivOfDirectSecondPartDef.

Referenced by G4AdjointComptonModel::AdjointCrossSection(), G4AdjointBremsstrahlungModel::AdjointCrossSection(), G4AdjointhIonisationModel::AdjointCrossSection(), AdjointCrossSection(), G4AdjointComptonModel::RapidSampleSecondaries(), G4AdjointBremsstrahlungModel::RapidSampleSecondaries(), G4AdjointhIonisationModel::RapidSampleSecondaries(), and G4AdjointBremsstrahlungModel::SampleSecondaries().

684 { if(couple != currentCouple) {
685  currentCouple = const_cast<G4MaterialCutsCouple*> (couple);
686  currentMaterial = const_cast<G4Material*> (couple->GetMaterial());
687  currentCoupleIndex = couple->GetIndex();
689  size_t idx=56;
690  currentTcutForDirectSecond =0.00000000001;
695  if (idx <56){
696  const std::vector<G4double>* aVec = G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(idx);
698  }
699  }
700 
701 
702  }
703 }
static G4AdjointGamma * AdjointGamma()
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
size_t GetIndex() const
Definition: G4Material.hh:260
static G4AdjointElectron * AdjointElectron()
G4Material * currentMaterial
G4MaterialCutsCouple * currentCouple
static G4ProductionCutsTable * GetProductionCutsTable()
G4ParticleDefinition * theAdjEquivOfDirectSecondPartDef
G4double currentTcutForDirectSecond
static G4AdjointPositron * AdjointPositron()
const G4Material * GetMaterial() const
void G4VEmAdjointModel::DefineDirectEMModel ( G4VEmModel aModel)
inline

Definition at line 193 of file G4VEmAdjointModel.hh.

References theDirectEMModel.

Referenced by G4AdjointPhotoElectricModel::SetTheDirectPEEffectModel().

193 {theDirectEMModel = aModel;}
G4VEmModel * theDirectEMModel
G4double G4VEmAdjointModel::DiffCrossSectionFunction1 ( G4double  kinEnergyProj)
protected

Definition at line 201 of file G4VEmAdjointModel.cc.

References ASelectedNucleus, CS_biasing_factor, DiffCrossSectionPerAtomPrimToSecond(), DiffCrossSectionPerVolumePrimToSecond(), kinEnergyProdForIntegration, SelectedMaterial, UseMatrixPerElement, and ZSelectedNucleus.

Referenced by ComputeAdjointCrossSectionVectorPerAtomForSecond(), and ComputeAdjointCrossSectionVectorPerVolumeForSecond().

201  {
202 
203 
204  G4double bias_factor = CS_biasing_factor*kinEnergyProdForIntegration/kinEnergyProj;
205 
206 
207  if (UseMatrixPerElement ) {
209  }
210  else {
212  }
213 }
virtual G4double DiffCrossSectionPerVolumePrimToSecond(const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyProd)
virtual G4double DiffCrossSectionPerAtomPrimToSecond(G4double kinEnergyProj, G4double kinEnergyProd, G4double Z, G4double A=0.)
G4Material * SelectedMaterial
double G4double
Definition: G4Types.hh:76
G4double kinEnergyProdForIntegration
G4double G4VEmAdjointModel::DiffCrossSectionFunction2 ( G4double  kinEnergyProj)
protected

Definition at line 217 of file G4VEmAdjointModel.cc.

References ASelectedNucleus, CS_biasing_factor, DiffCrossSectionPerAtomPrimToScatPrim(), DiffCrossSectionPerVolumePrimToScatPrim(), kinEnergyScatProjForIntegration, SelectedMaterial, UseMatrixPerElement, and ZSelectedNucleus.

Referenced by ComputeAdjointCrossSectionVectorPerAtomForScatProj(), and ComputeAdjointCrossSectionVectorPerVolumeForScatProj().

217  {
218 
219  G4double bias_factor = CS_biasing_factor*kinEnergyScatProjForIntegration/kinEnergyProj;
220  if (UseMatrixPerElement ) {
222  }
223  else {
225 
226  }
227 
228 }
G4double kinEnergyScatProjForIntegration
virtual G4double DiffCrossSectionPerVolumePrimToScatPrim(const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyScatProj)
G4Material * SelectedMaterial
virtual G4double DiffCrossSectionPerAtomPrimToScatPrim(G4double kinEnergyProj, G4double kinEnergyScatProj, G4double Z, G4double A=0.)
double G4double
Definition: G4Types.hh:76
G4double G4VEmAdjointModel::DiffCrossSectionPerAtomPrimToScatPrim ( G4double  kinEnergyProj,
G4double  kinEnergyScatProj,
G4double  Z,
G4double  A = 0. 
)
virtual

Reimplemented in G4AdjointComptonModel.

Definition at line 144 of file G4VEmAdjointModel.cc.

References DiffCrossSectionPerAtomPrimToSecond().

Referenced by DiffCrossSectionFunction2(), and SampleAdjSecEnergyFromDiffCrossSectionPerAtom().

149 { G4double kinEnergyProd = kinEnergyProj - kinEnergyScatProj;
150  G4double dSigmadEprod;
151  if (kinEnergyProd <=0) dSigmadEprod=0;
152  else dSigmadEprod=DiffCrossSectionPerAtomPrimToSecond(kinEnergyProj,kinEnergyProd,Z,A);
153  return dSigmadEprod;
154 
155 }
virtual G4double DiffCrossSectionPerAtomPrimToSecond(G4double kinEnergyProj, G4double kinEnergyProd, G4double Z, G4double A=0.)
double G4double
Definition: G4Types.hh:76
G4double G4VEmAdjointModel::DiffCrossSectionPerAtomPrimToSecond ( G4double  kinEnergyProj,
G4double  kinEnergyProd,
G4double  Z,
G4double  A = 0. 
)
virtual

Reimplemented in G4AdjointhIonisationModel, G4AdjointIonIonisationModel, G4AdjointComptonModel, and G4AdjointeIonisationModel.

Definition at line 112 of file G4VEmAdjointModel.cc.

References G4VEmModel::ComputeCrossSectionPerAtom(), GetSecondAdjEnergyMaxForProdToProjCase(), GetSecondAdjEnergyMinForProdToProjCase(), theDirectEMModel, and theDirectPrimaryPartDef.

Referenced by DiffCrossSectionFunction1(), DiffCrossSectionPerAtomPrimToScatPrim(), and SampleAdjSecEnergyFromDiffCrossSectionPerAtom().

117 {
118  G4double dSigmadEprod=0;
119  G4double Emax_proj = GetSecondAdjEnergyMaxForProdToProjCase(kinEnergyProd);
120  G4double Emin_proj = GetSecondAdjEnergyMinForProdToProjCase(kinEnergyProd);
121 
122 
123  if (kinEnergyProj>Emin_proj && kinEnergyProj<=Emax_proj){ //the produced particle should have a kinetic energy smaller than the projectile
124 
125  /*G4double Tmax=kinEnergyProj;
126  if (second_part_of_same_type) Tmax = kinEnergyProj/2.;*/
127 
128  G4double E1=kinEnergyProd;
129  G4double E2=kinEnergyProd*1.000001;
130  G4double dE=(E2-E1);
133 
134  dSigmadEprod=(sigma1-sigma2)/dE;
135  }
136  return dSigmadEprod;
137 
138 
139 
140 }
virtual G4double GetSecondAdjEnergyMaxForProdToProjCase(G4double PrimAdjEnergy)
G4ParticleDefinition * theDirectPrimaryPartDef
G4VEmModel * theDirectEMModel
virtual G4double GetSecondAdjEnergyMinForProdToProjCase(G4double PrimAdjEnergy)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:300
double G4double
Definition: G4Types.hh:76
G4double G4VEmAdjointModel::DiffCrossSectionPerVolumeFunctionForIntegrationOverEkinProj ( G4double  EkinProd)
protected

Definition at line 232 of file G4VEmAdjointModel.cc.

References DiffCrossSectionPerVolumePrimToSecond(), kinEnergyProjForIntegration, and SelectedMaterial.

233 {
235 }
virtual G4double DiffCrossSectionPerVolumePrimToSecond(const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyProd)
G4Material * SelectedMaterial
G4double kinEnergyProjForIntegration
G4double G4VEmAdjointModel::DiffCrossSectionPerVolumePrimToScatPrim ( const G4Material aMaterial,
G4double  kinEnergyProj,
G4double  kinEnergyScatProj 
)
virtual

Definition at line 188 of file G4VEmAdjointModel.cc.

References DiffCrossSectionPerVolumePrimToSecond().

Referenced by DiffCrossSectionFunction2(), and G4AdjointeIonisationModel::SampleSecondaries().

192 { G4double kinEnergyProd = kinEnergyProj - kinEnergyScatProj;
193  G4double dSigmadEprod;
194  if (kinEnergyProd <=0) dSigmadEprod=0;
195  else dSigmadEprod=DiffCrossSectionPerVolumePrimToSecond(aMaterial,kinEnergyProj,kinEnergyProd);
196  return dSigmadEprod;
197 
198 }
virtual G4double DiffCrossSectionPerVolumePrimToSecond(const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyProd)
double G4double
Definition: G4Types.hh:76
G4double G4VEmAdjointModel::DiffCrossSectionPerVolumePrimToSecond ( const G4Material aMaterial,
G4double  kinEnergyProj,
G4double  kinEnergyProd 
)
virtual

Reimplemented in G4AdjointBremsstrahlungModel.

Definition at line 160 of file G4VEmAdjointModel.cc.

References G4VEmModel::CrossSectionPerVolume(), GetSecondAdjEnergyMaxForProdToProjCase(), GetSecondAdjEnergyMinForProdToProjCase(), theDirectEMModel, and theDirectPrimaryPartDef.

Referenced by DiffCrossSectionFunction1(), DiffCrossSectionPerVolumeFunctionForIntegrationOverEkinProj(), DiffCrossSectionPerVolumePrimToScatPrim(), and G4AdjointeIonisationModel::SampleSecondaries().

164 {
165  G4double dSigmadEprod=0;
166  G4double Emax_proj = GetSecondAdjEnergyMaxForProdToProjCase(kinEnergyProd);
167  G4double Emin_proj = GetSecondAdjEnergyMinForProdToProjCase(kinEnergyProd);
168 
169 
170  if (kinEnergyProj>Emin_proj && kinEnergyProj<=Emax_proj){
171  /*G4double Tmax=kinEnergyProj;
172  if (second_part_of_same_type) Tmax = kinEnergyProj/2.;*/
173  G4double E1=kinEnergyProd;
174  G4double E2=kinEnergyProd*1.0001;
175  G4double dE=(E2-E1);
176  G4double sigma1=theDirectEMModel->CrossSectionPerVolume(aMaterial,theDirectPrimaryPartDef,kinEnergyProj,E1,1.e20);
177  G4double sigma2=theDirectEMModel->CrossSectionPerVolume(aMaterial,theDirectPrimaryPartDef,kinEnergyProj,E2,1.e20);
178  dSigmadEprod=(sigma1-sigma2)/dE;
179  }
180  return dSigmadEprod;
181 
182 
183 
184 }
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:245
virtual G4double GetSecondAdjEnergyMaxForProdToProjCase(G4double PrimAdjEnergy)
G4ParticleDefinition * theDirectPrimaryPartDef
G4VEmModel * theDirectEMModel
virtual G4double GetSecondAdjEnergyMinForProdToProjCase(G4double PrimAdjEnergy)
double G4double
Definition: G4Types.hh:76
G4double G4VEmAdjointModel::GetAdjointCrossSection ( const G4MaterialCutsCouple aCouple,
G4double  primEnergy,
G4bool  IsScatProjToProjCase 
)
virtual

Reimplemented in G4AdjointBremsstrahlungModel, G4AdjointComptonModel, and G4AdjointPhotoElectricModel.

Definition at line 83 of file G4VEmAdjointModel.cc.

References AdjointCrossSection().

Referenced by G4AdjointBremsstrahlungModel::GetAdjointCrossSection(), and G4VAdjointReverseReaction::GetMeanFreePath().

86 {
87  return AdjointCrossSection(aCouple, primEnergy,
88  IsScatProjToProjCase);
89 
90  /*
91  //To continue
92  DefineCurrentMaterial(aCouple);
93  preStepEnergy=primEnergy;
94  if (IsScatProjToProjCase){
95  G4double ekin=primEnergy*mass_ratio_projectile;
96  lastCS = G4AdjointCSManager::GetAdjointCSManager()->GetAdjointSigma(ekin, model_index,true, aCouple);
97  lastAdjointCSForScatProjToProjCase = lastCS;
98  //G4cout<<ekin<<std::endl;
99  }
100  else {
101  G4double ekin=primEnergy*mass_ratio_product;
102  lastCS = G4AdjointCSManager::GetAdjointCSManager()->GetAdjointSigma(ekin, model_index,false, aCouple);
103  lastAdjointCSForProdToProjCase = lastCS;
104  //G4cout<<ekin<<std::endl;
105  }
106  return lastCS;
107  */
108 }
virtual G4double AdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)
G4ParticleDefinition* G4VEmAdjointModel::GetAdjointEquivalentOfDirectPrimaryParticleDefinition ( )
inline

Definition at line 181 of file G4VEmAdjointModel.hh.

References theAdjEquivOfDirectPrimPartDef.

G4ParticleDefinition * theAdjEquivOfDirectPrimPartDef
G4ParticleDefinition* G4VEmAdjointModel::GetAdjointEquivalentOfDirectSecondaryParticleDefinition ( )
inline

Definition at line 183 of file G4VEmAdjointModel.hh.

References theAdjEquivOfDirectSecondPartDef.

G4ParticleDefinition * theAdjEquivOfDirectSecondPartDef
G4bool G4VEmAdjointModel::GetApplyCutInRange ( )
inline

Definition at line 214 of file G4VEmAdjointModel.hh.

References ApplyCutInRange.

Referenced by G4AdjointCSManager::ComputeAdjointCS().

214 { return ApplyCutInRange;}
G4double G4VEmAdjointModel::GetHighEnergyLimit ( )
inline

Definition at line 185 of file G4VEmAdjointModel.hh.

References HighEnergyLimit.

185 {return HighEnergyLimit;}
G4double G4VEmAdjointModel::GetLowEnergyLimit ( )
inline
G4String G4VEmAdjointModel::GetName ( void  )
inline

Definition at line 216 of file G4VEmAdjointModel.hh.

References name.

Referenced by G4AdjointCSManager::BuildCrossSectionMatrices().

216 { return name;}
G4double G4VEmAdjointModel::GetSecondAdjEnergyMaxForProdToProjCase ( G4double  PrimAdjEnergy)
virtual
G4double G4VEmAdjointModel::GetSecondAdjEnergyMaxForScatProjToProjCase ( G4double  PrimAdjEnergy)
virtual
G4double G4VEmAdjointModel::GetSecondAdjEnergyMinForProdToProjCase ( G4double  PrimAdjEnergy)
virtual
G4double G4VEmAdjointModel::GetSecondAdjEnergyMinForScatProjToProjCase ( G4double  PrimAdjEnergy,
G4double  Tcut = 0 
)
virtual
G4bool G4VEmAdjointModel::GetSecondPartOfSameType ( )
inline

Definition at line 203 of file G4VEmAdjointModel.hh.

References second_part_of_same_type.

203 {return second_part_of_same_type;}
G4bool G4VEmAdjointModel::GetUseMatrix ( )
inline
G4bool G4VEmAdjointModel::GetUseMatrixPerElement ( )
inline
G4bool G4VEmAdjointModel::GetUseOnlyOneMatrixForAllElements ( )
inline
G4double G4VEmAdjointModel::SampleAdjSecEnergyFromCSMatrix ( size_t  MatrixIndex,
G4double  prim_energy,
G4bool  IsScatProjToProjCase 
)
protected

Definition at line 454 of file G4VEmAdjointModel.cc.

References ApplyCutInRange, currentTcutForDirectSecond, G4AdjointInterpolator::FindPositionForLogVector(), G4cout, G4endl, G4UniformRand, G4AdjointCSMatrix::GetData(), G4AdjointInterpolator::GetInstance(), G4AdjointCSMatrix::GetLogPrimEnergyVector(), GetSecondAdjEnergyMaxForProdToProjCase(), GetSecondAdjEnergyMaxForScatProjToProjCase(), GetSecondAdjEnergyMinForProdToProjCase(), GetSecondAdjEnergyMinForScatProjToProjCase(), G4AdjointInterpolator::Interpolate(), G4AdjointInterpolator::InterpolateForLogVector(), G4AdjointCSMatrix::IsScatProjToProjCase(), G4AdjointInterpolator::LinearInterpolation(), G4INCL::Math::max(), G4INCL::Math::min(), and second_part_of_same_type.

Referenced by SampleAdjSecEnergyFromCSMatrix(), G4AdjointComptonModel::SampleSecondaries(), G4AdjointeIonisationModel::SampleSecondaries(), G4AdjointBremsstrahlungModel::SampleSecondaries(), G4AdjointIonIonisationModel::SampleSecondaries(), and G4AdjointhIonisationModel::SampleSecondaries().

455 {
456 
457 
458  G4AdjointCSMatrix* theMatrix= (*pOnCSMatrixForProdToProjBackwardScattering)[MatrixIndex];
459  if (IsScatProjToProjCase) theMatrix= (*pOnCSMatrixForScatProjToProjBackwardScattering)[MatrixIndex];
460  std::vector< double>* theLogPrimEnergyVector = theMatrix->GetLogPrimEnergyVector();
461 
462  if (theLogPrimEnergyVector->size() ==0){
463  G4cout<<"No data are contained in the given AdjointCSMatrix!"<<G4endl;
464  G4cout<<"The sampling procedure will be stopped."<<G4endl;
465  return 0.;
466 
467  }
468 
470  G4double aLogPrimEnergy = std::log(aPrimEnergy);
471  size_t ind =theInterpolator->FindPositionForLogVector(aLogPrimEnergy,*theLogPrimEnergyVector);
472 
473 
474  G4double aLogPrimEnergy1,aLogPrimEnergy2;
475  G4double aLogCS1,aLogCS2;
476  G4double log01,log02;
477  std::vector< double>* aLogSecondEnergyVector1 =0;
478  std::vector< double>* aLogSecondEnergyVector2 =0;
479  std::vector< double>* aLogProbVector1=0;
480  std::vector< double>* aLogProbVector2=0;
481  std::vector< size_t>* aLogProbVectorIndex1=0;
482  std::vector< size_t>* aLogProbVectorIndex2=0;
483 
484  theMatrix->GetData(ind, aLogPrimEnergy1,aLogCS1,log01, aLogSecondEnergyVector1,aLogProbVector1,aLogProbVectorIndex1);
485  theMatrix->GetData(ind+1, aLogPrimEnergy2,aLogCS2,log02, aLogSecondEnergyVector2,aLogProbVector2,aLogProbVectorIndex2);
486 
487  G4double rand_var = G4UniformRand();
488  G4double log_rand_var= std::log(rand_var);
489  G4double log_Tcut =std::log(currentTcutForDirectSecond);
490  G4double Esec=0;
491  G4double log_dE1,log_dE2;
492  G4double log_rand_var1,log_rand_var2;
493  G4double log_E1,log_E2;
494  log_rand_var1=log_rand_var;
495  log_rand_var2=log_rand_var;
496 
497  G4double Emin=0.;
498  G4double Emax=0.;
499  if (theMatrix->IsScatProjToProjCase()){ //case where Tcut plays a role
502  G4double dE=0;
503  if (Emin < Emax ){
504  if (ApplyCutInRange) {
505  if (second_part_of_same_type && currentTcutForDirectSecond>aPrimEnergy) return aPrimEnergy;
506 
507  log_rand_var1=log_rand_var+theInterpolator->InterpolateForLogVector(log_Tcut,*aLogSecondEnergyVector1,*aLogProbVector1);
508  log_rand_var2=log_rand_var+theInterpolator->InterpolateForLogVector(log_Tcut,*aLogSecondEnergyVector2,*aLogProbVector2);
509 
510  }
511  log_dE1 = theInterpolator->Interpolate(log_rand_var1,*aLogProbVector1,*aLogSecondEnergyVector1,"Lin");
512  log_dE2 = theInterpolator->Interpolate(log_rand_var2,*aLogProbVector2,*aLogSecondEnergyVector2,"Lin");
513  dE=std::exp(theInterpolator->LinearInterpolation(aLogPrimEnergy,aLogPrimEnergy1,aLogPrimEnergy2,log_dE1,log_dE2));
514  }
515 
516  Esec = aPrimEnergy +dE;
517  Esec=std::max(Esec,Emin);
518  Esec=std::min(Esec,Emax);
519 
520  }
521  else { //Tcut condition is already full-filled
522 
523  log_E1 = theInterpolator->Interpolate(log_rand_var,*aLogProbVector1,*aLogSecondEnergyVector1,"Lin");
524  log_E2 = theInterpolator->Interpolate(log_rand_var,*aLogProbVector2,*aLogSecondEnergyVector2,"Lin");
525 
526  Esec = std::exp(theInterpolator->LinearInterpolation(aLogPrimEnergy,aLogPrimEnergy1,aLogPrimEnergy2,log_E1,log_E2));
527  Emin=GetSecondAdjEnergyMinForProdToProjCase(aPrimEnergy);
528  Emax=GetSecondAdjEnergyMaxForProdToProjCase(aPrimEnergy);
529  Esec=std::max(Esec,Emin);
530  Esec=std::min(Esec,Emax);
531 
532  }
533 
534  return Esec;
535 
536 
537 
538 
539 
540 }
virtual G4double GetSecondAdjEnergyMaxForProdToProjCase(G4double PrimAdjEnergy)
G4bool GetData(unsigned int i, G4double &aPrimEnergy, G4double &aCS, G4double &log0, std::vector< double > *&aLogSecondEnergyVector, std::vector< double > *&aLogProbVector, std::vector< size_t > *&aLogProbVectorIndex)
virtual G4double GetSecondAdjEnergyMaxForScatProjToProjCase(G4double PrimAdjEnergy)
G4double InterpolateForLogVector(G4double &x, std::vector< G4double > &x_vec, std::vector< G4double > &y_vec)
G4double Interpolate(G4double &x, std::vector< G4double > &x_vec, std::vector< G4double > &y_vec, G4String InterPolMethod="Log")
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
virtual G4double GetSecondAdjEnergyMinForProdToProjCase(G4double PrimAdjEnergy)
std::vector< double > * GetLogPrimEnergyVector()
size_t FindPositionForLogVector(G4double &x, std::vector< G4double > &x_vec)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
virtual G4double GetSecondAdjEnergyMinForScatProjToProjCase(G4double PrimAdjEnergy, G4double Tcut=0)
G4double LinearInterpolation(G4double &x, G4double &x1, G4double &x2, G4double &y1, G4double &y2)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static G4AdjointInterpolator * GetInstance()
#define G4endl
Definition: G4ios.hh:61
G4double currentTcutForDirectSecond
double G4double
Definition: G4Types.hh:76
G4double G4VEmAdjointModel::SampleAdjSecEnergyFromCSMatrix ( G4double  prim_energy,
G4bool  IsScatProjToProjCase 
)
protected

Definition at line 544 of file G4VEmAdjointModel.cc.

References indexOfUsedCrossSectionMatrix, SampleAdjSecEnergyFromCSMatrix(), and SelectCSMatrix().

545 { SelectCSMatrix(IsScatProjToProjCase);
546  return SampleAdjSecEnergyFromCSMatrix(indexOfUsedCrossSectionMatrix, aPrimEnergy, IsScatProjToProjCase);
547 }
G4double SampleAdjSecEnergyFromCSMatrix(size_t MatrixIndex, G4double prim_energy, G4bool IsScatProjToProjCase)
size_t indexOfUsedCrossSectionMatrix
void SelectCSMatrix(G4bool IsScatProjToProjCase)
G4double G4VEmAdjointModel::SampleAdjSecEnergyFromDiffCrossSectionPerAtom ( G4double  prim_energy,
G4bool  IsScatProjToProjCase 
)
protectedvirtual

Definition at line 576 of file G4VEmAdjointModel.cc.

References currentTcutForDirectSecond, DiffCrossSectionPerAtomPrimToScatPrim(), DiffCrossSectionPerAtomPrimToSecond(), G4UniformRand, GetSecondAdjEnergyMaxForProdToProjCase(), GetSecondAdjEnergyMaxForScatProjToProjCase(), GetSecondAdjEnergyMinForProdToProjCase(), and test::x.

577 {
578  // here we try to use the rejection method
579  //-----------------------------------------
580 
581  G4double E=0;
582  G4double x,xmin,greject,q;
583  if ( IsScatProjToProjCase){
585  G4double Emin= prim_energy+currentTcutForDirectSecond;
586  xmin=Emin/Emax;
587  G4double grejmax = DiffCrossSectionPerAtomPrimToScatPrim(Emin,prim_energy,1)*prim_energy;
588 
589  do {
590  q = G4UniformRand();
591  x = 1./(q*(1./xmin -1.) +1.);
592  E=x*Emax;
593  greject = DiffCrossSectionPerAtomPrimToScatPrim( E,prim_energy ,1)*prim_energy;
594 
595  }
596 
597  while( greject < G4UniformRand()*grejmax );
598 
599  }
600  else {
603  xmin=Emin/Emax;
604  G4double grejmax = DiffCrossSectionPerAtomPrimToSecond(Emin,prim_energy,1);
605  do {
606  q = G4UniformRand();
607  x = std::pow(xmin, q);
608  E=x*Emax;
609  greject = DiffCrossSectionPerAtomPrimToSecond( E,prim_energy ,1);
610 
611  }
612 
613  while( greject < G4UniformRand()*grejmax );
614 
615 
616 
617  }
618 
619  return E;
620 }
virtual G4double GetSecondAdjEnergyMaxForProdToProjCase(G4double PrimAdjEnergy)
virtual G4double DiffCrossSectionPerAtomPrimToSecond(G4double kinEnergyProj, G4double kinEnergyProd, G4double Z, G4double A=0.)
virtual G4double GetSecondAdjEnergyMaxForScatProjToProjCase(G4double PrimAdjEnergy)
#define G4UniformRand()
Definition: Randomize.hh:87
virtual G4double GetSecondAdjEnergyMinForProdToProjCase(G4double PrimAdjEnergy)
virtual G4double DiffCrossSectionPerAtomPrimToScatPrim(G4double kinEnergyProj, G4double kinEnergyScatProj, G4double Z, G4double A=0.)
G4double currentTcutForDirectSecond
double G4double
Definition: G4Types.hh:76
virtual void G4VEmAdjointModel::SampleSecondaries ( const G4Track aTrack,
G4bool  IsScatProjToProjCase,
G4ParticleChange fParticleChange 
)
pure virtual
void G4VEmAdjointModel::SelectCSMatrix ( G4bool  IsScatProjToProjCase)
protected

Definition at line 550 of file G4VEmAdjointModel.cc.

References CS_Vs_ElementForProdToProjCase, CS_Vs_ElementForScatProjToProjCase, currentMaterial, currentMaterialIndex, G4UniformRand, G4Material::GetElement(), G4Element::GetIndex(), indexOfUsedCrossSectionMatrix, lastAdjointCSForProdToProjCase, lastAdjointCSForScatProjToProjCase, lastCS, UseMatrixPerElement, and UseOnlyOneMatrixForAllElements.

Referenced by SampleAdjSecEnergyFromCSMatrix().

551 {
554  else if (!UseOnlyOneMatrixForAllElements) { //Select Material
555  std::vector<G4double>* CS_Vs_Element = &CS_Vs_ElementForScatProjToProjCase;
557  if ( !IsScatProjToProjCase) {
558  CS_Vs_Element = &CS_Vs_ElementForProdToProjCase;
560  }
561  G4double rand_var= G4UniformRand();
562  G4double SumCS=0.;
563  size_t ind=0;
564  for (size_t i=0;i<CS_Vs_Element->size();i++){
565  SumCS+=(*CS_Vs_Element)[i];
566  if (rand_var<=SumCS/lastCS){
567  ind=i;
568  break;
569  }
570  }
572  }
573 }
G4double lastAdjointCSForScatProjToProjCase
G4Material * currentMaterial
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:200
size_t indexOfUsedCrossSectionMatrix
#define G4UniformRand()
Definition: Randomize.hh:87
std::vector< G4double > CS_Vs_ElementForProdToProjCase
std::vector< G4double > CS_Vs_ElementForScatProjToProjCase
size_t GetIndex() const
Definition: G4Element.hh:181
G4bool UseOnlyOneMatrixForAllElements
G4double lastAdjointCSForProdToProjCase
double G4double
Definition: G4Types.hh:76
void G4VEmAdjointModel::SetAdjointEquivalentOfDirectPrimaryParticleDefinition ( G4ParticleDefinition aPart)

Definition at line 719 of file G4VEmAdjointModel.cc.

References G4Electron::Electron(), G4Gamma::Gamma(), G4ParticleDefinition::GetParticleName(), theAdjEquivOfDirectPrimPartDef, and theDirectPrimaryPartDef.

720 {
724  if (theAdjEquivOfDirectPrimPartDef->GetParticleName() =="adj_gamma")
726 }
G4ParticleDefinition * theDirectPrimaryPartDef
G4ParticleDefinition * theAdjEquivOfDirectPrimPartDef
const G4String & GetParticleName() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4Electron * Electron()
Definition: G4Electron.cc:94
void G4VEmAdjointModel::SetAdjointEquivalentOfDirectSecondaryParticleDefinition ( G4ParticleDefinition aPart)
inline

Definition at line 197 of file G4VEmAdjointModel.hh.

References theAdjEquivOfDirectSecondPartDef.

197  {
199  }
G4ParticleDefinition * theAdjEquivOfDirectSecondPartDef
void G4VEmAdjointModel::SetApplyCutInRange ( G4bool  aBool)
inline
virtual void G4VEmAdjointModel::SetCSBiasingFactor ( G4double  aVal)
inlinevirtual

Definition at line 217 of file G4VEmAdjointModel.hh.

References CS_biasing_factor.

Referenced by G4AdjointPhysicsList::ConstructEM().

217 {CS_biasing_factor = aVal;}
void G4VEmAdjointModel::SetCSMatrices ( std::vector< G4AdjointCSMatrix * > *  Vec1CSMatrix,
std::vector< G4AdjointCSMatrix * > *  Vec2CSMatrix 
)
inline

Definition at line 174 of file G4VEmAdjointModel.hh.

References pOnCSMatrixForProdToProjBackwardScattering, and pOnCSMatrixForScatProjToProjBackwardScattering.

Referenced by G4AdjointCSManager::BuildCrossSectionMatrices().

174  {
177 
178 
179  };
std::vector< G4AdjointCSMatrix * > * pOnCSMatrixForProdToProjBackwardScattering
std::vector< G4AdjointCSMatrix * > * pOnCSMatrixForScatProjToProjBackwardScattering
void G4VEmAdjointModel::SetHighEnergyLimit ( G4double  aVal)

Definition at line 706 of file G4VEmAdjointModel.cc.

References HighEnergyLimit, G4VEmModel::SetHighEnergyLimit(), and theDirectEMModel.

Referenced by G4AdjointPhysicsList::ConstructEM().

707 { HighEnergyLimit=aVal;
709 }
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
G4VEmModel * theDirectEMModel
void G4VEmAdjointModel::SetLowEnergyLimit ( G4double  aVal)

Definition at line 712 of file G4VEmAdjointModel.cc.

References LowEnergyLimit, G4VEmModel::SetLowEnergyLimit(), and theDirectEMModel.

Referenced by G4AdjointPhysicsList::ConstructEM().

713 {
714  LowEnergyLimit=aVal;
716 }
G4VEmModel * theDirectEMModel
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
void G4VEmAdjointModel::SetSecondPartOfSameType ( G4bool  aBool)
inline
void G4VEmAdjointModel::SetUseMatrix ( G4bool  aBool)
inline
void G4VEmAdjointModel::SetUseMatrixPerElement ( G4bool  aBool)
inline
void G4VEmAdjointModel::SetUseOnlyOneMatrixForAllElements ( G4bool  aBool)
inline

Field Documentation

G4bool G4VEmAdjointModel::ApplyCutInRange
protected
G4int G4VEmAdjointModel::ASelectedNucleus
protected
G4double G4VEmAdjointModel::CS_biasing_factor
protected
std::vector<G4double> G4VEmAdjointModel::CS_Vs_ElementForProdToProjCase
protected

Definition at line 284 of file G4VEmAdjointModel.hh.

Referenced by AdjointCrossSection(), and SelectCSMatrix().

std::vector<G4double> G4VEmAdjointModel::CS_Vs_ElementForScatProjToProjCase
protected

Definition at line 283 of file G4VEmAdjointModel.hh.

Referenced by AdjointCrossSection(), and SelectCSMatrix().

G4MaterialCutsCouple* G4VEmAdjointModel::currentCouple
protected
size_t G4VEmAdjointModel::currentCoupleIndex
protected

Definition at line 311 of file G4VEmAdjointModel.hh.

Referenced by DefineCurrentMaterial().

G4Material* G4VEmAdjointModel::currentMaterial
protected
size_t G4VEmAdjointModel::currentMaterialIndex
protected

Definition at line 310 of file G4VEmAdjointModel.hh.

Referenced by DefineCurrentMaterial(), and SelectCSMatrix().

G4double G4VEmAdjointModel::currentTcutForDirectPrim
protected

Definition at line 312 of file G4VEmAdjointModel.hh.

G4double G4VEmAdjointModel::currentTcutForDirectSecond
protected
G4double G4VEmAdjointModel::HighEnergyLimit
protected
size_t G4VEmAdjointModel::indexOfUsedCrossSectionMatrix
protected

Definition at line 347 of file G4VEmAdjointModel.hh.

Referenced by SampleAdjSecEnergyFromCSMatrix(), and SelectCSMatrix().

G4double G4VEmAdjointModel::kinEnergyProdForIntegration
protected
G4double G4VEmAdjointModel::kinEnergyProjForIntegration
protected
G4double G4VEmAdjointModel::kinEnergyScatProjForIntegration
protected
G4double G4VEmAdjointModel::lastAdjointCSForProdToProjCase
protected
G4double G4VEmAdjointModel::lastAdjointCSForScatProjToProjCase
protected
G4double G4VEmAdjointModel::lastCS
protected
G4double G4VEmAdjointModel::LowEnergyLimit
protected

Definition at line 329 of file G4VEmAdjointModel.hh.

Referenced by GetLowEnergyLimit(), and SetLowEnergyLimit().

G4double G4VEmAdjointModel::mass_ratio_product
protected

Definition at line 321 of file G4VEmAdjointModel.hh.

Referenced by G4VEmAdjointModel().

G4double G4VEmAdjointModel::mass_ratio_projectile
protected

Definition at line 322 of file G4VEmAdjointModel.hh.

Referenced by G4VEmAdjointModel().

size_t G4VEmAdjointModel::model_index
protected

Definition at line 349 of file G4VEmAdjointModel.hh.

Referenced by G4VEmAdjointModel().

const G4String G4VEmAdjointModel::name
protected

Definition at line 264 of file G4VEmAdjointModel.hh.

Referenced by GetName().

std::vector< G4AdjointCSMatrix* >* G4VEmAdjointModel::pOnCSMatrixForProdToProjBackwardScattering
protected

Definition at line 281 of file G4VEmAdjointModel.hh.

Referenced by SetCSMatrices().

std::vector< G4AdjointCSMatrix* >* G4VEmAdjointModel::pOnCSMatrixForScatProjToProjBackwardScattering
protected

Definition at line 282 of file G4VEmAdjointModel.hh.

Referenced by SetCSMatrices().

G4VParticleChange* G4VEmAdjointModel::pParticleChange
protected

Definition at line 256 of file G4VEmAdjointModel.hh.

G4double G4VEmAdjointModel::preStepEnergy
protected

Definition at line 304 of file G4VEmAdjointModel.hh.

Referenced by AdjointCrossSection(), and CorrectPostStepWeight().

G4bool G4VEmAdjointModel::second_part_of_same_type
protected
G4Material* G4VEmAdjointModel::SelectedMaterial
protected
G4ParticleDefinition* G4VEmAdjointModel::theAdjEquivOfDirectPrimPartDef
protected
G4ParticleDefinition* G4VEmAdjointModel::theAdjEquivOfDirectSecondPartDef
protected
G4VEmModel* G4VEmAdjointModel::theDirectEMModel
protected
G4ParticleDefinition* G4VEmAdjointModel::theDirectPrimaryPartDef
protected
G4bool G4VEmAdjointModel::UseMatrix
protected
G4bool G4VEmAdjointModel::UseMatrixPerElement
protected
G4bool G4VEmAdjointModel::UseOnlyOneMatrixForAllElements
protected
G4int G4VEmAdjointModel::ZSelectedNucleus
protected

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