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

#include <G4DNAChampionElasticModel.hh>

Inheritance diagram for G4DNAChampionElasticModel:
G4VEmModel

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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *)
 
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
 

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 *)
 

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 40 of file G4DNAChampionElasticModel.cc.

References python.hepunit::eV, fParticleChangeForGamma, G4cout, G4endl, python.hepunit::MeV, G4VEmModel::SetHighEnergyLimit(), and G4VEmModel::SetLowEnergyLimit().

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;
49  SetLowEnergyLimit(lowEnergyLimit);
50  SetHighEnergyLimit(highEnergyLimit);
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
63  << "Energy range: "
64  << lowEnergyLimit / eV << " eV - "
65  << highEnergyLimit / MeV << " MeV"
66  << G4endl;
67  }
69  fpMolWaterDensity = 0;
70 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
#define G4endl
Definition: G4ios.hh:61
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
G4DNAChampionElasticModel::~G4DNAChampionElasticModel ( )
virtual

Definition at line 74 of file G4DNAChampionElasticModel.cc.

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

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 203 of file G4DNAChampionElasticModel.cc.

References python.hepunit::cm, DBL_MAX, python.hepunit::eV, FatalException, G4DNACrossSectionDataSet::FindValue(), G4cout, G4endl, G4Exception(), G4Material::GetIndex(), and G4ParticleDefinition::GetParticleName().

208 {
209  if (verboseLevel > 3)
210  G4cout << "Calling CrossSectionPerVolume() of G4DNAChampionElasticModel" << G4endl;
211 
212  // Calculate total cross section for model
213 
214  G4double sigma=0;
215 
216  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
217 
218  if(waterDensity!= 0.0)
219 // if (material == nistwater || material->GetBaseMaterial() == nistwater)
220  {
221  const G4String& particleName = p->GetParticleName();
222 
223  if (ekin < highEnergyLimit)
224  {
225  //SI : XS must not be zero otherwise sampling of secondaries method ignored
226  if (ekin < killBelowEnergy) return DBL_MAX;
227  //
228 
229  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
230  pos = tableData.find(particleName);
231 
232  if (pos != tableData.end())
233  {
234  G4DNACrossSectionDataSet* table = pos->second;
235  if (table != 0)
236  {
237  sigma = table->FindValue(ekin);
238  }
239  }
240  else
241  {
242  G4Exception("G4DNAChampionElasticModel::ComputeCrossSectionPerVolume","em0002",
243  FatalException,"Model not applicable to particle type.");
244  }
245  }
246 
247  if (verboseLevel > 2)
248  {
249  G4cout << "__________________________________" << G4endl;
250  G4cout << "°°° G4DNAChampionElasticModel - XS INFO START" << G4endl;
251  G4cout << "°°° Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
252  G4cout << "°°° Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
253  G4cout << "°°° Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
254  // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
255  G4cout << "°°° G4DNAChampionElasticModel - XS INFO END" << G4endl;
256  }
257 
258  }
259 
260  return sigma*waterDensity;
261 // return sigma*material->GetAtomicNumDensityVector()[1];
262 }
size_t GetIndex() const
Definition: G4Material.hh:260
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
virtual G4double FindValue(G4double e, G4int componentId=0) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
G4double G4DNAChampionElasticModel::GetKillBelowThreshold ( )
inline

Definition at line 66 of file G4DNAChampionElasticModel.hh.

66 { return killBelowEnergy; }
void G4DNAChampionElasticModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 93 of file G4DNAChampionElasticModel.cc.

References python.hepunit::cm, G4InuclParticleNames::electron, G4Electron::ElectronDefinition(), python.hepunit::eV, FatalException, fParticleChangeForGamma, G4cout, G4endl, G4Exception(), G4Material::GetMaterial(), G4DNAMolecularMaterial::GetNumMolPerVolTableFor(), G4VEmModel::GetParticleChangeForGamma(), G4ParticleDefinition::GetParticleName(), G4VEmModel::HighEnergyLimit(), G4DNAMolecularMaterial::Instance(), G4DNACrossSectionDataSet::LoadData(), G4VEmModel::LowEnergyLimit(), python.hepunit::MeV, G4VEmModel::SetHighEnergyLimit(), and G4VEmModel::SetLowEnergyLimit().

