KaliVeda
Toolkit for HIC analysis
KVedaLossMaterial.cpp
1 //Created by KVClassFactory on Wed Feb 2 16:13:08 2011
2 //Author: frankland,,,,
3 
4 #include "KVedaLossMaterial.h"
5 #include <TMath.h>
6 #include "KVIonRangeTable.h"
7 #include <TEnv.h>
8 #include <TF1.h>
9 #include "KVedaLossInverseRangeFunction.h"
10 #include "KVedaLoss.h"
11 
13 
16 
17 
20 
22  : KVIonRangeTableMaterial(), fInvRange(ZMAX_VEDALOSS, 1),
23  fEmin(ZMAX_VEDALOSS), fEmax(ZMAX_VEDALOSS), fCoeff(ZMAX_VEDALOSS, std::vector<Double_t>(14))
24 {
25  // Default constructor
26  for (int i = 0; i < ZMAX_VEDALOSS; i++) {
27  fEmin[i] = 0.0;
28  fEmax[i] = 500.0;
29  }
31 }
32 
33 
34 
37 
38 KVedaLossMaterial::KVedaLossMaterial(const KVIonRangeTable* t, const Char_t* name, const Char_t* type, const Char_t* state,
39  Double_t density, Double_t Z, Double_t A, Double_t)
40  : KVIonRangeTableMaterial(t, name, type, state, density, Z, A), fInvRange(ZMAX_VEDALOSS, 1),
41  fEmin(ZMAX_VEDALOSS), fEmax(ZMAX_VEDALOSS), fCoeff(ZMAX_VEDALOSS, std::vector<Double_t>(14))
42 {
43  // create new material
44  fgTable = static_cast<KVedaLoss*>(const_cast<KVIonRangeTable*>(t));
45  for (int i = 0; i < ZMAX_VEDALOSS; i++) {
46  fEmin[i] = 0.0;
47  fEmax[i] = 500.0;
48  }
50 }
51 
52 
53 
56 
58 {
59  // Destructor
60 }
61 
62 
63 
79 
81 {
82  // Read Z- & A-dependent range parameters for material
83  //
84  // For each material we create 4 TF1 objects:
85  // - KVedaLossMaterial:[type]:Range - gives range in g/cm**2 as a function of particle energy
86  // - KVedaLossMaterial:[type]:StoppingPower - gives dE/dx in MeV/(g/cm**2) as a function of particle energy
87  // - KVedaLossMaterial:[type]:EnergyLoss - gives dE as a function of particle energy
88  // - KVedaLossMaterial:[type]:ResidualEnergy - gives energy after material (0 if particle stops)
89  //
90  // The TF1::fNpx parameter for these functions is defined by the environment variables
91  //
92  // - KVedaLoss.Range.Npx: 20 /* also used for StoppingPower */
93  // - KVedaLoss.EnergyLoss.Npx: 50
94  // - KVedaLoss.ResidualEnergy.Npx: 20
95  //
96 
97  char line[132];
98 
99  //look for energy limits to calculation validity
100  if (!fgets(line, 132, fp)) {
101  Warning("ReadRangeTable", "Problem reading energy limits in range table file for %s (%s)",
102  GetName(), GetType());
103  return kFALSE;
104  }
105  else {
106  if (!strncmp(line, "COMPOUND", 8)) {
107  // material is compound. read composition for TGeoMixture.
108  if (fgets(line, 132, fp)) {}
109  int nel = atoi(line); // read number of elements
110  for (int el = 0; el < nel; el++) {
111  if (fgets(line, 132, fp)) {}
112  int z, a, w;
113  sscanf(line, "%d %d %d", &z, &a, &w);
114  AddCompoundElement(z, a, w);
115  }
116  }
117  else if (!strncmp(line, "MIXTURE", 7)) {
118  // material is mixture. read composition for TGeoMixture.
119  if (fgets(line, 132, fp)) {}
120  int nel = atoi(line); // read number of elements
121  for (int el = 0; el < nel; el++) {
122  if (fgets(line, 132, fp)) {}
123  int z, a, nat;
124  float w;
125  sscanf(line, "%d %d %d %f", &z, &a, &nat, &w);
126  AddMixtureElement(z, a, nat, w);
127  }
128  }
129  }
130  if (!fgets(line, 132, fp)) {
131  Warning("ReadRangeTable", "Problem reading energy limits in range table file for %s (%s)",
132  GetName(), GetType());
133  return kFALSE;
134  }
135  else {
136  while (line[0] == 'Z') {
137  Int_t z1, z2;
138  Float_t e1, e2;
139  if (sscanf(line, "Z = %d,%d %f < E/A < %f MeV", &z1,
140  &z2, &e1, &e2) != 4) {
141  Error("ReadRangeTable", "Problem reading energy limits in range table file for %s (%s)",
142  GetName(), GetType());
143  return kFALSE;
144  }
145  line[0] = '\0';
146  if (!fgets(line, 132, fp)) break;
147  for (int i = z1; i <= z2; i++) {
148  fEmin[i - 1] = e1;
149  fEmax[i - 1] = e2;
150  }
151  }
152  if (line[0] == '\0') {
153  Error("ReadRangeTable", "Problem with range table file for %s (%s)",
154  GetName(), GetType());
155  return kFALSE;
156  }
157  }
158 
159  // get require Npx value from (user-defined) environment variables
160  Int_t my_npx = gEnv->GetValue("KVedaLoss.Range.Npx", 100);
161 
162  fRange = new TF1(Form("KVedaLossMaterial:%s:Range", GetType()), this, &KVedaLossMaterial::RangeFunc,
163  0., 1.e+03, 0, "KVedaLossMaterial", "RangeFunc");
164  fRange->SetNpx(my_npx);
165 
166  fStopping = new TF1(Form("KVedaLossMaterial:%s:StoppingPower", GetType()), this, &KVedaLossMaterial::StoppingFunc,
167  0., 1.e+03, 0, "KVedaLossMaterial", "StoppingFunc");
168  fStopping->SetNpx(my_npx);
169 
170  my_npx = gEnv->GetValue("KVedaLoss.EnergyLoss.Npx", 100);
171  fDeltaE = new TF1(Form("KVedaLossMaterial:%s:EnergyLoss", GetType()), this, &KVedaLossMaterial::DeltaEFunc,
172  0., 1.e+03, 0, "KVedaLossMaterial", "DeltaEFunc");
173  fDeltaE->SetNpx(my_npx);
174 
175  my_npx = gEnv->GetValue("KVedaLoss.ResidualEnergy.Npx", 100);
176  fEres = new TF1(Form("KVedaLossMaterial:%s:ResidualEnergy", GetType()), this, &KVedaLossMaterial::EResFunc,
177  0., 1.e+03, 0, "KVedaLossMaterial", "EResFunc");
178  fEres->SetNpx(my_npx);
179 
180  for (int count = 0; count < ZMAX_VEDALOSS; count++) {
181 
182  if (sscanf(line, "%lf %lf %lf %lf %lf %lf %lf %lf",
183  &fCoeff[count][0], &fCoeff[count][1],
184  &fCoeff[count][2], &fCoeff[count][3],
185  &fCoeff[count][4], &fCoeff[count][5],
186  &fCoeff[count][6], &fCoeff[count][7])
187  != 8) {
188  Error("ReadRangeTable", "problem reading coeffs 0-7 in range table for %s (%s)", GetName(), GetType());
189  return kFALSE;
190  }
191  if (!fgets(line, 132, fp)) {
192  Warning("ReadRangeTable", "file too short ??? %s (%s)", GetName(), GetType());
193  return kFALSE;
194  }
195  else {
196  if (sscanf(line, "%lf %lf %lf %lf %lf %lf",
197  &fCoeff[count][8], &fCoeff[count][9],
198  &fCoeff[count][10], &fCoeff[count][11],
199  &fCoeff[count][12], &fCoeff[count][13])
200  != 6) {
201  Error("ReadRangeTable", "problem reading coeffs 8-13 in range table for %s (%s)", GetName(), GetType());
202  return kFALSE;
203  }
204  }
205  if (fNoLimits) {
206  // if we ignore nominal validity limits on incident energy, we must still use energy limits
207  // such that all range functions increase monotonically in the energy interval
208  GetRangeFunction(fCoeff[count][0], fCoeff[count][1])->SetRange(GetEminValid(fCoeff[count][0], fCoeff[count][1]), VERY_BIG_ENERGY);
209  Double_t emax = fRange->GetMaximumX() - 1;
210  emax /= fCoeff[count][1];
211  Double_t original_emax = fEmax[count];
212  // the new emax is only accepted if it is > than the nominal emax (400 or 250 AMeV),
213  // and at most 1 GeV/nucleon
214  fEmax[count] = TMath::Min(TMath::Max(original_emax, emax), 1000.);
215  // we may further reduce the upper limit to correspond to the minimum of stopping,
216  // if one exists
217  GetStoppingFunction(fCoeff[count][0], fCoeff[count][1])->SetRange(GetEminValid(fCoeff[count][0], fCoeff[count][1]), GetEmaxValid(fCoeff[count][0], fCoeff[count][1]));
218  emax = fStopping->GetMinimumX();
219  emax /= fCoeff[count][1];
220  // again, the new emax is only accepted if it is > than the nominal emax (400 or 250 AMeV),
221  // and at most 1 GeV/nucleon
222  fEmax[count] = TMath::Min(TMath::Max(original_emax, emax), 1000.);
223  //if(fEmax[count]!=original_emax) Info("ReadRangeTable", "Max. incident E for Z=%d ===> E/A = %f", count+1, fEmax[count]);
224  }
225  if (fgets(line, 132, fp)) {}
226  }
227 
228  return kTRUE;
229 }
230 
231 
232 
238 
240 {
241  // Function parameterising the energy loss of charged particles in this material.
242  //
243  // \param E[0] incident energy given in MeV.
244  // \returns energy loss calculated in MeV.
245 
246  return (E[0] - EResFunc(E, nullptr));
247 }
248 
249 
250 
256 
258 {
259  // Function parameterising the residual energy of charged particles in this material.
260  //
261  // \param E[0] incident energy given in MeV.
262  // \returns residual energy calculated in MeV.
263 
264  Double_t R0 = RangeFunc(E, nullptr);
265  if (R0 < thickness) {// if range < thickness, particle stops: Eres=0
266  fRangeOfLastDE = R0;
267  return 0.0;
268  }
269  else
271 
272  // calculate energy after absorber - invert range function to find Eres corresponding to (R0 - thickness)
273  R0 -= thickness;
274 
275  // invert range function to get energy after absorber
277  return static_cast<KVedaLossInverseRangeFunction*>(fInvRange[RF_Z])->GetEnergyPerNucleon(R0, riso) * RF_A;
278  }
279  return fRange->GetX(R0);
280 }
281 
282 
283 
284 
289 
291 {
292  // Return function giving range (in \f$g/cm^2\f$) as a function of energy (in MeV) for
293  // charged particles \f$Z,A\f$ in this material.
294  // \param isoAmat If required, the isotopic mass of the material can be given.
295 
296  RF_Z = Z;
297  RF_A = A;
298  // get parameters for this Z
299  par = &fCoeff[Z - 1];
300  // set up polynomial
301  Double_t x1 = TMath::Log(0.1);
302  Double_t x2 = TMath::Log(0.2);
303  ran = 0.0;
304  for (int j = 2; j < 7; j++)
305  ran += (*par)[j + 1] * TMath::Power(x2, (Double_t)(j - 1));
306  ran += (*par)[2];
307  Double_t y2 = ran;
308  ran = 0.0;
309  for (int jj = 2; jj < 7; jj++)
310  ran += (*par)[jj + 1] * TMath::Power(x1, (Double_t)(jj - 1));
311  ran += (*par)[2];
312  Double_t y1 = ran;
313  adm = (y2 - y1) / (x2 - x1);
314  adn = (y1 - adm * x1);
315  riso = RF_A / (*par)[1];
316  if (isoAmat > 0.0) riso *= (isoAmat / fAmat);
317 
318  fRange->SetRange(0., GetEmaxValid(Z, A));
319 
320  /*
321  * New inversion of range-energy curve (if KVedaLoss::fgNewRangeInversion=kTRUE)
322  */
324  if (!fInvRange[Z]) {
326  }
327  }
328 
329  return fRange;
330 }
331 
332 
333 
338 
340 {
341  // Return function giving stopping power (in \f$MeV/(g/cm^2)\f$) as a function of energy (in MeV) for
342  // charged particles \f$Z,A\f$ in this material.
343  // \param isoAmat If required, the isotopic mass of the material can be given.
344 
345  RF_Z = Z;
346  RF_A = A;
347  // get parameters for this Z
348  par = &fCoeff[Z - 1];
349  // set up polynomial
350  Double_t x1 = TMath::Log(0.1);
351  Double_t x2 = TMath::Log(0.2);
352  ran = 0.0;
353  for (int j = 2; j < 7; j++)
354  ran += (*par)[j + 1] * TMath::Power(x2, (Double_t)(j - 1));
355  ran += (*par)[2];
356  Double_t y2 = ran;
357  ran = 0.0;
358  for (int jj = 2; jj < 7; jj++)
359  ran += (*par)[jj + 1] * TMath::Power(x1, (Double_t)(jj - 1));
360  ran += (*par)[2];
361  Double_t y1 = ran;
362  adm = (y2 - y1) / (x2 - x1);
363  adn = (y1 - adm * x1);
364  riso = RF_A / (*par)[1];
365  if (isoAmat > 0.0) riso *= (isoAmat / fAmat);
366  fStopping->SetRange(0., GetEmaxValid(Z, A));
367  return fStopping;
368 }
369 
370 
371 
376 
378 {
379  // Function parameterising the range of charged particles in this material.
380  // \param E[0] energy is given in MeV.
381  // \returns range calculated in units of \f$g/cm^2\f$
382 
383  eps = E[0] / RF_A;
384  dleps = TMath::Log(eps);
385  if (eps < 0.1)
386  ran = adm * dleps + adn;
387  else {
388  DLEP = dleps;
389  ran = (*par)[2] + (*par)[3] * DLEP;
390  ran += (*par)[4] * (DLEP *= dleps);
391  ran += (*par)[5] * (DLEP *= dleps);
392  ran += (*par)[6] * (DLEP *= dleps);
393  ran += (*par)[7] * (DLEP *= dleps);
394  }
395 
396  // range in g/cm**2
397  return riso * TMath::Exp(ran) * KVUnits::mg;
398 }
399 
400 
401 
406 
408 {
409  // Function parameterising the stopping power of charged particles in this material.
410  // \param E[0] energy is given in MeV.
411  // \returns stopping power calculated in units of \f$Mev/(g/cm^2)\f$
412 
413  eps = E[0] / RF_A;
414  dleps = TMath::Log(eps);
415  if (eps < 0.1) {
416  ran = adm * dleps + adn;
417  drande = E[0] / (riso * TMath::Exp(ran) * KVUnits::mg) / adm;
418  return drande;
419  }
420  DLEP = dleps;
421  ran = (*par)[2] + (*par)[3] * DLEP;
422  drande = (*par)[3];
423  for (int i = 4; i < 8; i++) {
424  drande += (i - 2) * (*par)[i] * DLEP;
425  ran += (*par)[i] * (DLEP *= dleps);
426  }
427  // range in g/cm**2
428  return E[0] / (riso * TMath::Exp(ran) * KVUnits::mg) / drande;
429 }
430 
431 
432 
438 
440 {
441  // Return function giving energy loss (in MeV) as a function of incident energy (in MeV) for
442  // charged particles (Z,A) traversing (or not) material
443  // \param e thickness (in \f$g/cm^2\f$) of this material.
444  // \param isoAmat If required, the isotopic mass of the material can be given.
445 
446  GetRangeFunction(Z, A, isoAmat);
447  thickness = e;
448  fDeltaE->SetRange(0., GetEmaxValid(Z, A));
449  return fDeltaE;
450 }
451 
452 
453 
459 
461 {
462  // Return function giving residual energy (in MeV) as a function of incident energy (in MeV) for
463  // charged particles (Z,A) traversing (or not) material
464  // \param e thickness (in \f$g/cm^2\f$) of this material.
465  // \param isoAmat If required, the isotopic mass of the material can be given.
466 
467  GetRangeFunction(Z, A, isoAmat);
468  thickness = e;
469  fEres->SetRange(0., GetEmaxValid(Z, A));
470  return fEres;
471 }
472 
473 
474 
478 
480 {
481  // \returns range (in \f$g/cm^2\f$) of ion (Z,A) with energy E (MeV) in material.
482  // \param isoAmat change default (isotopic) mass of material,
483 
484  if (Z == 0) return 0.0; //only charged particles
485  TF1* f = GetRangeFunction(Z, A, isoAmat);
486  return f->Eval(E);
487 }
488 
489 
490 
496 
498 {
499  // \returns energy lost (in MeV) by ion (Z,A)
500  // \param E incident energy (MeV)
501  // \param e thickness (in \f$g/cm^2\f$)
502  // \param isoAmat change default (isotopic) mass of material,
503 
504  if (Z == 0) return 0.0; //only charged particles
505  TF1* f = GetDeltaEFunction(e, Z, A, isoAmat);
506  return f->Eval(E);
507 }
508 
509 
510 
516 
518  Double_t isoAmat)
519 {
520  // \returns energy lost (in MeV) by ion (Z,A)
521  // \param E incident energy (MeV)
522  // \param e thickness (in \f$g/cm^2\f$)
523  // \param isoAmat change default (isotopic) mass of material,
524 
525  if (Z == 0) return 0.0; //only charged particles
526  TF1* f = GetEResFunction(e, Z, A, isoAmat);
527  return f->Eval(E);
528 }
529 
530 
531 
537 
539 {
540  // Calculate incident energy (in MeV) for ion (Z,A) for which the range is equal to the
541  // given thickness e (in \f$g/cm^2\f$). At this energy the residual energy of the ion is (just) zero,
542  // for all energies above this energy the residual energy is > 0.
543  // \param isoAmat change default (isotopic) mass of material,
544 
545  if (Z == 0) return 0.0; //only charged particles
546  GetRangeFunction(Z, A, isoAmat);
548  return static_cast<KVedaLossInverseRangeFunction*>(fInvRange[Z])->GetEnergyPerNucleon(e, riso) * A;
549  }
550  return fRange->GetX(e);
551 }
552 
553 
554 
558 
560 {
561  // Calculates incident energy (in MeV) of an ion (Z,A) with residual energy Eres (MeV) after thickness e (in \f$g/cm^2\f$).
562  // \param isoAmat change default (isotopic) mass of material,
563 
564  if (Z == 0) return 0.0; //only charged particles
565  GetRangeFunction(Z, A, isoAmat);
566  Double_t R0 = fRange->Eval(Eres) + e;
568  return static_cast<KVedaLossInverseRangeFunction*>(fInvRange[Z])->GetEnergyPerNucleon(R0, riso) * A;
569  }
570  return fRange->GetX(R0);
571 }
572 
573 
574 
581 
582 void KVedaLossMaterial::GetParameters(Int_t Zion, Int_t& Aion, std::vector<Double_t>& rangepar)
583 {
584  // For the given ion atomic number, give the reference mass used and the six
585  // parameters for the range function fit
586  // \param[in] Zion ion atomic number
587  // \param[out] Aion reference mass used by table
588  // \param[out] rangepar vector containing the six parameters for the range function fit
589 
590  rangepar = std::vector<Double_t>(fCoeff[Zion - 1].begin() + 2, fCoeff[Zion - 1].end());
591  Aion = fCoeff[Zion - 1][1];
592 }
593 
594 
595 
int Int_t
#define f(i)
#define e(i)
bool Bool_t
char Char_t
float Float_t
constexpr Bool_t kFALSE
double Double_t
constexpr Bool_t kTRUE
R__EXTERN TEnv * gEnv
winID w
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
Option_t Option_t TPoint TPoint const char y1
char name[80]
char * Form(const char *fmt,...)
virtual const Char_t * GetType() const
Definition: KVBase.h:177
Material for use in energy loss & range calculations.
void AddCompoundElement(Int_t Z, Int_t A, Int_t Natoms)
TF1 * fDeltaE
function parameterising energy loss in material
void AddMixtureElement(Int_t Z, Int_t A, Int_t Natoms, Double_t Proportion)
Double_t fAmat
effective mass number of material
TF1 * fEres
function parameterising residual energy after crossing material
TF1 * fRange
function parameterising range of charged particles in material
TF1 * fStopping
function parameterising stopping power of charged particles in material
Double_t fRangeOfLastDE
range corresponding to last calculated DE
Abstract base class for calculation of range & energy loss of charged particles in matter.
Dedicated optimised inversion of range-energy function for KVedaLoss.
Description of material in the KVedaLoss range table.
Double_t EResFunc(Double_t *, Double_t *)
TObjArray fInvRange
KVedaLossInverseRangeFunction objects.
Float_t GetEminValid(Int_t Z, Int_t A) const
virtual Double_t GetEResOfIon(Int_t Z, Int_t A, Double_t E, Double_t e, Double_t isoAmat=0.)
static KVedaLoss * fgTable
Double_t thickness
in g/cm**2
virtual TF1 * GetDeltaEFunction(Double_t e, Int_t Z, Int_t A, Double_t isoAmat=0)
virtual Double_t GetDeltaEOfIon(Int_t Z, Int_t A, Double_t E, Double_t e, Double_t isoAmat=0.)
Double_t StoppingFunc(Double_t *, Double_t *)
Bool_t ReadRangeTable(FILE *fp)
Double_t RF_Z
internal variables used by RangeFunc/DeltaEFunc
std::vector< Double_t > fEmin
Z-dependent minimum energy/nucleon for calculation to be valid.
virtual TF1 * GetRangeFunction(Int_t Z, Int_t A, Double_t isoAmat=0)
virtual ~KVedaLossMaterial()
Destructor.
Double_t RangeFunc(Double_t *, Double_t *)
std::vector< Double_t > fEmax
Z-dependent maximum energy/nucleon for calculation to be valid.
virtual Double_t GetRangeOfIon(Int_t Z, Int_t A, Double_t E, Double_t isoAmat=0.)
static Bool_t fNoLimits
if kTRUE, ignore max E limit for validity of calculation
virtual Double_t GetPunchThroughEnergy(Int_t Z, Int_t A, Double_t e, Double_t isoAmat=0.)
void GetParameters(Int_t Zion, Int_t &Aion, std::vector< Double_t > &rangepar)
std::vector< Double_t > * par
Double_t DeltaEFunc(Double_t *, Double_t *)
KVedaLossMaterial()
Default constructor.
virtual TF1 * GetStoppingFunction(Int_t Z, Int_t A, Double_t isoAmat=0)
Float_t GetEmaxValid(Int_t Z, Int_t A) const
std::vector< std::vector< Double_t > > fCoeff
parameters for range tables
virtual Double_t GetEIncFromEResOfIon(Int_t Z, Int_t A, Double_t Eres, Double_t e, Double_t isoAmat=0.)
virtual TF1 * GetEResFunction(Double_t e, Int_t Z, Int_t A, Double_t isoAmat=0)
C++ implementation of VEDALOSS stopping power calculation.
Definition: KVedaLoss.h:63
static Bool_t IsUseNewRangeInversion()
Definition: KVedaLoss.h:92
virtual void SetOwner(Bool_t enable=kTRUE)
virtual const char * GetValue(const char *name, const char *dflt) const
virtual Double_t GetMinimumX(Double_t xmin=0, Double_t xmax=0, Double_t epsilon=1.E-10, Int_t maxiter=100, Bool_t logx=false) const
virtual void SetRange(Double_t xmin, Double_t xmax)
virtual void SetNpx(Int_t npx=100)
virtual Double_t GetX(Double_t y, Double_t xmin=0, Double_t xmax=0, Double_t epsilon=1.E-10, Int_t maxiter=100, Bool_t logx=false) const
virtual Double_t Eval(Double_t x, Double_t y=0, Double_t z=0, Double_t t=0) const
virtual Double_t GetMaximumX(Double_t xmin=0, Double_t xmax=0, Double_t epsilon=1.E-10, Int_t maxiter=100, Bool_t logx=false) const
const char * GetName() const override
virtual void Warning(const char *method, const char *msgfmt,...) const
virtual void Error(const char *method, const char *msgfmt,...) const
TLine * line
const long double mg
Definition: KVUnits.h:74
end
Double_t Min(Double_t a, Double_t b)
Double_t Exp(Double_t x)
Double_t Power(Double_t x, Double_t y)
constexpr Double_t E()
Double_t Log(Double_t x)
Double_t Max(Double_t a, Double_t b)
#define R0(v, w, x, y, z, i)
TArc a
ClassImp(TPyArg)