Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DNABornExcitationModel1 Class Reference

#include <G4DNABornExcitationModel1.hh>

Inheritance diagram for G4DNABornExcitationModel1:
Collaboration diagram for G4DNABornExcitationModel1:

Public Member Functions

 G4DNABornExcitationModel1 (const G4ParticleDefinition *p=0, const G4String &nam="DNABornExcitationModel")
 
virtual ~G4DNABornExcitationModel1 ()
 
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 G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SelectStationary (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 * > *)
 
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 43 of file G4DNABornExcitationModel1.hh.

Constructor & Destructor Documentation

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

Definition at line 41 of file G4DNABornExcitationModel1.cc.

42  :
43  G4VEmModel(nam), isInitialised(false), fTableData(0)
44 {
45  fpMolWaterDensity = 0;
46  fHighEnergy = 0;
47  fLowEnergy = 0;
48  fParticleDefinition = 0;
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  if (verboseLevel > 0)
59  {
60  G4cout << "Born excitation model is constructed " << G4endl;
61  }
63 
64  // Selection of stationary mode
65 
66  statCode = false;
67 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
G4ParticleChangeForGamma * fParticleChangeForGamma
G4DNABornExcitationModel1::~G4DNABornExcitationModel1 ( )
virtual

Definition at line 71 of file G4DNABornExcitationModel1.cc.

72 {
73  // Cross section
74  if (fTableData)
75  delete fTableData;
76 }

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 138 of file G4DNABornExcitationModel1.cc.

143 {
144  if (verboseLevel > 3)
145  {
146  G4cout << "Calling CrossSectionPerVolume() of G4DNABornExcitationModel1"
147  << G4endl;
148  }
149 
150  if(particleDefinition != fParticleDefinition) return 0;
151 
152  // Calculate total cross section for model
153 
154  G4double sigma=0;
155 
156  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
157 
158  if(waterDensity!= 0.0)
159  {
160  if (ekin >= fLowEnergy && ekin < fHighEnergy)
161  {
162  sigma = fTableData->FindValue(ekin);
163  }
164 
165  if (verboseLevel > 2)
166  {
167  G4cout << "__________________________________" << G4endl;
168  G4cout << "G4DNABornExcitationModel1 - XS INFO START" << G4endl;
169  G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
170  G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
171  G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
172  G4cout << "G4DNABornExcitationModel1 - XS INFO END" << G4endl;
173  }
174  } // if (waterMaterial)
175 
176  return sigma*waterDensity;
177 }
size_t GetIndex() const
Definition: G4Material.hh:262
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 G4DNABornExcitationModel1::GetPartialCrossSection ( const G4Material ,
G4int  level,
const G4ParticleDefinition particle,
G4double  kineticEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 218 of file G4DNABornExcitationModel1.cc.

222 {
223  if (fParticleDefinition != particle)
224  {
225  G4Exception("G4DNABornExcitationModel1::GetPartialCrossSection",
226  "bornParticleType",
228  "Model initialized for another particle type.");
229  }
230 
231  return fTableData->GetComponent(level)->FindValue(kineticEnergy);
232 }
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
virtual const G4VEMDataSet * GetComponent(G4int componentId) const
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 G4DNABornExcitationModel1::Initialise ( const G4ParticleDefinition particle,
const G4DataVector = *(new G4DataVector()) 
)
virtual

Implements G4VEmModel.

Definition at line 80 of file G4DNABornExcitationModel1.cc.

82 {
83 
84  if (verboseLevel > 3)
85  {
86  G4cout << "Calling G4DNABornExcitationModel1::Initialise()" << G4endl;
87  }
88 
89  if(fParticleDefinition != 0 && fParticleDefinition != particle)
90  {
91  G4Exception("G4DNABornExcitationModel1::Initialise","em0001",
92  FatalException,"Model already initialized for another particle type.");
93  }
94 
95  fParticleDefinition = particle;
96 
97  if(particle->GetParticleName() == "e-")
98  {
99  fTableFile = "dna/sigma_excitation_e_born";
100  fLowEnergy = 9*eV;
101  fHighEnergy = 1*MeV;
102  }
103  else if(particle->GetParticleName() == "proton")
104  {
105  fTableFile = "dna/sigma_excitation_p_born";
106  fLowEnergy = 500. * keV;
107  fHighEnergy = 100. * MeV;
108  }
109 
110  SetLowEnergyLimit(fLowEnergy);
111  SetHighEnergyLimit(fHighEnergy);
112 
113  G4double scaleFactor = (1.e-22 / 3.343) * m*m;
114  fTableData = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
115  fTableData->LoadData(fTableFile);
116 
117  if( verboseLevel>0 )
118  {
119  G4cout << "Born excitation model is initialized " << G4endl
120  << "Energy range: "
121  << LowEnergyLimit() / eV << " eV - "
122  << HighEnergyLimit() / keV << " keV for "
123  << particle->GetParticleName()
124  << G4endl;
125  }
126 
127  // Initialize water density pointer
129 
130  if (isInitialised)
131  { return;}
133  isInitialised = true;
134 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
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
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
static constexpr double m
Definition: G4SIunits.hh:129
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()
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4ParticleChangeForGamma * fParticleChangeForGamma
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
static constexpr double keV
Definition: G4SIunits.hh:216
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 181 of file G4DNABornExcitationModel1.cc.

186 {
187 
188  if (verboseLevel > 3)
189  {
190  G4cout << "Calling SampleSecondaries() of G4DNABornExcitationModel1"
191  << G4endl;
192  }
193 
194  G4double k = aDynamicParticle->GetKineticEnergy();
195 
196  G4int level = RandomSelect(k);
197  G4double excitationEnergy = waterStructure.ExcitationEnergy(level);
198  G4double newEnergy = k - excitationEnergy;
199 
200  if (newEnergy > 0)
201  {
203 
204  if (!statCode) fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
206 
208  }
209 
210  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
212  level,
213  theIncomingTrack);
214 }
G4double GetKineticEnergy() const
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
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
G4ParticleChangeForGamma * fParticleChangeForGamma
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

void G4DNABornExcitationModel1::SelectStationary ( G4bool  input)
inline

Definition at line 105 of file G4DNABornExcitationModel1.hh.

106 {
107  statCode = input;
108 }

Member Data Documentation

G4ParticleChangeForGamma* G4DNABornExcitationModel1::fParticleChangeForGamma
protected

Definition at line 75 of file G4DNABornExcitationModel1.hh.


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