Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 level, const G4ParticleDefinition *, G4double kineticEnergy)
 
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 * > *)
 
virtual 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=nullptr)
 
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

G4ParticleChangeForGammafParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
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.

Constructor & Destructor Documentation

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

Definition at line 44 of file G4DNAChampionElasticModel.cc.

44  :
45  G4VEmModel(nam), isInitialised(false)
46 {
47  SetLowEnergyLimit(7.4 * eV);
48  SetHighEnergyLimit(1. * MeV);
49 
50  verboseLevel = 0;
51  // Verbosity scale:
52  // 0 = nothing
53  // 1 = warning for energy non-conservation
54  // 2 = details of energy budget
55  // 3 = calculation of cross sections, file openings, sampling of atoms
56  // 4 = entering in methods
57 
58 #ifdef CHAMPION_VERBOSE
59  if (verboseLevel > 0)
60  {
61  G4cout << "Champion Elastic model is constructed "
62  << G4endl
63  << "Energy range: "
64  << LowEnergyLimit() / eV << " eV - "
65  << HighEnergyLimit() / MeV << " MeV"
66  << G4endl;
67  }
68 #endif
69 
71  fpMolWaterDensity = 0;
72  fpData = 0;
73 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
G4GLOB_DLL std::ostream G4cout
static constexpr double eV
Definition: G4SIunits.hh:215
G4ParticleChangeForGamma * fParticleChangeForGamma
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731

Here is the call graph for this function:

G4DNAChampionElasticModel::~G4DNAChampionElasticModel ( )
virtual

Definition at line 77 of file G4DNAChampionElasticModel.cc.

78 {
79  // For total cross section
80  if(fpData) delete fpData;
81 
82  // For final state
83  eVecm.clear();
84 }

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 230 of file G4DNAChampionElasticModel.cc.

239 {
240 #ifdef CHAMPION_VERBOSE
241  if (verboseLevel > 3)
242  {
243  G4cout << "Calling CrossSectionPerVolume() of G4DNAChampionElasticModel"
244  << G4endl;
245  }
246 #endif
247 
248  // Calculate total cross section for model
249 
250  G4double sigma = 0.;
251  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
252 
253  if(waterDensity!= 0.0)
254  {
255  if (ekin < HighEnergyLimit() && ekin >= LowEnergyLimit())
256  {
257  //SI : XS must not be zero otherwise sampling of secondaries method ignored
258  //
259  sigma = fpData->FindValue(ekin);
260  }
261 
262 #ifdef CHAMPION_VERBOSE
263  if (verboseLevel > 2)
264  {
265  G4cout << "__________________________________" << G4endl;
266  G4cout << "=== G4DNAChampionElasticModel - XS INFO START" << G4endl;
267  G4cout << "=== Kinetic energy(eV)=" << ekin/eV << " particle : " << p->GetParticleName() << G4endl;
268  G4cout << "=== Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
269  G4cout << "=== Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
270  G4cout << "=== G4DNAChampionElasticModel - XS INFO END" << G4endl;
271  }
272 #endif
273  }
274 
275  return sigma*waterDensity;
276 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
size_t GetIndex() const
Definition: G4Material.hh:262
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double eV
Definition: G4SIunits.hh:215
virtual G4double FindValue(G4double e, G4int componentId=0) const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double G4DNAChampionElasticModel::GetKillBelowThreshold ( )
inline

Definition at line 68 of file G4DNAChampionElasticModel.hh.

69  {
71  errMsg << "The method G4DNAChampionElasticModel::"
72  "GetKillBelowThreshold is deprecated";
73 
74  G4Exception("G4DNAChampionElasticModel::GetKillBelowThreshold",
75  "deprecated",
77  errMsg);
78  return 0.;
79  }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 88 of file G4DNAChampionElasticModel.cc.

90 {
91 #ifdef CHAMPION_VERBOSE
92  if (verboseLevel > 3)
93  {
94  G4cout << "Calling G4DNAChampionElasticModel::Initialise()" << G4endl;
95  }
96 #endif
97 
98  if(particle->GetParticleName() != "e-")
99  {
100  G4Exception("G4DNAChampionElasticModel::Initialise",
101  "em0002",
103  "Model not applicable to particle type.");
104  }
105 
106  // Energy limits
107 
108  if (LowEnergyLimit() < 7.4*eV)
109  {
110  G4cout << "G4DNAChampionElasticModel: low energy limit increased from "
111  << LowEnergyLimit()/eV << " eV to " << 7.4 << " eV"
112  << G4endl;
113  SetLowEnergyLimit(7.4*eV);
114  }
115 
116  if (HighEnergyLimit() > 1.*MeV)
117  {
118  G4cout << "G4DNAChampionElasticModel: high energy limit decreased from "
119  << HighEnergyLimit()/MeV << " MeV to " << 1. << " MeV"
120  << G4endl;
122  }
123 
124  if (isInitialised) { return; }
125 
126  // *** ELECTRON
127  // For total cross section
128  // Reading of data files
129 
130  G4double scaleFactor = 1e-16*cm*cm;
131 
132  G4String fileElectron("dna/sigma_elastic_e_champion");
133 
135  eV,
136  scaleFactor );
137  fpData->LoadData(fileElectron);
138 
139  // For final state
140 
141  char *path = getenv("G4LEDATA");
142 
143  if (!path)
144  {
145  G4Exception("G4ChampionElasticModel::Initialise",
146  "em0006",
148  "G4LEDATA environment variable not set.");
149  return;
150  }
151 
152  std::ostringstream eFullFileName;
153  eFullFileName << path << "/dna/sigmadiff_cumulated_elastic_e_champion_hp.dat";
154  std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
155 
156  if (!eDiffCrossSection)
157  {
158  G4ExceptionDescription errMsg;
159  errMsg << "Missing data file:/dna/sigmadiff_cumulated_elastic_e_champion_hp.dat; "
160  << "please use G4EMLOW6.36 and above.";
161 
162  G4Exception("G4DNAChampionElasticModel::Initialise",
163  "em0003",
165  errMsg);
166  }
167 
168  // March 25th, 2014 - Vaclav Stepan, Sebastien Incerti
169  // Added clear for MT
170 
171  eTdummyVec.clear();
172  eVecm.clear();
173  eDiffCrossSectionData.clear();
174 
175  //
176 
177  eTdummyVec.push_back(0.);
178 
179  while(!eDiffCrossSection.eof())
180  {
181  double tDummy;
182  double eDummy;
183  eDiffCrossSection >> tDummy >> eDummy;
184 
185  // SI : mandatory eVecm initialization
186 
187  if (tDummy != eTdummyVec.back())
188  {
189  eTdummyVec.push_back(tDummy);
190  eVecm[tDummy].push_back(0.);
191  }
192 
193  eDiffCrossSection >> eDiffCrossSectionData[tDummy][eDummy];
194 
195  if (eDummy != eVecm[tDummy].back()) eVecm[tDummy].push_back(eDummy);
196  }
197 
198  // End final state
199 #ifdef CHAMPION_VERBOSE
200  if (verboseLevel>0)
201  {
202  if (verboseLevel > 2)
203  {
204  G4cout << "Loaded cross section files for Champion Elastic model" << G4endl;
205  }
206 
207  G4cout << "Champion Elastic model is initialized " << G4endl
208  << "Energy range: "
209  << LowEnergyLimit() / eV << " eV - "
210  << HighEnergyLimit() / MeV << " MeV"
211  << G4endl;
212  }
213 #endif
214 
215  // Initialize water density pointer
217 
218  fpMolWaterDensity = G4DNAMolecularMaterial::Instance()->
219  GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
220 
222  isInitialised = true;
223 
224 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:602
virtual G4bool LoadData(const G4String &argFileName)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
G4GLOB_DLL std::ostream G4cout
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double eV
Definition: G4SIunits.hh:215
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4DNAMolecularMaterial * Instance()
G4ParticleChangeForGamma * fParticleChangeForGamma
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 280 of file G4DNAChampionElasticModel.cc.

285 {
286 
287 #ifdef CHAMPION_VERBOSE
288  if (verboseLevel > 3)
289  {
290  G4cout << "Calling SampleSecondaries() of G4DNAChampionElasticModel" << G4endl;
291  }
292 #endif
293 
294  G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
295 
296 // if (electronEnergy0 < HighEnergyLimit()) // necessaire ?
297  {
298  G4double cosTheta = RandomizeCosTheta(electronEnergy0);
299 
300  G4double phi = 2. * pi * G4UniformRand();
301 
302  G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
303  G4ThreeVector xVers = zVers.orthogonal();
304  G4ThreeVector yVers = zVers.cross(xVers);
305 
306  G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
307  G4double yDir = xDir;
308  xDir *= std::cos(phi);
309  yDir *= std::sin(phi);
310 
311  G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
312 
314 
316  // necessaire ?
317  }
318 }
G4double GetKineticEnergy() const
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector orthogonal() const
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
static constexpr double pi
Definition: G4SIunits.hh:75
Hep3Vector cross(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

void G4DNAChampionElasticModel::SetKillBelowThreshold ( G4double  threshold)

Definition at line 477 of file G4DNAChampionElasticModel.cc.

478 {
479  G4ExceptionDescription errMsg;
480  errMsg << "The method G4DNAChampionElasticModel::SetKillBelowThreshold is deprecated";
481 
482  G4Exception("G4DNAChampionElasticModel::SetKillBelowThreshold",
483  "deprecated",
484  JustWarning,
485  errMsg);
486 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Here is the call graph for this function:

Member Data Documentation

G4ParticleChangeForGamma* G4DNAChampionElasticModel::fParticleChangeForGamma
protected

Definition at line 96 of file G4DNAChampionElasticModel.hh.


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