Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TFormulaPrimitive_v5.cxx
Go to the documentation of this file.
1// @(#)root/hist:$Id$
2// Author: Marian Ivanov, 2005
3
4/*************************************************************************
5* Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6* All rights reserved. *
7* *
8* For the licensing terms see $ROOTSYS/LICENSE. *
9* For the list of contributors see $ROOTSYS/README/CREDITS. *
10*************************************************************************/
11
13
14#include "TMath.h"
15#include "TNamed.h"
16#include "TObjArray.h"
17#include "TVirtualMutex.h"
18
19#include <cmath>
20
21#ifdef WIN32
22#pragma optimize("",off)
23#endif
24
26
27
29
30namespace ROOT {
31 namespace v5 {
32
34
35/** \class TFormulaPrimitive TFormulaPrimitive.h "inc/v5/TFormulaPrimitive.h"
36 \ingroup Hist
37The Formula Primitive class
38
39Helper class for TFormula to speed up TFormula evaluation
40TFormula can use all functions registered in the list of TFormulaPrimitives
41User can add new function to the list of primitives
42if FormulaPrimitive with given name is already defined new primitive is ignored
43
44Example:
45
46~~~ {.cpp}
47 TFormulaPrimitive::AddFormula(new TFormulaPrimitive("Pow2","Pow2",TFastFun::Pow2));
48 TF1 f1("f1","Pow2(x)");
49~~~
50
51 - TFormulaPrimitive is used to get direct acces to the function pointers
52 - GenFunc - pointers to the static function
53 - TFunc - pointers to the data member functions
54
55The following sufixes are currently used, to describe function arguments:
56
57 - G - generic layout - pointer to double (arguments), pointer to double (parameters)
58 - 10 - double
59 - 110 - double, double
60 - 1110 - double, double, double
61*/
62
63//______________________________________________________________________________
64// TFormula primitive
65//
67#ifdef R__COMPLETE_MEM_TERMINATION
68namespace {
69 class TFormulaPrimitiveCleanup {
70 TObjArray **fListOfFunctions;
71 public:
72 TFormulaPrimitiveCleanup(TObjArray **functions) : fListOfFunctions(functions) {}
73 ~TFormulaPrimitiveCleanup() {
74 delete *fListOfFunctions;
75 }
76 };
77
78}
79#endif
80
81////////////////////////////////////////////////////////////////////////////////
82/// Default constructor.
83
84TFormulaPrimitive::TFormulaPrimitive() : fFuncG(nullptr), fType(0), fNArguments(0), fNParameters(0), fIsStatic(true) {}
85
86////////////////////////////////////////////////////////////////////////////////
87/// Constructor.
88
89TFormulaPrimitive::TFormulaPrimitive(const char *name,const char *formula,
90 GenFunc0 fpointer) : TNamed(name,formula),
91 fFunc0(fpointer),
92 fType(0),fNArguments(0),fNParameters(0),fIsStatic(kTRUE)
93{
94}
95
96////////////////////////////////////////////////////////////////////////////////
97/// Constructor.
98
99TFormulaPrimitive::TFormulaPrimitive(const char *name,const char *formula,
100 GenFunc10 fpointer) : TNamed(name,formula),
101 fFunc10(fpointer),
102 fType(10),fNArguments(1),fNParameters(0),fIsStatic(kTRUE)
103{
104}
105
106////////////////////////////////////////////////////////////////////////////////
107/// Constructor.
108
109TFormulaPrimitive::TFormulaPrimitive(const char *name,const char *formula,
110 GenFunc110 fpointer) : TNamed(name,formula),
111 fFunc110(fpointer),
112 fType(110),fNArguments(2),fNParameters(0),fIsStatic(kTRUE)
113{
114}
115
116////////////////////////////////////////////////////////////////////////////////
117/// Constructor.
118
119TFormulaPrimitive::TFormulaPrimitive(const char *name,const char *formula,
120 GenFunc1110 fpointer) : TNamed(name,formula),
121 fFunc1110(fpointer),
122 fType(1110),fNArguments(3),fNParameters(0),fIsStatic(kTRUE)
123{
124}
125
126////////////////////////////////////////////////////////////////////////////////
127/// Constructor.
128
129TFormulaPrimitive::TFormulaPrimitive(const char *name,const char *formula,
130 GenFuncG fpointer,Int_t npar) : TNamed(name,formula),
131 fFuncG(fpointer),
132 fType(-1),fNArguments(2),fNParameters(npar),fIsStatic(kTRUE)
133{
134}
135
136////////////////////////////////////////////////////////////////////////////////
137/// Constructor.
138
139TFormulaPrimitive::TFormulaPrimitive(const char *name,const char *formula,
140 TFuncG fpointer) : TNamed(name,formula),
141 fTFuncG(fpointer),
142 fType(0),fNArguments(0),fNParameters(0),fIsStatic(kFALSE)
143{
144}
145
146////////////////////////////////////////////////////////////////////////////////
147/// Constructor.
148
149TFormulaPrimitive::TFormulaPrimitive(const char *name,const char *formula,
150 TFunc0 fpointer) : TNamed(name,formula),
151 fTFunc0(fpointer),
152 fType(0),fNArguments(0),fNParameters(0),fIsStatic(kFALSE)
153{
154}
155
156////////////////////////////////////////////////////////////////////////////////
157/// Constructor.
158
159TFormulaPrimitive::TFormulaPrimitive(const char *name,const char *formula,
160 TFunc10 fpointer) : TNamed(name,formula),
161 fTFunc10(fpointer),
162 fType(-10),fNArguments(1),fNParameters(0),fIsStatic(kFALSE)
163{
164}
165
166////////////////////////////////////////////////////////////////////////////////
167/// Constructor.
168
169TFormulaPrimitive::TFormulaPrimitive(const char *name,const char *formula,
170 TFunc110 fpointer) : TNamed(name,formula),
171 fTFunc110(fpointer),
172 fType(-110),fNArguments(2),fNParameters(0),fIsStatic(kFALSE)
173{
174}
175
176////////////////////////////////////////////////////////////////////////////////
177/// Constructor.
178
179TFormulaPrimitive::TFormulaPrimitive(const char *name,const char *formula,
180 TFunc1110 fpointer) :TNamed(name,formula),
181 fTFunc1110(fpointer),
182 fType(-1110),fNArguments(3),fNParameters(0),fIsStatic(kFALSE)
183{
184}
185
186////////////////////////////////////////////////////////////////////////////////
187/// Add formula to the list of primitive formulas.
188/// If primitive formula already defined do nothing.
189
191{
193 if (fgListOfFunction == nullptr) BuildBasicFormulas();
194 if (FindFormula(formula->GetName(),formula->fNArguments)){
195 delete formula;
196 return 0;
197 }
198 fgListOfFunction->AddLast(formula);
199 return 1;
200}
201
202////////////////////////////////////////////////////////////////////////////////
203/// Eval primitive function at point x.
204
206{
207 if (fIsStatic == kFALSE) return 0;
208
209 if (fType==0) return fFunc0();
210 if (fType==10) {
211 return fFunc10(x[0]);
212 }
213 if (fType==110) {
214 return fFunc110(x[0],x[1]);
215 }
216 if (fType==1110) {
217 return fFunc1110(x[0],x[1],x[2]);
218 }
219 return 0;
220}
221
222////////////////////////////////////////////////////////////////////////////////
223/// Eval member function of object o at point x.
224
226{
227 if (fIsStatic == kTRUE) return 0;
228 if (fType== 0) return (*o.*fTFunc0)();
229 if (fType==-10) return (*o.*fTFunc10)(*x);
230 if (fType==-110) return (*o.*fTFunc110)(x[0],x[1]);
231 if (fType==-1110) return (*o.*fTFunc1110)(x[0],x[1],x[2]);
232 return 0;
233}
234
235////////////////////////////////////////////////////////////////////////////////
236/// Eval primitive parametric function.
237
239{
240 return fFuncG(x,param);
241}
242
243#define RTFastFun__POLY(var) \
244{ \
245 Double_t res= param[var-1]+param[var]*x[0]; \
246 for (Int_t j=var-1 ;j>0;j--) res = param[j-1]+x[0]*res; \
247 return res; \
248}
249
250namespace TFastFun {
251 //
252 // Namespace with basic primitive functions registered by TFormulaPrimitive
253 // all function registered by TFormulaPrimitive can be used in TFormula
254 //
259 inline Double_t FPoln(Double_t *x, Double_t *param, Int_t npar);
260 Double_t FPol0(Double_t * /*x*/, Double_t *param){ return param[0];}
261 Double_t FPol1(Double_t *x, Double_t *param){ return param[0]+param[1]*x[0];}
262 Double_t FPol2(Double_t *x, Double_t *param){ return param[0]+x[0]*(param[1]+param[2]*x[0]);}
263 Double_t FPol3(Double_t *x, Double_t *param){ return param[0]+x[0]*(param[1]+x[0]*(param[2]+param[3]*x[0]));}
271 //
272 //
285 Double_t Sqrt(Double_t x) {return x>0?sqrt(x):0;}
286 //
287 Double_t Sign(Double_t x){return (x<0)? -1:1;}
290 //logical
291 Double_t XandY(Double_t x, Double_t y){ return (x*y>0.1);}
292 Double_t XorY(Double_t x, Double_t y) { return (x+y>0.1);}
299 Double_t XNot(Double_t x){ return (x<0.1);}
300};
301
302////////////////////////////////////////////////////////////////////////////////
303/// Find the formula in the list of formulas.
304
306{
308 if (!fgListOfFunction) {
310 }
311 Int_t nobjects = fgListOfFunction->GetEntries();
312 for (Int_t i = 0; i < nobjects; ++i) {
314 if (formula && 0==strcmp(name, formula->GetName())) return formula;
315 }
316 return nullptr;
317}
318
319////////////////////////////////////////////////////////////////////////////////
320/// Find the formula in the list of formulas.
321
323{
325 if (!fgListOfFunction) {
327 }
328 Int_t nobjects = fgListOfFunction->GetEntries();
329 for (Int_t i = 0; i < nobjects; ++i) {
331 if (prim) {
332 bool match = ( ((UInt_t)prim->fNArguments) == nargs );
333 if (match && 0==strcmp(name, prim->GetName())) return prim;
334 }
335 }
336 return nullptr;
337}
338
339////////////////////////////////////////////////////////////////////////////////
340/// Find the formula in the list of formulas.
341
343{
344 // let's count the argument(s)
345 if (args) {
346 Int_t nargs = 0;
347 if (args[0]!=')') {
348 nargs = 1;
349 int nest = 0;
350 for(UInt_t c = 0; c < strlen(args); ++c ) {
351 switch (args[c]) {
352 case '(': ++nest; break;
353 case ')': --nest; break;
354 case '<': ++nest; break;
355 case '>': --nest; break;
356 case ',': nargs += (nest==0); break;
357 }
358 }
359 }
360 return FindFormula(name,nargs);
361 } else {
362 return FindFormula(name);
363 }
364 return nullptr;
365}
366
367////////////////////////////////////////////////////////////////////////////////
368/// FPoln.
369
371{
372 Double_t res = 0; Double_t temp=1;
373 for (Int_t j=npar ;j>=0;j--) {
374 res += temp*param[j];
375 temp *= *x;
376 }
377 return res;
378}
379
380////////////////////////////////////////////////////////////////////////////////
381/// Gauss.
382
384{
385 if (sigma == 0) return 1.e30;
386 Double_t arg = (x-mean)/sigma;
387 return TMath::Exp(-0.5*arg*arg);
388}
389
390////////////////////////////////////////////////////////////////////////////////
391/// Normalize gauss.
392
394{
395 if (sigma == 0) return 0;
396 Double_t arg = (x-mean)/sigma;
397 return TMath::Exp(-0.5*arg*arg)/(2.50662827463100024*sigma); //sqrt(2*Pi)=2.50662827463100024
398}
399
400////////////////////////////////////////////////////////////////////////////////
401/// Built-in functions.
402
404{
406 if (fgListOfFunction==nullptr) {
407 fgListOfFunction = new TObjArray(1000);
409 }
410#ifdef R__COMPLETE_MEM_TERMINATION
411 static TFormulaPrimitiveCleanup gCleanup(&fgListOfFunction);
412#endif
413
414 //
415 // logical
416 //
417 AddFormula(new TFormulaPrimitive("XandY","XandY",TFastFun::XandY));
418 AddFormula(new TFormulaPrimitive("XorY","XorY",TFastFun::XorY));
419 AddFormula(new TFormulaPrimitive("XNot","XNot",TFastFun::XNot));
421 AddFormula(new TFormulaPrimitive("XleY","XleY",TFastFun::XleY));
423 AddFormula(new TFormulaPrimitive("XgeY","XgeY",TFastFun::XgeY));
425 AddFormula(new TFormulaPrimitive("XneY","XneY",TFastFun::XneY));
426 // addition + multiplication
427 AddFormula(new TFormulaPrimitive("PlusXY","PlusXY",TFastFun::PlusXY));
428 AddFormula(new TFormulaPrimitive("MinusXY","MinusXY",TFastFun::MinusXY));
429 AddFormula(new TFormulaPrimitive("MultXY","MultXY",TFastFun::MultXY));
430 AddFormula(new TFormulaPrimitive("DivXY","DivXY",TFastFun::DivXY));
431 AddFormula(new TFormulaPrimitive("XpYpZ","XpYpZ",TFastFun::XpYpZ));
432 AddFormula(new TFormulaPrimitive("XxYxZ","XxYxZ",TFastFun::XxYxZ));
433 AddFormula(new TFormulaPrimitive("XxYpZ","XxYpZ",TFastFun::XxYpZ));
434 AddFormula(new TFormulaPrimitive("XpYxZ","XpYxZ",TFastFun::XpYxZ));
435 //
436 //
437 AddFormula(new TFormulaPrimitive("Gaus","Gaus",TFastFun::Gaus));
438 AddFormula(new TFormulaPrimitive("Gausn","Gausn",TFastFun::Gausn));
439 AddFormula(new TFormulaPrimitive("Landau","Landau",TFastFun::Landau));
440 AddFormula(new TFormulaPrimitive("Landaun","Landaun",TFastFun::Landaun));
441 //
442 //
443 // polynoms
444 //
445 //
456 AddFormula(new TFormulaPrimitive("Pol10","Pol10",(GenFuncG)TFastFun::FPol10,11));
457 //
458 // pows
459 AddFormula(new TFormulaPrimitive("Pow2","Pow2",TFastFun::Pow2));
460 AddFormula(new TFormulaPrimitive("Pow3","Pow3",TFastFun::Pow3));
461 AddFormula(new TFormulaPrimitive("Pow4","Pow4",TFastFun::Pow4));
462 AddFormula(new TFormulaPrimitive("Pow5","Pow5",TFastFun::Pow5));
463 //
464 //
465 AddFormula(new TFormulaPrimitive("TMath::Cos","TMath::Cos",cos)); // 10
466 AddFormula(new TFormulaPrimitive("cos","cos",cos)); // 10
467 AddFormula(new TFormulaPrimitive("TMath::Sin","TMath::Sin",sin)); // 11
468 AddFormula(new TFormulaPrimitive("sin","sin",sin)); // 11
469 AddFormula(new TFormulaPrimitive("TMath::Tan","TMath::Tan",tan)); // 12
470 AddFormula(new TFormulaPrimitive("tan","tan",tan)); // 12
471 AddFormula(new TFormulaPrimitive("TMath::ACos","TMath::ACos",acos)); // 13
472 AddFormula(new TFormulaPrimitive("acos","acos",acos)); // 13
473 AddFormula(new TFormulaPrimitive("TMath::ASin","TMath::ASin",asin)); // 14
474 AddFormula(new TFormulaPrimitive("asin","asin",asin)); // 14
475 AddFormula(new TFormulaPrimitive("TMath::ATan","TMath::ATan",atan)); // 15
476 AddFormula(new TFormulaPrimitive("atan","atan",atan)); // 15
477 AddFormula(new TFormulaPrimitive("TMath::ATan2","TMath::ATan2",atan2)); // 16
478 AddFormula(new TFormulaPrimitive("atan2","atan2",atan2)); // 16
479 // kpow = 20, ksq = 21, ksqrt = 22,
480 AddFormula(new TFormulaPrimitive("pow","pow",TMath::Power)); //20
481 AddFormula(new TFormulaPrimitive("sq","sq",TFastFun::Pow2)); //21
482 AddFormula(new TFormulaPrimitive("sqrt","sqrt",TFastFun::Sqrt)); //22
483 // kmin = 24, kmax = 25,
484 AddFormula(new TFormulaPrimitive("min","min",(GenFunc110)TMath::Min)); //24
485 AddFormula(new TFormulaPrimitive("max","max",(GenFunc110)TMath::Max)); //25
486 // klog = 30, kexp = 31, klog10 = 32,
487 AddFormula(new TFormulaPrimitive("log","log",TMath::Log)); //30
488 AddFormula(new TFormulaPrimitive("exp","exp",TMath::Exp)); //31
489 AddFormula(new TFormulaPrimitive("log10","log10",TMath::Log10)); //32
490 //
491 // cosh 70 acosh 73
492 // sinh 71 asinh 74
493 // tanh 72 atanh 75
494 //
495 AddFormula(new TFormulaPrimitive("TMath::CosH","TMath::Cosh",cosh)); // 70
496 AddFormula(new TFormulaPrimitive("cosh","cosh",cosh)); // 70
497 AddFormula(new TFormulaPrimitive("TMath::SinH","TMath::SinH",sinh)); // 71
498 AddFormula(new TFormulaPrimitive("sinh","sinh",sinh)); // 71
499 AddFormula(new TFormulaPrimitive("TMath::TanH","TMath::Tanh",tanh)); // 72
500 AddFormula(new TFormulaPrimitive("tanh","tanh",tanh)); // 72
501 AddFormula(new TFormulaPrimitive("TMath::ACosH","TMath::ACosh",TMath::ACosH)); // 73
502 AddFormula(new TFormulaPrimitive("acosh","acosH",TMath::ACosH)); // 73
503 AddFormula(new TFormulaPrimitive("TMath::ASinH","TMath::ASinh",TMath::ASinH)); // 74
504 AddFormula(new TFormulaPrimitive("acosh","acosH",TMath::ASinH)); // 74
505 AddFormula(new TFormulaPrimitive("TMath::ATanH","TMath::ATanh",TMath::ATanH)); // 75
506 AddFormula(new TFormulaPrimitive("atanh","atanh",TMath::ATanH)); // 75
507 //
508 AddFormula(new TFormulaPrimitive("TMath::Abs","TMath::Abs",TMath::Abs));
509 AddFormula(new TFormulaPrimitive("TMath::BreitWigner","TMath::BreitWigner",TMath::BreitWigner));
510
511 //Disable direct access to TMath::Landau for now because of the default parameter.
512 //AddFormula(new TFormulaPrimitive("TMath::Landau","TMath::Landau",(TFormulaPrimitive::GenFunc1110)TMath::Landau));
513
515 return 1;
516}
517
518 } // end namespace v5
519
520} // end namespace ROOT
#define c(i)
Definition RSha256.hxx:101
int Int_t
Definition RtypesCore.h:45
unsigned int UInt_t
Definition RtypesCore.h:46
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
double Double_t
Definition RtypesCore.h:59
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
#define ClassImp(name)
Definition Rtypes.h:382
static TVirtualMutex * gTFormulaPrimativeListMutex
#define RTFastFun__POLY(var)
char name[80]
Definition TGX11.cxx:110
#define R__LOCKGUARD2(mutex)
The Formula Primitive class.
GenFunc1110 fFunc1110
pointer to the function
static Int_t AddFormula(TFormulaPrimitive *formula)
Add formula to the list of primitive formulas.
Double_t(* GenFuncG)(const Double_t *, const Double_t *)
TFormulaPrimitive()
Default constructor.
TFunc110 fTFunc110
pointer to member function
GenFunc10 fFunc10
pointer to the function
static TFormulaPrimitive * FindFormula(const char *name)
Find the formula in the list of formulas.
TFunc10 fTFunc10
pointer to member function
GenFunc110 fFunc110
pointer to the function
TFunc0 fTFunc0
pointer to the TFormula generic function
static Int_t BuildBasicFormulas()
list of global primitive formulas
GenFunc0 fFunc0
pointer to the TFormula generic function
Double_t Eval(Double_t *x)
Eval primitive function at point x.
Double_t(* GenFunc110)(Double_t, Double_t)
TFunc1110 fTFunc1110
pointer to member function
static TObjArray * fgListOfFunction
virtual void SetOwner(Bool_t enable=kTRUE)
Set whether this collection is the owner (enable==true) of its content.
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
An array of TObjects.
Definition TObjArray.h:31
Int_t GetEntries() const override
Return the number of objects in array (i.e.
TObject * At(Int_t idx) const override
Definition TObjArray.h:164
void AddLast(TObject *obj) override
Add object in the next empty slot in the array.
Mother of all ROOT objects.
Definition TObject.h:41
This class implements a mutex interface.
const Double_t sigma
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
Double_t Landau(Double_t x, Double_t mean, Double_t sigma)
Double_t Pow3(Double_t x)
Double_t XorY(Double_t x, Double_t y)
Double_t Abs(Double_t x)
Double_t XpYpZ(Double_t x, Double_t y, Double_t z)
Double_t Pow4(Double_t x)
Double_t FPol7(Double_t *x, Double_t *param)
Double_t XneY(Double_t x, Double_t y)
Double_t Pow2(Double_t x)
Double_t FPol0(Double_t *, Double_t *param)
Double_t XlY(Double_t x, Double_t y)
Double_t XxYpZ(Double_t x, Double_t y, Double_t z)
Double_t FPoln(Double_t *x, Double_t *param, Int_t npar)
FPoln.
Double_t XNot(Double_t x)
Double_t XxYxZ(Double_t x, Double_t y, Double_t z)
Double_t FPol1(Double_t *x, Double_t *param)
Double_t Nint(Double_t x)
Double_t XgeY(Double_t x, Double_t y)
Double_t XpYxZ(Double_t x, Double_t y, Double_t z)
Double_t Gausn(Double_t x, Double_t mean, Double_t sigma)
Normalize gauss.
Double_t PlusXY(Double_t x, Double_t y)
Double_t Gaus(Double_t x, Double_t mean, Double_t sigma)
Gauss.
Double_t MinusXY(Double_t x, Double_t y)
Double_t Landaun(Double_t x, Double_t mean, Double_t sigma)
Double_t FPol4(Double_t *x, Double_t *param)
Double_t XandY(Double_t x, Double_t y)
Double_t XeY(Double_t x, Double_t y)
Double_t XleY(Double_t x, Double_t y)
Double_t Sqrt(Double_t x)
Double_t FPol9(Double_t *x, Double_t *param)
Double_t FPol6(Double_t *x, Double_t *param)
Double_t FPol10(Double_t *x, Double_t *param)
Double_t DivXY(Double_t x, Double_t y)
Double_t FPol3(Double_t *x, Double_t *param)
Double_t XgY(Double_t x, Double_t y)
Double_t Pow5(Double_t x)
Double_t FPol2(Double_t *x, Double_t *param)
Double_t FPol8(Double_t *x, Double_t *param)
Double_t Sign(Double_t x)
Double_t FPol5(Double_t *x, Double_t *param)
Double_t MultXY(Double_t x, Double_t y)
void TMath_GenerInterface()
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...
Int_t Nint(T x)
Round to nearest integer. Rounds half integers to the nearest even integer.
Definition TMath.h:697
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:250
Double_t Exp(Double_t x)
Returns the base-e exponential function of x, which is e raised to the power x.
Definition TMath.h:713
Double_t ASinH(Double_t)
Returns the area hyperbolic sine of x.
Definition TMath.cxx:67
Double_t BreitWigner(Double_t x, Double_t mean=0, Double_t gamma=1)
Calculates a Breit Wigner function with mean and gamma.
Definition TMath.cxx:442
Double_t Landau(Double_t x, Double_t mpv=0, Double_t sigma=1, Bool_t norm=kFALSE)
The LANDAU function.
Definition TMath.cxx:492
Double_t ACosH(Double_t)
Returns the nonnegative area hyperbolic cosine of x.
Definition TMath.cxx:81
Double_t Log(Double_t x)
Returns the natural logarithm of x.
Definition TMath.h:760
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Definition TMath.h:725
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:198
Double_t ATanH(Double_t)
Returns the area hyperbolic tangent of x.
Definition TMath.cxx:95
Double_t Log10(Double_t x)
Returns the common (base-10) logarithm of x.
Definition TMath.h:766
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123