Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4DNABornExcitationModel2 Class Reference

#include <G4DNABornExcitationModel2.hh>

Inheritance diagram for G4DNABornExcitationModel2:
Collaboration diagram for G4DNABornExcitationModel2:

Public Member Functions

 G4DNABornExcitationModel2 (const G4ParticleDefinition *p=0, const G4String &nam="DNABornExcitationModel")
 
virtual ~G4DNABornExcitationModel2 ()
 
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 SetFluctuationFlag (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
 
G4bool lossFlucFlag
 

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 45 of file G4DNABornExcitationModel2.hh.

Constructor & Destructor Documentation

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

Definition at line 44 of file G4DNABornExcitationModel2.cc.

45  :
46  G4VEmModel(nam), isInitialised(false), fTableData(0)
47 {
48  fpMolWaterDensity = 0;
49  fHighEnergy = 0;
50  fLowEnergy = 0;
51  fParticleDefinition = 0;
52 
53  verboseLevel = 0;
54  // Verbosity scale:
55  // 0 = nothing
56  // 1 = warning for energy non-conservation
57  // 2 = details of energy budget
58  // 3 = calculation of cross sections, file openings, sampling of atoms
59  // 4 = entering in methods
60 
61  if (verboseLevel > 0)
62  {
63  G4cout << "Born excitation model is constructed " << G4endl;
64  }
66  fLastBinCallForFinalXS = 0;
67  fTotalXS = 0;
68  fTableData = 0;
69 
70  // Selection of stationary mode
71 
72  statCode = false;
73 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
G4ParticleChangeForGamma * fParticleChangeForGamma
G4DNABornExcitationModel2::~G4DNABornExcitationModel2 ( )
virtual

Definition at line 77 of file G4DNABornExcitationModel2.cc.

78 {
79  // Cross section
80  if (fTableData)
81  delete fTableData;
82 }

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 188 of file G4DNABornExcitationModel2.cc.

193 {
194  if (verboseLevel > 3)
195  {
196  G4cout << "Calling CrossSectionPerVolume() of G4DNABornExcitationModel2"
197  << G4endl;
198  }
199 
200  if(particleDefinition != fParticleDefinition) return 0;
201 
202  // Calculate total cross section for model
203 
204  G4double sigma=0;
205 
206  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
207 
208  if(waterDensity!= 0.0)
209  {
210  if (ekin >= fLowEnergy && ekin < fHighEnergy)
211  {
212  sigma = fTotalXS->Value(ekin, fLastBinCallForFinalXS);
213 
214 // for(size_t i = 0; i < 5; ++i)
215 // sigma += (*fTableData)[i]->Value(ekin);
216 
217  if(sigma == 0)
218  {
219  G4cerr << "PROBLEM SIGMA = 0 at " << G4BestUnit(ekin, "Energy")<< G4endl;
220  }
221  }
222 
223  if (verboseLevel > 2)
224  {
225  G4cout << "__________________________________" << G4endl;
226  G4cout << "G4DNABornExcitationModel2 - XS INFO START" << G4endl;
227  G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
228  G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
229  G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
230  G4cout << "G4DNABornExcitationModel2 - XS INFO END" << G4endl;
231  }
232  } // if (waterMaterial)
233  else
234  {
235  G4cerr << "PROBLEM NO WATER MATERIAL ?" << G4endl;
236  }
237 
238  return sigma*waterDensity;
239 }
size_t GetIndex() const
Definition: G4Material.hh:262
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4GLOB_DLL std::ostream G4cout
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double eV
Definition: G4SIunits.hh:215
G4double Value(G4double theEnergy, size_t &lastidx) const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4GLOB_DLL std::ostream G4cerr

Here is the call graph for this function:

G4double G4DNABornExcitationModel2::GetPartialCrossSection ( const G4Material ,
G4int  level,
const G4ParticleDefinition particle,
G4double  kineticEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 280 of file G4DNABornExcitationModel2.cc.

284 {
285  if (fParticleDefinition != particle)
286  {
287  G4Exception("G4DNABornExcitationModel2::GetPartialCrossSection",
288  "bornParticleType",
290  "Model initialized for another particle type.");
291  }
292 
293  return (*fTableData)(level)->Value(kineticEnergy);
294 }
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:377

Here is the call graph for this function:

void G4DNABornExcitationModel2::Initialise ( const G4ParticleDefinition particle,
const G4DataVector = *(new G4DataVector()) 
)
virtual

Implements G4VEmModel.

Definition at line 86 of file G4DNABornExcitationModel2.cc.

88 {
89 
90  if (verboseLevel > 3)
91  {
92  G4cout << "Calling G4DNABornExcitationModel2::Initialise()" << G4endl;
93  }
94 
95  if(fParticleDefinition != 0 && fParticleDefinition != particle)
96  {
97  G4Exception("G4DNABornExcitationModel2::Initialise","em0001",
98  FatalException,"Model already initialized for another particle type.");
99  }
100 
101  fParticleDefinition = particle;
102 
103  std::ostringstream fullFileName;
104  char *path = getenv("G4LEDATA");
105 
106  if(G4String(path) == "")
107  {
108  G4Exception("G4DNABornExcitationModel2::Initialise","G4LEDATA-CHECK",
109  FatalException, "G4LEDATA not defined in environment variables");
110  }
111 
112  fullFileName << path;
113 
114  if(particle->GetParticleName() == "e-")
115  {
116  fullFileName << "/dna/bornExcitation-e.dat";
117  fLowEnergy = 9*eV;
118  fHighEnergy = 1*MeV;
119  }
120  else if(particle->GetParticleName() == "proton")
121  {
122  fullFileName << "/dna/bornExcitation-p.dat";
123  fLowEnergy = 500. * keV;
124  fHighEnergy = 100. * MeV;
125  }
126 
127  SetLowEnergyLimit(fLowEnergy);
128  SetHighEnergyLimit(fHighEnergy);
129 
130 // G4double scaleFactor = (1.e-22 / 3.343) * m*m;
131 
132  fTableData = new G4PhysicsTable();
133  fTableData->RetrievePhysicsTable(fullFileName.str().c_str(), true);
134  for(size_t level = 0; level<fTableData->size(); ++level)
135  {
136 // (*fTableData)(level)->ScaleVector(1,scaleFactor);
137  (*fTableData)(level)->SetSpline(true);
138  }
139 
140  size_t finalBin_i = 2000;
141  G4double E_min = fLowEnergy;
142  G4double E_max = fHighEnergy;
143  fTotalXS = new G4PhysicsLogVector(E_min, E_max, finalBin_i);
144  fTotalXS->SetSpline(true);
146  G4double finalXS;
147 
148  for(size_t energy_i = 0; energy_i < finalBin_i; ++energy_i)
149  {
150  energy = fTotalXS->Energy(energy_i);
151  finalXS = 0;
152 
153  for(size_t level = 0; level<fTableData->size(); ++level)
154  {
155  finalXS += (*fTableData)(level)->Value(energy);
156  }
157  fTotalXS->PutValue(energy_i, finalXS);
158 // G4cout << "energy = " << energy << " " << fTotalXS->Value(energy)
159 // << " " << energy_i << " " << finalXS << G4endl;
160  }
161 
162 // for(energy = LowEnergyLimit() ; energy < HighEnergyLimit() ; energy += 1*pow(10,log10(energy)))
163 // {
164 // G4cout << "energy = " << energy << " " << fTotalXS->Value(energy) << G4endl;
165 // }
166 
167  if( verboseLevel>0 )
168  {
169  G4cout << "Born excitation model is initialized " << G4endl
170  << "Energy range: "
171  << LowEnergyLimit() / eV << " eV - "
172  << HighEnergyLimit() / keV << " keV for "
173  << particle->GetParticleName()
174  << G4endl;
175  }
176 
177  // Initialize water density pointer
179 
180  if (isInitialised)
181  { return;}
183  isInitialised = true;
184 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:643
G4bool RetrievePhysicsTable(const G4String &filename, G4bool ascii=false)
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:636
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:602
const G4String & GetParticleName() const
void SetSpline(G4bool)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
G4GLOB_DLL std::ostream G4cout
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
void PutValue(size_t index, G4double theValue)
static constexpr double eV
Definition: G4SIunits.hh:215
G4double Energy(size_t index) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4DNAMolecularMaterial * Instance()
G4double energy(const ThreeVector &p, const G4double m)
#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:739
static constexpr double keV
Definition: G4SIunits.hh:216
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:377
G4ParticleChangeForGamma * fParticleChangeForGamma

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 243 of file G4DNABornExcitationModel2.cc.

248 {
249 
250  if (verboseLevel > 3)
251  {
252  G4cout << "Calling SampleSecondaries() of G4DNABornExcitationModel2"
253  << G4endl;
254  }
255 
256  G4double k = aDynamicParticle->GetKineticEnergy();
257 
258  G4int level = RandomSelect(k);
259  G4double excitationEnergy = waterStructure.ExcitationEnergy(level);
260  G4double newEnergy = k - excitationEnergy;
261 
262  if (newEnergy > 0)
263  {
265 
266  if (!statCode) fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
268 
270  }
271 
272  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
274  level,
275  theIncomingTrack);
276 }
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
double G4double
Definition: G4Types.hh:76
G4ParticleChangeForGamma * fParticleChangeForGamma

Here is the call graph for this function:

void G4DNABornExcitationModel2::SelectStationary ( G4bool  input)
inline

Definition at line 109 of file G4DNABornExcitationModel2.hh.

110 {
111  statCode = input;
112 }

Member Data Documentation

G4ParticleChangeForGamma* G4DNABornExcitationModel2::fParticleChangeForGamma
protected

Definition at line 77 of file G4DNABornExcitationModel2.hh.


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