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

#include <G4DNADingfelderChargeDecreaseModel.hh>

Inheritance diagram for G4DNADingfelderChargeDecreaseModel:
G4VEmModel

Public Member Functions

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

Constructor & Destructor Documentation

G4DNADingfelderChargeDecreaseModel::G4DNADingfelderChargeDecreaseModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNADingfelderChargeDecreaseModel" 
)

Definition at line 40 of file G4DNADingfelderChargeDecreaseModel.cc.

References fParticleChangeForGamma, G4cout, and G4endl.

42  :G4VEmModel(nam),isInitialised(false)
43 {
44  // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
45  fpMolWaterDensity = 0;
46  numberOfPartialCrossSections[0] = 0;
47  numberOfPartialCrossSections[1] = 0;
48  numberOfPartialCrossSections[2] = 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 << "Dingfelder charge decrease model is constructed " << G4endl;
61  }
63 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
G4DNADingfelderChargeDecreaseModel::~G4DNADingfelderChargeDecreaseModel ( )
virtual

Definition at line 67 of file G4DNADingfelderChargeDecreaseModel.cc.

68 {}

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 201 of file G4DNADingfelderChargeDecreaseModel.cc.

References python.hepunit::cm, python.hepunit::eV, G4cout, G4endl, G4Material::GetIndex(), G4DNAGenericIonsManager::GetIon(), G4ParticleDefinition::GetParticleName(), G4DNAGenericIonsManager::Instance(), and G4Proton::ProtonDefinition().

206 {
207  if (verboseLevel > 3)
208  G4cout << "Calling CrossSectionPerVolume() of G4DNADingfelderChargeDecreaseModel" << G4endl;
209 
210  // Calculate total cross section for model
211 
212  G4DNAGenericIonsManager *instance;
214 
215  if (
216  particleDefinition != G4Proton::ProtonDefinition()
217  &&
218  particleDefinition != instance->GetIon("alpha++")
219  &&
220  particleDefinition != instance->GetIon("alpha+")
221  )
222 
223  return 0;
224 
225  G4double lowLim = 0;
226  G4double highLim = 0;
227  G4double crossSection = 0.;
228 
229  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
230 
231  if(waterDensity!= 0.0)
232  // if (material == nistwater || material->GetBaseMaterial() == nistwater)
233  {
234  const G4String& particleName = particleDefinition->GetParticleName();
235 
236  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
237  pos1 = lowEnergyLimit.find(particleName);
238 
239  if (pos1 != lowEnergyLimit.end())
240  {
241  lowLim = pos1->second;
242  }
243 
244  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
245  pos2 = highEnergyLimit.find(particleName);
246 
247  if (pos2 != highEnergyLimit.end())
248  {
249  highLim = pos2->second;
250  }
251 
252  if (k >= lowLim && k < highLim)
253  {
254  crossSection = Sum(k,particleDefinition);
255  }
256 
257  if (verboseLevel > 2)
258  {
259  G4cout << "_______________________________________" << G4endl;
260  G4cout << "°°° G4DNADingfelderChargeDecreaeModel" << G4endl;
261  G4cout << "---> Kinetic energy(eV)=" << k/eV << "particle :" << particleName << G4endl;
262  G4cout << " - Cross section per water molecule (cm^2)=" << crossSection/cm/cm << G4endl;
263  G4cout << " - Cross section per water molecule (cm^-1)=" << crossSection*
264  waterDensity/(1./cm) << G4endl;
265 // material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
266  }
267  }
268 
269  return crossSection*waterDensity;
270 // return crossSection*material->GetAtomicNumDensityVector()[1];
271 
272 }
size_t GetIndex() const
Definition: G4Material.hh:260
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
G4GLOB_DLL std::ostream G4cout
static G4DNAGenericIonsManager * Instance(void)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4ParticleDefinition * GetIon(const G4String &name)
void G4DNADingfelderChargeDecreaseModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 72 of file G4DNADingfelderChargeDecreaseModel.cc.

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

