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

#include <G4DNABornExcitationModel.hh>

Inheritance diagram for G4DNABornExcitationModel:
G4VEmModel

Public Member Functions

 G4DNABornExcitationModel (const G4ParticleDefinition *p=0, const G4String &nam="DNABornExcitationModel")
 
virtual ~G4DNABornExcitationModel ()
 
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)
 
- 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 43 of file G4DNABornExcitationModel.hh.

Constructor & Destructor Documentation

G4DNABornExcitationModel::G4DNABornExcitationModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNABornExcitationModel" 
)

Definition at line 40 of file G4DNABornExcitationModel.cc.

References fParticleChangeForGamma, G4cout, and G4endl.

42  :G4VEmModel(nam),isInitialised(false)
43 {
44  // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
45  fpMolWaterDensity = 0;
46 
47  verboseLevel= 0;
48  // Verbosity scale:
49  // 0 = nothing
50  // 1 = warning for energy non-conservation
51  // 2 = details of energy budget
52  // 3 = calculation of cross sections, file openings, sampling of atoms
53  // 4 = entering in methods
54 
55  if( verboseLevel>0 )
56  {
57  G4cout << "Born excitation model is constructed " << G4endl;
58  }
60 }
G4ParticleChangeForGamma * fParticleChangeForGamma
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
G4DNABornExcitationModel::~G4DNABornExcitationModel ( )
virtual

Definition at line 64 of file G4DNABornExcitationModel.cc.

65 {
66  // Cross section
67 
68  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
69  for (pos = tableData.begin(); pos != tableData.end(); ++pos)
70  {
71  G4DNACrossSectionDataSet* table = pos->second;
72  delete table;
73  }
74 
75 }

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 163 of file G4DNABornExcitationModel.cc.

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

168 {
169  if (verboseLevel > 3)
170  G4cout << "Calling CrossSectionPerVolume() of G4DNABornExcitationModel" << G4endl;
171 
172  if (
173  particleDefinition != G4Proton::ProtonDefinition()
174  &&
175  particleDefinition != G4Electron::ElectronDefinition()
176  )
177 
178  return 0;
179 
180  // Calculate total cross section for model
181 
182  G4double lowLim = 0;
183  G4double highLim = 0;
184  G4double sigma=0;
185 
186  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
187 
188  if(waterDensity!= 0.0)
189  // if (material == nistwater || material->GetBaseMaterial() == nistwater)
190  {
191  const G4String& particleName = particleDefinition->GetParticleName();
192 
193  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
194  pos1 = lowEnergyLimit.find(particleName);
195  if (pos1 != lowEnergyLimit.end())
196  {
197  lowLim = pos1->second;
198  }
199 
200  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
201  pos2 = highEnergyLimit.find(particleName);
202  if (pos2 != highEnergyLimit.end())
203  {
204  highLim = pos2->second;
205  }
206 
207  if (ekin >= lowLim && ekin < highLim)
208  {
209  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
210  pos = tableData.find(particleName);
211 
212  if (pos != tableData.end())
213  {
214  G4DNACrossSectionDataSet* table = pos->second;
215  if (table != 0)
216  {
217  sigma = table->FindValue(ekin);
218  }
219  }
220  else
221  {
222  G4Exception("G4DNABornExcitationModel::CrossSectionPerVolume","em0002",
223  FatalException,"Model not applicable to particle type.");
224  }
225  }
226 
227  if (verboseLevel > 2)
228  {
229  G4cout << "__________________________________" << G4endl;
230  G4cout << "°°° G4DNABornExcitationModel - XS INFO START" << G4endl;
231  G4cout << "°°° Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
232  G4cout << "°°° Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
233  G4cout << "°°° Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
234  // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
235  G4cout << "°°° G4DNABornExcitationModel - XS INFO END" << G4endl;
236  }
237 
238  } // if (waterMaterial)
239 
240  return sigma*waterDensity;
241  // return sigma*material->GetAtomicNumDensityVector()[1];
242 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
size_t GetIndex() const
Definition: G4Material.hh:260
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
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
void G4DNABornExcitationModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector = *(new G4DataVector()) 
)
virtual

