Geant4-11
G4BetheBlochModel.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// -------------------------------------------------------------------
27//
28// GEANT4 Class header file
29//
30//
31// File name: G4BetheBlochModel
32//
33// Author: Vladimir Ivanchenko on base of Laszlo Urban code
34//
35// Creation date: 03.01.2002
36//
37// Modifications:
38//
39// 04-12-02 Fix problem of G4DynamicParticle constructor (V.Ivanchenko)
40// 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
41// 27-01-03 Make models region aware (V.Ivanchenko)
42// 13-02-03 Add name (V.Ivanchenko)
43// 24-03-05 Add G4EmCorrections (V.Ivanchenko)
44// 11-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
45// 11-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
46// 12-02-06 move G4LossTableManager::Instance()->EmCorrections()
47// in constructor (mma)
48// 12-08-08 Added methods GetParticleCharge, GetChargeSquareRatio,
49// CorrectionsAlongStep needed for ions(V.Ivanchenko)
50//
51// -------------------------------------------------------------------
52//
53
54
55//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57
58#include "G4BetheBlochModel.hh"
59#include "Randomize.hh"
61#include "G4SystemOfUnits.hh"
62#include "G4NistManager.hh"
63#include "G4Electron.hh"
64#include "G4LossTableManager.hh"
65#include "G4EmCorrections.hh"
66#include "G4EmParameters.hh"
69#include "G4Log.hh"
70#include "G4DeltaAngle.hh"
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73
74using namespace std;
75
77 const G4String& nam)
78 : G4VEmModel(nam),
79 twoln10(2.0*G4Log(10.0)),
80 fAlphaTlimit(1*CLHEP::GeV),
81 fProtonTlimit(10*CLHEP::GeV)
82{
87}
88
89//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
90
92{}
93
94//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
95
97 const G4DataVector&)
98{
99 if(p != particle) { SetupParameters(p); }
100
101 //G4cout << "G4BetheBlochModel::Initialise for " << p->GetParticleName()
102 // << " isIon= " << isIon
103 // << G4endl;
104
105 // always false before the run
106 SetDeexcitationFlag(false);
107
108 // initialisation once
109 if(nullptr == fParticleChange) {
111 if(IsMaster() && G4EmParameters::Instance()->UseICRU90Data() &&
112 (pname == "proton" || pname == "GenericIon" || pname == "alpha")) {
115 }
117 pname == "GenericIon") { isIon = true; }
118 if(pname == "alpha") { isAlpha = true; }
119
121 if(UseAngularGeneratorFlag() && nullptr == GetAngularDistribution()) {
123 }
124 }
125}
126
127//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
128
130 const G4Material* mat,
131 G4double kineticEnergy)
132{
133 // this method is called only for ions, so no check if it is an ion
134 return
135 (!isAlpha) ? corr->EffectiveChargeSquareRatio(p,mat,kineticEnergy) : 1.0;
136}
137
138//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
139
141 const G4Material* mat,
142 G4double kineticEnergy)
143{
144 // this method is called only for ions, so no check if it is an ion
145 return corr->GetParticleCharge(p,mat,kineticEnergy);
146}
147
148//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
149
151{
152 particle = p;
156 if(!isIon && q > 1.1) { isIon = true; }
157 chargeSquare = q*q;
159 static const G4double aMag = 1./(0.5*eplus*CLHEP::hbar_Planck*CLHEP::c_squared);
161 magMoment2 = magmom*magmom - 1.0;
162 formfact = 0.0;
163 tlimit = DBL_MAX;
164 if(particle->GetLeptonNumber() == 0) {
165 G4double x = 0.8426*CLHEP::GeV;
166 if(spin == 0.0 && mass < CLHEP::GeV) { x = 0.736*CLHEP::GeV; }
167 else if (mass > CLHEP::GeV) {
168 G4int iz = G4lrint(std::abs(q));
169 if(iz > 1) { x /= nist->GetA27(iz); }
170 }
172 tlimit = 2.0/formfact;
173 }
174}
175
176//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
177
179 const G4MaterialCutsCouple* couple)
180{
182}
183
184//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
185
188 G4double kineticEnergy,
189 G4double cutEnergy,
190 G4double maxKinEnergy)
191{
192 G4double cross = 0.0;
193 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
194 G4double maxEnergy = std::min(tmax, maxKinEnergy);
195 if(cutEnergy < maxEnergy) {
196
197 G4double totEnergy = kineticEnergy + mass;
198 G4double energy2 = totEnergy*totEnergy;
199 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
200
201 cross = (maxEnergy - cutEnergy)/(cutEnergy*maxEnergy)
202 - beta2*G4Log(maxEnergy/cutEnergy)/tmax;
203
204 // +term for spin=1/2 particle
205 if( 0.0 < spin ) { cross += 0.5*(maxEnergy - cutEnergy)/energy2; }
206
207 cross *= CLHEP::twopi_mc2_rcl2*chargeSquare/beta2;
208 }
209
210 // G4cout << "BB: e= " << kineticEnergy << " tmin= " << cutEnergy
211 // << " tmax= " << tmax << " cross= " << cross << G4endl;
212
213 return cross;
214}
215
216//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
217
219 const G4ParticleDefinition* p,
220 G4double kinEnergy,
222 G4double cutEnergy,
223 G4double maxEnergy)
224{
225 return Z*ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,maxEnergy);
226}
227
228//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
229
231 const G4Material* mat,
232 const G4ParticleDefinition* p,
233 G4double kinEnergy,
234 G4double cutEnergy,
235 G4double maxEnergy)
236{
237 G4double sigma = mat->GetElectronDensity()
238 *ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,maxEnergy);
239 if(isAlpha) {
240 sigma *= corr->EffectiveChargeSquareRatio(p,mat,kinEnergy)/chargeSquare;
241 }
242 return sigma;
243}
244
245//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
246
248 const G4ParticleDefinition* p,
249 G4double kineticEnergy,
250 G4double cut)
251{
252 const G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
253 // projectile formfactor limit energy loss
254 const G4double cutEnergy = std::min(std::min(cut,tmax), tlimit);
255
256 G4double tau = kineticEnergy/mass;
257 G4double gam = tau + 1.0;
258 G4double bg2 = tau * (tau+2.0);
259 G4double beta2 = bg2/(gam*gam);
260 G4double xc = cutEnergy/tmax;
261
262 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
263 G4double eexc2 = eexc*eexc;
264
265 G4double eDensity = material->GetElectronDensity();
266
267 // added ICRU90 stopping data for limited list of materials
268 /*
269 G4cout << "### DEDX ICRI90:" << (nullptr != fICRU90)
270 << " Ekin=" << kineticEnergy
271 << " " << p->GetParticleName()
272 << " q2=" << chargeSquare
273 << " inside " << material->GetName() << G4endl;
274 */
275 if(nullptr != fICRU90 && kineticEnergy < fProtonTlimit) {
278 baseMaterial = material->GetBaseMaterial()
279 ? material->GetBaseMaterial() : material;
281 }
282 if(iICRU90 >= 0) {
283 G4double dedx = 0.0;
284 // only for alpha
285 if(isAlpha) {
286 if(kineticEnergy <= fAlphaTlimit) {
287 dedx = fICRU90->GetElectronicDEDXforAlpha(iICRU90, kineticEnergy);
288 } else {
289 const G4double e = kineticEnergy*CLHEP::proton_mass_c2/mass;
291 }
292 } else {
293 dedx = fICRU90->GetElectronicDEDXforProton(iICRU90, kineticEnergy)
295 }
296 dedx *= material->GetDensity();
297 if(cutEnergy < tmax) {
298 dedx += (G4Log(xc) + (1.0 - xc)*beta2)*CLHEP::twopi_mc2_rcl2
299 *(eDensity*chargeSquare/beta2);
300 }
301 //G4cout << " iICRU90=" << iICRU90 << " dedx=" << dedx << G4endl;
302 if(dedx > 0.0) { return dedx; }
303 }
304 }
305 // general Bethe-Bloch formula
306 G4double dedx = G4Log(2.0*CLHEP::electron_mass_c2*bg2*cutEnergy/eexc2)
307 - (1.0 + xc)*beta2;
308
309 if(0.0 < spin) {
310 G4double del = 0.5*cutEnergy/(kineticEnergy + mass);
311 dedx += del*del;
312 }
313
314 // density correction
315 G4double x = G4Log(bg2)/twoln10;
316 dedx -= material->GetIonisation()->DensityCorrection(x);
317
318 // shell correction
319 dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
320
321 // now compute the total ionization loss
322 dedx *= CLHEP::twopi_mc2_rcl2*chargeSquare*eDensity/beta2;
323
324 //High order correction different for hadrons and ions
325 if(isIon) {
326 dedx += corr->IonBarkasCorrection(p,material,kineticEnergy);
327 } else {
328 dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy);
329 }
330
331 dedx = std::max(dedx, 0.0);
332 /*
333 G4cout << "E(MeV)= " << kineticEnergy/CLHEP::MeV << " dedx= " << dedx
334 << " " << material->GetName() << G4endl;
335 */
336 return dedx;
337}
338
339//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
340
342 const G4DynamicParticle* dp,
343 const G4double& length,
344 G4double& eloss)
345{
346 // no correction at the last step
347 const G4double preKinEnergy = dp->GetKineticEnergy();
348 if(eloss >= preKinEnergy) { return; }
349
350 const G4ParticleDefinition* p = dp->GetDefinition();
351 if(p != particle) { SetupParameters(p); }
352 if(isIon) {
353 const G4Material* mat = couple->GetMaterial();
354 G4double e = std::max(preKinEnergy - eloss*0.5, preKinEnergy*0.75);
355
356 const G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,e);
358 const G4double qfactor =
360 const G4double highOrder = length*corr->IonHighOrderCorrections(p,couple,e);
361 const G4double elossnew = eloss*qfactor + highOrder;
362 /*
363 G4cout << "G4BetheBlochModel::CorrectionsAlongStep: e= " << preKinEnergy
364 << " eloss=" << eloss << " elossnew=" << elossnew
365 << " qfactor= " << qfactor << " highOrder= " << highOrder << " (ratio="
366 << highOrder/eloss << ")"
367 " q2= " << q2 << " chargeSquare= " << chargeSquare << G4endl;
368 */
369 eloss = elossnew;
370 }
371}
372
373//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
374
375void G4BetheBlochModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
376 const G4MaterialCutsCouple* couple,
377 const G4DynamicParticle* dp,
378 G4double minKinEnergy,
379 G4double maxEnergy)
380{
381 G4double kineticEnergy = dp->GetKineticEnergy();
382 const G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy);
383 const G4double maxKinEnergy = std::min(maxEnergy,tmax);
384 if(minKinEnergy >= maxKinEnergy) { return; }
385
386 //G4cout << "G4BetheBlochModel::SampleSecondaries Emin= " << minKinEnergy
387 // << " Emax= " << maxKinEnergy << G4endl;
388
389 const G4double totEnergy = kineticEnergy + mass;
390 const G4double etot2 = totEnergy*totEnergy;
391 const G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
392
393 G4double deltaKinEnergy, f;
394 G4double f1 = 0.0;
395 G4double fmax = 1.0;
396 if( 0.0 < spin ) { fmax += 0.5*maxKinEnergy*maxKinEnergy/etot2; }
397
398 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
399 G4double rndm[2];
400
401 // sampling without nuclear size effect
402 do {
403 rndmEngineMod->flatArray(2, rndm);
404 deltaKinEnergy = minKinEnergy*maxKinEnergy
405 /(minKinEnergy*(1.0 - rndm[0]) + maxKinEnergy*rndm[0]);
406
407 f = 1.0 - beta2*deltaKinEnergy/tmax;
408 if( 0.0 < spin ) {
409 f1 = 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
410 f += f1;
411 }
412
413 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
414 } while( fmax*rndm[1] > f);
415
416 // projectile formfactor - suppresion of high energy
417 // delta-electron production at high energy
418
419 G4double x = formfact*deltaKinEnergy;
420 if(x > 1.e-6) {
421
422 G4double x1 = 1.0 + x;
423 G4double grej = 1.0/(x1*x1);
424 if( 0.0 < spin ) {
425 G4double x2 = 0.5*electron_mass_c2*deltaKinEnergy/(mass*mass);
426 grej *= (1.0 + magMoment2*(x2 - f1/f)/(1.0 + x2));
427 }
428 if(grej > 1.1) {
429 G4cout << "### G4BetheBlochModel WARNING: grej= " << grej
430 << " " << dp->GetDefinition()->GetParticleName()
431 << " Ekin(MeV)= " << kineticEnergy
432 << " delEkin(MeV)= " << deltaKinEnergy
433 << G4endl;
434 }
435 if(rndmEngineMod->flat() > grej) { return; }
436 }
437
438 G4ThreeVector deltaDirection;
439
441 const G4Material* mat = couple->GetMaterial();
442 deltaDirection =
443 GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy,
445 mat);
446 } else {
447
448 G4double deltaMomentum =
449 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
450 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
451 (deltaMomentum * dp->GetTotalMomentum());
452 cost = std::min(cost, 1.0);
453 const G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
454 const G4double phi = twopi*rndmEngineMod->flat();
455
456 deltaDirection.set(sint*std::cos(phi),sint*std::sin(phi), cost) ;
457 deltaDirection.rotateUz(dp->GetMomentumDirection());
458 }
459 /*
460 G4cout << "### G4BetheBlochModel "
461 << dp->GetDefinition()->GetParticleName()
462 << " Ekin(MeV)= " << kineticEnergy
463 << " delEkin(MeV)= " << deltaKinEnergy
464 << " tmin(MeV)= " << minKinEnergy
465 << " tmax(MeV)= " << maxKinEnergy
466 << " dir= " << dp->GetMomentumDirection()
467 << " dirDelta= " << deltaDirection
468 << G4endl;
469 */
470 // create G4DynamicParticle object for delta ray
471 G4DynamicParticle* delta =
472 new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
473
474 vdp->push_back(delta);
475
476 // Change kinematics of primary particle
477 kineticEnergy -= deltaKinEnergy;
478 G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
479 finalP = finalP.unit();
480
483}
484
485//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
486
488 G4double kinEnergy)
489{
490 // here particle type is checked for the case,
491 // when this model is shared between particles
492 if(pd != particle) { SetupParameters(pd); }
493 G4double tau = kinEnergy/mass;
494 return 2.0*CLHEP::electron_mass_c2*tau*(tau + 2.) /
495 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
496}
497
498//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4double G4Log(G4double x)
Definition: G4Log.hh:226
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double eplus
Definition: G4SIunits.hh:184
static constexpr double GeV
Definition: G4SIunits.hh:203
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Hep3Vector unit() const
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
virtual double flat()=0
virtual void flatArray(const int size, double *vect)=0
void CorrectionsAlongStep(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, const G4double &length, G4double &eloss) override
G4EmCorrections * corr
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
const G4ParticleDefinition * theElectron
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
G4double GetChargeSquareRatio() const
G4NistManager * nist
~G4BetheBlochModel() override
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
const G4ParticleDefinition * particle
G4ParticleChangeForLoss * fParticleChange
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
const G4Material * currentMaterial
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
G4BetheBlochModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="BetheBloch")
G4ICRU90StoppingData * fICRU90
void SetupParameters(const G4ParticleDefinition *p)
const G4Material * baseMaterial
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double GetTotalMomentum() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double IonHighOrderCorrections(const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double kineticEnergy)
G4double IonBarkasCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double HighOrderCorrections(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy, G4double cutEnergy)
G4double ShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double EffectiveChargeCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
static G4EmParameters * Instance()
G4double GetElectronicDEDXforProton(const G4Material *, G4double kinEnergy) const
G4int GetIndex(const G4Material *) const
G4double GetElectronicDEDXforAlpha(const G4Material *, G4double scaledKinEnergy) const
G4double GetMeanExcitationEnergy() const
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
const G4Material * GetMaterial() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:222
G4double GetElectronDensity() const
Definition: G4Material.hh:213
G4ICRU90StoppingData * GetICRU90StoppingData()
G4double GetA27(G4int Z) const
static G4NistManager * Instance()
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4double GetPDGMagneticMoment() const
G4double GetPDGCharge() const
const G4String & GetParticleName() const
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
virtual void SetParticleAndCharge(const G4ParticleDefinition *, G4double q2)
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:614
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:621
G4double inveplus
Definition: G4VEmModel.hh:431
G4bool IsMaster() const
Definition: G4VEmModel.hh:746
G4int SelectRandomAtomNumber(const G4Material *)
Definition: G4VEmModel.cc:297
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:774
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:823
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:628
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:718
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:108
Definition: DoubConv.h:17
static constexpr double eplus
static constexpr double electron_mass_c2
static constexpr double proton_mass_c2
static constexpr double GeV
static constexpr double hbar_Planck
static constexpr double MeV
static constexpr double c_squared
static constexpr double twopi_mc2_rcl2
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
string material
Definition: eplot.py:19
string pname
Definition: eplot.py:33
float electron_mass_c2
Definition: hepunit.py:273
int G4lrint(double ad)
Definition: templates.hh:134
#define DBL_MAX
Definition: templates.hh:62