74 {
75 
76  if (verboseLevel > 3)
77  G4cout << "Calling G4DNADingfelderChargeDecreaseModel::Initialise()" << G4endl;
78 
79  // Energy limits
80 
81  G4DNAGenericIonsManager *instance;
84  G4ParticleDefinition* alphaPlusPlusDef = instance->GetIon("alpha++");
85  G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
86 
88  G4String alphaPlusPlus;
89  G4String alphaPlus;
90 
91  // LIMITS
92 
93  proton = protonDef->GetParticleName();
94  lowEnergyLimit[proton] = 100. * eV;
95  highEnergyLimit[proton] = 100. * MeV;
96 
97  alphaPlusPlus = alphaPlusPlusDef->GetParticleName();
98  lowEnergyLimit[alphaPlusPlus] = 1. * keV;
99  highEnergyLimit[alphaPlusPlus] = 400. * MeV;
100 
101  alphaPlus = alphaPlusDef->GetParticleName();
102  lowEnergyLimit[alphaPlus] = 1. * keV;
103  highEnergyLimit[alphaPlus] = 400. * MeV;
104 
105  //
106 
107  if (particle==protonDef)
108  {
109  SetLowEnergyLimit(lowEnergyLimit[proton]);
110  SetHighEnergyLimit(highEnergyLimit[proton]);
111  }
112 
113  if (particle==alphaPlusPlusDef)
114  {
115  SetLowEnergyLimit(lowEnergyLimit[alphaPlusPlus]);
116  SetHighEnergyLimit(highEnergyLimit[alphaPlusPlus]);
117  }
118 
119  if (particle==alphaPlusDef)
120  {
121  SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
122  SetHighEnergyLimit(highEnergyLimit[alphaPlus]);
123  }
124 
125  // Final state
126 
127  //PROTON
128  f0[0][0]=1.;
129  a0[0][0]=-0.180;
130  a1[0][0]=-3.600;
131  b0[0][0]=-18.22;
132  b1[0][0]=-1.997;
133  c0[0][0]=0.215;
134  d0[0][0]=3.550;
135  x0[0][0]=3.450;
136  x1[0][0]=5.251;
137 
138  numberOfPartialCrossSections[0] = 1;
139 
140  //ALPHA++
141  f0[0][1]=1.; a0[0][1]=0.95;
142  a1[0][1]=-2.75;
143  b0[0][1]=-23.00;
144  c0[0][1]=0.215;
145  d0[0][1]=2.95;
146  x0[0][1]=3.50;
147 
148  f0[1][1]=1.;
149  a0[1][1]=0.95;
150  a1[1][1]=-2.75;
151  b0[1][1]=-23.73;
152  c0[1][1]=0.250;
153  d0[1][1]=3.55;
154  x0[1][1]=3.72;
155 
156  x1[0][1]=-1.;
157  b1[0][1]=-1.;
158 
159  x1[1][1]=-1.;
160  b1[1][1]=-1.;
161 
162  numberOfPartialCrossSections[1] = 2;
163 
164  // ALPHA+
165  f0[0][2]=1.;
166  a0[0][2]=0.65;
167  a1[0][2]=-2.75;
168  b0[0][2]=-21.81;
169  c0[0][2]=0.232;
170  d0[0][2]=2.95;
171  x0[0][2]=3.53;
172 
173  x1[0][2]=-1.;
174  b1[0][2]=-1.;
175 
176  numberOfPartialCrossSections[2] = 1;
177 
178  //
179 
180  if( verboseLevel>0 )
181  {
182  G4cout << "Dingfelder charge decrease model is initialized " << G4endl
183  << "Energy range: "
184  << LowEnergyLimit() / keV << " keV - "
185  << HighEnergyLimit() / MeV << " MeV for "
186  << particle->GetParticleName()
187  << G4endl;
188  }
189 
190  // Initialize water density pointer
192 
193  if (isInitialised) { return; }
195  isInitialised = true;
196 
197 }
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
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 G4DNAGenericIonsManager * Instance(void)
static G4DNAMolecularMaterial * Instance()
#define G4endl
Definition: G4ios.hh:61
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
G4ParticleDefinition * GetIon(const G4String &name)
void G4DNADingfelderChargeDecreaseModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicParticle,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 276 of file G4DNADingfelderChargeDecreaseModel.cc.

References python.hepunit::electron_mass_c2, FatalException, fParticleChangeForGamma, fStopAndKill, G4cout, G4endl, G4Exception(), G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4ParticleDefinition::GetPDGMass(), n, G4VParticleChange::ProposeLocalEnergyDeposit(), G4VParticleChange::ProposeTrackStatus(), G4Proton::Proton(), and python.hepunit::proton_mass_c2.

281 {
282  if (verboseLevel > 3)
283  G4cout << "Calling SampleSecondaries() of G4DNADingfelderChargeDecreaseModel" << G4endl;
284 
285  G4double inK = aDynamicParticle->GetKineticEnergy();
286 
287  G4ParticleDefinition* definition = aDynamicParticle->GetDefinition();
288 
289  G4double particleMass = definition->GetPDGMass();
290 
291  G4int finalStateIndex = RandomSelect(inK,definition);
292 
293  G4int n = NumberOfFinalStates(definition, finalStateIndex);
294  G4double waterBindingEnergy = WaterBindingEnergyConstant(definition, finalStateIndex);
295  G4double outgoingParticleBindingEnergy = OutgoingParticleBindingEnergyConstant(definition, finalStateIndex);
296 
297  G4double outK = 0.;
298  if (definition==G4Proton::Proton())
299  outK = inK - n*(inK*electron_mass_c2/proton_mass_c2) - waterBindingEnergy + outgoingParticleBindingEnergy;
300  else
301  outK = inK - n*(inK*electron_mass_c2/particleMass) - waterBindingEnergy + outgoingParticleBindingEnergy;
302 
303  if (outK<0)
304  {
305  G4Exception("G4DNADingfelderChargeDecreaseModel::SampleSecondaries","em0004",
306  FatalException,"Final kinetic energy is negative.");
307  }
308 
311 
312  G4DynamicParticle* dp = new G4DynamicParticle (OutgoingParticleDefinition(definition, finalStateIndex),
313  aDynamicParticle->GetMomentumDirection(),
314  outK) ;
315  fvect->push_back(dp);
316 
317 }
G4double GetKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
static G4Proton * Proton()
Definition: G4Proton.cc:93
float proton_mass_c2
Definition: hepunit.py:275
float electron_mass_c2
Definition: hepunit.py:274
const G4int n
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double GetPDGMass() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)

Field Documentation

G4ParticleChangeForGamma* G4DNADingfelderChargeDecreaseModel::fParticleChangeForGamma
protected

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