Implements G4VEmModel.

Definition at line 79 of file G4DNABornExcitationModel.cc.

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

81 {
82 
83  if (verboseLevel > 3)
84  G4cout << "Calling G4DNABornExcitationModel::Initialise()" << G4endl;
85 
86  G4String fileElectron("dna/sigma_excitation_e_born");
87  G4String fileProton("dna/sigma_excitation_p_born");
88 
91 
94 
95  G4double scaleFactor = (1.e-22 / 3.343) * m*m;
96 
97  // *** ELECTRON
98 
99  electron = electronDef->GetParticleName();
100 
101  tableFile[electron] = fileElectron;
102 
103  lowEnergyLimit[electron] = 9. * eV;
104  highEnergyLimit[electron] = 1. * MeV;
105 
106  // Cross section
107 
109  tableE->LoadData(fileElectron);
110 
111  tableData[electron] = tableE;
112 
113  // *** PROTON
114 
115  proton = protonDef->GetParticleName();
116 
117  tableFile[proton] = fileProton;
118 
119  lowEnergyLimit[proton] = 500. * keV;
120  highEnergyLimit[proton] = 100. * MeV;
121 
122  // Cross section
123 
125  tableP->LoadData(fileProton);
126 
127  tableData[proton] = tableP;
128 
129  //
130 
131  if (particle==electronDef)
132  {
133  SetLowEnergyLimit(lowEnergyLimit[electron]);
134  SetHighEnergyLimit(highEnergyLimit[electron]);
135  }
136 
137  if (particle==protonDef)
138  {
139  SetLowEnergyLimit(lowEnergyLimit[proton]);
140  SetHighEnergyLimit(highEnergyLimit[proton]);
141  }
142 
143  if( verboseLevel>0 )
144  {
145  G4cout << "Born excitation model is initialized " << G4endl
146  << "Energy range: "
147  << LowEnergyLimit() / eV << " eV - "
148  << HighEnergyLimit() / keV << " keV for "
149  << particle->GetParticleName()
150  << G4endl;
151  }
152 
153  // Initialize water density pointer
155 
156  if (isInitialised) { return; }
158  isInitialised = true;
159 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
G4ParticleChangeForGamma * fParticleChangeForGamma
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
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
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
static G4DNAMolecularMaterial * Instance()
#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 G4DNABornExcitationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicParticle,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 246 of file G4DNABornExcitationModel.cc.

References G4DNAChemistryManager::CreateWaterMolecule(), eExcitedMolecule, G4DNAWaterExcitationStructure::ExcitationEnergy(), fParticleChangeForGamma, G4cout, G4endl, G4ParticleChangeForGamma::GetCurrentTrack(), G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4ParticleDefinition::GetParticleName(), G4DNAChemistryManager::Instance(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), and G4ParticleChangeForGamma::SetProposedKineticEnergy().

251 {
252 
253  if (verboseLevel > 3)
254  G4cout << "Calling SampleSecondaries() of G4DNABornExcitationModel" << G4endl;
255 
256  G4double k = aDynamicParticle->GetKineticEnergy();
257 
258  const G4String& particleName = aDynamicParticle->GetDefinition()->GetParticleName();
259 
260  G4int level = RandomSelect(k,particleName);
261  G4double excitationEnergy = waterStructure.ExcitationEnergy(level);
262  G4double newEnergy = k - excitationEnergy;
263 
264  if (newEnergy > 0)
265  {
269  }
270 
271  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
273  level,
274  theIncomingTrack);
275 }
G4ParticleChangeForGamma * fParticleChangeForGamma
G4double GetKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4String & GetParticleName() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Field Documentation

G4ParticleChangeForGamma* G4DNABornExcitationModel::fParticleChangeForGamma
protected

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