Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56
57#include "G4BetheBlochModel.hh"
58#include "Randomize.hh"
60#include "G4SystemOfUnits.hh"
61#include "G4NistManager.hh"
62#include "G4Electron.hh"
63#include "G4LossTableManager.hh"
64#include "G4EmCorrections.hh"
65#include "G4EmParameters.hh"
68#include "G4Log.hh"
69#include "G4DeltaAngle.hh"
70#include <vector>
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73
75 const G4String& nam)
76 : G4VEmModel(nam),
77 twoln10(2.0*G4Log(10.0)),
78 fAlphaTlimit(1*CLHEP::GeV),
79 fProtonTlimit(10*CLHEP::GeV)
80{
81 theElectron = G4Electron::Electron();
84 SetLowEnergyLimit(2.0*CLHEP::MeV);
85}
86
87//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
88
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
92
94 const G4DataVector&)
95{
96 if(p != particle) { SetupParameters(p); }
97
98 // always false before the run
100
101 // initialisation once
102 if(nullptr == fParticleChange) {
103 const G4String& pname = particle->GetParticleName();
104 if(G4EmParameters::Instance()->UseICRU90Data() &&
105 (pname == "proton" || pname == "GenericIon" || pname == "alpha")) {
106 fICRU90 = nist->GetICRU90StoppingData();
107 }
108 if(particle->GetPDGCharge() > CLHEP::eplus ||
109 pname == "GenericIon") { isIon = true; }
110 if(pname == "alpha") { isAlpha = true; }
111
112 fParticleChange = GetParticleChangeForLoss();
113 if(UseAngularGeneratorFlag() && nullptr == GetAngularDistribution()) {
115 }
116 }
117 // initialisation for each new run
118 if(IsMaster() && nullptr != fICRU90) {
119 fICRU90->Initialise();
120 }
121}
122
123//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
124
126 const G4Material* mat,
127 G4double kinEnergy)
128{
129 // this method is called only for ions, so no check if it is an ion
130 if(isAlpha) { return 1.0; }
131 chargeSquare = corr->EffectiveChargeSquareRatio(p, mat, kinEnergy);
132 return chargeSquare;
133}
134
135//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
136
138 const G4Material* mat,
139 G4double kineticEnergy)
140{
141 // this method is called only for ions, so no check if it is an ion
142 return corr->GetParticleCharge(p, mat, kineticEnergy);
143}
144
145//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
146
147void G4BetheBlochModel::SetupParameters(const G4ParticleDefinition* p)
148{
149 particle = p;
150 mass = particle->GetPDGMass();
151 spin = particle->GetPDGSpin();
152 G4double q = particle->GetPDGCharge()*inveplus;
153 isIon = (!isAlpha && q > 1.1);
154 chargeSquare = q*q;
155 ratio = electron_mass_c2/mass;
156 constexpr G4double aMag = 1./(0.5*eplus*CLHEP::hbar_Planck*CLHEP::c_squared);
157 G4double magmom = particle->GetPDGMagneticMoment()*mass*aMag;
158 magMoment2 = magmom*magmom - 1.0;
159 formfact = 0.0;
160 tlimit = DBL_MAX;
161 if(particle->GetLeptonNumber() == 0) {
162 G4double x = 0.8426*CLHEP::GeV;
163 if(spin == 0.0 && mass < CLHEP::GeV) { x = 0.736*CLHEP::GeV; }
164 else if (mass > CLHEP::GeV) {
165 G4int iz = G4lrint(std::abs(q));
166 if(iz > 1) { x /= nist->GetA27(iz); }
167 }
168 formfact = 2.0*CLHEP::electron_mass_c2/(x*x);
169 tlimit = 2.0/formfact;
170 }
171}
172
173//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
174
180
181//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
182
185 G4double kineticEnergy,
186 G4double cut,
187 G4double maxKinEnergy)
188{
189 G4double cross = 0.0;
190 const G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
191 const G4double cutEnergy = std::min(std::min(cut,tmax), tlimit);
192 const G4double maxEnergy = std::min(tmax, maxKinEnergy);
193 if(cutEnergy < maxEnergy) {
194
195 G4double totEnergy = kineticEnergy + mass;
196 G4double energy2 = totEnergy*totEnergy;
197 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
198
199 cross = (maxEnergy - cutEnergy)/(cutEnergy*maxEnergy)
200 - beta2*G4Log(maxEnergy/cutEnergy)/tmax;
201
202 // +term for spin=1/2 particle
203 if( 0.0 < spin ) { cross += 0.5*(maxEnergy - cutEnergy)/energy2; }
204
205 cross *= CLHEP::twopi_mc2_rcl2*chargeSquare/beta2;
206 }
207
208 // G4cout << "BB: e= " << kineticEnergy << " tmin= " << cutEnergy
209 // << " tmax= " << tmax << " cross= " << cross << G4endl;
210
211 return cross;
212}
213
214//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
215
217 const G4ParticleDefinition* p,
218 G4double kinEnergy,
220 G4double cutEnergy,
221 G4double maxEnergy)
222{
223 return Z*ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,maxEnergy);
224}
225
226//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
227
229 const G4Material* mat,
230 const G4ParticleDefinition* p,
231 G4double kinEnergy,
232 G4double cutEnergy,
233 G4double maxEnergy)
234{
235 G4double sigma = mat->GetElectronDensity()
236 *ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,maxEnergy);
237 if(isAlpha) {
238 sigma *= corr->EffectiveChargeSquareRatio(p,mat,kinEnergy)/chargeSquare;
239 }
240 return sigma;
241}
242
243//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
244
246 const G4ParticleDefinition* p,
247 G4double kineticEnergy,
248 G4double cut)
249{
250 const G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
251 // projectile formfactor limit energy loss
252 const G4double cutEnergy = std::min(std::min(cut,tmax), tlimit);
253
254 G4double tau = kineticEnergy/mass;
255 G4double gam = tau + 1.0;
256 G4double bg2 = tau * (tau+2.0);
257 G4double beta2 = bg2/(gam*gam);
258 G4double xc = cutEnergy/tmax;
259
260 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
261 G4double eexc2 = eexc*eexc;
262
263 G4double eDensity = material->GetElectronDensity();
264
265 // added ICRU90 stopping data for limited list of materials
266 /*
267 G4cout << "### DEDX ICRI90:" << (nullptr != fICRU90)
268 << " Ekin=" << kineticEnergy
269 << " " << p->GetParticleName()
270 << " q2=" << chargeSquare
271 << " inside " << material->GetName() << G4endl;
272 */
273 if(nullptr != fICRU90 && kineticEnergy < fProtonTlimit) {
274 if(material != currentMaterial) {
275 currentMaterial = material;
276 baseMaterial = material->GetBaseMaterial()
277 ? material->GetBaseMaterial() : material;
278 iICRU90 = fICRU90->GetIndex(baseMaterial);
279 }
280 if(iICRU90 >= 0) {
281 G4double dedx = 0.0;
282 // only for alpha
283 if(isAlpha) {
284 if(kineticEnergy <= fAlphaTlimit) {
285 dedx = fICRU90->GetElectronicDEDXforAlpha(iICRU90, kineticEnergy);
286 } else {
287 const G4double e = kineticEnergy*CLHEP::proton_mass_c2/mass;
288 dedx = fICRU90->GetElectronicDEDXforProton(iICRU90, e)*chargeSquare;
289 }
290 } else {
291 dedx = fICRU90->GetElectronicDEDXforProton(iICRU90, kineticEnergy)
292 *chargeSquare;
293 }
294 dedx *= material->GetDensity();
295 if(cutEnergy < tmax) {
296 dedx += (G4Log(xc) + (1.0 - xc)*beta2)*CLHEP::twopi_mc2_rcl2
297 *(eDensity*chargeSquare/beta2);
298 }
299 //G4cout << " iICRU90=" << iICRU90 << " dedx=" << dedx << G4endl;
300 if(dedx > 0.0) { return dedx; }
301 }
302 }
303 // general Bethe-Bloch formula
304 G4double dedx = G4Log(2.0*CLHEP::electron_mass_c2*bg2*cutEnergy/eexc2)
305 - (1.0 + xc)*beta2;
306
307 if(0.0 < spin) {
308 G4double del = 0.5*cutEnergy/(kineticEnergy + mass);
309 dedx += del*del;
310 }
311
312 // density correction
313 G4double x = G4Log(bg2)/twoln10;
314 dedx -= material->GetIonisation()->DensityCorrection(x);
315
316 // shell correction
317 dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
318
319 // now compute the total ionization loss
320 dedx *= CLHEP::twopi_mc2_rcl2*chargeSquare*eDensity/beta2;
321
322 //High order correction different for hadrons and ions
323 if(isIon) {
324 dedx += corr->IonBarkasCorrection(p,material,kineticEnergy);
325 } else {
326 dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy);
327 }
328
329 dedx = std::max(dedx, 0.0);
330 /*
331 G4cout << "E(MeV)= " << kineticEnergy/CLHEP::MeV << " dedx= " << dedx
332 << " " << material->GetName() << G4endl;
333 */
334 return dedx;
335}
336
337//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
338
340 const G4DynamicParticle* dp,
341 const G4double& /*length*/,
342 G4double& eloss)
343{
344 // no correction for alpha
345 if(isAlpha) { return; }
346
347 // no correction at the last step or at small step
348 const G4double preKinEnergy = dp->GetKineticEnergy();
349 if(eloss >= preKinEnergy || eloss < preKinEnergy*0.05) { return; }
350
351 // corrections for all charged particles with Q > 1
352 const G4ParticleDefinition* p = dp->GetDefinition();
353 if(p != particle) { SetupParameters(p); }
354 if(!isIon) { return; }
355
356 // effective energy and charge at a step
357 const G4double e = std::max(preKinEnergy - eloss*0.5, preKinEnergy*0.5);
358 const G4Material* mat = couple->GetMaterial();
359 const G4double q20 = corr->EffectiveChargeSquareRatio(p, mat, preKinEnergy);
360 const G4double q2 = corr->EffectiveChargeSquareRatio(p, mat, e);
361 const G4double qfactor = q2/q20;
362
363 /*
364 G4cout << "G4BetheBlochModel::CorrectionsAlongStep: Epre(MeV)="
365 << preKinEnergy << " Eeff(MeV)=" << e
366 << " eloss=" << eloss << " elossnew=" << eloss*qfactor
367 << " qfactor=" << qfactor << " Qpre=" << q20
368 << p->GetParticleName() <<G4endl;
369 */
370 eloss *= qfactor;
371}
372
373//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
374
375void G4BetheBlochModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
376 const G4MaterialCutsCouple* couple,
377 const G4DynamicParticle* dp,
378 G4double cut,
379 G4double maxEnergy)
380{
381 G4double kinEnergy = dp->GetKineticEnergy();
382 const G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(), kinEnergy);
383 const G4double minKinEnergy = std::min(cut, tmax);
384 const G4double maxKinEnergy = std::min(maxEnergy, tmax);
385 if(minKinEnergy >= maxKinEnergy) { return; }
386
387 //G4cout << "G4BetheBlochModel::SampleSecondaries Emin= " << minKinEnergy
388 // << " Emax= " << maxKinEnergy << G4endl;
389
390 const G4double totEnergy = kinEnergy + mass;
391 const G4double etot2 = totEnergy*totEnergy;
392 const G4double beta2 = kinEnergy*(kinEnergy + 2.0*mass)/etot2;
393
394 G4double deltaKinEnergy, f;
395 G4double f1 = 0.0;
396 G4double fmax = 1.0;
397 if( 0.0 < spin ) { fmax += 0.5*maxKinEnergy*maxKinEnergy/etot2; }
398
399 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
400 G4double rndm[2];
401
402 // sampling without nuclear size effect
403 do {
404 rndmEngineMod->flatArray(2, rndm);
405 deltaKinEnergy = minKinEnergy*maxKinEnergy
406 /(minKinEnergy*(1.0 - rndm[0]) + maxKinEnergy*rndm[0]);
407
408 f = 1.0 - beta2*deltaKinEnergy/tmax;
409 if( 0.0 < spin ) {
410 f1 = 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
411 f += f1;
412 }
413
414 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
415 } while( fmax*rndm[1] > f);
416
417 // projectile formfactor - suppresion of high energy
418 // delta-electron production at high energy
419
420 G4double x = formfact*deltaKinEnergy;
421 if(x > 1.e-6) {
422
423 G4double x1 = 1.0 + x;
424 G4double grej = 1.0/(x1*x1);
425 if( 0.0 < spin ) {
426 G4double x2 = 0.5*electron_mass_c2*deltaKinEnergy/(mass*mass);
427 grej *= (1.0 + magMoment2*(x2 - f1/f)/(1.0 + x2));
428 }
429 if(grej > 1.1) {
430 G4cout << "### G4BetheBlochModel WARNING: grej= " << grej
431 << " " << dp->GetDefinition()->GetParticleName()
432 << " Ekin(MeV)= " << kinEnergy
433 << " delEkin(MeV)= " << deltaKinEnergy
434 << G4endl;
435 }
436 if(rndmEngineMod->flat() > grej) { return; }
437 }
438
439 G4ThreeVector deltaDirection;
440
442 const G4Material* mat = couple->GetMaterial();
443 deltaDirection =
444 GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy,
446 mat);
447 } else {
448
449 G4double deltaMomentum =
450 std::sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
451 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
452 (deltaMomentum * dp->GetTotalMomentum());
453 cost = std::min(cost, 1.0);
454 const G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
455 const G4double phi = twopi*rndmEngineMod->flat();
456
457 deltaDirection.set(sint*std::cos(phi),sint*std::sin(phi), cost) ;
458 deltaDirection.rotateUz(dp->GetMomentumDirection());
459 }
460 /*
461 G4cout << "### G4BetheBlochModel "
462 << dp->GetDefinition()->GetParticleName()
463 << " Ekin(MeV)= " << kinEnergy
464 << " delEkin(MeV)= " << deltaKinEnergy
465 << " tmin(MeV)= " << minKinEnergy
466 << " tmax(MeV)= " << maxKinEnergy
467 << " dir= " << dp->GetMomentumDirection()
468 << " dirDelta= " << deltaDirection
469 << G4endl;
470 */
471 // create G4DynamicParticle object for delta ray
472 auto delta = new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
473
474 vdp->push_back(delta);
475
476 // Change kinematics of primary particle
477 kinEnergy -= deltaKinEnergy;
478 G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
479 finalP = finalP.unit();
480
481 fParticleChange->SetProposedKineticEnergy(kinEnergy);
482 fParticleChange->SetProposedMomentumDirection(finalP);
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:227
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
Hep3Vector unit() const
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
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
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)
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
~G4BetheBlochModel() override
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
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")
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double GetTotalMomentum() const
static G4Electron * Electron()
Definition G4Electron.cc:91
G4double HighOrderCorrections(const G4ParticleDefinition *, const G4Material *, const G4double kineticEnergy, const G4double cutEnergy)
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, const G4double kineticEnergy)
G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, const G4double kineticEnergy)
G4double ShellCorrection(const G4ParticleDefinition *, const G4Material *, const G4double kineticEnergy)
G4double IonBarkasCorrection(const G4ParticleDefinition *, const G4Material *, const 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 DensityCorrection(G4double x) const
G4double GetMeanExcitationEnergy() const
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
const G4Material * GetMaterial() const
G4double GetDensity() const
const G4Material * GetBaseMaterial() const
G4IonisParamMat * GetIonisation() const
G4double GetElectronDensity() const
G4ICRU90StoppingData * GetICRU90StoppingData()
G4double GetA27(G4int Z) const
static G4NistManager * Instance()
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4double GetPDGMagneticMoment() const
const G4String & GetParticleName() const
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4VEmAngularDistribution * GetAngularDistribution()
G4double inveplus
G4int SelectRandomAtomNumber(const G4Material *) const
G4bool IsMaster() const
void SetLowEnergyLimit(G4double)
void SetDeexcitationFlag(G4bool val)
void SetAngularDistribution(G4VEmAngularDistribution *)
G4bool UseAngularGeneratorFlag() const
G4ParticleChangeForLoss * GetParticleChangeForLoss()
int G4lrint(double ad)
Definition templates.hh:134
#define DBL_MAX
Definition templates.hh:62