Logo ROOT   6.07/09
Reference Guide
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 
12 #include "v5/TFormulaPrimitive.h"
13 
14 #include "TMath.h"
15 #include "TNamed.h"
16 #include "TObjArray.h"
17 #include "TVirtualMutex.h"
18 
19 #include <math.h>
20 
21 #ifdef WIN32
22 #pragma optimize("",off)
23 #endif
24 
26 
27 
29 
30 namespace ROOT {
31  namespace v5 {
32 
33  void TMath_GenerInterface();
34 
35 /** \class TFormulaPrimitive TFormulaPrimitive.h "inc/v5/TFormulaPrimitive.h"
36  \ingroup Hist
37 The Formula Primitive class
38 
39 Helper class for TFormula to speed up TFormula evaluation
40 TFormula can use all functions registered in the list of TFormulaPrimitives
41 User can add new function to the list of primitives
42 if FormulaPrimitive with given name is already defined new primitive is ignored
43 
44 Example:
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 
55 The 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 //
66 TObjArray * TFormulaPrimitive::fgListOfFunction = 0;
67 #ifdef R__COMPLETE_MEM_TERMINATION
68 namespace {
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 
84 TFormulaPrimitive::TFormulaPrimitive() : TNamed(),
85  fFuncG(0),
86  fType(0),fNArguments(0),fNParameters(0),fIsStatic(kTRUE)
87 {
88 }
89 
90 ////////////////////////////////////////////////////////////////////////////////
91 /// Constructor.
92 
93 TFormulaPrimitive::TFormulaPrimitive(const char *name,const char *formula,
94  GenFunc0 fpointer) : TNamed(name,formula),
95  fFunc0(fpointer),
97 {
98 }
99 
100 ////////////////////////////////////////////////////////////////////////////////
101 /// Constructor.
102 
103 TFormulaPrimitive::TFormulaPrimitive(const char *name,const char *formula,
104  GenFunc10 fpointer) : TNamed(name,formula),
105  fFunc10(fpointer),
107 {
108 }
109 
110 ////////////////////////////////////////////////////////////////////////////////
111 /// Constructor.
112 
113 TFormulaPrimitive::TFormulaPrimitive(const char *name,const char *formula,
114  GenFunc110 fpointer) : TNamed(name,formula),
115  fFunc110(fpointer),
117 {
118 }
119 
120 ////////////////////////////////////////////////////////////////////////////////
121 /// Constructor.
122 
123 TFormulaPrimitive::TFormulaPrimitive(const char *name,const char *formula,
124  GenFunc1110 fpointer) : TNamed(name,formula),
125  fFunc1110(fpointer),
127 {
128 }
129 
130 ////////////////////////////////////////////////////////////////////////////////
131 /// Constructor.
132 
133 TFormulaPrimitive::TFormulaPrimitive(const char *name,const char *formula,
134  GenFuncG fpointer,Int_t npar) : TNamed(name,formula),
135  fFuncG(fpointer),
137 {
138 }
139 
140 ////////////////////////////////////////////////////////////////////////////////
141 /// Constructor.
142 
143 TFormulaPrimitive::TFormulaPrimitive(const char *name,const char *formula,
144  TFuncG fpointer) : TNamed(name,formula),
145  fTFuncG(fpointer),
147 {
148 }
149 
150 ////////////////////////////////////////////////////////////////////////////////
151 /// Constructor.
152 
153 TFormulaPrimitive::TFormulaPrimitive(const char *name,const char *formula,
154  TFunc0 fpointer) : TNamed(name,formula),
155  fTFunc0(fpointer),
157 {
158 }
159 
160 ////////////////////////////////////////////////////////////////////////////////
161 /// Constructor.
162 
163 TFormulaPrimitive::TFormulaPrimitive(const char *name,const char *formula,
164  TFunc10 fpointer) : TNamed(name,formula),
165  fTFunc10(fpointer),
167 {
168 }
169 
170 ////////////////////////////////////////////////////////////////////////////////
171 /// Constructor.
172 
173 TFormulaPrimitive::TFormulaPrimitive(const char *name,const char *formula,
174  TFunc110 fpointer) : TNamed(name,formula),
175  fTFunc110(fpointer),
177 {
178 }
179 
180 ////////////////////////////////////////////////////////////////////////////////
181 /// Constructor.
182 
183 TFormulaPrimitive::TFormulaPrimitive(const char *name,const char *formula,
184  TFunc1110 fpointer) :TNamed(name,formula),
185  fTFunc1110(fpointer),
187 {
188 }
189 
190 ////////////////////////////////////////////////////////////////////////////////
191 /// Add formula to the list of primitive formulas.
192 /// If primitive formula already defined do nothing.
193 
195 {
196  R__LOCKGUARD2(gTFormulaPrimativeListMutex);
198  if (FindFormula(formula->GetName(),formula->fNArguments)){
199  delete formula;
200  return 0;
201  }
202  fgListOfFunction->AddLast(formula);
203  return 1;
204 }
205 
206 ////////////////////////////////////////////////////////////////////////////////
207 /// Eval primitive function at point x.
208 
210 {
211  if (fIsStatic == kFALSE) return 0;
212 
213  if (fType==0) return fFunc0();
214  if (fType==10) {
215  return fFunc10(x[0]);
216  }
217  if (fType==110) {
218  return fFunc110(x[0],x[1]);
219  }
220  if (fType==1110) {
221  return fFunc1110(x[0],x[1],x[2]);
222  }
223  return 0;
224 }
225 
226 ////////////////////////////////////////////////////////////////////////////////
227 /// Eval member function of object o at point x.
228 
230 {
231  if (fIsStatic == kTRUE) return 0;
232  if (fType== 0) return (*o.*fTFunc0)();
233  if (fType==-10) return (*o.*fTFunc10)(*x);
234  if (fType==-110) return (*o.*fTFunc110)(x[0],x[1]);
235  if (fType==-1110) return (*o.*fTFunc1110)(x[0],x[1],x[2]);
236  return 0;
237 }
238 
239 ////////////////////////////////////////////////////////////////////////////////
240 /// Eval primitive parametric function.
241 
243 {
244  return fFuncG(x,param);
245 }
246 
247 #define RTFastFun__POLY(var) \
248 { \
249  Double_t res= param[var-1]+param[var]*x[0]; \
250  for (Int_t j=var-1 ;j>0;j--) res = param[j-1]+x[0]*res; \
251  return res; \
252 }
253 
254 namespace TFastFun {
255  //
256  // Namespace with basic primitive functions registered by TFormulaPrimitive
257  // all function registered by TFormulaPrimitive can be used in TFormula
258  //
259  Double_t Pow2(Double_t x){return x*x;}
260  Double_t Pow3(Double_t x){return x*x*x;}
261  Double_t Pow4(Double_t x){return x*x*x*x;}
262  Double_t Pow5(Double_t x){return x*x*x*x*x;}
263  inline Double_t FPoln(Double_t *x, Double_t *param, Int_t npar);
264  Double_t FPol0(Double_t * /*x*/, Double_t *param){ return param[0];}
265  Double_t FPol1(Double_t *x, Double_t *param){ return param[0]+param[1]*x[0];}
266  Double_t FPol2(Double_t *x, Double_t *param){ return param[0]+x[0]*(param[1]+param[2]*x[0]);}
267  Double_t FPol3(Double_t *x, Double_t *param){ return param[0]+x[0]*(param[1]+x[0]*(param[2]+param[3]*x[0]));}
275  //
276  //
280  Double_t DivXY(Double_t x, Double_t y){return TMath::Abs(y)>0 ? x/y:0;}
283  Double_t XxYpZ(Double_t x, Double_t y, Double_t z){ return x*(y+z);}
284  Double_t XpYxZ(Double_t x, Double_t y, Double_t z){ return x+(y*z);}
286  Double_t Gausn(Double_t x, Double_t mean, Double_t sigma);
287  Double_t Landau(Double_t x, Double_t mean, Double_t sigma){return TMath::Landau(x,mean,sigma,kFALSE);}
288  Double_t Landaun(Double_t x, Double_t mean, Double_t sigma){return TMath::Landau(x,mean,sigma,kTRUE);}
289  Double_t Sqrt(Double_t x) {return x>0?sqrt(x):0;}
290  //
291  Double_t Sign(Double_t x){return (x<0)? -1:1;}
294  //logical
295  Double_t XandY(Double_t x, Double_t y){ return (x*y>0.1);}
296  Double_t XorY(Double_t x, Double_t y) { return (x+y>0.1);}
297  Double_t XgY(Double_t x, Double_t y) {return (x>y);}
298  Double_t XgeY(Double_t x, Double_t y) {return (x>=y);}
299  Double_t XlY(Double_t x, Double_t y) {return (x<y);}
300  Double_t XleY(Double_t x, Double_t y) {return (x<=y);}
301  Double_t XeY(Double_t x,Double_t y) {return (x==y);}
302  Double_t XneY(Double_t x,Double_t y) {return (x!=y);}
303  Double_t XNot(Double_t x){ return (x<0.1);}
304 };
305 
306 ////////////////////////////////////////////////////////////////////////////////
307 /// Find the formula in the list of formulas.
308 
310 {
311  R__LOCKGUARD2(gTFormulaPrimativeListMutex);
312  if (!fgListOfFunction) {
314  }
315  Int_t nobjects = fgListOfFunction->GetEntries();
316  for (Int_t i = 0; i < nobjects; ++i) {
318  if (formula && 0==strcmp(name, formula->GetName())) return formula;
319  }
320  return 0;
321 }
322 
323 ////////////////////////////////////////////////////////////////////////////////
324 /// Find the formula in the list of formulas.
325 
327 {
328  R__LOCKGUARD2(gTFormulaPrimativeListMutex);
329  if (!fgListOfFunction) {
331  }
332  Int_t nobjects = fgListOfFunction->GetEntries();
333  for (Int_t i = 0; i < nobjects; ++i) {
335  if (prim) {
336  bool match = ( ((UInt_t)prim->fNArguments) == nargs );
337  if (match && 0==strcmp(name, prim->GetName())) return prim;
338  }
339  }
340  return 0;
341 }
342 
343 ////////////////////////////////////////////////////////////////////////////////
344 /// Find the formula in the list of formulas.
345 
347 {
348  // let's count the argument(s)
349  if (args) {
350  Int_t nargs = 0;
351  if (args[0]!=')') {
352  nargs = 1;
353  int nest = 0;
354  for(UInt_t c = 0; c < strlen(args); ++c ) {
355  switch (args[c]) {
356  case '(': ++nest; break;
357  case ')': --nest; break;
358  case '<': ++nest; break;
359  case '>': --nest; break;
360  case ',': nargs += (nest==0); break;
361  }
362  }
363  }
364  return FindFormula(name,nargs);
365  } else {
366  return FindFormula(name);
367  }
368  return 0;
369 }
370 
371 ////////////////////////////////////////////////////////////////////////////////
372 /// FPoln.
373 
375 {
376  Double_t res = 0; Double_t temp=1;
377  for (Int_t j=npar ;j>=0;j--) {
378  res += temp*param[j];
379  temp *= *x;
380  }
381  return res;
382 }
383 
384 ////////////////////////////////////////////////////////////////////////////////
385 /// Gauss.
386 
388 {
389  if (sigma == 0) return 1.e30;
390  Double_t arg = (x-mean)/sigma;
391  return TMath::Exp(-0.5*arg*arg);
392 }
393 
394 ////////////////////////////////////////////////////////////////////////////////
395 /// Normalize gauss.
396 
398 {
399  if (sigma == 0) return 0;
400  Double_t arg = (x-mean)/sigma;
401  return TMath::Exp(-0.5*arg*arg)/(2.50662827463100024*sigma); //sqrt(2*Pi)=2.50662827463100024
402 }
403 
404 ////////////////////////////////////////////////////////////////////////////////
405 /// Built-in functions.
406 
408 {
409  R__LOCKGUARD2(gTFormulaPrimativeListMutex);
410  if (fgListOfFunction==0) {
411  fgListOfFunction = new TObjArray(1000);
413  }
414 #ifdef R__COMPLETE_MEM_TERMINATION
415  static TFormulaPrimitiveCleanup gCleanup(&fgListOfFunction);
416 #endif
417 
418  //
419  // logical
420  //
421  AddFormula(new TFormulaPrimitive("XandY","XandY",TFastFun::XandY));
422  AddFormula(new TFormulaPrimitive("XorY","XorY",TFastFun::XorY));
423  AddFormula(new TFormulaPrimitive("XNot","XNot",TFastFun::XNot));
424  AddFormula(new TFormulaPrimitive("XlY","XlY",TFastFun::XlY));
425  AddFormula(new TFormulaPrimitive("XleY","XleY",TFastFun::XleY));
426  AddFormula(new TFormulaPrimitive("XgY","XgY",TFastFun::XgY));
427  AddFormula(new TFormulaPrimitive("XgeY","XgeY",TFastFun::XgeY));
428  AddFormula(new TFormulaPrimitive("XeY","XeY",TFastFun::XeY));
429  AddFormula(new TFormulaPrimitive("XneY","XneY",TFastFun::XneY));
430  // addition + multiplication
431  AddFormula(new TFormulaPrimitive("PlusXY","PlusXY",TFastFun::PlusXY));
432  AddFormula(new TFormulaPrimitive("MinusXY","MinusXY",TFastFun::MinusXY));
433  AddFormula(new TFormulaPrimitive("MultXY","MultXY",TFastFun::MultXY));
434  AddFormula(new TFormulaPrimitive("DivXY","DivXY",TFastFun::DivXY));
435  AddFormula(new TFormulaPrimitive("XpYpZ","XpYpZ",TFastFun::XpYpZ));
436  AddFormula(new TFormulaPrimitive("XxYxZ","XxYxZ",TFastFun::XxYxZ));
437  AddFormula(new TFormulaPrimitive("XxYpZ","XxYpZ",TFastFun::XxYpZ));
438  AddFormula(new TFormulaPrimitive("XpYxZ","XpYxZ",TFastFun::XpYxZ));
439  //
440  //
441  AddFormula(new TFormulaPrimitive("Gaus","Gaus",TFastFun::Gaus));
442  AddFormula(new TFormulaPrimitive("Gausn","Gausn",TFastFun::Gausn));
443  AddFormula(new TFormulaPrimitive("Landau","Landau",TFastFun::Landau));
444  AddFormula(new TFormulaPrimitive("Landaun","Landaun",TFastFun::Landaun));
445  //
446  //
447  // polynoms
448  //
449  //
450  AddFormula(new TFormulaPrimitive("Pol0","Pol0",(GenFuncG)TFastFun::FPol0,1));
451  AddFormula(new TFormulaPrimitive("Pol1","Pol1",(GenFuncG)TFastFun::FPol1,2));
452  AddFormula(new TFormulaPrimitive("Pol2","Pol2",(GenFuncG)TFastFun::FPol2,3));
453  AddFormula(new TFormulaPrimitive("Pol3","Pol3",(GenFuncG)TFastFun::FPol3,4));
454  AddFormula(new TFormulaPrimitive("Pol4","Pol4",(GenFuncG)TFastFun::FPol4,5));
455  AddFormula(new TFormulaPrimitive("Pol5","Pol5",(GenFuncG)TFastFun::FPol5,6));
456  AddFormula(new TFormulaPrimitive("Pol6","Pol6",(GenFuncG)TFastFun::FPol6,7));
457  AddFormula(new TFormulaPrimitive("Pol7","Pol7",(GenFuncG)TFastFun::FPol7,8));
458  AddFormula(new TFormulaPrimitive("Pol8","Pol8",(GenFuncG)TFastFun::FPol8,9));
459  AddFormula(new TFormulaPrimitive("Pol9","Pol9",(GenFuncG)TFastFun::FPol9,10));
460  AddFormula(new TFormulaPrimitive("Pol10","Pol10",(GenFuncG)TFastFun::FPol10,11));
461  //
462  // pows
463  AddFormula(new TFormulaPrimitive("Pow2","Pow2",TFastFun::Pow2));
464  AddFormula(new TFormulaPrimitive("Pow3","Pow3",TFastFun::Pow3));
465  AddFormula(new TFormulaPrimitive("Pow4","Pow4",TFastFun::Pow4));
466  AddFormula(new TFormulaPrimitive("Pow5","Pow5",TFastFun::Pow5));
467  //
468  //
469  AddFormula(new TFormulaPrimitive("TMath::Cos","TMath::Cos",cos)); // 10
470  AddFormula(new TFormulaPrimitive("cos","cos",cos)); // 10
471  AddFormula(new TFormulaPrimitive("TMath::Sin","TMath::Sin",sin)); // 11
472  AddFormula(new TFormulaPrimitive("sin","sin",sin)); // 11
473  AddFormula(new TFormulaPrimitive("TMath::Tan","TMath::Tan",tan)); // 12
474  AddFormula(new TFormulaPrimitive("tan","tan",tan)); // 12
475  AddFormula(new TFormulaPrimitive("TMath::ACos","TMath::ACos",acos)); // 13
476  AddFormula(new TFormulaPrimitive("acos","acos",acos)); // 13
477  AddFormula(new TFormulaPrimitive("TMath::ASin","TMath::ASin",asin)); // 14
478  AddFormula(new TFormulaPrimitive("asin","asin",asin)); // 14
479  AddFormula(new TFormulaPrimitive("TMath::ATan","TMath::ATan",atan)); // 15
480  AddFormula(new TFormulaPrimitive("atan","atan",atan)); // 15
481  AddFormula(new TFormulaPrimitive("TMath::ATan2","TMath::ATan2",atan2)); // 16
482  AddFormula(new TFormulaPrimitive("atan2","atan2",atan2)); // 16
483  // kpow = 20, ksq = 21, ksqrt = 22,
484  AddFormula(new TFormulaPrimitive("pow","pow",TMath::Power)); //20
485  AddFormula(new TFormulaPrimitive("sq","sq",TFastFun::Pow2)); //21
486  AddFormula(new TFormulaPrimitive("sqrt","sqrt",TFastFun::Sqrt)); //22
487  // kmin = 24, kmax = 25,
488  AddFormula(new TFormulaPrimitive("min","min",(GenFunc110)TMath::Min)); //24
489  AddFormula(new TFormulaPrimitive("max","max",(GenFunc110)TMath::Max)); //25
490  // klog = 30, kexp = 31, klog10 = 32,
491  AddFormula(new TFormulaPrimitive("log","log",TMath::Log)); //30
492  AddFormula(new TFormulaPrimitive("exp","exp",TMath::Exp)); //31
493  AddFormula(new TFormulaPrimitive("log10","log10",TMath::Log10)); //32
494  //
495  // cosh 70 acosh 73
496  // sinh 71 asinh 74
497  // tanh 72 atanh 75
498  //
499  AddFormula(new TFormulaPrimitive("TMath::CosH","TMath::Cosh",cosh)); // 70
500  AddFormula(new TFormulaPrimitive("cosh","cosh",cosh)); // 70
501  AddFormula(new TFormulaPrimitive("TMath::SinH","TMath::SinH",sinh)); // 71
502  AddFormula(new TFormulaPrimitive("sinh","sinh",sinh)); // 71
503  AddFormula(new TFormulaPrimitive("TMath::TanH","TMath::Tanh",tanh)); // 72
504  AddFormula(new TFormulaPrimitive("tanh","tanh",tanh)); // 72
505  AddFormula(new TFormulaPrimitive("TMath::ACosH","TMath::ACosh",TMath::ACosH)); // 73
506  AddFormula(new TFormulaPrimitive("acosh","acosH",TMath::ACosH)); // 73
507  AddFormula(new TFormulaPrimitive("TMath::ASinH","TMath::ASinh",TMath::ASinH)); // 74
508  AddFormula(new TFormulaPrimitive("acosh","acosH",TMath::ASinH)); // 74
509  AddFormula(new TFormulaPrimitive("TMath::ATanH","TMath::ATanh",TMath::ATanH)); // 75
510  AddFormula(new TFormulaPrimitive("atanh","atanh",TMath::ATanH)); // 75
511  //
512  AddFormula(new TFormulaPrimitive("TMath::Abs","TMath::Abs",TMath::Abs));
513  AddFormula(new TFormulaPrimitive("TMath::BreitWigner","TMath::BreitWigner",TMath::BreitWigner));
514 
515  //Disable direct access to TMath::Landau for now because of the default parameter.
516  //AddFormula(new TFormulaPrimitive("TMath::Landau","TMath::Landau",(TFormulaPrimitive::GenFunc1110)TMath::Landau));
517 
519  return 1;
520 }
521 
522  } // end namespace v5
523 
524 } // end namespace ROOT
Double_t ACosH(Double_t)
Definition: TMath.cxx:80
Double_t MinusXY(Double_t x, Double_t y)
Double_t Landau(Double_t x, Double_t mpv=0, Double_t sigma=1, Bool_t norm=kFALSE)
The LANDAU function.
Definition: TMath.cxx:472
An array of TObjects.
Definition: TObjArray.h:39
Double_t BreitWigner(Double_t x, Double_t mean=0, Double_t gamma=1)
Calculate a Breit Wigner function with mean and gamma.
Definition: TMath.cxx:441
Double_t(TObject::* TFunc10)(Double_t) const
double tanh(double)
Double_t XgeY(Double_t x, Double_t y)
Double_t PlusXY(Double_t x, Double_t y)
This namespace contains pre-defined functions to be used in conjuction with TExecutor::Map and TExecu...
Definition: StringConv.hxx:21
Double_t Log(Double_t x)
Definition: TMath.h:526
Double_t XgY(Double_t x, Double_t y)
Double_t(TObject::* TFunc1110)(Double_t, Double_t, Double_t) const
Double_t Landaun(Double_t x, Double_t mean, Double_t sigma)
Double_t FPol8(Double_t *x, Double_t *param)
Double_t XneY(Double_t x, Double_t y)
Double_t XpYpZ(Double_t x, Double_t y, Double_t z)
return c
Double_t Landau(Double_t x, Double_t mean, Double_t sigma)
Double_t FPol9(Double_t *x, Double_t *param)
virtual void SetOwner(Bool_t enable=kTRUE)
Set whether this collection is the owner (enable==true) of its content.
GenFunc10 fFunc10
pointer to the function
Double_t FPol5(Double_t *x, Double_t *param)
Double_t FPol3(Double_t *x, Double_t *param)
This class implements a mutex interface.
Definition: TVirtualMutex.h:34
Double_t XxYxZ(Double_t x, Double_t y, Double_t z)
Double_t Gaus(Double_t x, Double_t mean, Double_t sigma)
Gauss.
Double_t(* GenFuncG)(const Double_t *, const Double_t *)
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:170
int Int_t
Definition: RtypesCore.h:41
Double_t FPoln(Double_t *x, Double_t *param, Int_t npar)
FPoln.
TFunc110 fTFunc110
pointer to member function
const Bool_t kFALSE
Definition: Rtypes.h:92
Double_t MultXY(Double_t x, Double_t y)
Double_t(* GenFunc1110)(Double_t, Double_t, Double_t)
TFunc0 fTFunc0
pointer to the TFormula generic function
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
double cos(double)
Double_t XxYpZ(Double_t x, Double_t y, Double_t z)
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
GenFunc0 fFunc0
pointer to the TFormula generic function
Double_t(TObject::* TFunc0)() const
Double_t Pow5(Double_t x)
The Formula Primitive class.
double sqrt(double)
TFunc10 fTFunc10
pointer to member function
Double_t FPol10(Double_t *x, Double_t *param)
Double_t x[n]
Definition: legend1.C:17
double acos(double)
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:33
double sinh(double)
Double_t FPol1(Double_t *x, Double_t *param)
Double_t ASinH(Double_t)
Definition: TMath.cxx:67
Double_t Log10(Double_t x)
Definition: TMath.h:529
Double_t XlY(Double_t x, Double_t y)
const Double_t sigma
Double_t FPol6(Double_t *x, Double_t *param)
TFunc1110 fTFunc1110
pointer to member function
Double_t(* GenFunc110)(Double_t, Double_t)
double sin(double)
Double_t Eval(Double_t *x)
Eval primitive function at point x.
Double_t Sqrt(Double_t x)
Double_t XNot(Double_t x)
Double_t Pow4(Double_t x)
double cosh(double)
Double_t XleY(Double_t x, Double_t y)
Double_t FPol2(Double_t *x, Double_t *param)
unsigned int UInt_t
Definition: RtypesCore.h:42
Double_t Nint(Double_t x)
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
double asin(double)
Double_t XandY(Double_t x, Double_t y)
#define RTFastFun__POLY(var)
#define R__LOCKGUARD2(mutex)
PyObject * fType
Double_t(* GenFunc10)(Double_t)
Double_t Exp(Double_t x)
Definition: TMath.h:495
#define ClassImp(name)
Definition: Rtypes.h:279
double Double_t
Definition: RtypesCore.h:55
Double_t Gausn(Double_t x, Double_t mean, Double_t sigma)
Normalize gauss.
double atan2(double, double)
Double_t y[n]
Definition: legend1.C:17
double atan(double)
Double_t Pow2(Double_t x)
Double_t FPol0(Double_t *, Double_t *param)
Int_t GetEntries() const
Return the number of objects in array (i.e.
Definition: TObjArray.cxx:494
static Int_t AddFormula(TFormulaPrimitive *formula)
Add formula to the list of primitive formulas.
static Int_t BuildBasicFormulas()
list of global primitive formulas
Mother of all ROOT objects.
Definition: TObject.h:44
you should not use this method at all Int_t Int_t z
Definition: TRolke.cxx:630
Double_t(TObject::* TFunc110)(Double_t, Double_t) const
TFormulaPrimitive()
Default constructor.
Double_t Abs(Double_t x)
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
double tan(double)
static TVirtualMutex * gTFormulaPrimativeListMutex
Double_t Pow3(Double_t x)
Double_t FPol4(Double_t *x, Double_t *param)
Double_t Sign(Double_t x)
Double_t FPol7(Double_t *x, Double_t *param)
virtual void AddLast(TObject *obj)
Add object in the next empty slot in the array.
Definition: TObjArray.cxx:169
GenFunc110 fFunc110
pointer to the function
Double_t(TObject::* TFuncG)(const Double_t *, const Double_t *) const
static TObjArray * fgListOfFunction
Double_t DivXY(Double_t x, Double_t y)
TObject * At(Int_t idx) const
Definition: TObjArray.h:167
const Bool_t kTRUE
Definition: Rtypes.h:91
Double_t XpYxZ(Double_t x, Double_t y, Double_t z)
Double_t XeY(Double_t x, Double_t y)
Int_t Nint(T x)
Definition: TMath.h:480
Double_t ATanH(Double_t)
Definition: TMath.cxx:93
Double_t XorY(Double_t x, Double_t y)
char name[80]
Definition: TGX11.cxx:109
static TFormulaPrimitive * FindFormula(const char *name)
Find the formula in the list of formulas.
TFuncG fTFuncG
pointer to the function
void TMath_GenerInterface()
GenFunc1110 fFunc1110
pointer to the function