Geant4  10.02.p03
G4DNAChampionElasticModel Class Reference

#include <G4DNAChampionElasticModel.hh>

Inheritance diagram for G4DNAChampionElasticModel:
Collaboration diagram for G4DNAChampionElasticModel:

Public Member Functions

 G4DNAChampionElasticModel (const G4ParticleDefinition *p=0, const G4String &nam="DNAChampionElasticModel")
 
virtual ~G4DNAChampionElasticModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const 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)
 
void SetKillBelowThreshold (G4double threshold)
 
G4double GetKillBelowThreshold ()
 
- 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 GetPartialCrossSection (const G4Material *, G4int, const G4ParticleDefinition *, G4double)
 
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 Theta (G4ParticleDefinition *aParticleDefinition, G4double k, G4double integrDiff)
 
G4double LinLinInterpolate (G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
 
G4double LinLogInterpolate (G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
 
G4double LogLogInterpolate (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)
 
G4double RandomizeCosTheta (G4double k)
 
G4DNAChampionElasticModeloperator= (const G4DNAChampionElasticModel &right)
 
 G4DNAChampionElasticModel (const G4DNAChampionElasticModel &)
 

Private Attributes

const std::vector< G4double > * fpMolWaterDensity
 
G4double killBelowEnergy
 
G4double lowEnergyLimit
 
G4double highEnergyLimit
 
G4bool isInitialised
 
G4int verboseLevel
 
MapFile tableFile
 
MapData tableData
 
TriDimensionMap eDiffCrossSectionData
 
std::vector< double > eTdummyVec
 
VecMap eVecm
 

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 41 of file G4DNAChampionElasticModel.hh.

Member Typedef Documentation

◆ MapData

Definition at line 91 of file G4DNAChampionElasticModel.hh.

◆ MapFile

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

Definition at line 88 of file G4DNAChampionElasticModel.hh.

◆ TriDimensionMap

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

Definition at line 133 of file G4DNAChampionElasticModel.hh.

◆ VecMap

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

Definition at line 138 of file G4DNAChampionElasticModel.hh.

Constructor & Destructor Documentation

◆ G4DNAChampionElasticModel() [1/2]

G4DNAChampionElasticModel::G4DNAChampionElasticModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNAChampionElasticModel" 
)

Definition at line 40 of file G4DNAChampionElasticModel.cc.

41  :
42  G4VEmModel(nam), isInitialised(false)
43 {
44 // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
45 
46  killBelowEnergy = 7.4 * eV;
47  lowEnergyLimit = 0 * eV;
48  highEnergyLimit = 1. * MeV;
51 
52  verboseLevel = 0;
53  // Verbosity scale:
54  // 0 = nothing
55  // 1 = warning for energy non-conservation
56  // 2 = details of energy budget
57  // 3 = calculation of cross sections, file openings, sampling of atoms
58  // 4 = entering in methods
59 
60  if (verboseLevel > 0)
61  {
62  G4cout << "Champion Elastic model is constructed " << G4endl<< "Energy range: "
63  << lowEnergyLimit / eV << " eV - "
64  << highEnergyLimit / MeV << " MeV"
65  << G4endl;
66  }
69 }
static const double MeV
Definition: G4SIunits.hh:211
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
G4GLOB_DLL std::ostream G4cout
static const double eV
Definition: G4SIunits.hh:212
G4ParticleChangeForGamma * fParticleChangeForGamma
#define G4endl
Definition: G4ios.hh:61
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
const std::vector< G4double > * fpMolWaterDensity
Here is the call graph for this function:

◆ ~G4DNAChampionElasticModel()

G4DNAChampionElasticModel::~G4DNAChampionElasticModel ( )
virtual

Definition at line 73 of file G4DNAChampionElasticModel.cc.

74 {
75  // For total cross section
76 
77  std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> >::iterator pos;
78  for (pos = tableData.begin(); pos != tableData.end(); ++pos)
79  {
80  G4DNACrossSectionDataSet* table = pos->second;
81  delete table;
82  }
83 
84  // For final state
85 
86  eVecm.clear();
87 
88 }
static const G4double pos

◆ G4DNAChampionElasticModel() [2/2]

G4DNAChampionElasticModel::G4DNAChampionElasticModel ( const G4DNAChampionElasticModel )
private

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 215 of file G4DNAChampionElasticModel.cc.

220 {
221  if (verboseLevel > 3)
222  G4cout << "Calling CrossSectionPerVolume() of G4DNAChampionElasticModel"
223  << G4endl;
224 
225  // Calculate total cross section for model
226 
227  G4double sigma=0;
228 
229  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
230 
231  if(waterDensity!= 0.0)
232 // if (material == nistwater || material->GetBaseMaterial() == nistwater)
233  {
234  const G4String& particleName = p->GetParticleName();
235 
236  if (ekin < highEnergyLimit)
237  {
238  //SI : XS must not be zero otherwise sampling of secondaries method ignored
239  if (ekin < killBelowEnergy) return DBL_MAX;
240  //
241 
242  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
243  pos = tableData.find(particleName);
244 
245  if (pos != tableData.end())
246  {
247  G4DNACrossSectionDataSet* table = pos->second;
248  if (table != 0)
249  {
250  sigma = table->FindValue(ekin);
251  }
252  }
253  else
254  {
255  G4Exception("G4DNAChampionElasticModel::ComputeCrossSectionPerVolume","em0002",
256  FatalException,"Model not applicable to particle type.");
257  }
258  }
259 
260  if (verboseLevel > 2)
261  {
262  G4cout << "__________________________________" << G4endl;
263  G4cout << "=== G4DNAChampionElasticModel - XS INFO START" << G4endl;
264  G4cout << "=== Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
265  G4cout << "=== Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
266  G4cout << "=== Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
267  // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
268  G4cout << "=== G4DNAChampionElasticModel - XS INFO END" << G4endl;
269  }
270 
271  }
272 
273  return sigma*waterDensity;
274 // return sigma*material->GetAtomicNumDensityVector()[1];
275 }
static const double cm
Definition: G4SIunits.hh:118
size_t GetIndex() const
Definition: G4Material.hh:262
const G4String & GetParticleName() const
virtual G4double FindValue(G4double e, G4int componentId=0) const
G4GLOB_DLL std::ostream G4cout
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
#define DBL_MAX
Definition: templates.hh:83
static const G4double pos
Here is the call graph for this function:

◆ GetKillBelowThreshold()

G4double G4DNAChampionElasticModel::GetKillBelowThreshold ( )
inline

Definition at line 67 of file G4DNAChampionElasticModel.hh.

68  {
69  return killBelowEnergy;
70  }

◆ Initialise()

void G4DNAChampionElasticModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 92 of file G4DNAChampionElasticModel.cc.

94 {
95 
96  if (verboseLevel > 3)
97  G4cout << "Calling G4DNAChampionElasticModel::Initialise()" << G4endl;
98 
99  // Energy limits
100 
102  {
103  G4cout << "G4DNAChampionElasticModel: low energy limit increased from " <<
104  LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
106  }
107 
109  {
110  G4cout << "G4DNAChampionElasticModel: high energy limit decreased from " <<
111  HighEnergyLimit()/MeV << " MeV to " << highEnergyLimit/MeV << " MeV" << G4endl;
113  }
114 
115  // Reading of data files
116 
117  G4double scaleFactor = 1e-16*cm*cm;
118 
119  G4String fileElectron("dna/sigma_elastic_e_champion");
120 
123 
124  // *** ELECTRON
125 
126  // For total cross section
127 
128  electron = electronDef->GetParticleName();
129 
130  tableFile[electron] = fileElectron;
131 
132  G4DNACrossSectionDataSet* tableE =
133  new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
134  tableE->LoadData(fileElectron);
135  tableData[electron] = tableE;
136 
137  // For final state
138 
139  char *path = getenv("G4LEDATA");
140 
141  if (!path)
142  {
143  G4Exception("G4ChampionElasticModel::Initialise","em0006",
144  FatalException,"G4LEDATA environment variable not set.");
145  return;
146  }
147 
148  std::ostringstream eFullFileName;
149  eFullFileName << path << "/dna/sigmadiff_cumulated_elastic_e_champion_hp.dat";
150  std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
151 
152  if (!eDiffCrossSection)
153  G4Exception("G4DNAChampionElasticModel::Initialise","em0003",
155  "Missing data file:/dna/sigmadiff_cumulated_elastic_e_champion_hp.dat; please use G4EMLOW6.36 and above.");
156 
157  // March 25th, 2014 - Vaclav Stepan, Sebastien Incerti
158  // Added clear for MT
159 
160  eTdummyVec.clear();
161  eVecm.clear();
162  eDiffCrossSectionData.clear();
163 
164  //
165 
166  eTdummyVec.push_back(0.);
167 
168  while(!eDiffCrossSection.eof())
169  {
170  double tDummy;
171  double eDummy;
172  eDiffCrossSection>>tDummy>>eDummy;
173 
174  // SI : mandatory eVecm initialization
175 
176  if (tDummy != eTdummyVec.back())
177  {
178  eTdummyVec.push_back(tDummy);
179  eVecm[tDummy].push_back(0.);
180  }
181 
182  eDiffCrossSection>>eDiffCrossSectionData[tDummy][eDummy];
183 
184  if (eDummy != eVecm[tDummy].back()) eVecm[tDummy].push_back(eDummy);
185 
186  }
187 
188  // End final state
189 
190  if (verboseLevel > 2)
191  G4cout << "Loaded cross section files for Champion Elastic model" << G4endl;
192 
193  if( verboseLevel>0 )
194  {
195  G4cout << "Champion Elastic model is initialized " << G4endl
196  << "Energy range: "
197  << LowEnergyLimit() / eV << " eV - "
198  << HighEnergyLimit() / MeV << " MeV"
199  << G4endl;
200  }
201 
202  // Initialize water density pointer
205 
206  if (isInitialised)
207  { return;}
209  isInitialised = true;
210 
211 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
static const double cm
Definition: G4SIunits.hh:118
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static const double MeV
Definition: G4SIunits.hh:211
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
virtual G4bool LoadData(const G4String &argFileName)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
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()
G4ParticleChangeForGamma * fParticleChangeForGamma
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
const std::vector< G4double > * fpMolWaterDensity
Here is the call graph for this function:

◆ LinLinInterpolate()

G4double G4DNAChampionElasticModel::LinLinInterpolate ( G4double  e1,
G4double  e2,
G4double  e,
G4double  xs1,
G4double  xs2 
)
private

Definition at line 395 of file G4DNAChampionElasticModel.cc.

400 {
401  G4double d1 = xs1;
402  G4double d2 = xs2;
403  G4double value = (d1 + (d2 - d1) * (e - e1) / (e2 - e1));
404  return value;
405 }
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:

◆ LinLogInterpolate()

G4double G4DNAChampionElasticModel::LinLogInterpolate ( G4double  e1,
G4double  e2,
G4double  e,
G4double  xs1,
G4double  xs2 
)
private

Definition at line 381 of file G4DNAChampionElasticModel.cc.

386 {
387  G4double d1 = std::log(xs1);
388  G4double d2 = std::log(xs2);
389  G4double value = std::exp(d1 + (d2 - d1) * (e - e1) / (e2 - e1));
390  return value;
391 }
static const G4double d1
static const G4double e2
static const G4double e1
double G4double
Definition: G4Types.hh:76
static const G4double d2

◆ LogLogInterpolate()

G4double G4DNAChampionElasticModel::LogLogInterpolate ( G4double  e1,
G4double  e2,
G4double  e,
G4double  xs1,
G4double  xs2 
)
private

Definition at line 409 of file G4DNAChampionElasticModel.cc.

414 {
415  G4double a = (std::log10(xs2) - std::log10(xs1))
416  / (std::log10(e2) - std::log10(e1));
417  G4double b = std::log10(xs2) - a * std::log10(e2);
418  G4double sigma = a * std::log10(e) + b;
419  G4double value = (std::pow(10., sigma));
420  return value;
421 }
static const G4double e2
static const G4double e1
double G4double
Definition: G4Types.hh:76

◆ operator=()

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

◆ QuadInterpolator()

G4double G4DNAChampionElasticModel::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 425 of file G4DNAChampionElasticModel.cc.

437 {
438  // Log-Log
439  /*
440  G4double interpolatedvalue1 = LogLogInterpolate(e11, e12, e, xs11, xs12);
441  G4double interpolatedvalue2 = LogLogInterpolate(e21, e22, e, xs21, xs22);
442  G4double value = LogLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
443 
444 
445  // Lin-Log
446  G4double interpolatedvalue1 = LinLogInterpolate(e11, e12, e, xs11, xs12);
447  G4double interpolatedvalue2 = LinLogInterpolate(e21, e22, e, xs21, xs22);
448  G4double value = LinLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
449  */
450 
451  // Lin-Lin
452  G4double interpolatedvalue1 = LinLinInterpolate(e11, e12, e, xs11, xs12);
453  G4double interpolatedvalue2 = LinLinInterpolate(e21, e22, e, xs21, xs22);
454  G4double value = LinLinInterpolate(t1, t2, t, interpolatedvalue1,
455  interpolatedvalue2);
456 
457  return value;
458 }
TTree * t1
Definition: plottest35.C:26
TTree * t2
Definition: plottest35.C:36
double G4double
Definition: G4Types.hh:76
G4double LinLinInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ RandomizeCosTheta()

G4double G4DNAChampionElasticModel::RandomizeCosTheta ( G4double  k)
private

Definition at line 462 of file G4DNAChampionElasticModel.cc.

463 {
464 
465  G4double integrdiff = 0;
466  G4double uniformRand = G4UniformRand();
467  integrdiff = uniformRand;
468 
469  G4double theta = 0.;
470  G4double cosTheta = 0.;
471  theta = Theta(G4Electron::ElectronDefinition(), k / eV, integrdiff);
472 
473  cosTheta = std::cos(theta * pi / 180);
474 
475  return cosTheta;
476 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
#define G4UniformRand()
Definition: Randomize.hh:97
G4double Theta(G4ParticleDefinition *aParticleDefinition, G4double k, G4double integrDiff)
static const double eV
Definition: G4SIunits.hh:212
static const double pi
Definition: G4SIunits.hh:74
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SampleSecondaries()

void G4DNAChampionElasticModel::SampleSecondaries ( std::vector< G4DynamicParticle *> *  ,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicElectron,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 279 of file G4DNAChampionElasticModel.cc.

284 {
285 
286  if (verboseLevel > 3)
287  G4cout << "Calling SampleSecondaries() of G4DNAChampionElasticModel" << G4endl;
288 
289  G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
290 
291  if (electronEnergy0 < killBelowEnergy)
292  {
293  fParticleChangeForGamma->SetProposedKineticEnergy(0.);
294  fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
295  fParticleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0);
296  return;
297  }
298 
299  if (electronEnergy0>= killBelowEnergy && electronEnergy0 < highEnergyLimit)
300  {
301 
302  G4double cosTheta = RandomizeCosTheta(electronEnergy0);
303 
304  G4double phi = 2. * pi * G4UniformRand();
305 
306  G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
307  G4ThreeVector xVers = zVers.orthogonal();
308  G4ThreeVector yVers = zVers.cross(xVers);
309 
310  G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
311  G4double yDir = xDir;
312  xDir *= std::cos(phi);
313  yDir *= std::sin(phi);
314 
315  G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
316 
317  fParticleChangeForGamma->ProposeMomentumDirection(zPrimeVers.unit());
318 
319  fParticleChangeForGamma->SetProposedKineticEnergy(electronEnergy0);
320  }
321 
322 }
G4double GetKineticEnergy() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
Hep3Vector cross(const Hep3Vector &) const
Hep3Vector orthogonal() const
static const double pi
Definition: G4SIunits.hh:74
const G4ThreeVector & GetMomentumDirection() const
G4ParticleChangeForGamma * fParticleChangeForGamma
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ SetKillBelowThreshold()

void G4DNAChampionElasticModel::SetKillBelowThreshold ( G4double  threshold)

Definition at line 480 of file G4DNAChampionElasticModel.cc.

481 {
482  killBelowEnergy = threshold;
483 
484  if (threshold < 7.4 * eV)
485  {
486  G4cout<< "*** WARNING : the G4DNAChampionElasticModel class is not "
487  "activated below 7.4 eV !" << G4endl;
488  G4cout<< "*** NOTE: if you are using G4EmDNAChemistry, do not worry, this is"
489  " the expected behavior" << G4endl;
490  }
491 }
G4GLOB_DLL std::ostream G4cout
static const double eV
Definition: G4SIunits.hh:212
#define G4endl
Definition: G4ios.hh:61

◆ Theta()

G4double G4DNAChampionElasticModel::Theta ( G4ParticleDefinition aParticleDefinition,
G4double  k,
G4double  integrDiff 
)
private

Definition at line 326 of file G4DNAChampionElasticModel.cc.

329 {
330  G4double theta = 0.;
331  G4double valueT1 = 0;
332  G4double valueT2 = 0;
333  G4double valueE21 = 0;
334  G4double valueE22 = 0;
335  G4double valueE12 = 0;
336  G4double valueE11 = 0;
337  G4double xs11 = 0;
338  G4double xs12 = 0;
339  G4double xs21 = 0;
340  G4double xs22 = 0;
341 
342  if (particleDefinition == G4Electron::ElectronDefinition())
343  {
344  std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),
345  eTdummyVec.end(), k);
346  std::vector<double>::iterator t1 = t2 - 1;
347 
348  std::vector<double>::iterator e12 = std::upper_bound(eVecm[(*t1)].begin(),
349  eVecm[(*t1)].end(),
350  integrDiff);
351  std::vector<double>::iterator e11 = e12 - 1;
352 
353  std::vector<double>::iterator e22 = std::upper_bound(eVecm[(*t2)].begin(),
354  eVecm[(*t2)].end(),
355  integrDiff);
356  std::vector<double>::iterator e21 = e22 - 1;
357 
358  valueT1 = *t1;
359  valueT2 = *t2;
360  valueE21 = *e21;
361  valueE22 = *e22;
362  valueE12 = *e12;
363  valueE11 = *e11;
364 
365  xs11 = eDiffCrossSectionData[valueT1][valueE11];
366  xs12 = eDiffCrossSectionData[valueT1][valueE12];
367  xs21 = eDiffCrossSectionData[valueT2][valueE21];
368  xs22 = eDiffCrossSectionData[valueT2][valueE22];
369  }
370 
371  if (xs11 == 0 && xs12 == 0 && xs21 == 0 && xs22 == 0) return (0.);
372 
373  theta = QuadInterpolator(valueE11, valueE12, valueE21, valueE22, xs11, xs12,
374  xs21, xs22, valueT1, valueT2, k, integrDiff);
375 
376  return theta;
377 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
TTree * t1
Definition: plottest35.C:26
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)
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:

Member Data Documentation

◆ eDiffCrossSectionData

TriDimensionMap G4DNAChampionElasticModel::eDiffCrossSectionData
private

Definition at line 135 of file G4DNAChampionElasticModel.hh.

◆ eTdummyVec

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

Definition at line 136 of file G4DNAChampionElasticModel.hh.

◆ eVecm

VecMap G4DNAChampionElasticModel::eVecm
private

Definition at line 139 of file G4DNAChampionElasticModel.hh.

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAChampionElasticModel::fParticleChangeForGamma
protected

Definition at line 74 of file G4DNAChampionElasticModel.hh.

◆ fpMolWaterDensity

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

Definition at line 78 of file G4DNAChampionElasticModel.hh.

◆ highEnergyLimit

G4double G4DNAChampionElasticModel::highEnergyLimit
private

Definition at line 82 of file G4DNAChampionElasticModel.hh.

◆ isInitialised

G4bool G4DNAChampionElasticModel::isInitialised
private

Definition at line 83 of file G4DNAChampionElasticModel.hh.

◆ killBelowEnergy

G4double G4DNAChampionElasticModel::killBelowEnergy
private

Definition at line 80 of file G4DNAChampionElasticModel.hh.

◆ lowEnergyLimit

G4double G4DNAChampionElasticModel::lowEnergyLimit
private

Definition at line 81 of file G4DNAChampionElasticModel.hh.

◆ tableData

MapData G4DNAChampionElasticModel::tableData
private

Definition at line 92 of file G4DNAChampionElasticModel.hh.

◆ tableFile

MapFile G4DNAChampionElasticModel::tableFile
private

Definition at line 89 of file G4DNAChampionElasticModel.hh.

◆ verboseLevel

G4int G4DNAChampionElasticModel::verboseLevel
private

Definition at line 84 of file G4DNAChampionElasticModel.hh.


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