Geant4  10.02.p03
G4DNABornIonisationModel1 Class Reference

#include <G4DNABornIonisationModel1.hh>

Inheritance diagram for G4DNABornIonisationModel1:
Collaboration diagram for G4DNABornIonisationModel1:

Public Member Functions

 G4DNABornIonisationModel1 (const G4ParticleDefinition *p=0, const G4String &nam="DNABornIonisationModel")
 
virtual ~G4DNABornIonisationModel1 ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector()))
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle *> *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int, const G4ParticleDefinition *, G4double)
 
double DifferentialCrossSection (G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
 
G4double TransferedEnergy (G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell, G4double random)
 
void SelectFasterComputation (G4bool input)
 
void SelectStationary (G4bool input)
 
void SelectSPScaling (G4bool input)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector *> *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 

Protected Attributes

G4ParticleChangeForGamma * fParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChange * pParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Private Types

typedef std::map< G4String, G4String, std::less< G4String > > MapFile
 
typedef std::map< G4String, G4DNACrossSectionDataSet *, std::less< G4String > > MapData
 
typedef std::map< double, std::map< double, double > > TriDimensionMap
 
typedef std::map< double, std::vector< double > > VecMap
 

Private Member Functions

G4double RandomizeEjectedElectronEnergy (G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell)
 
G4double RandomizeEjectedElectronEnergyFromCumulatedDcs (G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell)
 
G4double Interpolate (G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
 
G4double QuadInterpolator (G4double e11, G4double e12, G4double e21, G4double e22, G4double x11, G4double x12, G4double x21, G4double x22, G4double t1, G4double t2, G4double t, G4double e)
 
G4int RandomSelect (G4double energy, const G4String &particle)
 
G4DNABornIonisationModel1operator= (const G4DNABornIonisationModel1 &right)
 
 G4DNABornIonisationModel1 (const G4DNABornIonisationModel1 &)
 

Private Attributes

G4bool fasterCode
 
G4bool statCode
 
G4bool spScaling
 
const std::vector< G4double > * fpMolWaterDensity
 
G4VAtomDeexcitationfAtomDeexcitation
 
std::map< G4String, G4double, std::less< G4String > > lowEnergyLimit
 
std::map< G4String, G4double, std::less< G4String > > highEnergyLimit
 
G4bool isInitialised
 
G4int verboseLevel
 
MapFile tableFile
 
MapData tableData
 
G4DNAWaterIonisationStructure waterStructure
 
TriDimensionMap eDiffCrossSectionData [6]
 
TriDimensionMap eNrjTransfData [6]
 
TriDimensionMap pDiffCrossSectionData [6]
 
TriDimensionMap pNrjTransfData [6]
 
std::vector< double > eTdummyVec
 
std::vector< double > pTdummyVec
 
VecMap eVecm
 
VecMap pVecm
 
VecMap eProbaShellMap [6]
 
VecMap pProbaShellMap [6]
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLoss * GetParticleChangeForLoss ()
 
G4ParticleChangeForGamma * GetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 48 of file G4DNABornIonisationModel1.hh.

Member Typedef Documentation

◆ MapData

Definition at line 119 of file G4DNABornIonisationModel1.hh.

◆ MapFile

typedef std::map<G4String,G4String,std::less<G4String> > G4DNABornIonisationModel1::MapFile
private

Definition at line 116 of file G4DNABornIonisationModel1.hh.

◆ TriDimensionMap

typedef std::map<double, std::map<double, double> > G4DNABornIonisationModel1::TriDimensionMap
private

Definition at line 145 of file G4DNABornIonisationModel1.hh.

◆ VecMap

typedef std::map<double, std::vector<double> > G4DNABornIonisationModel1::VecMap
private

Definition at line 156 of file G4DNABornIonisationModel1.hh.

Constructor & Destructor Documentation

◆ G4DNABornIonisationModel1() [1/2]

G4DNABornIonisationModel1::G4DNABornIonisationModel1 ( const G4ParticleDefinition p = 0,
const G4String nam = "DNABornIonisationModel" 
)

Definition at line 45 of file G4DNABornIonisationModel1.cc.

46  :
47  G4VEmModel(nam), isInitialised(false)
48 {
49  verboseLevel = 0;
50  // Verbosity scale:
51  // 0 = nothing
52  // 1 = warning for energy non-conservation
53  // 2 = details of energy budget
54  // 3 = calculation of cross sections, file openings, sampling of atoms
55  // 4 = entering in methods
56 
57  if (verboseLevel > 0)
58  {
59  G4cout << "Born ionisation model is constructed " << G4endl;
60  }
61 
62  //Mark this model as "applicable" for atomic deexcitation
63  SetDeexcitationFlag(true);
67 
68  // define default angular generator
70 
71  // Selection of computation method
72 
73  fasterCode = false;
74 
75  // Selection of stationary mode
76 
77  statCode = false;
78 
79  // Selection of SP scaling
80 
81  spScaling = true;
82 }
const std::vector< G4double > * fpMolWaterDensity
G4VAtomDeexcitation * fAtomDeexcitation
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
G4ParticleChangeForGamma * fParticleChangeForGamma
G4GLOB_DLL std::ostream G4cout
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:624
#define G4endl
Definition: G4ios.hh:61
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:781
Here is the call graph for this function:

◆ ~G4DNABornIonisationModel1()

G4DNABornIonisationModel1::~G4DNABornIonisationModel1 ( )
virtual

Definition at line 86 of file G4DNABornIonisationModel1.cc.

87 {
88  // Cross section
89 
90  std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> >::iterator pos;
91  for (pos = tableData.begin(); pos != tableData.end(); ++pos)
92  {
93  G4DNACrossSectionDataSet* table = pos->second;
94  delete table;
95  }
96 
97  // Final state
98 
99  eVecm.clear();
100  pVecm.clear();
101 }
static const G4double pos

◆ G4DNABornIonisationModel1() [2/2]

G4DNABornIonisationModel1::G4DNABornIonisationModel1 ( const G4DNABornIonisationModel1 )
private

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4DNABornIonisationModel1::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  ekin,
G4double  emin,
G4double  emax 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 314 of file G4DNABornIonisationModel1.cc.

319 {
320  if (verboseLevel > 3)
321  {
322  G4cout << "Calling CrossSectionPerVolume() of G4DNABornIonisationModel1"
323  << G4endl;
324 
325  }
326 
327  if (
328  particleDefinition != G4Proton::ProtonDefinition()
329  &&
330  particleDefinition != G4Electron::ElectronDefinition()
331  )
332 
333  return 0;
334 
335  // Calculate total cross section for model
336 
337  G4double lowLim = 0;
338  G4double highLim = 0;
339  G4double sigma=0;
340 
341  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
342 
343  if(waterDensity!= 0.0)
344  // if (material == nistwater || material->GetBaseMaterial() == nistwater)
345  {
346  const G4String& particleName = particleDefinition->GetParticleName();
347 
348  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
349  pos1 = lowEnergyLimit.find(particleName);
350  if (pos1 != lowEnergyLimit.end())
351  {
352  lowLim = pos1->second;
353  }
354 
355  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
356  pos2 = highEnergyLimit.find(particleName);
357  if (pos2 != highEnergyLimit.end())
358  {
359  highLim = pos2->second;
360  }
361 
362  if (ekin >= lowLim && ekin < highLim)
363  {
364  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
365  pos = tableData.find(particleName);
366 
367  if (pos != tableData.end())
368  {
369  G4DNACrossSectionDataSet* table = pos->second;
370  if (table != 0)
371  {
372  sigma = table->FindValue(ekin);
373 
374  // ICRU49 electronic SP scaling - ZF, SI
375 
376  if (particleDefinition == G4Proton::ProtonDefinition() && ekin < 70*MeV && spScaling)
377  {
378  G4double A = 1.39241700556072800000E-009 ;
379  G4double B = -8.52610412942622630000E-002 ;
380  sigma = sigma * std::exp(A*(ekin/eV)+B);
381  }
382  //
383 
384  }
385  }
386  else
387  {
388  G4Exception("G4DNABornIonisationModel1::CrossSectionPerVolume","em0002",
389  FatalException,"Model not applicable to particle type.");
390  }
391  }
392 
393  if (verboseLevel > 2)
394  {
395  G4cout << "__________________________________" << G4endl;
396  G4cout << "G4DNABornIonisationModel1 - XS INFO START" << G4endl;
397  G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
398  G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
399  G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
400  G4cout << "G4DNABornIonisationModel1 - XS INFO END" << G4endl;
401  }
402  } // if (waterMaterial)
403 
404  return sigma*waterDensity;
405 }
static const double cm
Definition: G4SIunits.hh:118
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static const double MeV
Definition: G4SIunits.hh:211
std::map< G4String, G4double, std::less< G4String > > lowEnergyLimit
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
size_t GetIndex() const
Definition: G4Material.hh:262
virtual G4double FindValue(G4double e, G4int componentId=0) const
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
std::map< G4String, G4double, std::less< G4String > > highEnergyLimit
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const double eV
Definition: G4SIunits.hh:212
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static const G4double pos
Here is the call graph for this function:

◆ DifferentialCrossSection()

double G4DNABornIonisationModel1::DifferentialCrossSection ( G4ParticleDefinition aParticleDefinition,
G4double  k,
G4double  energyTransfer,
G4int  shell 
)

Definition at line 724 of file G4DNABornIonisationModel1.cc.

728 {
729  G4double sigma = 0.;
730 
731  if (energyTransfer >= waterStructure.IonisationEnergy(ionizationLevelIndex))
732  {
733  G4double valueT1 = 0;
734  G4double valueT2 = 0;
735  G4double valueE21 = 0;
736  G4double valueE22 = 0;
737  G4double valueE12 = 0;
738  G4double valueE11 = 0;
739 
740  G4double xs11 = 0;
741  G4double xs12 = 0;
742  G4double xs21 = 0;
743  G4double xs22 = 0;
744 
745  if (particleDefinition == G4Electron::ElectronDefinition())
746  {
747  // k should be in eV and energy transfer eV also
748 
749  std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),
750  eTdummyVec.end(),
751  k);
752 
753  std::vector<double>::iterator t1 = t2 - 1;
754 
755  // SI : the following condition avoids situations where energyTransfer >last vector element
756  if (energyTransfer <= eVecm[(*t1)].back()
757  && energyTransfer <= eVecm[(*t2)].back())
758  {
759  std::vector<double>::iterator e12 =
760  std::upper_bound(eVecm[(*t1)].begin(),
761  eVecm[(*t1)].end(),
762  energyTransfer);
763  std::vector<double>::iterator e11 = e12 - 1;
764 
765  std::vector<double>::iterator e22 =
766  std::upper_bound(eVecm[(*t2)].begin(),
767  eVecm[(*t2)].end(),
768  energyTransfer);
769  std::vector<double>::iterator e21 = e22 - 1;
770 
771  valueT1 = *t1;
772  valueT2 = *t2;
773  valueE21 = *e21;
774  valueE22 = *e22;
775  valueE12 = *e12;
776  valueE11 = *e11;
777 
778  xs11 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
779  xs12 = eDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
780  xs21 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
781  xs22 = eDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
782 
783  }
784 
785  }
786 
787  if (particleDefinition == G4Proton::ProtonDefinition())
788  {
789  // k should be in eV and energy transfer eV also
790  std::vector<double>::iterator t2 = std::upper_bound(pTdummyVec.begin(),
791  pTdummyVec.end(),
792  k);
793  std::vector<double>::iterator t1 = t2 - 1;
794 
795  std::vector<double>::iterator e12 = std::upper_bound(pVecm[(*t1)].begin(),
796  pVecm[(*t1)].end(),
797  energyTransfer);
798  std::vector<double>::iterator e11 = e12 - 1;
799 
800  std::vector<double>::iterator e22 = std::upper_bound(pVecm[(*t2)].begin(),
801  pVecm[(*t2)].end(),
802  energyTransfer);
803  std::vector<double>::iterator e21 = e22 - 1;
804 
805  valueT1 = *t1;
806  valueT2 = *t2;
807  valueE21 = *e21;
808  valueE22 = *e22;
809  valueE12 = *e12;
810  valueE11 = *e11;
811 
812  xs11 = pDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE11];
813  xs12 = pDiffCrossSectionData[ionizationLevelIndex][valueT1][valueE12];
814  xs21 = pDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE21];
815  xs22 = pDiffCrossSectionData[ionizationLevelIndex][valueT2][valueE22];
816 
817  }
818 
819  G4double xsProduct = xs11 * xs12 * xs21 * xs22;
820  if (xsProduct != 0.)
821  {
822  sigma = QuadInterpolator(valueE11,
823  valueE12,
824  valueE21,
825  valueE22,
826  xs11,
827  xs12,
828  xs21,
829  xs22,
830  valueT1,
831  valueT2,
832  k,
833  energyTransfer);
834  }
835  }
836 
837  return sigma;
838 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
TTree * t1
Definition: plottest35.C:26
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
G4double QuadInterpolator(G4double e11, G4double e12, G4double e21, G4double e22, G4double x11, G4double x12, G4double x21, G4double x22, G4double t1, G4double t2, G4double t, G4double e)
G4DNAWaterIonisationStructure waterStructure
TTree * t2
Definition: plottest35.C:36
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetPartialCrossSection()