95 {
96 
97  if (verboseLevel > 3)
98  G4cout << "Calling G4DNAChampionElasticModel::Initialise()" << G4endl;
99 
100  // Energy limits
101 
102  if (LowEnergyLimit() < lowEnergyLimit)
103  {
104  G4cout << "G4DNAChampionElasticModel: low energy limit increased from " <<
105  LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
106  SetLowEnergyLimit(lowEnergyLimit);
107  }
108 
109  if (HighEnergyLimit() > highEnergyLimit)
110  {
111  G4cout << "G4DNAChampionElasticModel: high energy limit decreased from " <<
112  HighEnergyLimit()/MeV << " MeV to " << highEnergyLimit/MeV << " MeV" << G4endl;
113  SetHighEnergyLimit(highEnergyLimit);
114  }
115 
116  // Reading of data files
117 
118  G4double scaleFactor = 1e-16*cm*cm;
119 
120  G4String fileElectron("dna/sigma_elastic_e_champion");
121 
124 
125  // *** ELECTRON
126 
127  // For total cross section
128 
129  electron = electronDef->GetParticleName();
130 
131  tableFile[electron] = fileElectron;
132 
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.dat";
150  std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
151 
152  if (!eDiffCrossSection)
153  G4Exception("G4DNAChampionElasticModel::Initialise","em0003",
154  FatalException,"Missing data file:/dna/sigmadiff_cumulated_elastic_e_champion.dat");
155 
156  eTdummyVec.push_back(0.);
157 
158  while(!eDiffCrossSection.eof())
159  {
160  double tDummy;
161  double eDummy;
162  eDiffCrossSection>>tDummy>>eDummy;
163 
164  // SI : mandatory eVecm initialization
165 
166  if (tDummy != eTdummyVec.back())
167  {
168  eTdummyVec.push_back(tDummy);
169  eVecm[tDummy].push_back(0.);
170  }
171 
172  eDiffCrossSection>>eDiffCrossSectionData[tDummy][eDummy];
173 
174  if (eDummy != eVecm[tDummy].back()) eVecm[tDummy].push_back(eDummy);
175 
176  }
177 
178  // End final state
179 
180  if (verboseLevel > 2)
181  G4cout << "Loaded cross section files for Champion Elastic model" << G4endl;
182 
183  if( verboseLevel>0 )
184  {
185  G4cout << "Champion Elastic model is initialized " << G4endl
186  << "Energy range: "
187  << LowEnergyLimit() / eV << " eV - "
188  << HighEnergyLimit() / MeV << " MeV"
189  << G4endl;
190  }
191 
192  // Initialize water density pointer
194 
195  if (isInitialised) { return; }
197  isInitialised = true;
198 
199 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:578
virtual G4bool LoadData(const G4String &argFileName)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
G4GLOB_DLL std::ostream G4cout
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
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
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
void G4DNAChampionElasticModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicElectron,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 266 of file G4DNAChampionElasticModel.cc.

References CLHEP::Hep3Vector::cross(), fParticleChangeForGamma, fStopAndKill, G4cout, G4endl, G4UniformRand, G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), CLHEP::Hep3Vector::orthogonal(), python.hepunit::pi, G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), G4VParticleChange::ProposeTrackStatus(), G4ParticleChangeForGamma::SetProposedKineticEnergy(), and CLHEP::Hep3Vector::unit().

271 {
272 
273  if (verboseLevel > 3)
274  G4cout << "Calling SampleSecondaries() of G4DNAChampionElasticModel" << G4endl;
275 
276  G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
277 
278  if (electronEnergy0 < killBelowEnergy)
279  {
283  return ;
284  }
285 
286  if (electronEnergy0>= killBelowEnergy && electronEnergy0 < highEnergyLimit)
287  {
288 
289  G4double cosTheta = RandomizeCosTheta(electronEnergy0);
290 
291  G4double phi = 2. * pi * G4UniformRand();
292 
293  G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
294  G4ThreeVector xVers = zVers.orthogonal();
295  G4ThreeVector yVers = zVers.cross(xVers);
296 
297  G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
298  G4double yDir = xDir;
299  xDir *= std::cos(phi);
300  yDir *= std::sin(phi);
301 
302  G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
303 
305 
307  }
308 
309 }
G4double GetKineticEnergy() const
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector orthogonal() const
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
Hep3Vector cross(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void G4DNAChampionElasticModel::SetKillBelowThreshold ( G4double  threshold)
inline

Definition at line 133 of file G4DNAChampionElasticModel.hh.

References G4Exception(), and JustWarning.

134 {
135 
136 // SI - commented on 19/06/2013
137 /*
138  killBelowEnergy = threshold;
139 
140  if (threshold < 1*eV)
141  G4Exception ("*** WARNING : the G4DNAChampionElasticModel class is not validated below 1 eV !","",JustWarning,"") ;
142 
143  if (threshold < 0.025*eV) threshold = 0.025*eV;
144 */
145 
146  G4Exception ("*** WARNING : G4DNAChampionElasticModel::SetKillBelowThreshold INACTIVE for now","",JustWarning,"") ;
147 
148 }
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Field Documentation

G4ParticleChangeForGamma* G4DNAChampionElasticModel::fParticleChangeForGamma
protected

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