G4double G4DNABornIonisationModel1::GetPartialCrossSection ( const G4Material ,
G4int  level,
const G4ParticleDefinition particle,
G4double  kineticEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 930 of file G4DNABornIonisationModel1.cc.

934 {
935  std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> >::iterator pos;
936  pos = tableData.find(particle->GetParticleName());
937 
938  if (pos != tableData.end())
939  {
940  G4DNACrossSectionDataSet* table = pos->second;
941  return table->GetComponent(level)->FindValue(kineticEnergy);
942  }
943 
944  return 0;
945 }
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
const G4String & GetParticleName() const
virtual const G4VEMDataSet * GetComponent(G4int componentId) const
static const G4double pos
Here is the call graph for this function:

◆ Initialise()

void G4DNABornIonisationModel1::Initialise ( const G4ParticleDefinition particle,
const G4DataVector = *(new G4DataVector()) 
)
virtual

Implements G4VEmModel.

Definition at line 105 of file G4DNABornIonisationModel1.cc.

107 {
108 
109  if (verboseLevel > 3)
110  {
111  G4cout << "Calling G4DNABornIonisationModel1::Initialise()" << G4endl;
112  }
113 
114  // Energy limits
115 
116  G4String fileElectron("dna/sigma_ionisation_e_born");
117  G4String fileProton("dna/sigma_ionisation_p_born");
118 
121 
124 
125  G4double scaleFactor = (1.e-22 / 3.343) * m*m;
126 
127  char *path = getenv("G4LEDATA");
128 
129  // *** ELECTRON
130 
131  electron = electronDef->GetParticleName();
132 
133  tableFile[electron] = fileElectron;
134 
135  lowEnergyLimit[electron] = 11. * eV;
136  highEnergyLimit[electron] = 1. * MeV;
137 
138  // Cross section
139 
141  tableE->LoadData(fileElectron);
142 
143  tableData[electron] = tableE;
144 
145  // Final state
146 
147  std::ostringstream eFullFileName;
148 
149  if (fasterCode) eFullFileName << path << "/dna/sigmadiff_cumulated_ionisation_e_born_hp.dat";
150  if (!fasterCode) eFullFileName << path << "/dna/sigmadiff_ionisation_e_born.dat";
151 
152  std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
153 
154  if (!eDiffCrossSection)
155  {
156  if (fasterCode) G4Exception("G4DNABornIonisationModel1::Initialise","em0003",
157  FatalException,"Missing data file:/dna/sigmadiff_cumulated_ionisation_e_born_hp.dat");
158 
159  if (!fasterCode) G4Exception("G4DNABornIonisationModel1::Initialise","em0003",
160  FatalException,"Missing data file:/dna/sigmadiff_ionisation_e_born.dat");
161  }
162 
163  //
164 
165  // Clear the arrays for re-initialization case (MT mode)
166  // March 25th, 2014 - Vaclav Stepan, Sebastien Incerti
167 
168  eTdummyVec.clear();
169  pTdummyVec.clear();
170 
171  eVecm.clear();
172  pVecm.clear();
173 
174  for (int j=0; j<5; j++)
175  {
176  eProbaShellMap[j].clear();
177  pProbaShellMap[j].clear();
178 
179  eDiffCrossSectionData[j].clear();
180  pDiffCrossSectionData[j].clear();
181 
182  eNrjTransfData[j].clear();
183  pNrjTransfData[j].clear();
184  }
185  //
186 
187  eTdummyVec.push_back(0.);
188  while(!eDiffCrossSection.eof())
189  {
190  double tDummy;
191  double eDummy;
192  eDiffCrossSection>>tDummy>>eDummy;
193  if (tDummy != eTdummyVec.back()) eTdummyVec.push_back(tDummy);
194 
195  double tmp;
196  for (int j=0; j<5; j++)
197  {
198  eDiffCrossSection>> tmp;
199 
200  eDiffCrossSectionData[j][tDummy][eDummy] = tmp;
201 
202  if (fasterCode)
203  {
204  eNrjTransfData[j][tDummy][eDiffCrossSectionData[j][tDummy][eDummy]]=eDummy;
205  eProbaShellMap[j][tDummy].push_back(eDiffCrossSectionData[j][tDummy][eDummy]);
206  }
207 
208  // SI - only if eof is not reached
209  if (!eDiffCrossSection.eof() && !fasterCode) eDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor;
210 
211  if (!fasterCode) eVecm[tDummy].push_back(eDummy);
212 
213  }
214  }
215 
216  // *** PROTON
217 
218  proton = protonDef->GetParticleName();
219 
220  tableFile[proton] = fileProton;
221 
222  lowEnergyLimit[proton] = 500. * keV;
223  highEnergyLimit[proton] = 100. * MeV;
224 
225  // Cross section
226 
228  tableP->LoadData(fileProton);
229 
230  tableData[proton] = tableP;
231 
232  // Final state
233 
234  std::ostringstream pFullFileName;
235 
236  if (fasterCode) pFullFileName << path << "/dna/sigmadiff_cumulated_ionisation_p_born_hp.dat";
237 
238  if (!fasterCode) pFullFileName << path << "/dna/sigmadiff_ionisation_p_born.dat";
239 
240  std::ifstream pDiffCrossSection(pFullFileName.str().c_str());
241 
242  if (!pDiffCrossSection)
243  {
244  if (fasterCode) G4Exception("G4DNABornIonisationModel1::Initialise","em0003",
245  FatalException,"Missing data file:/dna/sigmadiff_cumulated_ionisation_p_born_hp.dat");
246 
247  if (!fasterCode) G4Exception("G4DNABornIonisationModel1::Initialise","em0003",
248  FatalException,"Missing data file:/dna/sigmadiff_ionisation_p_born.dat");
249  }
250 
251  pTdummyVec.push_back(0.);
252  while(!pDiffCrossSection.eof())
253  {
254  double tDummy;
255  double eDummy;
256  pDiffCrossSection>>tDummy>>eDummy;
257  if (tDummy != pTdummyVec.back()) pTdummyVec.push_back(tDummy);
258  for (int j=0; j<5; j++)
259  {
260  pDiffCrossSection>>pDiffCrossSectionData[j][tDummy][eDummy];
261 
262  if (fasterCode)
263  {
264  pNrjTransfData[j][tDummy][pDiffCrossSectionData[j][tDummy][eDummy]]=eDummy;
265  pProbaShellMap[j][tDummy].push_back(pDiffCrossSectionData[j][tDummy][eDummy]);
266  }
267 
268  // SI - only if eof is not reached !
269  if (!pDiffCrossSection.eof() && !fasterCode) pDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor;
270 
271  if (!fasterCode) pVecm[tDummy].push_back(eDummy);
272  }
273  }
274 
275  //
276 
277  if (particle==electronDef)
278  {
281  }
282 
283  if (particle==protonDef)
284  {
287  }
288 
289  if( verboseLevel>0 )
290  {
291  G4cout << "Born ionisation model is initialized " << G4endl
292  << "Energy range: "
293  << LowEnergyLimit() / eV << " eV - "
294  << HighEnergyLimit() / keV << " keV for "
295  << particle->GetParticleName()
296  << G4endl;
297  }
298 
299  // Initialize water density pointer
301  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
302 
303  //
305 
306  if (isInitialised)
307  { return;}
309  isInitialised = true;
310 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static const double MeV
Definition: G4SIunits.hh:211
std::map< G4String, G4double, std::less< G4String > > lowEnergyLimit
Float_t tmp
static G4LossTableManager * Instance()
const std::vector< G4double > * fpMolWaterDensity
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
virtual G4bool LoadData(const G4String &argFileName)
G4VAtomDeexcitation * fAtomDeexcitation
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
G4ParticleChangeForGamma * fParticleChangeForGamma
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
std::map< G4String, G4double, std::less< G4String > > highEnergyLimit
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const double eV
Definition: G4SIunits.hh:212
static G4DNAMolecularMaterial * Instance()
#define G4endl
Definition: G4ios.hh:61
static const double m
Definition: G4SIunits.hh:128
G4VAtomDeexcitation * AtomDeexcitation()
static const double keV
Definition: G4SIunits.hh:213
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
Here is the call graph for this function:

◆ Interpolate()

G4double G4DNABornIonisationModel1::Interpolate ( G4double  e1,
G4double  e2,
G4double  e,
G4double  xs1,
G4double  xs2 
)
private

Definition at line 842 of file G4DNABornIonisationModel1.cc.

847 {
848  G4double value = 0.;
849 
850  // Log-log interpolation by default
851 
852  if (e1 != 0 && e2 != 0 && (std::log10(e2) - std::log10(e1)) != 0
853  && !fasterCode)
854  {
855  G4double a = (std::log10(xs2) - std::log10(xs1))
856  / (std::log10(e2) - std::log10(e1));
857  G4double b = std::log10(xs2) - a * std::log10(e2);
858  G4double sigma = a * std::log10(e) + b;
859  value = (std::pow(10., sigma));
860  }
861 
862  // Switch to lin-lin interpolation
863  /*
864  if ((e2-e1)!=0)
865  {
866  G4double d1 = xs1;
867  G4double d2 = xs2;
868  value = (d1 + (d2 - d1)*(e - e1)/ (e2 - e1));
869  }
870  */
871 
872  // Switch to log-lin interpolation for faster code
873  if ((e2 - e1) != 0 && xs1 != 0 && xs2 != 0 && fasterCode)
874  {
875  G4double d1 = std::log10(xs1);
876  G4double d2 = std::log10(xs2);
877  value = std::pow(10., (d1 + (d2 - d1) * (e - e1) / (e2 - e1)));
878  }
879 
880  // Switch to lin-lin interpolation for faster code
881  // in case one of xs1 or xs2 (=cum proba) value is zero
882 
883  if ((e2 - e1) != 0 && (xs1 == 0 || xs2 == 0) && fasterCode)
884  {
885  G4double d1 = xs1;
886  G4double d2 = xs2;
887  value = (d1 + (d2 - d1) * (e - e1) / (e2 - e1));
888  }
889 
890  /*
891  G4cout
892  << e1 << " "
893  << e2 << " "
894  << e << " "
895  << xs1 << " "
896  << xs2 << " "
897  << value
898  << G4endl;
899  */
900 
901  return value;
902 }
static const G4double d1
static const G4double e2
static const G4double e1
double G4double
Definition: G4Types.hh:76
static const G4double d2
Here is the caller graph for this function:

◆ operator=()

G4DNABornIonisationModel1& G4DNABornIonisationModel1::operator= ( const G4DNABornIonisationModel1 right)
private

◆ QuadInterpolator()

G4double G4DNABornIonisationModel1::QuadInterpolator ( G4double  e11,
G4double  e12,
G4double  e21,
G4double  e22,
G4double  x11,
G4double  x12,
G4double  x21,
G4double  x22,
G4double  t1,
G4double  t2,
G4double  t,
G4double  e 
)
private

Definition at line 906 of file G4DNABornIonisationModel1.cc.

918 {
919  G4double interpolatedvalue1 = Interpolate(e11, e12, e, xs11, xs12);
920  G4double interpolatedvalue2 = Interpolate(e21, e22, e, xs21, xs22);
921  G4double value = Interpolate(t1,
922  t2,
923  t,
924  interpolatedvalue1,
925  interpolatedvalue2);
926 
927  return value;
928 }
TTree * t1
Definition: plottest35.C:26
G4double Interpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
TTree * t2
Definition: plottest35.C:36
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ RandomizeEjectedElectronEnergy()

G4double G4DNABornIonisationModel1::RandomizeEjectedElectronEnergy ( G4ParticleDefinition aParticleDefinition,
G4double  incomingParticleEnergy,
G4int  shell 
)
private

Definition at line 586 of file G4DNABornIonisationModel1.cc.

589 {
590  // G4cout << "*** SLOW computation for " << " " << particleDefinition->GetParticleName() << G4endl;
591 
592  if (particleDefinition == G4Electron::ElectronDefinition())
593  {
594  G4double maximumEnergyTransfer = 0.;
595  if ((k + waterStructure.IonisationEnergy(shell)) / 2. > k)
596  maximumEnergyTransfer = k;
597  else
598  maximumEnergyTransfer = (k + waterStructure.IonisationEnergy(shell)) / 2.;
599 
600  // SI : original method
601  /*
602  G4double crossSectionMaximum = 0.;
603  for(G4double value=waterStructure.IonisationEnergy(shell); value<=maximumEnergyTransfer; value+=0.1*eV)
604  {
605  G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
606  if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
607  }
608  */
609 
610  // SI : alternative method
611  G4double crossSectionMaximum = 0.;
612 
613  G4double minEnergy = waterStructure.IonisationEnergy(shell);
614  G4double maxEnergy = maximumEnergyTransfer;
615  G4int nEnergySteps = 50;
616 
617  G4double value(minEnergy);
618  G4double stpEnergy(std::pow(maxEnergy / value,
619  1. / static_cast<G4double>(nEnergySteps - 1)));
620  G4int step(nEnergySteps);
621  while (step > 0)
622  {
623  step--;
624  G4double differentialCrossSection =
625  DifferentialCrossSection(particleDefinition,
626  k / eV,
627  value / eV,
628  shell);
629  if (differentialCrossSection >= crossSectionMaximum)
630  crossSectionMaximum = differentialCrossSection;
631  value *= stpEnergy;
632  }
633  //
634 
635  G4double secondaryElectronKineticEnergy = 0.;
636  do
637  {
638  secondaryElectronKineticEnergy = G4UniformRand()* (maximumEnergyTransfer-waterStructure.IonisationEnergy(shell));
639  }while(G4UniformRand()*crossSectionMaximum >
640  DifferentialCrossSection(particleDefinition, k/eV,
641  (secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
642 
643  return secondaryElectronKineticEnergy;
644 
645  }
646 
647  else if (particleDefinition == G4Proton::ProtonDefinition())
648  {
649  G4double maximumKineticEnergyTransfer = 4.
651 
652  G4double crossSectionMaximum = 0.;
653  for (G4double value = waterStructure.IonisationEnergy(shell);
654  value <= 4. * waterStructure.IonisationEnergy(shell); value += 0.1 * eV)
655  {
656  G4double differentialCrossSection =
657  DifferentialCrossSection(particleDefinition,
658  k / eV,
659  value / eV,
660  shell);
661  if (differentialCrossSection >= crossSectionMaximum)
662  crossSectionMaximum = differentialCrossSection;
663  }
664 
665  G4double secondaryElectronKineticEnergy = 0.;
666  do
667  {
668  secondaryElectronKineticEnergy = G4UniformRand()* maximumKineticEnergyTransfer;
669  }while(G4UniformRand()*crossSectionMaximum >=
670  DifferentialCrossSection(particleDefinition, k/eV,
671  (secondaryElectronKineticEnergy+waterStructure.IonisationEnergy(shell))/eV,shell));
672 
673  return secondaryElectronKineticEnergy;
674  }
675 
676  return 0;
677 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:97
float proton_mass_c2
Definition: hepunit.py:275
float electron_mass_c2
Definition: hepunit.py:274
G4DNAWaterIonisationStructure waterStructure
static const double eV
Definition: G4SIunits.hh:212
double DifferentialCrossSection(G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ RandomizeEjectedElectronEnergyFromCumulatedDcs()

G4double G4DNABornIonisationModel1::RandomizeEjectedElectronEnergyFromCumulatedDcs ( G4ParticleDefinition aParticleDefinition,
G4double  incomingParticleEnergy,
G4int  shell 
)
private

Definition at line 1008 of file G4DNABornIonisationModel1.cc.

1011 {
1012  //G4cout << "*** FAST computation for " << " " << particleDefinition->GetParticleName() << G4endl;
1013 
1014  G4double secondaryElectronKineticEnergy = 0.;
1015 
1016  G4double random = G4UniformRand();
1017 
1018  secondaryElectronKineticEnergy = TransferedEnergy(particleDefinition,
1019  k / eV,
1020  shell,
1021  random) * eV
1022  - waterStructure.IonisationEnergy(shell);
1023 
1024  //G4cout << RandomTransferedEnergy(particleDefinition, k/eV, shell) << G4endl;
1025  // SI - 01/04/2014
1026  if (secondaryElectronKineticEnergy < 0.)
1027  return 0.;
1028  //
1029 
1030  return secondaryElectronKineticEnergy;
1031 }
G4double TransferedEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell, G4double random)
#define G4UniformRand()
Definition: Randomize.hh:97
G4DNAWaterIonisationStructure waterStructure
static const double eV
Definition: G4SIunits.hh:212
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ RandomSelect()

G4int G4DNABornIonisationModel1::RandomSelect ( G4double  energy,
const G4String particle 
)
private

Definition at line 949 of file G4DNABornIonisationModel1.cc.

951 {
952  G4int level = 0;
953 
954  std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> >::iterator pos;
955  pos = tableData.find(particle);
956 
957  if (pos != tableData.end())
958  {
959  G4DNACrossSectionDataSet* table = pos->second;
960 
961  if (table != 0)
962  {
963  G4double* valuesBuffer = new G4double[table->NumberOfComponents()];
964  const size_t n(table->NumberOfComponents());
965  size_t i(n);
966  G4double value = 0.;
967 
968  while (i > 0)
969  {
970  i--;
971  valuesBuffer[i] = table->GetComponent(i)->FindValue(k);
972  value += valuesBuffer[i];
973  }
974 
975  value *= G4UniformRand();
976 
977  i = n;
978 
979  while (i > 0)
980  {
981  i--;
982 
983  if (valuesBuffer[i] > value)
984  {
985  delete[] valuesBuffer;
986  return i;
987  }
988  value -= valuesBuffer[i];
989  }
990 
991  if (valuesBuffer)
992  delete[] valuesBuffer;
993 
994  }
995  } else
996  {
997  G4Exception("G4DNABornIonisationModel1::RandomSelect",
998  "em0002",
1000  "Model not applicable to particle type.");
1001  }
1002 
1003  return level;
1004 }
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
int G4int
Definition: G4Types.hh:78
Char_t n[5]
#define G4UniformRand()
Definition: Randomize.hh:97
virtual const G4VEMDataSet * GetComponent(G4int componentId) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
double G4double
Definition: G4Types.hh:76
static const G4double pos
virtual size_t NumberOfComponents(void) const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SampleSecondaries()

void G4DNABornIonisationModel1::SampleSecondaries ( std::vector< G4DynamicParticle *> *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle particle,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 409 of file G4DNABornIonisationModel1.cc.

414 {
415 
416  if (verboseLevel > 3)
417  {
418  G4cout << "Calling SampleSecondaries() of G4DNABornIonisationModel1"
419  << G4endl;
420  }
421 
422  G4double lowLim = 0;
423  G4double highLim = 0;
424 
425  G4double k = particle->GetKineticEnergy();
426 
427  const G4String& particleName = particle->GetDefinition()->GetParticleName();
428 
429  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
430  pos1 = lowEnergyLimit.find(particleName);
431 
432  if (pos1 != lowEnergyLimit.end())
433  {
434  lowLim = pos1->second;
435  }
436 
437  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
438  pos2 = highEnergyLimit.find(particleName);
439 
440  if (pos2 != highEnergyLimit.end())
441  {
442  highLim = pos2->second;
443  }
444 
445  if (k >= lowLim && k < highLim)
446  {
447  G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
448  G4double particleMass = particle->GetDefinition()->GetPDGMass();
449  G4double totalEnergy = k + particleMass;
450  G4double pSquare = k * (totalEnergy + particleMass);
451  G4double totalMomentum = std::sqrt(pSquare);
452 
453  G4int ionizationShell = 0;
454 
455  if (!fasterCode) ionizationShell = RandomSelect(k,particleName);
456 
457  // SI: The following protection is necessary to avoid infinite loops :
458  // sigmadiff_ionisation_e_born.dat has non zero partial xs at 18 eV for shell 3 (ionizationShell ==2)
459  // sigmadiff_cumulated_ionisation_e_born.dat has zero cumulated partial xs at 18 eV for shell 3 (ionizationShell ==2)
460  // this is due to the fact that the max allowed transfered energy is (18+10.79)/2=17.025 eV and only transfered energies
461  // strictly above this value have non zero partial xs in sigmadiff_ionisation_e_born.dat (starting at trans = 17.12 eV)
462 
463  if (fasterCode)
464  do
465  {
466  ionizationShell = RandomSelect(k,particleName);
467  }while (k<19*eV && ionizationShell==2 && particle->GetDefinition()==G4Electron::ElectronDefinition());
468 
469  // AM: sample deexcitation
470  // here we assume that H_{2}O electronic levels are the same as Oxygen.
471  // this can be considered true with a rough 10% error in energy on K-shell,
472 
473  G4int secNumberInit = 0;// need to know at a certain point the energy of secondaries
474  G4int secNumberFinal = 0;// So I'll make the diference and then sum the energies
475 
476  G4double bindingEnergy = 0;
477  bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
478 
479  //SI: additional protection if tcs interpolation method is modified
480  if (k<bindingEnergy) return;
481  //
482 
483  G4int Z = 8;
485  {
487 
488  if (ionizationShell <5 && ionizationShell >1)
489  {
490  as = G4AtomicShellEnumerator(4-ionizationShell);
491  }
492  else if (ionizationShell <2)
493  {
494  as = G4AtomicShellEnumerator(3);
495  }
496 
497  // FOR DEBUG ONLY
498  // if (ionizationShell == 4) {
499  //
500  // G4cout << "Z: " << Z << " as: " << as
501  // << " ionizationShell: " << ionizationShell << " bindingEnergy: "<< bindingEnergy/eV << G4endl;
502  // G4cout << "Press <Enter> key to continue..." << G4endl;
503  // G4cin.ignore();
504  // }
505 
506  const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
507  secNumberInit = fvect->size();
508  fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
509  secNumberFinal = fvect->size();
510  }
511 
512  G4double secondaryKinetic=-1000*eV;
513 
514  if (fasterCode == false)
515  {
516  secondaryKinetic = RandomizeEjectedElectronEnergy(particle->GetDefinition(),k,ionizationShell);
517  }
518  // SI - 01/04/2014
519  else
520  {
521  secondaryKinetic = RandomizeEjectedElectronEnergyFromCumulatedDcs(particle->GetDefinition(),k,ionizationShell);
522  }
523  //
524 
525  G4ThreeVector deltaDirection =
526  GetAngularDistribution()->SampleDirectionForShell(particle, secondaryKinetic,
527  Z, ionizationShell,
528  couple->GetMaterial());
529 
530  if (particle->GetDefinition() == G4Electron::ElectronDefinition())
531  {
532  G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
533 
534  G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
535  G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
536  G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
537  G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
538  finalPx /= finalMomentum;
539  finalPy /= finalMomentum;
540  finalPz /= finalMomentum;
541 
542  G4ThreeVector direction;
543  direction.set(finalPx,finalPy,finalPz);
544 
545  fParticleChangeForGamma->ProposeMomentumDirection(direction.unit());
546  }
547 
548  else fParticleChangeForGamma->ProposeMomentumDirection(primaryDirection);
549 
550  // note that secondaryKinetic is the energy of the delta ray, not of all secondaries.
551  G4double scatteredEnergy = k-bindingEnergy-secondaryKinetic;
552  G4double deexSecEnergy = 0;
553  for (G4int j=secNumberInit; j < secNumberFinal; j++)
554  {
555  deexSecEnergy = deexSecEnergy + (*fvect)[j]->GetKineticEnergy();
556  }
557 
558  if (!statCode)
559  {
560  fParticleChangeForGamma->SetProposedKineticEnergy(scatteredEnergy);
561  fParticleChangeForGamma->ProposeLocalEnergyDeposit(k-scatteredEnergy-secondaryKinetic-deexSecEnergy);
562  }
563  else
564  {
565  fParticleChangeForGamma->SetProposedKineticEnergy(k);
566  fParticleChangeForGamma->ProposeLocalEnergyDeposit(k-scatteredEnergy);
567  }
568 
569  // SI - 01/04/2014
570  if (secondaryKinetic>0)
571  {
572  G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic);
573  fvect->push_back(dp);
574  }
575  //
576 
577  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
579  ionizationShell,
580  theIncomingTrack);
581  }
582 }
void set(double x, double y, double z)
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
std::map< G4String, G4double, std::less< G4String > > lowEnergyLimit
const G4Material * GetMaterial() const
G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell)
G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell)
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:617
G4VAtomDeexcitation * fAtomDeexcitation
int G4int
Definition: G4Types.hh:78
virtual G4ThreeVector & SampleDirectionForShell(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
G4ParticleChangeForGamma * fParticleChangeForGamma
G4double GetKineticEnergy() const
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
Float_t Z
Hep3Vector unit() const
std::map< G4String, G4double, std::less< G4String > > highEnergyLimit
float electron_mass_c2
Definition: hepunit.py:274
G4DNAWaterIonisationStructure waterStructure
double x() const
static G4DNAChemistryManager * Instance()
static const double eV
Definition: G4SIunits.hh:212
double y() const
const G4ThreeVector & GetMomentumDirection() const
void GenerateParticles(std::vector< G4DynamicParticle *> *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
double z() const
G4int RandomSelect(G4double energy, const G4String &particle)
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4AtomicShellEnumerator
Here is the call graph for this function:

◆ SelectFasterComputation()

void G4DNABornIonisationModel1::SelectFasterComputation ( G4bool  input)
inline

Definition at line 175 of file G4DNABornIonisationModel1.hh.

176 {
177  fasterCode = input;
178 }

◆ SelectSPScaling()

void G4DNABornIonisationModel1::SelectSPScaling ( G4bool  input)
inline

Definition at line 189 of file G4DNABornIonisationModel1.hh.

190 {
191  spScaling = input;
192 }

◆ SelectStationary()

void G4DNABornIonisationModel1::SelectStationary ( G4bool  input)
inline

Definition at line 182 of file G4DNABornIonisationModel1.hh.

183 {
184  statCode = input;
185 }

◆ TransferedEnergy()

G4double G4DNABornIonisationModel1::TransferedEnergy ( G4ParticleDefinition aParticleDefinition,
G4double  incomingParticleEnergy,
G4int  shell,
G4double  random 
)

Definition at line 1035 of file G4DNABornIonisationModel1.cc.

1039 {
1040  G4double nrj = 0.;
1041 
1042  G4double valueK1 = 0;
1043  G4double valueK2 = 0;
1044  G4double valuePROB21 = 0;
1045  G4double valuePROB22 = 0;
1046  G4double valuePROB12 = 0;
1047  G4double valuePROB11 = 0;
1048 
1049  G4double nrjTransf11 = 0;
1050  G4double nrjTransf12 = 0;
1051  G4double nrjTransf21 = 0;
1052  G4double nrjTransf22 = 0;
1053 
1054  if (particleDefinition == G4Electron::ElectronDefinition())
1055  {
1056  // k should be in eV
1057  std::vector<double>::iterator k2 = std::upper_bound(eTdummyVec.begin(),
1058  eTdummyVec.end(),
1059  k);
1060  std::vector<double>::iterator k1 = k2 - 1;
1061 
1062  /*
1063  G4cout << "----> k=" << k
1064  << " " << *k1
1065  << " " << *k2
1066  << " " << random
1067  << " " << ionizationLevelIndex
1068  << " " << eProbaShellMap[ionizationLevelIndex][(*k1)].back()
1069  << " " << eProbaShellMap[ionizationLevelIndex][(*k2)].back()
1070  << G4endl;
1071  */
1072 
1073  // SI : the following condition avoids situations where random >last vector element
1074  if (random <= eProbaShellMap[ionizationLevelIndex][(*k1)].back()
1075  && random <= eProbaShellMap[ionizationLevelIndex][(*k2)].back())
1076  {
1077  std::vector<double>::iterator prob12 =
1078  std::upper_bound(eProbaShellMap[ionizationLevelIndex][(*k1)].begin(),
1079  eProbaShellMap[ionizationLevelIndex][(*k1)].end(),
1080  random);
1081 
1082  std::vector<double>::iterator prob11 = prob12 - 1;
1083 
1084  std::vector<double>::iterator prob22 =
1085  std::upper_bound(eProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
1086  eProbaShellMap[ionizationLevelIndex][(*k2)].end(),
1087  random);
1088 
1089  std::vector<double>::iterator prob21 = prob22 - 1;
1090 
1091  valueK1 = *k1;
1092  valueK2 = *k2;
1093  valuePROB21 = *prob21;
1094  valuePROB22 = *prob22;
1095  valuePROB12 = *prob12;
1096  valuePROB11 = *prob11;
1097 
1098  /*
1099  G4cout << " " << random << " " << valuePROB11 << " "
1100  << valuePROB12 << " " << valuePROB21 << " " << valuePROB22 << G4endl;
1101  */
1102 
1103  nrjTransf11 = eNrjTransfData[ionizationLevelIndex][valueK1][valuePROB11];
1104  nrjTransf12 = eNrjTransfData[ionizationLevelIndex][valueK1][valuePROB12];
1105  nrjTransf21 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
1106  nrjTransf22 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
1107 
1108  /*
1109  G4cout << " " << ionizationLevelIndex << " "
1110  << random << " " <<valueK1 << " " << valueK2 << G4endl;
1111 
1112  G4cout << " " << random << " " << nrjTransf11 << " "
1113  << nrjTransf12 << " " << nrjTransf21 << " " <<nrjTransf22 << G4endl;
1114  */
1115 
1116  }
1117  // Avoids cases where cum xs is zero for k1 and is not for k2 (with always k1<k2)
1118  if (random > eProbaShellMap[ionizationLevelIndex][(*k1)].back())
1119  {
1120  std::vector<double>::iterator prob22 =
1121  std::upper_bound(eProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
1122  eProbaShellMap[ionizationLevelIndex][(*k2)].end(),
1123  random);
1124 
1125  std::vector<double>::iterator prob21 = prob22 - 1;
1126 
1127  valueK1 = *k1;
1128  valueK2 = *k2;
1129  valuePROB21 = *prob21;
1130  valuePROB22 = *prob22;
1131 
1132  //G4cout << " " << random << " " << valuePROB21 << " " << valuePROB22 << G4endl;
1133 
1134  nrjTransf21 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
1135  nrjTransf22 = eNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
1136 
1137  G4double interpolatedvalue2 = Interpolate(valuePROB21,
1138  valuePROB22,
1139  random,
1140  nrjTransf21,
1141  nrjTransf22);
1142 
1143  // zeros are explicitely set
1144 
1145  G4double value = Interpolate(valueK1, valueK2, k, 0., interpolatedvalue2);
1146 
1147  /*
1148  G4cout << " " << ionizationLevelIndex << " "
1149  << random << " " <<valueK1 << " " << valueK2 << G4endl;
1150 
1151  G4cout << " " << random << " " << nrjTransf11 << " "
1152  << nrjTransf12 << " " << nrjTransf21 << " " <<nrjTransf22 << G4endl;
1153 
1154  G4cout << "ici" << " " << value << G4endl;
1155  */
1156 
1157  return value;
1158  }
1159  }
1160  //
1161  else if (particleDefinition == G4Proton::ProtonDefinition())
1162  {
1163  // k should be in eV
1164 
1165  std::vector<double>::iterator k2 = std::upper_bound(pTdummyVec.begin(),
1166  pTdummyVec.end(),
1167  k);
1168 
1169  std::vector<double>::iterator k1 = k2 - 1;
1170 
1171  /*
1172  G4cout << "----> k=" << k
1173  << " " << *k1
1174  << " " << *k2
1175  << " " << random
1176  << " " << ionizationLevelIndex
1177  << " " << pProbaShellMap[ionizationLevelIndex][(*k1)].back()
1178  << " " << pProbaShellMap[ionizationLevelIndex][(*k2)].back()
1179  << G4endl;
1180  */
1181 
1182  // SI : the following condition avoids situations where random > last vector element,
1183  // for eg. when the last element is zero
1184  if (random <= pProbaShellMap[ionizationLevelIndex][(*k1)].back()
1185  && random <= pProbaShellMap[ionizationLevelIndex][(*k2)].back())
1186  {
1187  std::vector<double>::iterator prob12 =
1188  std::upper_bound(pProbaShellMap[ionizationLevelIndex][(*k1)].begin(),
1189  pProbaShellMap[ionizationLevelIndex][(*k1)].end(),
1190  random);
1191 
1192  std::vector<double>::iterator prob11 = prob12 - 1;
1193 
1194  std::vector<double>::iterator prob22 =
1195  std::upper_bound(pProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
1196  pProbaShellMap[ionizationLevelIndex][(*k2)].end(),
1197  random);
1198 
1199  std::vector<double>::iterator prob21 = prob22 - 1;
1200 
1201  valueK1 = *k1;
1202  valueK2 = *k2;
1203  valuePROB21 = *prob21;
1204  valuePROB22 = *prob22;
1205  valuePROB12 = *prob12;
1206  valuePROB11 = *prob11;
1207 
1208  /*
1209  G4cout << " " << random << " " << valuePROB11 << " "
1210  << valuePROB12 << " " << valuePROB21 << " " << valuePROB22 << G4endl;
1211  */
1212 
1213  nrjTransf11 = pNrjTransfData[ionizationLevelIndex][valueK1][valuePROB11];
1214  nrjTransf12 = pNrjTransfData[ionizationLevelIndex][valueK1][valuePROB12];
1215  nrjTransf21 = pNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
1216  nrjTransf22 = pNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
1217 
1218  /*
1219  G4cout << " " << ionizationLevelIndex << " "
1220  << random << " " <<valueK1 << " " << valueK2 << G4endl;
1221 
1222  G4cout << " " << random << " " << nrjTransf11 << " "
1223  << nrjTransf12 << " " << nrjTransf21 << " " <<nrjTransf22 << G4endl;
1224  */
1225  }
1226 
1227  // Avoids cases where cum xs is zero for k1 and is not for k2 (with always k1<k2)
1228 
1229  if (random > pProbaShellMap[ionizationLevelIndex][(*k1)].back())
1230  {
1231  std::vector<double>::iterator prob22 =
1232  std::upper_bound(pProbaShellMap[ionizationLevelIndex][(*k2)].begin(),
1233  pProbaShellMap[ionizationLevelIndex][(*k2)].end(),
1234  random);
1235 
1236  std::vector<double>::iterator prob21 = prob22 - 1;
1237 
1238  valueK1 = *k1;
1239  valueK2 = *k2;
1240  valuePROB21 = *prob21;
1241  valuePROB22 = *prob22;
1242 
1243  //G4cout << " " << random << " " << valuePROB21 << " " << valuePROB22 << G4endl;
1244 
1245  nrjTransf21 = pNrjTransfData[ionizationLevelIndex][valueK2][valuePROB21];
1246  nrjTransf22 = pNrjTransfData[ionizationLevelIndex][valueK2][valuePROB22];
1247 
1248  G4double interpolatedvalue2 = Interpolate(valuePROB21,
1249  valuePROB22,
1250  random,
1251  nrjTransf21,
1252  nrjTransf22);
1253 
1254  // zeros are explicitely set
1255 
1256  G4double value = Interpolate(valueK1, valueK2, k, 0., interpolatedvalue2);
1257 
1258  /*
1259  G4cout << " " << ionizationLevelIndex << " "
1260  << random << " " <<valueK1 << " " << valueK2 << G4endl;
1261 
1262  G4cout << " " << random << " " << nrjTransf11 << " "
1263  << nrjTransf12 << " " << nrjTransf21 << " " <<nrjTransf22 << G4endl;
1264 
1265  G4cout << "ici" << " " << value << G4endl;
1266  */
1267 
1268  return value;
1269  }
1270  }
1271  // End electron and proton cases
1272 
1273  G4double nrjTransfProduct = nrjTransf11 * nrjTransf12 * nrjTransf21
1274  * nrjTransf22;
1275  //G4cout << "nrjTransfProduct=" << nrjTransfProduct << G4endl;
1276 
1277  if (nrjTransfProduct != 0.)
1278  {
1279  nrj = QuadInterpolator(valuePROB11,
1280  valuePROB12,
1281  valuePROB21,
1282  valuePROB22,
1283  nrjTransf11,
1284  nrjTransf12,
1285  nrjTransf21,
1286  nrjTransf22,
1287  valueK1,
1288  valueK2,
1289  k,
1290  random);
1291  }
1292  //G4cout << nrj << endl;
1293 
1294  return nrj;
1295 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
G4double QuadInterpolator(G4double e11, G4double e12, G4double e21, G4double e22, G4double x11, G4double x12, G4double x21, G4double x22, G4double t1, G4double t2, G4double t, G4double e)
G4double Interpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ eDiffCrossSectionData

TriDimensionMap G4DNABornIonisationModel1::eDiffCrossSectionData[6]
private

Definition at line 147 of file G4DNABornIonisationModel1.hh.

◆ eNrjTransfData

TriDimensionMap G4DNABornIonisationModel1::eNrjTransfData[6]
private

Definition at line 148 of file G4DNABornIonisationModel1.hh.

◆ eProbaShellMap

VecMap G4DNABornIonisationModel1::eProbaShellMap[6]
private

Definition at line 161 of file G4DNABornIonisationModel1.hh.

◆ eTdummyVec

std::vector<double> G4DNABornIonisationModel1::eTdummyVec
private

Definition at line 153 of file G4DNABornIonisationModel1.hh.

◆ eVecm

VecMap G4DNABornIonisationModel1::eVecm
private

Definition at line 158 of file G4DNABornIonisationModel1.hh.

◆ fasterCode

G4bool G4DNABornIonisationModel1::fasterCode
private

Definition at line 94 of file G4DNABornIonisationModel1.hh.

◆ fAtomDeexcitation

G4VAtomDeexcitation* G4DNABornIonisationModel1::fAtomDeexcitation
private

Definition at line 102 of file G4DNABornIonisationModel1.hh.

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNABornIonisationModel1::fParticleChangeForGamma
protected

Definition at line 90 of file G4DNABornIonisationModel1.hh.

◆ fpMolWaterDensity

const std::vector<G4double>* G4DNABornIonisationModel1::fpMolWaterDensity
private

Definition at line 99 of file G4DNABornIonisationModel1.hh.

◆ highEnergyLimit

std::map<G4String,G4double,std::less<G4String> > G4DNABornIonisationModel1::highEnergyLimit
private

Definition at line 105 of file G4DNABornIonisationModel1.hh.

◆ isInitialised

G4bool G4DNABornIonisationModel1::isInitialised
private

Definition at line 111 of file G4DNABornIonisationModel1.hh.

◆ lowEnergyLimit

std::map<G4String,G4double,std::less<G4String> > G4DNABornIonisationModel1::lowEnergyLimit
private

Definition at line 104 of file G4DNABornIonisationModel1.hh.

◆ pDiffCrossSectionData

TriDimensionMap G4DNABornIonisationModel1::pDiffCrossSectionData[6]
private

Definition at line 150 of file G4DNABornIonisationModel1.hh.

◆ pNrjTransfData

TriDimensionMap G4DNABornIonisationModel1::pNrjTransfData[6]
private

Definition at line 151 of file G4DNABornIonisationModel1.hh.

◆ pProbaShellMap

VecMap G4DNABornIonisationModel1::pProbaShellMap[6]
private

Definition at line 162 of file G4DNABornIonisationModel1.hh.

◆ pTdummyVec

std::vector<double> G4DNABornIonisationModel1::pTdummyVec
private

Definition at line 154 of file G4DNABornIonisationModel1.hh.

◆ pVecm

VecMap G4DNABornIonisationModel1::pVecm
private

Definition at line 159 of file G4DNABornIonisationModel1.hh.

◆ spScaling

G4bool G4DNABornIonisationModel1::spScaling
private

Definition at line 96 of file G4DNABornIonisationModel1.hh.

◆ statCode

G4bool G4DNABornIonisationModel1::statCode
private

Definition at line 95 of file G4DNABornIonisationModel1.hh.

◆ tableData

MapData G4DNABornIonisationModel1::tableData
private

Definition at line 120 of file G4DNABornIonisationModel1.hh.

◆ tableFile

MapFile G4DNABornIonisationModel1::tableFile
private

Definition at line 117 of file G4DNABornIonisationModel1.hh.

◆ verboseLevel

G4int G4DNABornIonisationModel1::verboseLevel
private

Definition at line 112 of file G4DNABornIonisationModel1.hh.

◆ waterStructure

G4DNAWaterIonisationStructure G4DNABornIonisationModel1::waterStructure
private

Definition at line 124 of file G4DNABornIonisationModel1.hh.


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