ROOT  6.06/09
Reference Guide
TFormula.cxx
Go to the documentation of this file.
1 // @(#)root/hist:$Id$
2 // Author: Maciej Zimnoch 30/09/2013
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2013, 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 #if __cplusplus >= 201103L
13 #define ROOT_CPLUSPLUS11 1
14 #endif
15 
16 #include "TROOT.h"
17 #include "TClass.h"
18 #include "TMethod.h"
19 #include "TMath.h"
20 #include "TF1.h"
21 #include "TMethodCall.h"
22 #include <TBenchmark.h>
23 #include "TError.h"
24 #include "TInterpreter.h"
25 #include "TFormula.h"
26 #include <cassert>
27 #include <iostream>
28 #include <unordered_map>
29 #include <functional>
30 
31 using namespace std;
32 
33 // #define __STDC_LIMIT_MACROS
34 // #define __STDC_CONSTANT_MACROS
35 
36 // #include "cling/Interpreter/Interpreter.h"
37 // #include "cling/Interpreter/Value.h"
38 // #include "cling/Interpreter/StoredValueRef.h"
39 
40 
41 #ifdef WIN32
42 #pragma optimize("",off)
43 #endif
44 #include "v5/TFormula.h"
45 
47 /** \class TFormula TFormula.h "inc/TFormula.h"
48  \ingroup Hist
49 The F O R M U L A class
50 
51 <p>This is a new version of the TFormula class based on Cling.
52 This class is not 100% backward compatible with the old TFOrmula class, which is still available in ROOT as
53 <code>ROOT::v5::TFormula</code>. Some of the TFormula member funtions available in version 5, such as
54 <code>Analyze</code> and <code>AnalyzeFunction</code> are not available in the new TFormula.
55 On the other hand formula expressions which were valid in version 5 are still valid in TFormula version 6</p>
56 
57 <p>This class has been implemented during Google Summer of Code 2013 by Maciej Zimnoch.</p>
58 
59 <h3>Example of valid expressions:</h3>
60 
61 <ul>
62 <li><code>sin(x)/x</code></li>
63 <li><code>[0]*sin(x) + [1]*exp(-[2]*x)</code></li>
64 <li><code>x + y**2</code></li>
65 <li><code>x^2 + y^2</code></li>
66 <li><code>[0]*pow([1],4)</code></li>
67 <li><code>2*pi*sqrt(x/y)</code></li>
68 <li><code>gaus(0)*expo(3) + ypol3(5)*x</code></li>
69 <li><code>gausn(0)*expo(3) + ypol3(5)*x</code></li>
70 </ul>
71 
72 <p>In the last example above:</p>
73 
74 <ul>
75 <li><code>gaus(0)</code> is a substitute for <code>[0]*exp(-0.5*((x-[1])/[2])**2)</code>
76  and (0) means start numbering parameters at 0</li>
77 <li><code>gausn(0)</code> is a substitute for <code>[0]*exp(-0.5*((x-[1])/[2])**2)/(sqrt(2*pi)*[2]))</code>
78  and (0) means start numbering parameters at 0</li>
79 <li><code>expo(3)</code> is a substitute for <code>exp([3]+[4]*x</code></li>
80 <li><code>pol3(5)</code> is a substitute for <code>par[5]+par[6]*x+par[7]*x**2+par[8]*x**3</code>
81  (<code>PolN</code> stands for Polynomial of degree N)</li>
82 </ul>
83 
84 <p><code>TMath</code> functions can be part of the expression, eg:</p>
85 
86 <ul>
87 <li><code>TMath::Landau(x)*sin(x)</code></li>
88 <li><code>TMath::Erf(x)</code></li>
89 </ul>
90 
91 <p>Formula may contain constans, eg:</p>
92 
93 <ul>
94 <li><code>sqrt2</code></li>
95 <li><code>e</code></li>
96 <li><code>pi</code></li>
97 <li><code>ln10</code></li>
98 <li><code>infinity</code></li>
99 </ul>
100 
101 <p>and more.</p>
102 
103 <p>Comparisons operators are also supported <code>(&amp;&amp;, ||, ==, &lt;=, &gt;=, !)</code></p>
104 
105 <p>Examples:</p>
106 
107 <pre><code> `sin(x*(x&lt;0.5 || x&gt;1))`
108 </code></pre>
109 
110 <p>If the result of a comparison is TRUE, the result is 1, otherwise 0.</p>
111 
112 <p>Already predefined names can be given. For example, if the formula</p>
113 
114 <pre><code> `TFormula old("old",sin(x*(x&lt;0.5 || x&gt;1)))`
115 </code></pre>
116 
117 <p>one can assign a name to the formula. By default the name of the object = title = formula itself.</p>
118 
119 <pre><code> `TFormula new("new","x*old")`
120 </code></pre>
121 
122 <p>is equivalent to:</p>
123 
124 <pre><code> `TFormula new("new","x*sin(x*(x&lt;0.5 || x&gt;1))")`
125 </code></pre>
126 
127 <p>The class supports unlimited numer of variables and parameters.
128  By default the names which can be used for the variables are <code>x,y,z,t</code> or
129  <code>x[0],x[1],x[2],x[3],....x[N]</code> for N-dimensionals formula.</p>
130 
131 <p>This class is not anymore the base class for the function classes <code>TF1</code>, but it has now
132 adata member of TF1 which can be access via <code>TF1::GetFormula</code>. </p>
133 
134 
135 \class TFormulaFunction
136  Helper class for TFormula
137 \class TFormulaVariable
138  Another helper class for TFormula
139 \class TFormulaParamOrder
140  Functor defining the parameter order
141 
142 */
143 
144 // prefix used for function name passed to Cling
145 static const TString gNamePrefix = "TFormula__";
146 
147 // static map of function pointers and expressions
148 //static std::unordered_map<std::string, TInterpreter::CallFuncIFacePtr_t::Generic_t> gClingFunctions = std::unordered_map<TString, TInterpreter::CallFuncIFacePtr_t::Generic_t>();
149 static std::unordered_map<std::string, void *> gClingFunctions = std::unordered_map<std::string, void * >();
150 
151 Bool_t TFormula::IsOperator(const char c)
152 {
153  // operator ":" must be handled separatly
154  char ops[] = { '+','^', '-','/','*','<','>','|','&','!','=','?'};
155  Int_t opsLen = sizeof(ops)/sizeof(char);
156  for(Int_t i = 0; i < opsLen; ++i)
157  if(ops[i] == c)
158  return true;
159  return false;
160 }
161 
163 {
164  char brackets[] = { ')','(','{','}'};
165  Int_t bracketsLen = sizeof(brackets)/sizeof(char);
166  for(Int_t i = 0; i < bracketsLen; ++i)
167  if(brackets[i] == c)
168  return true;
169  return false;
170 }
171 
173 {
174  return !IsBracket(c) && !IsOperator(c) && c != ',' && c != ' ';
175 }
176 
178 {
179  return name == "x" || name == "z" || name == "y" || name == "t";
180 }
181 
182 
184 {
185  // check if the character at position i is part of a scientific notation
186  if ( (formula[i] == 'e' || formula[i] == 'E') && (i > 0 && i < formula.Length()-1) ) {
187  // handle cases: 2e+3 2e-3 2e3 and 2.e+3
188  if ( (isdigit(formula[i-1]) || formula[i-1] == '.') && ( isdigit(formula[i+1]) || formula[i+1] == '+' || formula[i+1] == '-' ) )
189  return true;
190  }
191  return false;
192 }
193 
194 Bool_t TFormula::IsHexadecimal(const TString & formula, int i)
195 {
196  // check if the character at position i is part of a scientific notation
197  if ( (formula[i] == 'x' || formula[i] == 'X') && (i > 0 && i < formula.Length()-1) && formula[i-1] == '0') {
198  if (isdigit(formula[i+1]) )
199  return true;
200  static char hex_values[12] = { 'a','A', 'b','B','c','C','d','D','e','E','f','F'};
201  for (int jjj = 0; jjj < 12; ++jjj) {
202  if (formula[i+1] == hex_values[jjj])
203  return true;
204  }
205  }
206  // else
207  // return false;
208  // // handle cases: 2e+3 2e-3 2e3 and 2.e+3
209  // if ( (isdigit(formula[i-1]) || formula[i-1] == '.') && ( isdigit(formula[i+1]) || formula[i+1] == '+' || formula[i+1] == '-' ) )
210  // return true;
211  // }
212  return false;
213 }
214 
215 bool TFormulaParamOrder::operator() (const TString& a, const TString& b) const {
216  // implement comparison used to set parameter orders in TFormula
217  // want p2 to be before p10
218 
219  // strip first character in case you have (p0, p1, pN)
220  if ( a[0] == 'p' && a.Length() > 1) {
221  if ( b[0] == 'p' && b.Length() > 1) {
222  // strip first character
223  TString lhs = a(1,a.Length()-1);
224  TString rhs = b(1,b.Length()-1);
225  if (lhs.IsDigit() && rhs.IsDigit() )
226  return (lhs.Atoi() < rhs.Atoi() );
227  }
228  else {
229  return true; // assume a(a numeric name) is always before b (an alphanumeric name)
230  }
231  }
232  else {
233  if ( b[0] == 'p' && b.Length() > 1)
234  // now b is numeric and a is not so return false
235  return false;
236 
237  // case both names are numeric
238  if (a.IsDigit() && b.IsDigit() )
239  return (a.Atoi() < b.Atoi() );
240 
241  }
242 
243  return a < b;
244 }
245 
247 {
248  fName = "";
249  fTitle = "";
250  fClingInput = "";
251  fReadyToExecute = false;
252  fClingInitialized = false;
253  fAllParametersSetted = false;
254  fMethod = 0;
255  fNdim = 0;
256  fNpar = 0;
257  fNumber = 0;
258  fClingName = "";
259  fFormula = "";
260  fLambdaPtr = nullptr;
261 }
262 
263 ////////////////////////////////////////////////////////////////////////////////
264 
265 static bool IsReservedName(const char* name){
266  if (strlen(name)!=1) return false;
267  for (auto const & specialName : {"x","y","z","t"}){
268  if (strcmp(name,specialName)==0) return true;
269  }
270  return false;
271 }
272 
273 
275 {
276 
277  // N.B. a memory leak may happen if user set bit after constructing the object,
278  // Setting of bit should be done only internally
279  if (!TestBit(TFormula::kNotGlobal) && gROOT ) {
281  gROOT->GetListOfFunctions()->Remove(this);
282  }
283 
284  if(fMethod)
285  {
286  fMethod->Delete();
287  }
288  int nLinParts = fLinearParts.size();
289  if (nLinParts > 0) {
290  for (int i = 0; i < nLinParts; ++i) delete fLinearParts[i];
291  }
292 }
293 
294 TFormula::TFormula(const char *name, const char *formula, bool addToGlobList) :
295  TNamed(name,formula),
296  fClingInput(formula),fFormula(formula)
297 {
298  fReadyToExecute = false;
299  fClingInitialized = false;
300  fMethod = 0;
301  fNdim = 0;
302  fNpar = 0;
303  fNumber = 0;
304  fLambdaPtr = nullptr;
305 
306  FillDefaults();
307 
308 
309  if (addToGlobList && gROOT) {
310  TFormula *old = 0;
312  old = dynamic_cast<TFormula*> ( gROOT->GetListOfFunctions()->FindObject(name) );
313  if (old)
314  gROOT->GetListOfFunctions()->Remove(old);
315  if (IsReservedName(name))
316  Error("TFormula","The name %s is reserved as a TFormula variable name.\n",name);
317  else
318  gROOT->GetListOfFunctions()->Add(this);
319  }
320  SetBit(kNotGlobal,!addToGlobList);
321 
322  //fName = gNamePrefix + name; // is this needed
323 
324  // do not process null formulas.
325  if (!fFormula.IsNull() ) {
327 
329  }
330 
331 }
332 
333 /// constructor from a full compileable C++ expression
334 TFormula::TFormula(const char *name, const char *formula, int ndim, int npar, bool addToGlobList) :
335  TNamed(name,formula),
336  fClingInput(formula),fFormula(formula)
337 {
338  fReadyToExecute = false;
339  fClingInitialized = false;
340  fNpar = 0;
341  fLambdaPtr = nullptr;
342 
343 
344  fNdim = ndim;
345  for (int i = 0; i < npar; ++i) {
346  DoAddParameter(TString::Format("p%d",i), 0, false);
347  }
348  fAllParametersSetted = true;
349  assert (fNpar == npar);
350 
351  bool ret = InitLambdaExpression(formula);
352 
353  if (ret) {
354 
356 
357  fReadyToExecute = true;
358 
359  if (addToGlobList && gROOT) {
360  TFormula *old = 0;
362  old = dynamic_cast<TFormula*> ( gROOT->GetListOfFunctions()->FindObject(name) );
363  if (old)
364  gROOT->GetListOfFunctions()->Remove(old);
365  if (IsReservedName(name))
366  Error("TFormula","The name %s is reserved as a TFormula variable name.\n",name);
367  else
368  gROOT->GetListOfFunctions()->Add(this);
369  }
370  SetBit(kNotGlobal,!addToGlobList);
371  }
372  else
373  Error("TFormula","Syntax error in building the lambda expression %s", formula );
374 }
375 
376 
377 TFormula::TFormula(const TFormula &formula) :
378  TNamed(formula.GetName(),formula.GetTitle())
379 {
380  fReadyToExecute = false;
381  fClingInitialized = false;
382  fMethod = 0;
383  fNdim = formula.GetNdim();
384  fNpar = formula.GetNpar();
385  fNumber = formula.GetNumber();
386  fFormula = formula.GetExpFormula(); // returns fFormula in case of Lambda's
387  fLambdaPtr = nullptr;
388 
389  // case of function based on a C++ expression (lambda's) which is ready to be compiled
390  if (formula.fLambdaPtr && formula.TestBit(TFormula::kLambda)) {
391 
393  fParams = formula.fParams;
396 
397  bool ret = InitLambdaExpression(fFormula);
398 
399  if (ret) {
401  fReadyToExecute = true;
402  }
403  else
404  Error("TFormula","Syntax error in building the lambda expression %s", fFormula.Data() );
405 
406  }
407  else {
408 
409  FillDefaults();
410 
413  }
414 
415 
416  if (!TestBit(TFormula::kNotGlobal) && gROOT ) {
418  TFormula *old = (TFormula*)gROOT->GetListOfFunctions()->FindObject(formula.GetName());
419  if (old)
420  gROOT->GetListOfFunctions()->Remove(old);
421 
422  if (IsReservedName(formula.GetName())) {
423  Error("TFormula","The name %s is reserved as a TFormula variable name.\n",formula.GetName());
424  } else
425  gROOT->GetListOfFunctions()->Add(this);
426  }
427 
428 }
429 
431 {
432  //*-*
433  //*-* = Operator
434  //*-*
435 
436  if (this != &rhs) {
437  rhs.Copy(*this);
438  }
439  return *this;
440 }
441 
442 Bool_t TFormula::InitLambdaExpression(const char * formula) {
443 
444  std::string lambdaExpression = formula;
445 
446  // check if formula exist already in the map
447  {
449 
450  auto funcit = gClingFunctions.find(lambdaExpression);
451  if (funcit != gClingFunctions.end() ) {
452  fLambdaPtr = funcit->second;
453  fClingInitialized = true;
454  return true;
455  }
456  }
457 
458  // set the cling name using hash of the static formulae map
459  auto hasher = gClingFunctions.hash_function();
460  TString lambdaName = TString::Format("lambda__id%zu",(unsigned long) hasher(lambdaExpression) );
461 
462  //lambdaExpression = TString::Format("[&](double * x, double *){ return %s ;}",formula);
463  //TString lambdaName = TString::Format("mylambda_%s",GetName() );
464  TString lineExpr = TString::Format("std::function<double(double*,double*)> %s = %s ;",lambdaName.Data(), lambdaExpression.c_str() );
465  gInterpreter->ProcessLine(lineExpr);
466  fLambdaPtr = (void*) gInterpreter->ProcessLine(TString(lambdaName)+TString(";")); // add ; to avoid printing
467  if (fLambdaPtr != nullptr) {
469  gClingFunctions.insert ( std::make_pair ( lambdaExpression, fLambdaPtr) );
470  fClingInitialized = true;
471  return true;
472  }
473  fClingInitialized = false;
474  return false;
475 }
476 
477 
478 
479 Int_t TFormula::Compile(const char *expression)
480 {
481  // Compile the given expression with Cling
482  // backward compatibility method to be used in combination with the empty constructor
483  // if no expression is given , the current stored formula (retrieved with GetExpFormula()) or the title is used.
484  // return 0 if the formula compilation is successfull
485 
486 
487  TString formula = expression;
488  if (formula.IsNull() ) {
489  formula = fFormula;
490  if (formula.IsNull() ) formula = GetTitle();
491  }
492 
493  if (formula.IsNull() ) return -1;
494 
495  // do not re-process if it was done before
496  if (IsValid() && formula == fFormula ) return 0;
497 
498  // clear if a formula was already existing
499  if (!fFormula.IsNull() ) Clear();
500 
501  fFormula = formula;
502 
503  if (TestBit(TFormula::kLambda) ) {
504  bool ret = InitLambdaExpression(fFormula);
505  return (ret) ? 0 : 1;
506  }
507 
508  if (fVars.empty() ) FillDefaults();
509  // prepare the formula for Cling
510  //printf("compile: processing formula %s\n",fFormula.Data() );
512  // pass formula in CLing
513  bool ret = PrepareFormula(fFormula);
514 
515  return (ret) ? 0 : 1;
516 }
517 
519 {
520  TNamed::Copy(obj);
521  // need to copy also cling parameters
522  TFormula & fnew = dynamic_cast<TFormula&>(obj);
523 
526 
527  fnew.fFuncs = fFuncs;
528  fnew.fVars = fVars;
529  fnew.fParams = fParams;
530  fnew.fConsts = fConsts;
532  fnew.fFormula = fFormula;
533  fnew.fNdim = fNdim;
534  fnew.fNpar = fNpar;
535  fnew.fNumber = fNumber;
537  // copy Linear parts (it is a vector of TFormula pointers) needs to be copied one by one
538  // looping at all the elements
539  // delete first previous elements
540  int nLinParts = fnew.fLinearParts.size();
541  if (nLinParts > 0) {
542  for (int i = 0; i < nLinParts; ++i) delete fnew.fLinearParts[i];
543  fnew.fLinearParts.clear();
544  }
545  // old size that needs to be copied
546  nLinParts = fLinearParts.size();
547  if (nLinParts > 0) {
548  fnew.fLinearParts.reserve(nLinParts);
549  for (int i = 0; i < nLinParts; ++i) {
550  TFormula * linearNew = new TFormula();
551  TFormula * linearOld = (TFormula*) fLinearParts[i];
552  if (linearOld) {
553  linearOld->Copy(*linearNew);
554  fnew.fLinearParts.push_back(linearNew);
555  }
556  else
557  Warning("Copy","Function %s - expr %s has a dummy linear part %d",GetName(),GetExpFormula().Data(),i);
558  }
559  }
560 
561  fnew.fClingInput = fClingInput;
565  fnew.fClingName = fClingName;
566 
567  // case of function based on a C++ expression (lambda's) which is ready to be compiled
569 
570  bool ret = fnew.InitLambdaExpression(fnew.fFormula);
571  if (ret) {
573  fnew.fReadyToExecute = true;
574  }
575  else {
576  Error("TFormula","Syntax error in building the lambda expression %s", fFormula.Data() );
577  fnew.fReadyToExecute = false;
578  }
579  }
580  else if (fMethod) {
581  if (fnew.fMethod) delete fnew.fMethod;
582  // use copy-constructor of TMethodCall
583  TMethodCall *m = new TMethodCall(*fMethod);
584  fnew.fMethod = m;
585  }
586 
587  fnew.fFuncPtr = fFuncPtr;
588 
589 }
590 
592 {
593  // clear the formula setting expression to empty and reset the variables and parameters containers
594  fNdim = 0;
595  fNpar = 0;
596  fNumber = 0;
597  fFormula = "";
598  fClingName = "";
599 
600 
601  if(fMethod) fMethod->Delete();
602  fMethod = nullptr;
603 
604  fClingVariables.clear();
605  fClingParameters.clear();
606  fReadyToExecute = false;
607  fClingInitialized = false;
608  fAllParametersSetted = false;
609  fFuncs.clear();
610  fVars.clear();
611  fParams.clear();
612  fConsts.clear();
613  fFunctionsShortcuts.clear();
614 
615  // delete linear parts
616  int nLinParts = fLinearParts.size();
617  if (nLinParts > 0) {
618  for (int i = 0; i < nLinParts; ++i) delete fLinearParts[i];
619  }
620  fLinearParts.clear();
621 
622 }
624 {
625  //*-*
626  //*-* Sets TMethodCall to function inside Cling environment
627  //*-* TFormula uses it to execute function.
628  //*-* After call, TFormula should be ready to evaluate formula.
629  //*-*
630  if(!fMethod)
631  {
632  fMethod = new TMethodCall();
633 
634  Bool_t hasParameters = (fNpar > 0);
635  Bool_t hasVariables = (fNdim > 0);
636  TString prototypeArguments = "";
637  if(hasVariables)
638  {
639  prototypeArguments.Append("Double_t*");
640  }
641  if(hasVariables && hasParameters)
642  {
643  prototypeArguments.Append(",");
644  }
645  if(hasParameters)
646  {
647  prototypeArguments.Append("Double_t*");
648  }
649  // init method call using real function name (cling name) which is defined in ProcessFormula
650  fMethod->InitWithPrototype(fClingName,prototypeArguments);
651  if(!fMethod->IsValid())
652  {
653  Error("Eval","Can't find %s function prototype with arguments %s",fClingName.Data(),prototypeArguments.Data());
654  return false;
655  }
656 
657  // not needed anymore since we use the function pointer
658  // if(hasParameters)
659  // {
660  // Long_t args[2];
661  // args[0] = (Long_t)fClingVariables.data();
662  // args[1] = (Long_t)fClingParameters.data();
663  // fMethod->SetParamPtrs(args,2);
664  // }
665  // else
666  // {
667  // Long_t args[1];
668  // args[0] = (Long_t)fClingVariables.data();
669  // fMethod->SetParamPtrs(args,1);
670  // }
671 
672  CallFunc_t * callfunc = fMethod->GetCallFunc();
674  fFuncPtr = faceptr.fGeneric;
675  }
676  return true;
677 }
678 
680 {
681  //*-*
682  //*-* Inputs formula, transfered to C++ code into Cling
683  //*-*
685  {
688  }
689 }
691 {
692  //*-*
693  //*-* Fill structures with default variables, constants and function shortcuts
694  //*-*
695 //#ifdef ROOT_CPLUSPLUS11
696 
697  const TString defvars[] = { "x","y","z","t"};
698  const pair<TString,Double_t> defconsts[] = { {"pi",TMath::Pi()}, {"sqrt2",TMath::Sqrt2()},
699  {"infinity",TMath::Infinity()}, {"e",TMath::E()}, {"ln10",TMath::Ln10()},
700  {"loge",TMath::LogE()}, {"c",TMath::C()}, {"g",TMath::G()},
701  {"h",TMath::H()}, {"k",TMath::K()},{"sigma",TMath::Sigma()},
702  {"r",TMath::R()}, {"eg",TMath::EulerGamma()},{"true",1},{"false",0} };
703  // const pair<TString,Double_t> defconsts[] = { {"pi",TMath::Pi()}, {"sqrt2",TMath::Sqrt2()},
704  // {"infinity",TMath::Infinity()}, {"ln10",TMath::Ln10()},
705  // {"loge",TMath::LogE()}, {"true",1},{"false",0} };
706  const pair<TString,TString> funShortcuts[] =
707  { {"sin","TMath::Sin" },
708  {"cos","TMath::Cos" }, {"exp","TMath::Exp"}, {"log","TMath::Log"}, {"log10","TMath::Log10"},
709  {"tan","TMath::Tan"}, {"sinh","TMath::SinH"}, {"cosh","TMath::CosH"},
710  {"tanh","TMath::TanH"}, {"asin","TMath::ASin"}, {"acos","TMath::ACos"},
711  {"atan","TMath::ATan"}, {"atan2","TMath::ATan2"}, {"sqrt","TMath::Sqrt"},
712  {"ceil","TMath::Ceil"}, {"floor","TMath::Floor"}, {"pow","TMath::Power"},
713  {"binomial","TMath::Binomial"},{"abs","TMath::Abs"},
714  {"min","TMath::Min"},{"max","TMath::Max"},{"sign","TMath::Sign" },
715  {"sq","TMath::Sq"}
716  };
717 
718  std::vector<TString> defvars2(10);
719  for (int i = 0; i < 9; ++i)
720  defvars2[i] = TString::Format("x[%d]",i);
721 
722  for(auto var : defvars)
723  {
724  int pos = fVars.size();
725  fVars[var] = TFormulaVariable(var,0,pos);
726  fClingVariables.push_back(0);
727  }
728  // add also the variables definesd like x[0],x[1],x[2],...
729  // support up to x[9] - if needed extend that to higher value
730  // const int maxdim = 10;
731  // for (int i = 0; i < maxdim; ++i) {
732  // TString xvar = TString::Format("x[%d]",i);
733  // fVars[xvar] = TFormulaVariable(xvar,0,i);
734  // fClingVariables.push_back(0);
735  // }
736 
737  for(auto con : defconsts)
738  {
739  fConsts[con.first] = con.second;
740  }
741  for(auto fun : funShortcuts)
742  {
743  fFunctionsShortcuts[fun.first] = fun.second;
744  }
745 
746 /*** - old code tu support C++03
747 #else
748 
749  TString defvarsNames[] = {"x","y","z","t"};
750  Int_t defvarsLength = sizeof(defvarsNames)/sizeof(TString);
751 
752  TString defconstsNames[] = {"pi","sqrt2","infinity","e","ln10","loge","c","g","h","k","sigma","r","eg","true","false"};
753  Double_t defconstsValues[] = {TMath::Pi(),TMath::Sqrt2(),TMath::Infinity(),TMath::E(),TMath::Ln10(),TMath::LogE(),
754  TMath::C(),TMath::G(),TMath::H(),TMath::K(),TMath::Sigma(),TMath::R(),TMath::EulerGamma(), 1, 0};
755  Int_t defconstsLength = sizeof(defconstsNames)/sizeof(TString);
756 
757  TString funShortcutsNames[] = {"sin","cos","exp","log","tan","sinh","cosh","tanh","asin","acos","atan","atan2","sqrt",
758  "ceil","floor","pow","binomial","abs"};
759  TString funShortcutsExtendedNames[] = {"TMath::Sin","TMath::Cos","TMath::Exp","TMath::Log","TMath::Tan","TMath::SinH",
760  "TMath::CosH","TMath::TanH","TMath::ASin","TMath::ACos","TMath::ATan","TMath::ATan2",
761  "TMath::Sqrt","TMath::Ceil","TMath::Floor","TMath::Power","TMath::Binomial","TMath::Abs"};
762  Int_t funShortcutsLength = sizeof(funShortcutsNames)/sizeof(TString);
763 
764  for(Int_t i = 0; i < defvarsLength; ++i)
765  {
766  TString var = defvarsNames[i];
767  Double_t value = 0;
768  unsigned int pos = fVars.size();
769  fVars[var] = TFormulaVariable(var,value,pos);
770  fClingVariables.push_back(value);
771  }
772 
773  for(Int_t i = 0; i < defconstsLength; ++i)
774  {
775  fConsts[defconstsNames[i]] = defconstsValues[i];
776  }
777  for(Int_t i = 0; i < funShortcutsLength; ++i)
778  {
779  pair<TString,TString> fun(funShortcutsNames[i],funShortcutsExtendedNames[i]);
780  fFunctionsShortcuts[fun.first] = fun.second;
781  }
782 
783 #endif
784 ***/
785 
786 }
787 
789 {
790  //*-*
791  //*-* Handling polN
792  //*-* If before 'pol' exist any name, this name will be treated as variable used in polynomial
793  //*-* eg.
794  //*-* varpol2(5) will be replaced with: [5] + [6]*var + [7]*var^2
795  //*-* Empty name is treated like variable x.
796  //*-*
797  Int_t polPos = formula.Index("pol");
798  while(polPos != kNPOS)
799  {
800 
801  Bool_t defaultVariable = false;
802  TString variable;
803  Int_t openingBracketPos = formula.Index('(',polPos);
804  Bool_t defaultCounter = openingBracketPos == kNPOS;
805  Bool_t defaultDegree = true;
806  Int_t degree,counter;
807  TString sdegree;
808  if(!defaultCounter)
809  {
810  // veryfy first of opening parenthesis belongs to pol expression
811  // character between 'pol' and '(' must all be digits
812  sdegree = formula(polPos + 3,openingBracketPos - polPos - 3);
813  if (!sdegree.IsDigit() ) defaultCounter = true;
814  }
815  if (!defaultCounter) {
816  degree = sdegree.Atoi();
817  counter = TString(formula(openingBracketPos+1,formula.Index(')',polPos) - openingBracketPos)).Atoi();
818  }
819  else
820  {
821  Int_t temp = polPos+3;
822  while(temp < formula.Length() && isdigit(formula[temp]))
823  {
824  defaultDegree = false;
825  temp++;
826  }
827  degree = TString(formula(polPos+3,temp - polPos - 3)).Atoi();
828  counter = 0;
829  }
830 
831  TString replacement = TString::Format("[%d]",counter);
832  if(polPos - 1 < 0 || !IsFunctionNameChar(formula[polPos-1]) || formula[polPos-1] == ':' )
833  {
834  variable = "x";
835  defaultVariable = true;
836  }
837  else
838  {
839  Int_t tmp = polPos - 1;
840  while(tmp >= 0 && IsFunctionNameChar(formula[tmp]) && formula[tmp] != ':')
841  {
842  tmp--;
843  }
844  variable = formula(tmp + 1, polPos - (tmp+1));
845  }
846  Int_t param = counter + 1;
847  Int_t tmp = 1;
848  while(tmp <= degree)
849  {
850  replacement.Append(TString::Format("+[%d]*%s^%d",param,variable.Data(),tmp));
851  param++;
852  tmp++;
853  }
855  if(defaultCounter && !defaultDegree)
856  {
857  pattern = TString::Format("%spol%d",(defaultVariable ? "" : variable.Data()),degree);
858  }
859  else if(defaultCounter && defaultDegree)
860  {
861  pattern = TString::Format("%spol",(defaultVariable ? "" : variable.Data()));
862  }
863  else
864  {
865  pattern = TString::Format("%spol%d(%d)",(defaultVariable ? "" : variable.Data()),degree,counter);
866  }
867 
868  if (!formula.Contains(pattern)) {
869  Error("HandlePolN","Error handling polynomial function - expression is %s - trying to replace %s with %s ", formula.Data(), pattern.Data(), replacement.Data() );
870  break;
871  }
872  if (formula == pattern) {
873  // case of single polynomial
874  SetBit(kLinear,1);
875  fNumber = 300 + degree;
876  }
877  formula.ReplaceAll(pattern,replacement);
878  polPos = formula.Index("pol");
879  }
880 }
882 {
883  //*-*
884  //*-* Handling parametrized functions
885  //*-* Function can be normalized, and have different variable then x.
886  //*-* Variables should be placed in brackets after function name.
887  //*-* No brackets are treated like [x].
888  //*-* Normalized function has char 'n' after name, eg.
889  //*-* gausn[var](0) will be replaced with [0]*exp(-0.5*((var-[1])/[2])^2)/(sqrt(2*pi)*[2])
890  //*-*
891  //*-* Adding function is easy, just follow these rules:
892  //*-* - Key for function map is pair of name and dimension of function
893  //*-* - value of key is a pair function body and normalized function body
894  //*-* - {Vn} is a place where variable appear, n represents n-th variable from variable list.
895  //*-* Count starts from 0.
896  //*-* - [num] stands for parameter number.
897  //*-* If user pass to function argument 5, num will stand for (5 + num) parameter.
898  //*-*
899 
900  map< pair<TString,Int_t> ,pair<TString,TString> > functions;
901  functions.insert(make_pair(make_pair("gaus",1),make_pair("[0]*exp(-0.5*(({V0}-[1])/[2])*(({V0}-[1])/[2]))","[0]*exp(-0.5*(({V0}-[1])/[2])*(({V0}-[1])/[2]))/(sqrt(2*pi)*[2])")));
902  functions.insert(make_pair(make_pair("landau",1),make_pair("[0]*TMath::Landau({V0},[1],[2],false)","[0]*TMath::Landau({V0},[1],[2],true)")));
903  functions.insert(make_pair(make_pair("expo",1),make_pair("exp([0]+[1]*{V0})","")));
904  functions.insert(make_pair(make_pair("crystalball",1),make_pair("[0]*ROOT::Math::crystalball_function({V0},[3],[4],[2],[1])","[0]*ROOT::Math::crystalball_pdf({V0},[3],[4],[2],[1])")));
905  functions.insert(make_pair(make_pair("breitwigner",1),make_pair("[0]*ROOT::Math::breitwigner_pdf({V0},[2],[1])","[0]*ROOT::Math::breitwigner_pdf({V0},[2],[4],[1])")));
906  // chebyshev polynomial
907  functions.insert(make_pair(make_pair("cheb0" ,1),make_pair("ROOT::Math::Chebyshev0({V0},[0])","")));
908  functions.insert(make_pair(make_pair("cheb1" ,1),make_pair("ROOT::Math::Chebyshev1({V0},[0],[1])","")));
909  functions.insert(make_pair(make_pair("cheb2" ,1),make_pair("ROOT::Math::Chebyshev2({V0},[0],[1],[2])","")));
910  functions.insert(make_pair(make_pair("cheb3" ,1),make_pair("ROOT::Math::Chebyshev3({V0},[0],[1],[2],[3])","")));
911  functions.insert(make_pair(make_pair("cheb4" ,1),make_pair("ROOT::Math::Chebyshev4({V0},[0],[1],[2],[3],[4])","")));
912  functions.insert(make_pair(make_pair("cheb5" ,1),make_pair("ROOT::Math::Chebyshev5({V0},[0],[1],[2],[3],[4],[5])","")));
913  functions.insert(make_pair(make_pair("cheb6" ,1),make_pair("ROOT::Math::Chebyshev6({V0},[0],[1],[2],[3],[4],[5],[6])","")));
914  functions.insert(make_pair(make_pair("cheb7" ,1),make_pair("ROOT::Math::Chebyshev7({V0},[0],[1],[2],[3],[4],[5],[6],[7])","")));
915  functions.insert(make_pair(make_pair("cheb8" ,1),make_pair("ROOT::Math::Chebyshev8({V0},[0],[1],[2],[3],[4],[5],[6],[7],[8])","")));
916  functions.insert(make_pair(make_pair("cheb9" ,1),make_pair("ROOT::Math::Chebyshev9({V0},[0],[1],[2],[3],[4],[5],[6],[7],[8],[9])","")));
917  functions.insert(make_pair(make_pair("cheb10",1),make_pair("ROOT::Math::Chebyshev10({V0},[0],[1],[2],[3],[4],[5],[6],[7],[8],[9],[10])","")));
918  // 2-dimensional functions
919  functions.insert(make_pair(make_pair("gaus",2),make_pair("[0]*exp(-0.5*(({V0}-[1])/[2])^2 - 0.5*(({V1}-[3])/[4])^2)","")));
920  functions.insert(make_pair(make_pair("landau",2),make_pair("[0]*TMath::Landau({V0},[1],[2],false)*TMath::Landau({V1},[3],[4],false)","")));
921  functions.insert(make_pair(make_pair("expo",2),make_pair("exp([0]+[1]*{V0})","exp([0]+[1]*{V0}+[2]*{V1})")));
922  // gaussian with correlations
923  functions.insert(make_pair(make_pair("bigaus",2),make_pair("[0]*ROOT::Math::bigaussian_pdf({V0},{V1},[2],[4],[5],[1],[3])","[0]*ROOT::Math::bigaussian_pdf({V0},{V1},[2],[4],[5],[1],[3])")));
924 
925  map<TString,Int_t> functionsNumbers;
926  functionsNumbers["gaus"] = 100;
927  functionsNumbers["bigaus"] = 102;
928  functionsNumbers["landau"] = 400;
929  functionsNumbers["expo"] = 200;
930  functionsNumbers["crystalball"] = 500;
931 
932  // replace old names xygaus -> gaus[x,y]
933  formula.ReplaceAll("xygaus","gaus[x,y]");
934  formula.ReplaceAll("xylandau","landau[x,y]");
935  formula.ReplaceAll("xyexpo","expo[x,y]");
936  // at the moment pre-defined functions have no more than 3 dimensions
937  const char * defaultVariableNames[] = { "x","y","z"};
938 
939  for(map<pair<TString,Int_t>,pair<TString,TString> >::iterator it = functions.begin(); it != functions.end(); ++it)
940  {
941 
942  TString funName = it->first.first;
943  Int_t funDim = it->first.second;
944  Int_t funPos = formula.Index(funName);
945 
946 
947 
948  //std::cout << formula << " ---- " << funName << " " << funPos << std::endl;
949  while(funPos != kNPOS)
950  {
951 
952  // should also check that function is not something else (e.g. exponential - parse the expo)
953  Int_t lastFunPos = funPos + funName.Length();
954 
955  // check that first and last character is not alphanumeric
956  Int_t iposBefore = funPos - 1;
957  //std::cout << "looping on funpos is " << funPos << " formula is " << formula << " function " << funName << std::endl;
958  if (iposBefore >= 0) {
959  assert( iposBefore < formula.Length() );
960  if (isalpha(formula[iposBefore] ) ) {
961  //std::cout << "previous character for function " << funName << " is " << formula[iposBefore] << "- skip " << std::endl;
962  funPos = formula.Index(funName,lastFunPos);
963  continue;
964  }
965  }
966 
967  Bool_t isNormalized = false;
968  if (lastFunPos < formula.Length() ) {
969  // check if function is normalized by looking at "n" character after function name (e.g. gausn)
970  isNormalized = (formula[lastFunPos] == 'n');
971  if (isNormalized) lastFunPos += 1;
972  if (lastFunPos < formula.Length() ) {
973  char c = formula[lastFunPos];
974  // check if also last character is not alphanumeric or is not an operator and not a parenthesis ( or [.
975  // Parenthesis [] are used to express the variables
976  if ( isalnum(c ) || ( ! IsOperator(c ) && c != '(' && c != ')' && c != '[' && c != ']' ) ) {
977  //std::cout << "last character for function " << funName << " is " << c << " skip .." << std::endl;
978  funPos = formula.Index(funName,lastFunPos);
979  continue;
980  }
981  }
982  }
983 
984  if(isNormalized)
985  {
986  SetBit(kNormalized,1);
987  }
988  std::vector<TString> variables;
989  Int_t dim = 0;
990  TString varList = "";
991  Bool_t defaultVariables = false;
992 
993  // check if function has specified the [...] e.g. gaus[x,y]
994  Int_t openingBracketPos = funPos + funName.Length() + (isNormalized ? 1 : 0);
995  Int_t closingBracketPos = kNPOS;
996  if(openingBracketPos > formula.Length() || formula[openingBracketPos] != '[')
997  {
998  dim = funDim;
999  variables.resize(dim);
1000  for (Int_t idim = 0; idim < dim; ++idim)
1001  variables[idim] = defaultVariableNames[idim];
1002  defaultVariables = true;
1003  }
1004  else
1005  {
1006  // in case of [..] found, assume they specify all the variables. Use it to get function dimension
1007  closingBracketPos = formula.Index(']',openingBracketPos);
1008  varList = formula(openingBracketPos+1,closingBracketPos - openingBracketPos - 1);
1009  dim = varList.CountChar(',') + 1;
1010  variables.resize(dim);
1011  Int_t Nvar = 0;
1012  TString varName = "";
1013  for(Int_t i = 0 ; i < varList.Length(); ++i)
1014  {
1015  if(IsFunctionNameChar(varList[i]))
1016  {
1017  varName.Append(varList[i]);
1018  }
1019  if(varList[i] == ',')
1020  {
1021  variables[Nvar] = varName;
1022  varName = "";
1023  Nvar++;
1024  }
1025  }
1026  if(varName != "") // we will miss last variable
1027  {
1028  variables[Nvar] = varName;
1029  }
1030  }
1031  // chech if dimension obtained from [...] is compatible with existing pre-defined functions
1032  if(dim != funDim)
1033  {
1034  pair<TString,Int_t> key = make_pair(funName,dim);
1035  if(functions.find(key) == functions.end())
1036  {
1037  Int_t funDim = it->first.second;
1038  Error("PreProcessFormula","Dimension of function %s is detected to be of dimension %d and is not compatible with existing pre-defined function which has dim %d",
1039  funName.Data(),dim,funDim);
1040  return;
1041  }
1042  // skip the particular function found - we might find later on the corresponding pre-defined function
1043  funPos = formula.Index(funName,lastFunPos);
1044  continue;
1045  }
1046  // look now for the (..) brackets to get the parameter counter (e.g. gaus(0) + gaus(3) )
1047  // need to start for a position
1048  Int_t openingParenthesisPos = (closingBracketPos == kNPOS) ? openingBracketPos : closingBracketPos + 1;
1049  bool defaultCounter = (openingParenthesisPos > formula.Length() || formula[openingParenthesisPos] != '(');
1050 
1051  //Int_t openingParenthesisPos = formula.Index('(',funPos);
1052  //Bool_t defaultCounter = (openingParenthesisPos == kNPOS);
1053  Int_t counter;
1054  if(defaultCounter)
1055  {
1056  counter = 0;
1057  }
1058  else
1059  {
1060  counter = TString(formula(openingParenthesisPos+1,formula.Index(')',funPos) - openingParenthesisPos -1)).Atoi();
1061  }
1062  //std::cout << "openingParenthesisPos " << openingParenthesisPos << " counter is " << counter << std::endl;
1063 
1064  TString body = (isNormalized ? it->second.second : it->second.first);
1065  if(isNormalized && body == "")
1066  {
1067  Error("PreprocessFormula","%d dimension function %s has no normalized form.",it->first.second,funName.Data());
1068  break;
1069  }
1070  for(int i = 0 ; i < body.Length() ; ++i)
1071  {
1072  if(body[i] == '{')
1073  {
1074  // replace {Vn} with variable names
1075  i += 2; // skip '{' and 'V'
1076  Int_t num = TString(body(i,body.Index('}',i) - i)).Atoi();
1077  TString variable = variables[num];
1078  TString pattern = TString::Format("{V%d}",num);
1079  i -= 2; // restore original position
1080  body.Replace(i, pattern.Length(),variable,variable.Length());
1081  i += variable.Length()-1; // update i to reflect change in body string
1082  }
1083  else if(body[i] == '[')
1084  {
1085  // update parameter counters in case of many functions (e.g. gaus(0)+gaus(3) )
1086  Int_t tmp = i;
1087  while(tmp < body.Length() && body[tmp] != ']')
1088  {
1089  tmp++;
1090  }
1091  Int_t num = TString(body(i+1,tmp - 1 - i)).Atoi();
1092  num += counter;
1093  TString replacement = TString::Format("%d",num);
1094 
1095  body.Replace(i+1,tmp - 1 - i,replacement,replacement.Length());
1096  i += replacement.Length() + 1;
1097  }
1098  }
1099  TString pattern;
1100  if(defaultCounter && defaultVariables)
1101  {
1102  pattern = TString::Format("%s%s",
1103  funName.Data(),
1104  (isNormalized ? "n" : ""));
1105  }
1106  if(!defaultCounter && defaultVariables)
1107  {
1108  pattern = TString::Format("%s%s(%d)",
1109  funName.Data(),
1110  (isNormalized ? "n" : ""),
1111  counter);
1112  }
1113  if(defaultCounter && !defaultVariables)
1114  {
1115  pattern = TString::Format("%s%s[%s]",
1116  funName.Data(),
1117  (isNormalized ? "n":""),
1118  varList.Data());
1119  }
1120  if(!defaultCounter && !defaultVariables)
1121  {
1122  pattern = TString::Format("%s%s[%s](%d)",
1123  funName.Data(),
1124  (isNormalized ? "n" : ""),
1125  varList.Data(),
1126  counter);
1127  }
1128  TString replacement = body;
1129 
1130  // set the number (only in case a function exists without anything else
1131  if (fNumber == 0 && formula.Length() <= (pattern.Length()-funPos) +1 ) { // leave 1 extra
1132  fNumber = functionsNumbers[funName] + 10*(dim-1);
1133  }
1134 
1135  //std::cout << " replace " << pattern << " with " << replacement << std::endl;
1136 
1137  formula.Replace(funPos,pattern.Length(),replacement,replacement.Length());
1138 
1139  funPos = formula.Index(funName);
1140  }
1141  //std::cout << " End loop of " << funName << " formula is now " << formula << std::endl;
1142  }
1143 
1144 }
1146 {
1147  //*-*
1148  //*-* Handling exponentiation
1149  //*-* Can handle multiple carets, eg.
1150  //*-* 2^3^4 will be treated like 2^(3^4)
1151  //*-*
1152  Int_t caretPos = formula.Last('^');
1153  while(caretPos != kNPOS)
1154  {
1155 
1156  TString right,left;
1157  Int_t temp = caretPos;
1158  temp--;
1159  // get the expression in ( ) which has the operator^ applied
1160  if(formula[temp] == ')')
1161  {
1162  Int_t depth = 1;
1163  temp--;
1164  while(depth != 0 && temp > 0)
1165  {
1166  if(formula[temp] == ')')
1167  depth++;
1168  if(formula[temp] == '(')
1169  depth--;
1170  temp--;
1171  }
1172  if (depth == 0) temp++;
1173  }
1174  // this in case of someting like sin(x+2)^2
1175  do {
1176  temp--; // go down one
1177  // handle scientific notation cases (1.e-2 ^ 3 )
1178  if (temp>=2 && IsScientificNotation(formula, temp-1) ) temp-=3;
1179  }
1180  while(temp >= 0 && !IsOperator(formula[temp]) && !IsBracket(formula[temp]) );
1181 
1182  assert(temp+1 >= 0);
1183  Int_t leftPos = temp+1;
1184  left = formula(leftPos, caretPos - leftPos);
1185  //std::cout << "left to replace is " << left << std::endl;
1186 
1187  // look now at the expression after the ^ operator
1188  temp = caretPos;
1189  temp++;
1190  if (temp >= formula.Length() ) {
1191  Error("HandleExponentiation","Invalid position of operator ^");
1192  return;
1193  }
1194  if(formula[temp] == '(')
1195  {
1196  Int_t depth = 1;
1197  temp++;
1198  while(depth != 0 && temp < formula.Length())
1199  {
1200  if(formula[temp] == ')')
1201  depth--;
1202  if(formula[temp] == '(')
1203  depth++;
1204  temp++;
1205  }
1206  temp--;
1207  }
1208  else {
1209  // handle case first character is operator - or + continue
1210  if (formula[temp] == '-' || formula[temp] == '+' ) temp++;
1211  // handle cases x^-2 or x^+2
1212  // need to handle also cases x^sin(x+y)
1213  Int_t depth = 0;
1214  // stop right expression if is an operator or if is a ")" from a zero depth
1215  while(temp < formula.Length() && ( (depth > 0) || !IsOperator(formula[temp]) ) )
1216  {
1217  temp++;
1218  // handle scientific notation cases (1.e-2 ^ 3 )
1219  if (temp>=2 && IsScientificNotation(formula, temp) ) temp+=2;
1220  // for internal parenthesis
1221  if (temp < formula.Length() && formula[temp] == '(') depth++;
1222  if (temp < formula.Length() && formula[temp] == ')') {
1223  if (depth > 0)
1224  depth--;
1225  else
1226  break; // case of end of a previously started expression e.g. sin(x^2)
1227  }
1228  }
1229  }
1230  right = formula(caretPos + 1, (temp - 1) - caretPos );
1231  //std::cout << "right to replace is " << right << std::endl;
1232 
1233  TString pattern = TString::Format("%s^%s",left.Data(),right.Data());
1234  TString replacement = TString::Format("pow(%s,%s)",left.Data(),right.Data());
1235 
1236  //std::cout << "pattern : " << pattern << std::endl;
1237  //std::cout << "replacement : " << replacement << std::endl;
1238  formula.Replace(leftPos,pattern.Length(),replacement,replacement.Length());
1239 
1240  caretPos = formula.Last('^');
1241  }
1242 }
1243 
1244 // handle linear functions defined with the operator ++
1246 {
1247  // Handle Linear functions identified with "@" operator
1248  Int_t linPos = formula.Index("@");
1249  if (linPos == kNPOS ) return; // function is not linear
1250  Int_t NofLinParts = formula.CountChar((int)'@');
1251  assert(NofLinParts > 0);
1252  fLinearParts.reserve(NofLinParts + 1);
1253  Int_t Nlinear = 0;
1254  bool first = true;
1255  while(linPos != kNPOS)
1256  {
1257  SetBit(kLinear,1);
1258  // analyze left part only the first time
1259  Int_t temp = 0;
1260  TString left;
1261  if (first) {
1262  temp = linPos - 1;
1263  while(temp >= 0 && formula[temp] != '@')
1264  {
1265  temp--;
1266  }
1267  left = formula(temp+1,linPos - (temp +1));
1268  }
1269  temp = linPos + 1;
1270  while(temp < formula.Length() && formula[temp] != '@')
1271  {
1272  temp++;
1273  }
1274  TString right = formula(linPos+1,temp - (linPos+1));
1275 
1276  TString pattern = (first) ? TString::Format("%s@%s",left.Data(),right.Data()) : TString::Format("@%s",right.Data());
1277  TString replacement = (first) ? TString::Format("([%d]*(%s))+([%d]*(%s))",Nlinear,left.Data(),Nlinear+1,right.Data()) : TString::Format("+([%d]*(%s))",Nlinear,right.Data());
1278  Nlinear += (first) ? 2 : 1;
1279 
1280  formula.ReplaceAll(pattern,replacement);
1281  if (first) {
1282  TFormula *lin1 = new TFormula("__linear1",left,false);
1283  fLinearParts.push_back(lin1);
1284  }
1285  TFormula *lin2 = new TFormula("__linear2",right,false);
1286  fLinearParts.push_back(lin2);
1287 
1288  linPos = formula.Index("@");
1289  first = false;
1290  }
1291 }
1292 
1294 {
1295  //*-*
1296  //*-* Preprocessing of formula
1297  //*-* Replace all ** by ^, and removes spaces.
1298  //*-* Handle also parametrized functions like polN,gaus,expo,landau
1299  //*-* and exponentiation.
1300  //*-* Similar functionality should be added here.
1301  //*-*
1302  formula.ReplaceAll("**","^");
1303  formula.ReplaceAll("++","@"); // for linear functions
1304  formula.ReplaceAll(" ","");
1305  HandlePolN(formula);
1306  HandleParametrizedFunctions(formula);
1307  HandleExponentiation(formula);
1308  // "++" wil be dealt with Handle Linear
1309  HandleLinear(formula);
1310  // special case for "--" and "++"
1311  // ("++" needs to be written with whitespace that is removed before but then we re-add it again
1312  formula.ReplaceAll("--","- -");
1313  formula.ReplaceAll("++","+ +");
1314 }
1316 {
1317  // prepare the formula to be executed
1318  // normally is called with fFormula
1319 
1320  fFuncs.clear();
1321  fReadyToExecute = false;
1322  ExtractFunctors(formula);
1323 
1324  // update the expression with the new formula
1325  fFormula = formula;
1326  // save formula to parse variable and parameters for Cling
1327  fClingInput = formula;
1328  // replace all { and }
1329  fFormula.ReplaceAll("{","");
1330  fFormula.ReplaceAll("}","");
1331 
1332  //std::cout << "functors are extracted formula is " << std::endl;
1333  //std::cout << fFormula << std::endl << std::endl;
1334 
1335  fFuncs.sort();
1336  fFuncs.unique();
1337 
1338  // use inputFormula for Cling
1340 
1341  // for pre-defined functions (need after processing)
1342  if (fNumber != 0) SetPredefinedParamNames();
1343 
1345 }
1347 {
1348  //*-*
1349  //*-* Extracts functors from formula, and put them in fFuncs.
1350  //*-* Simple grammar:
1351  //*-* <function> := name(arg1,arg2...)
1352  //*-* <variable> := name
1353  //*-* <parameter> := [number]
1354  //*-* <name> := String containing lower and upper letters, numbers, underscores
1355  //*-* <number> := Integer number
1356  //*-* Operators are omitted.
1357  //*-*
1358  TString name = "";
1359  TString body = "";
1360  //printf("formula is : %s \n",formula.Data() );
1361  for(Int_t i = 0 ; i < formula.Length(); ++i )
1362  {
1363 
1364  //std::cout << "loop on character : " << i << " " << formula[i] << std::endl;
1365  // case of parameters
1366  if(formula[i] == '[')
1367  {
1368  Int_t tmp = i;
1369  i++;
1370  TString param = "";
1371  while(formula[i] != ']' && i < formula.Length())
1372  {
1373  param.Append(formula[i++]);
1374  }
1375  i++;
1376  //rename parameter name XX to pXX
1377  if (param.IsDigit() ) param.Insert(0,'p');
1378  // handle whitespace characters in parname
1379  param.ReplaceAll("\\s"," ");
1380  DoAddParameter(param,0,false);
1381  TString replacement = TString::Format("{[%s]}",param.Data());
1382  formula.Replace(tmp,i - tmp, replacement,replacement.Length());
1383  fFuncs.push_back(TFormulaFunction(param));
1384  //printf("found parameter %s \n",param.Data() );
1385  continue;
1386  }
1387  // case of strings
1388  if (formula[i] == '\"') {
1389  // look for next instance of "\"
1390  do {
1391  i++;
1392  } while(formula[i] != '\"');
1393  }
1394  // case of e or E for numbers in exponential notaton (e.g. 2.2e-3)
1395  if (IsScientificNotation(formula, i) )
1396  continue;
1397  // case of x for hexadecimal numbers
1398  if (IsHexadecimal(formula, i) ) {
1399  // find position of operator
1400  // do not check cases if character is not only a to f, but accept anything
1401  while ( !IsOperator(formula[i]) && i < formula.Length() ) {
1402  i++;
1403  }
1404  continue;
1405  }
1406 
1407 
1408  //std::cout << "investigating character : " << i << " " << formula[i] << " of formula " << formula << std::endl;
1409  // look for variable and function names. They start in C++ with alphanumeric characters
1410  if(isalpha(formula[i]) && !IsOperator(formula[i])) // not really needed to check if operator (if isalpha is not an operator)
1411  {
1412  //std::cout << "character : " << i << " " << formula[i] << " is not an operator and is alpha " << std::endl;
1413 
1414  while( IsFunctionNameChar(formula[i]) && i < formula.Length())
1415  {
1416  // need special case for separting operator ":" from scope operator "::"
1417  if (formula[i] == ':' && ( (i+1) < formula.Length() ) ) {
1418  if ( formula[i+1] == ':' ) {
1419  // case of :: (scopeOperator)
1420  name.Append("::");
1421  i+=2;
1422  continue;
1423  }
1424  else
1425  break;
1426  }
1427 
1428  name.Append(formula[i++]);
1429  }
1430  //printf(" build a name %s \n",name.Data() );
1431  if(formula[i] == '(')
1432  {
1433  i++;
1434  if(formula[i] == ')')
1435  {
1436  fFuncs.push_back(TFormulaFunction(name,body,0));
1437  name = body = "";
1438  continue;
1439  }
1440  Int_t depth = 1;
1441  Int_t args = 1; // we will miss first argument
1442  while(depth != 0 && i < formula.Length())
1443  {
1444  switch(formula[i])
1445  {
1446  case '(': depth++; break;
1447  case ')': depth--; break;
1448  case ',': if(depth == 1) args++; break;
1449  }
1450  if(depth != 0) // we don't want last ')' inside body
1451  {
1452  body.Append(formula[i++]);
1453  }
1454  }
1455  Int_t originalBodyLen = body.Length();
1456  ExtractFunctors(body);
1457  formula.Replace(i-originalBodyLen,originalBodyLen,body,body.Length());
1458  i += body.Length() - originalBodyLen;
1459  fFuncs.push_back(TFormulaFunction(name,body,args));
1460  }
1461  else
1462  {
1463 
1464  //std::cout << "check if character : " << i << " " << formula[i] << " from name " << name << " is a function " << std::endl;
1465 
1466  // check if function is provided by gROOT
1467  TObject *obj = 0;
1468  {
1470  obj = gROOT->GetListOfFunctions()->FindObject(name);
1471  }
1472  TFormula * f = dynamic_cast<TFormula*> (obj);
1473  if (!f) {
1474  // maybe object is a TF1
1475  TF1 * f1 = dynamic_cast<TF1*> (obj);
1476  if (f1) f = f1->GetFormula();
1477  }
1478  if (f) {
1479  TString replacementFormula = f->GetExpFormula();
1480  // analyze expression string
1481  //std::cout << "formula to replace for " << f->GetName() << " is " << replacementFormula << std::endl;
1482  PreProcessFormula(replacementFormula);
1483  // we need to define different parameters if we use the unnamed default parameters ([0])
1484  // I need to replace all the terms in the functor for backward compatibility of the case
1485  // f1("[0]*x") f2("[0]*x") f1+f2 - it is weird but it is better to support
1486  //std::cout << "current number of parameter is " << fNpar << std::endl;
1487  int nparOffset = 0;
1488  //if (fParams.find("0") != fParams.end() ) {
1489  // do in any case if parameters are existing
1490  std::vector<TString> newNames;
1491  if (fNpar > 0) {
1492  nparOffset = fNpar;
1493  newNames.resize(f->GetNpar() );
1494  // start from higher number to avoid overlap
1495  for (int jpar = f->GetNpar()-1; jpar >= 0; --jpar ) {
1496  // parameters name have a "p" added in front
1497  TString pj = TString(f->GetParName(jpar));
1498  if ( pj[0] == 'p' && TString(pj(1,pj.Length())).IsDigit() ) {
1499  TString oldName = TString::Format("[%s]",f->GetParName(jpar));
1500  TString newName = TString::Format("[p%d]",nparOffset+jpar);
1501  //std::cout << "replace - parameter " << f->GetParName(jpar) << " with " << newName << std::endl;
1502  replacementFormula.ReplaceAll(oldName,newName);
1503  newNames[jpar] = newName;
1504  }
1505  else
1506  newNames[jpar] = f->GetParName(jpar);
1507  }
1508  //std::cout << "after replacing params " << replacementFormula << std::endl;
1509  }
1510  ExtractFunctors(replacementFormula);
1511  //std::cout << "after re-extracting functors " << replacementFormula << std::endl;
1512 
1513  // set parameter value from replacement formula
1514  for (int jpar = 0; jpar < f->GetNpar(); ++jpar) {
1515  if (nparOffset> 0) {
1516  // parameter have an offset- so take this into accound
1517  assert((int) newNames.size() == f->GetNpar() );
1518  SetParameter(newNames[jpar], f->GetParameter(jpar) );
1519  }
1520  else
1521  // names are the same between current formula and replaced one
1522  SetParameter(f->GetParName(jpar), f->GetParameter(jpar) );
1523  }
1524  // need to add parenthesis at begin end end of replacementFormula
1525  replacementFormula.Insert(0,'(');
1526  replacementFormula.Insert(replacementFormula.Length(),')');
1527  formula.Replace(i-name.Length(),name.Length(), replacementFormula, replacementFormula.Length());
1528  // move forward the index i of the main loop
1529  i += replacementFormula.Length()-name.Length();
1530 
1531  // we have extracted all the functor for "fname"
1532  //std::cout << " i = " << i << " f[i] = " << formula[i] << " - " << formula << std::endl;
1533  name = "";
1534 
1535  continue;
1536  }
1537 
1538  // add now functor in
1539  TString replacement = TString::Format("{%s}",name.Data());
1540  formula.Replace(i-name.Length(),name.Length(),replacement,replacement.Length());
1541  i += 2;
1542  fFuncs.push_back(TFormulaFunction(name));
1543  }
1544  }
1545  name = body = "";
1546 
1547  }
1548 }
1550 {
1551  //*-*
1552  //*-* Iterates through funtors in fFuncs and performs the appropriate action.
1553  //*-* If functor has 0 arguments (has only name) can be:
1554  //*-* - variable
1555  //*-* * will be replaced with x[num], where x is an array containing value of this variable under num.
1556  //*-* - pre-defined formula
1557  //*-* * will be replaced with formulas body
1558  //*-* - constant
1559  //*-* * will be replaced with constant value
1560  //*-* - parameter
1561  //*-* * will be replaced with p[num], where p is an array containing value of this parameter under num.
1562  //*-* If has arguments it can be :
1563  //*-* - function shortcut, eg. sin
1564  //*-* * will be replaced with fullname of function, eg. sin -> TMath::Sin
1565  //*-* - function from cling environment, eg. TMath::BreitWigner(x,y,z)
1566  //*-* * first check if function exists, and has same number of arguments, then accept it and set as found.
1567  //*-* If all functors after iteration are matched with corresponding action,
1568  //*-* it inputs C++ code of formula into cling, and sets flag that formula is ready to evaluate.
1569  //*-*
1570 
1571  // std::cout << "Begin: formula is " << formula << " list of functors " << fFuncs.size() << std::endl;
1572 
1573  for(list<TFormulaFunction>::iterator funcsIt = fFuncs.begin(); funcsIt != fFuncs.end(); ++funcsIt)
1574  {
1575  TFormulaFunction & fun = *funcsIt;
1576 
1577  //std::cout << "fun is " << fun.GetName() << std::endl;
1578 
1579  if(fun.fFound)
1580  continue;
1581  if(fun.IsFuncCall())
1582  {
1583  map<TString,TString>::iterator it = fFunctionsShortcuts.find(fun.GetName());
1584  if(it != fFunctionsShortcuts.end())
1585  {
1586  TString shortcut = it->first;
1587  TString full = it->second;
1588  //std::cout << " functor " << fun.GetName() << " found - replace " << shortcut << " with " << full << " in " << formula << std::endl;
1589  // replace all functors
1590  Ssiz_t index = formula.Index(shortcut,0);
1591  while ( index != kNPOS) {
1592  // check that function is not in a namespace and is not in other characters
1593  //std::cout << "analyzing " << shortcut << " in " << formula << std::endl;
1594  Ssiz_t i2 = index + shortcut.Length();
1595  if ( (index > 0) && (isalpha( formula[index-1] ) || formula[index-1] == ':' )) {
1596  index = formula.Index(shortcut,i2);
1597  continue;
1598  }
1599  if (i2 < formula.Length() && formula[i2] != '(') {
1600  index = formula.Index(shortcut,i2);
1601  continue;
1602  }
1603  // now replace the string
1604  formula.Replace(index, shortcut.Length(), full);
1605  Ssiz_t inext = index + full.Length();
1606  index = formula.Index(shortcut,inext);
1607  fun.fFound = true;
1608  }
1609  }
1610  if(fun.fName.Contains("::")) // add support for nested namespaces
1611  {
1612  // look for last occurence of "::"
1613  std::string name(fun.fName);
1614  size_t index = name.rfind("::");
1615  assert(index != std::string::npos);
1616  TString className = fun.fName(0,fun.fName(0,index).Length());
1617  TString functionName = fun.fName(index + 2, fun.fName.Length());
1618 
1619  Bool_t silent = true;
1620  TClass *tclass = TClass::GetClass(className,silent);
1621  // std::cout << "looking for class " << className << std::endl;
1622  const TList *methodList = tclass->GetListOfAllPublicMethods();
1623  TIter next(methodList);
1624  TMethod *p;
1625  while ((p = (TMethod*) next()))
1626  {
1627  if (strcmp(p->GetName(),functionName.Data()) == 0 &&
1628  (fun.GetNargs() <= p->GetNargs() && fun.GetNargs() >= p->GetNargs() - p->GetNargsOpt() ) )
1629  {
1630  fun.fFound = true;
1631  break;
1632  }
1633  }
1634  }
1635  if(!fun.fFound)
1636  {
1637  // try to look into all the global functions in gROOT
1638  TFunction* f;
1639  {
1641  f = (TFunction*) gROOT->GetListOfGlobalFunctions(true)->FindObject(fun.fName);
1642  }
1643  // if found a function with matching arguments
1644  if (f && fun.GetNargs() <= f->GetNargs() && fun.GetNargs() >= f->GetNargs() - f->GetNargsOpt() )
1645  {
1646  fun.fFound = true;
1647  }
1648  }
1649 
1650  if(!fun.fFound)
1651  {
1652  // ignore not found functions
1653  if (gDebug)
1654  Info("TFormula","Could not find %s function with %d argument(s)",fun.GetName(),fun.GetNargs());
1655  fun.fFound = false;
1656  }
1657  }
1658  else
1659  {
1660  TFormula* old = 0;
1661  {
1663  old = (TFormula*)gROOT->GetListOfFunctions()->FindObject(gNamePrefix + fun.fName);
1664  }
1665  if(old)
1666  {
1667  // we should not go here (this analysis is done before in ExtractFunctors)
1668  assert(false);
1669  fun.fFound = true;
1670  TString pattern = TString::Format("{%s}",fun.GetName());
1671  TString replacement = old->GetExpFormula();
1672  PreProcessFormula(replacement);
1673  ExtractFunctors(replacement);
1674  formula.ReplaceAll(pattern,replacement);
1675  continue;
1676  }
1677  // looking for default variables defined in fVars
1678 
1679  map<TString,TFormulaVariable>::iterator varsIt = fVars.find(fun.GetName());
1680  if(varsIt!= fVars.end())
1681  {
1682 
1683  TString name = (*varsIt).second.GetName();
1684  Double_t value = (*varsIt).second.fValue;
1685 
1686 
1687  AddVariable(name,value); // this set the cling variable
1688  if(!fVars[name].fFound)
1689  {
1690 
1691 
1692  fVars[name].fFound = true;
1693  int varDim = (*varsIt).second.fArrayPos; // variable dimenions (0 for x, 1 for y, 2, for z)
1694  if (varDim >= fNdim) {
1695  fNdim = varDim+1;
1696 
1697  // we need to be sure that all other variables are added with position less
1698  for ( auto &v : fVars) {
1699  if (v.second.fArrayPos < varDim && !v.second.fFound ) {
1700  AddVariable(v.first, v.second.fValue);
1701  v.second.fFound = true;
1702  }
1703  }
1704  }
1705  }
1706  // remove the "{.. }" added around the variable
1707  TString pattern = TString::Format("{%s}",name.Data());
1708  TString replacement = TString::Format("x[%d]",(*varsIt).second.fArrayPos);
1709  formula.ReplaceAll(pattern,replacement);
1710 
1711  //std::cout << "Found an observable for " << fun.GetName() << std::endl;
1712 
1713  fun.fFound = true;
1714  continue;
1715  }
1716  // check for observables defined as x[0],x[1],....
1717  // maybe could use a regular expression here
1718  // only in case match with defined variables is not successfull
1719  TString funname = fun.GetName();
1720  if (funname.Contains("x[") && funname.Contains("]") ) {
1721  TString sdigit = funname(2,funname.Index("]") );
1722  int digit = sdigit.Atoi();
1723  if (digit >= fNdim) {
1724  fNdim = digit+1;
1725  // we need to add the variables in fVars all of them before x[n]
1726  for (int j = 0; j < fNdim; ++j) {
1727  TString vname = TString::Format("x[%d]",j);
1728  if (fVars.find(vname) == fVars.end() ) {
1729  fVars[vname] = TFormulaVariable(vname,0,j);
1730  fVars[vname].fFound = true;
1731  AddVariable(vname,0.);
1732  }
1733  }
1734  }
1735  //std::cout << "Found matching observable for " << funname << std::endl;
1736  fun.fFound = true;
1737  // remove the "{.. }" added around the variable
1738  TString pattern = TString::Format("{%s}",funname.Data());
1739  formula.ReplaceAll(pattern,funname);
1740  continue;
1741  }
1742  //}
1743 
1744  auto paramsIt = fParams.find(fun.GetName());
1745  if(paramsIt != fParams.end())
1746  {
1747  //TString name = (*paramsIt).second.GetName();
1748  TString pattern = TString::Format("{[%s]}",fun.GetName());
1749  //std::cout << "pattern is " << pattern << std::endl;
1750  if(formula.Index(pattern) != kNPOS)
1751  {
1752  //TString replacement = TString::Format("p[%d]",(*paramsIt).second.fArrayPos);
1753  TString replacement = TString::Format("p[%d]",(*paramsIt).second);
1754  //std::cout << "replace pattern " << pattern << " with " << replacement << std::endl;
1755  formula.ReplaceAll(pattern,replacement);
1756 
1757  }
1758  fun.fFound = true;
1759  continue;
1760  }
1761  else {
1762  //std::cout << "functor " << fun.GetName() << " is not a parameter " << std::endl;
1763  }
1764 
1765  // looking for constants (needs to be done after looking at the parameters)
1766  map<TString,Double_t>::iterator constIt = fConsts.find(fun.GetName());
1767  if(constIt != fConsts.end())
1768  {
1769  TString pattern = TString::Format("{%s}",fun.GetName());
1770  TString value = TString::Format("%lf",(*constIt).second);
1771  formula.ReplaceAll(pattern,value);
1772  fun.fFound = true;
1773  //std::cout << "constant with name " << fun.GetName() << " is found " << std::endl;
1774  continue;
1775  }
1776 
1777 
1778  fun.fFound = false;
1779  }
1780  }
1781  //std::cout << "End: formula is " << formula << std::endl;
1782 
1783  // ignore case of functors have been matched - try to pass it to Cling
1784  if(!fReadyToExecute)
1785  {
1786  fReadyToExecute = true;
1787  Bool_t hasVariables = (fNdim > 0);
1788  Bool_t hasParameters = (fNpar > 0);
1789  if(!hasParameters)
1790  {
1791  fAllParametersSetted = true;
1792  }
1793  // assume a function without variables is always 1-dimensional
1794  if (hasParameters && ! hasVariables) {
1795  fNdim = 1;
1796  AddVariable("x",0);
1797  hasVariables = true;
1798  }
1799  Bool_t hasBoth = hasVariables && hasParameters;
1800  Bool_t inputIntoCling = (formula.Length() > 0);
1801  if (inputIntoCling) {
1802 
1803  // save copy of inputFormula in a std::strig for the unordered map
1804  // and also formula is same as FClingInput typically and it will be modified
1805  std::string inputFormula = std::string(formula);
1806 
1807 
1808  // valid input formula - try to put into Cling
1809  TString argumentsPrototype =
1810  TString::Format("%s%s%s",(hasVariables ? "Double_t *x" : ""), (hasBoth ? "," : ""),
1811  (hasParameters ? "Double_t *p" : ""));
1812 
1813 
1814  // set the name for Cling using the hash_function
1815  fClingName = gNamePrefix;
1816 
1817  // check if formula exist already in the map
1819 
1820  auto funcit = gClingFunctions.find(inputFormula);
1821 
1822  if (funcit != gClingFunctions.end() ) {
1824  fClingInitialized = true;
1825  inputIntoCling = false;
1826  }
1827 
1828  // set the cling name using hash of the static formulae map
1829  auto hasher = gClingFunctions.hash_function();
1830  fClingName = TString::Format("%s__id%zu",gNamePrefix.Data(),(unsigned long) hasher(inputFormula) );
1831 
1832  fClingInput = TString::Format("Double_t %s(%s){ return %s ; }", fClingName.Data(),argumentsPrototype.Data(),inputFormula.c_str());
1833 
1834  // this is not needed (maybe can be re-added in case of recompilation of identical expressions
1835  // // check in case of a change if need to re-initialize
1836  // if (fClingInitialized) {
1837  // if (oldClingInput == fClingInput)
1838  // inputIntoCling = false;
1839  // else
1840  // fClingInitialized = false;
1841  // }
1842 
1843 
1844  if(inputIntoCling) {
1846  if (fClingInitialized) {
1847  // if Cling has been succesfully initialized
1848  // dave function ptr in the static map
1850  gClingFunctions.insert ( std::make_pair ( inputFormula, (void*) fFuncPtr) );
1851  }
1852 
1853  }
1854  else {
1855  fAllParametersSetted = true;
1856  fClingInitialized = true;
1857  }
1858  }
1859  }
1860 
1861 
1862  // IN case of a Cling Error check components wich are not found in Cling
1863  // check that all formula components arematched otherwise emit an error
1864  if (!fClingInitialized) {
1865  Bool_t allFunctorsMatched = true;
1866  for(list<TFormulaFunction>::iterator it = fFuncs.begin(); it != fFuncs.end(); it++)
1867  {
1868  if(!it->fFound)
1869  {
1870  allFunctorsMatched = false;
1871  if (it->GetNargs() == 0)
1872  Error("ProcessFormula","\"%s\" has not been matched in the formula expression",it->GetName() );
1873  else
1874  Error("ProcessFormula","Could not find %s function with %d argument(s)",it->GetName(),it->GetNargs());
1875  }
1876  }
1877  if (!allFunctorsMatched) {
1878  Error("ProcessFormula","Formula \"%s\" is invalid !", GetExpFormula().Data() );
1879  fReadyToExecute = false;
1880  }
1881  }
1882 
1883  // clean up un-used default variables in case formula is valid
1885  auto itvar = fVars.begin();
1886  do
1887  {
1888  if ( ! itvar->second.fFound ) {
1889  //std::cout << "Erase variable " << itvar->first << std::endl;
1890  itvar = fVars.erase(itvar);
1891  }
1892  else
1893  itvar++;
1894  }
1895  while( itvar != fVars.end() );
1896  }
1897 
1898 }
1900 
1901  // set parameter names only in case of pre-defined functions
1902  if (fNumber == 0) return;
1903 
1904  if (fNumber == 100) { // Gaussian
1905  SetParName(0,"Constant");
1906  SetParName(1,"Mean");
1907  SetParName(2,"Sigma");
1908  return;
1909  }
1910  if (fNumber == 110) {
1911  SetParName(0,"Constant");
1912  SetParName(1,"MeanX");
1913  SetParName(2,"SigmaX");
1914  SetParName(3,"MeanY");
1915  SetParName(4,"SigmaY");
1916  return;
1917  }
1918  if (fNumber == 112) { // bigaus
1919  SetParName(0,"Constant");
1920  SetParName(1,"MeanX");
1921  SetParName(2,"SigmaX");
1922  SetParName(3,"MeanY");
1923  SetParName(4,"SigmaY");
1924  SetParName(5,"Rho");
1925  return;
1926  }
1927  if (fNumber == 200) { // exponential
1928  SetParName(0,"Constant");
1929  SetParName(1,"Slope");
1930  return;
1931  }
1932  if (fNumber == 400) { // landau
1933  SetParName(0,"Constant");
1934  SetParName(1,"MPV");
1935  SetParName(2,"Sigma");
1936  return;
1937  }
1938  if (fNumber == 500) { // crystal-ball
1939  SetParName(0,"Constant");
1940  SetParName(1,"Mean");
1941  SetParName(2,"Sigma");
1942  SetParName(3,"Alpha");
1943  SetParName(4,"N");
1944  return;
1945  }
1946  if (fNumber == 600) { // breit-wigner
1947  SetParName(0,"Constant");
1948  SetParName(1,"Mean");
1949  SetParName(2,"Gamma");
1950  return;
1951  }
1952  // if formula is a polynomial (or chebyshev), set parameter names
1953  // not needed anymore (p0 is assigned by default)
1954  // if (fNumber == (300+fNpar-1) ) {
1955  // for (int i = 0; i < fNpar; i++) SetParName(i,TString::Format("p%d",i));
1956  // return;
1957  // }
1958 
1959  // // general case if parameters are digits (XX) change to pXX
1960  // auto paramMap = fParams; // need to copy the map because SetParName is going to modify it
1961  // for ( auto & p : paramMap) {
1962  // if (p.first.IsDigit() )
1963  // SetParName(p.second,TString::Format("p%s",p.first.Data()));
1964  // }
1965 
1966  return;
1967 }
1969 {
1970  // Return linear part.
1971 
1972  if (!fLinearParts.empty()) {
1973  int n = fLinearParts.size();
1974  if (i < 0 || i >= n ) {
1975  Error("GetLinearPart","Formula %s has only %d linear parts - requested %d",GetName(),n,i);
1976  return nullptr;
1977  }
1978  return fLinearParts[i];
1979  }
1980  return nullptr;
1981 }
1983 {
1984  //*-*
1985  //*-* Adds variable to known variables, and reprocess formula.
1986  //*-*
1987 
1988  if(fVars.find(name) != fVars.end() )
1989  {
1990  TFormulaVariable & var = fVars[name];
1991  var.fValue = value;
1992 
1993  // If the position is not defined in the Cling vectors, make space for it
1994  // but normally is variable is defined in fVars a slot should be also present in fClingVariables
1995  if(var.fArrayPos < 0)
1996  {
1997  var.fArrayPos = fVars.size();
1998  }
1999  if(var.fArrayPos >= (int)fClingVariables.size())
2000  {
2001  fClingVariables.resize(var.fArrayPos+1);
2002  }
2004  }
2005  else
2006  {
2007  TFormulaVariable var(name,value,fVars.size());
2008  fVars[name] = var;
2009  fClingVariables.push_back(value);
2010  if (!fFormula.IsNull() ) {
2011  //printf("process formula again - %s \n",fClingInput.Data() );
2013  }
2014  }
2015 
2016 }
2017 void TFormula::AddVariables(const TString *vars, const Int_t size)
2018 {
2019  //*-*
2020  //*-* Adds multiple variables.
2021  //*-* First argument is an array of pairs<TString,Double>, where
2022  //*-* first argument is name of variable,
2023  //*-* second argument represents value.
2024  //*-* size - number of variables passed in first argument
2025  //*-*
2026 
2027  Bool_t anyNewVar = false;
2028  for(Int_t i = 0 ; i < size; ++i)
2029  {
2030 
2031  const TString & vname = vars[i];
2032 
2033  TFormulaVariable &var = fVars[vname];
2034  if(var.fArrayPos < 0)
2035  {
2036 
2037  var.fName = vname;
2038  var.fArrayPos = fVars.size();
2039  anyNewVar = true;
2040  var.fValue = 0;
2041  if(var.fArrayPos >= (int)fClingVariables.capacity())
2042  {
2043  Int_t multiplier = 2;
2044  if(fFuncs.size() > 100)
2045  {
2046  multiplier = TMath::Floor(TMath::Log10(fFuncs.size()) * 10);
2047  }
2048  fClingVariables.reserve(multiplier * fClingVariables.capacity());
2049  }
2050  fClingVariables.push_back(0.0);
2051  }
2052  // else
2053  // {
2054  // var.fValue = v.second;
2055  // fClingVariables[var.fArrayPos] = v.second;
2056  // }
2057  }
2058  if(anyNewVar && !fFormula.IsNull())
2059  {
2061  }
2062 
2063 }
2064 
2065 void TFormula::SetName(const char* name)
2066 {
2067  // Set the name of the formula. We need to allow the list of function to
2068  // properly handle the hashes.
2069  if (IsReservedName(name)) {
2070  Error("SetName","The name \'%s\' is reserved as a TFormula variable name.\n"
2071  "\tThis function will not be renamed.",name);
2072  } else {
2073  // Here we need to remove and re-add to keep the hashes consistent with
2074  // the underlying names.
2075  auto listOfFunctions = gROOT->GetListOfFunctions();
2076  TObject* thisAsFunctionInList = nullptr;
2078  if (listOfFunctions){
2079  thisAsFunctionInList = listOfFunctions->FindObject(this);
2080  if (thisAsFunctionInList) listOfFunctions->Remove(thisAsFunctionInList);
2081  }
2082  TNamed::SetName(name);
2083  if (thisAsFunctionInList) listOfFunctions->Add(thisAsFunctionInList);
2084  }
2085 }
2086 
2087 void TFormula::SetVariables(const pair<TString,Double_t> *vars, const Int_t size)
2088 {
2089  //*-*
2090  //*-* Sets multiple variables.
2091  //*-* First argument is an array of pairs<TString,Double>, where
2092  //*-* first argument is name of variable,
2093  //*-* second argument represents value.
2094  //*-* size - number of variables passed in first argument
2095  //*-*
2096  for(Int_t i = 0; i < size; ++i)
2097  {
2098  pair<TString,Double_t> v = vars[i];
2099  if(fVars.find(v.first) != fVars.end())
2100  {
2101  fVars[v.first].fValue = v.second;
2102  fClingVariables[fVars[v.first].fArrayPos] = v.second;
2103  }
2104  else
2105  {
2106  Error("SetVariables","Variable %s is not defined.",v.first.Data());
2107  }
2108  }
2109 }
2110 
2112 {
2113  //*-*
2114  //*-* Returns variable value.
2115  //*-*
2116  TString sname(name);
2117  if(fVars.find(sname) == fVars.end())
2118  {
2119  Error("GetVariable","Variable %s is not defined.",sname.Data());
2120  return -1;
2121  }
2122  return fVars.find(sname)->second.fValue;
2123 }
2125 {
2126  //*-*
2127  //*-* Returns variable number (positon in array) given its name
2128  //*-*
2129  TString sname(name);
2130  if(fVars.find(sname) == fVars.end())
2131  {
2132  Error("GetVarNumber","Variable %s is not defined.",sname.Data());
2133  return -1;
2134  }
2135  return fVars.find(sname)->second.fArrayPos;
2136 }
2137 
2139 {
2140  //*-*
2141  //*-* Returns variable name given its position in the array
2142  //*-*
2143 
2144  if (ivar < 0 || ivar >= fNdim) return "";
2145 
2146  // need to loop on the map to find corresponding variable
2147  for ( auto & v : fVars) {
2148  if (v.second.fArrayPos == ivar) return v.first;
2149  }
2150  Error("GetVarName","Variable with index %d not found !!",ivar);
2151  //return TString::Format("x%d",ivar);
2152  return TString();
2153 }
2154 
2156 {
2157  //*-*
2158  //*-* Sets variable value.
2159  //*-*
2160  if(fVars.find(name) == fVars.end())
2161  {
2162  Error("SetVariable","Variable %s is not defined.",name.Data());
2163  return;
2164  }
2165  fVars[name].fValue = value;
2166  fClingVariables[fVars[name].fArrayPos] = value;
2167 }
2168 
2170 {
2171  //*-*
2172  //*-* Adds parameter to known parameters.
2173  //*-* User should use SetParameter, because parameters are added during initialization part,
2174  //*-* and after that adding new will be pointless.
2175  //*-*
2176 
2177  //std::cout << "adding parameter " << name << std::endl;
2178 
2179  // if parameter is already defined in fParams - just set the new value
2180  if(fParams.find(name) != fParams.end() )
2181  {
2182  int ipos = fParams[name];
2183  //TFormulaVariable & par = fParams[name];
2184  //par.fValue = value;
2185  if (ipos < 0)
2186  {
2187  ipos = fParams.size();
2188  fParams[name] = ipos;
2189  }
2190 //
2191  if(ipos >= (int)fClingParameters.size())
2192  {
2193  if(ipos >= (int)fClingParameters.capacity())
2194  fClingParameters.reserve( TMath::Max(int(fParams.size()), ipos+1));
2195  fClingParameters.insert(fClingParameters.end(),ipos+1-fClingParameters.size(),0.0);
2196  }
2197  fClingParameters[ipos] = value;
2198  }
2199  else
2200  {
2201  // new parameter defined
2202  fNpar++;
2203  //TFormulaVariable(name,value,fParams.size());
2204  int pos = fParams.size();
2205  //fParams.insert(std::make_pair<TString,TFormulaVariable>(name,TFormulaVariable(name,value,pos)));
2206  auto ret = fParams.insert(std::make_pair(name,pos));
2207  // map returns a std::pair<iterator, bool>
2208  // use map the order for defult position of parameters in the vector
2209  // (i.e use the alphabetic order)
2210  if (ret.second) {
2211  // a new element is inserted
2212  if (ret.first == fParams.begin() )
2213  pos = 0;
2214  else {
2215  auto previous = (ret.first);
2216  --previous;
2217  pos = previous->second + 1;
2218  }
2219 
2220 
2221  if (pos < (int)fClingParameters.size() )
2222  fClingParameters.insert(fClingParameters.begin()+pos,value);
2223  else {
2224  // this should not happen
2225  if (pos > (int)fClingParameters.size() )
2226  Warning("inserting parameter %s at pos %d when vector size is %d \n",name.Data(),pos,(int)fClingParameters.size() );
2227 
2228  if(pos >= (int)fClingParameters.capacity())
2229  fClingParameters.reserve( TMath::Max(int(fParams.size()), pos+1));
2230  fClingParameters.insert(fClingParameters.end(),pos+1-fClingParameters.size(),0.0);
2231  fClingParameters[pos] = value;
2232  }
2233 
2234  // need to adjust all other positions
2235  for ( auto it = ret.first; it != fParams.end(); ++it ) {
2236  it->second = pos;
2237  pos++;
2238  }
2239  // for (auto & p : fParams)
2240  // std::cout << "Parameter " << p.first << " position " << p.second << std::endl;
2241  // printf("inserted parameters size params %d size cling %d \n",fParams.size(), fClingParameters.size() );
2242  }
2243  if (processFormula) {
2244  // replace first in input parameter name with [name]
2245  fClingInput.ReplaceAll(name,TString::Format("[%s]",name.Data() ) );
2247  }
2248  }
2249 
2250 }
2251 Int_t TFormula::GetParNumber(const char * name) const {
2252  // return parameter index given a name (return -1 for not existing parameters)
2253  // non need to print an error
2254  auto it = fParams.find(name);
2255  if(it == fParams.end())
2256  {
2257  return -1;
2258  }
2259  return it->second;
2260 
2261 }
2262 
2264 {
2265  //*-*
2266  //*-* Returns parameter value given by string.
2267  //*-*
2268  int i = GetParNumber(name);
2269  if (i == -1) {
2270  Error("GetParameter","Parameter %s is not defined.",name);
2271  return TMath::QuietNaN();
2272  }
2273 
2274  return GetParameter( GetParNumber(name) );
2275 }
2277 {
2278  //*-*
2279  //*-* Return parameter value given by integer.
2280  //*-*
2281  //*-*
2282  //TString name = TString::Format("%d",param);
2283  if(param >=0 && param < (int) fClingParameters.size())
2284  return fClingParameters[param];
2285  Error("GetParameter","wrong index used - use GetParameter(name)");
2286  return TMath::QuietNaN();
2287 }
2288 const char * TFormula::GetParName(Int_t ipar) const
2289 {
2290  //*-*
2291  //*-* Return parameter name given by integer.
2292  //*-*
2293  if (ipar < 0 || ipar >= fNpar) return "";
2294 
2295  // need to loop on the map to find corresponding parameter
2296  for ( auto & p : fParams) {
2297  if (p.second == ipar) return p.first.Data();
2298  }
2299  Error("GetParName","Parameter with index %d not found !!",ipar);
2300  //return TString::Format("p%d",ipar);
2301  return TString();
2302 }
2304 {
2305  if(!fClingParameters.empty())
2306  return const_cast<Double_t*>(&fClingParameters[0]);
2307  return 0;
2308 }
2309 
2311 {
2312  for(Int_t i = 0; i < fNpar; ++i)
2313  {
2314  if (Int_t(fClingParameters.size()) > i)
2315  params[i] = fClingParameters[i];
2316  else
2317  params[i] = -1;
2318  }
2319 }
2320 
2322 {
2323  //*-*
2324  //*-* Sets parameter value.
2325  //*-*
2326 
2327  SetParameter( GetParNumber(name), value);
2328 
2329  // do we need this ???
2330 #ifdef OLDPARAMS
2331  if(fParams.find(name) == fParams.end())
2332  {
2333  Error("SetParameter","Parameter %s is not defined.",name.Data());
2334  return;
2335  }
2336  fParams[name].fValue = value;
2337  fParams[name].fFound = true;
2338  fClingParameters[fParams[name].fArrayPos] = value;
2339  fAllParametersSetted = true;
2340  for(map<TString,TFormulaVariable>::iterator it = fParams.begin(); it != fParams.end(); it++)
2341  {
2342  if(!it->second.fFound)
2343  {
2344  fAllParametersSetted = false;
2345  break;
2346  }
2347  }
2348 #endif
2349 }
2350 #ifdef OLDPARAMS
2351 void TFormula::SetParameters(const pair<TString,Double_t> *params,const Int_t size)
2352 {
2353  //*-*
2354  //*-* Set multiple parameters.
2355  //*-* First argument is an array of pairs<TString,Double>, where
2356  //*-* first argument is name of parameter,
2357  //*-* second argument represents value.
2358  //*-* size - number of params passed in first argument
2359  //*-*
2360  for(Int_t i = 0 ; i < size ; ++i)
2361  {
2362  pair<TString,Double_t> p = params[i];
2363  if(fParams.find(p.first) == fParams.end())
2364  {
2365  Error("SetParameters","Parameter %s is not defined",p.first.Data());
2366  continue;
2367  }
2368  fParams[p.first].fValue = p.second;
2369  fParams[p.first].fFound = true;
2370  fClingParameters[fParams[p.first].fArrayPos] = p.second;
2371  }
2372  fAllParametersSetted = true;
2373  for(map<TString,TFormulaVariable>::iterator it = fParams.begin(); it != fParams.end(); it++)
2374  {
2375  if(!it->second.fFound)
2376  {
2377  fAllParametersSetted = false;
2378  break;
2379  }
2380  }
2381 }
2382 #endif
2383 
2384 void TFormula::DoSetParameters(const Double_t *params, Int_t size)
2385 {
2386  if(!params || size < 0 || size > fNpar) return;
2387  // reset vector of cling parameters
2388  if (size != (int) fClingParameters.size() ) {
2389  Warning("SetParameters","size is not same of cling parameter size %d - %d",size,int(fClingParameters.size()) );
2390  for(Int_t i = 0; i < size; ++i)
2391  {
2392  TString name = TString::Format("%d",i);
2393  SetParameter(name,params[i]);
2394  }
2395  return;
2396  }
2397  fAllParametersSetted = true;
2398  std::copy(params, params+size, fClingParameters.begin() );
2399 }
2400 
2402 {
2403  // set a vector of parameters value
2404  // Order in the vector is by default the aphabetic order given to the parameters
2405  // apart if the users has defined explicitly the parameter names
2406  DoSetParameters(params,fNpar);
2407 }
2409  Double_t p5,Double_t p6,Double_t p7,Double_t p8,
2410  Double_t p9,Double_t p10)
2411 {
2412  // Set a list of parameters.
2413  // The order is by default the aphabetic order given to the parameters
2414  // apart if the users has defined explicitly the parameter names
2415  if(fNpar >= 1) SetParameter(0,p0);
2416  if(fNpar >= 2) SetParameter(1,p1);
2417  if(fNpar >= 3) SetParameter(2,p2);
2418  if(fNpar >= 4) SetParameter(3,p3);
2419  if(fNpar >= 5) SetParameter(4,p4);
2420  if(fNpar >= 6) SetParameter(5,p5);
2421  if(fNpar >= 7) SetParameter(6,p6);
2422  if(fNpar >= 8) SetParameter(7,p7);
2423  if(fNpar >= 9) SetParameter(8,p8);
2424  if(fNpar >= 10) SetParameter(9,p9);
2425  if(fNpar >= 11) SetParameter(10,p10);
2426 }
2428 {
2429  // Set a parameter given a parameter index
2430  // The parameter index is by default the aphabetic order given to the parameters
2431  // apart if the users has defined explicitly the parameter names
2432  if (param < 0 || param >= fNpar) return;
2433  assert(int(fClingParameters.size()) == fNpar);
2434  fClingParameters[param] = value;
2435  // TString name = TString::Format("%d",param);
2436  // SetParameter(name,value);
2437 }
2438 void TFormula::SetParNames(const char *name0,const char *name1,const char *name2,const char *name3,
2439  const char *name4, const char *name5,const char *name6,const char *name7,
2440  const char *name8,const char *name9,const char *name10)
2441 {
2442  if(fNpar >= 1) SetParName(0,name0);
2443  if(fNpar >= 2) SetParName(1,name1);
2444  if(fNpar >= 3) SetParName(2,name2);
2445  if(fNpar >= 4) SetParName(3,name3);
2446  if(fNpar >= 5) SetParName(4,name4);
2447  if(fNpar >= 6) SetParName(5,name5);
2448  if(fNpar >= 7) SetParName(6,name6);
2449  if(fNpar >= 8) SetParName(7,name7);
2450  if(fNpar >= 9) SetParName(8,name8);
2451  if(fNpar >= 10) SetParName(9,name9);
2452  if(fNpar >= 11) SetParName(10,name10);
2453 }
2454 void TFormula::SetParName(Int_t ipar, const char * name)
2455 {
2456 
2457  if (ipar < 0 || ipar > fNpar) {
2458  Error("SetParName","Wrong Parameter index %d ",ipar);
2459  return;
2460  }
2461  TString oldName;
2462  // find parameter with given index
2463  for ( auto &it : fParams) {
2464  if (it.second == ipar) {
2465  oldName = it.first;
2466  fParams.erase(oldName);
2467  fParams.insert(std::make_pair(name, ipar) );
2468  break;
2469  }
2470  }
2471  if (oldName.IsNull() ) {
2472  Error("SetParName","Parameter %d is not existing.",ipar);
2473  return;
2474  }
2475 
2476  //replace also parameter name in formula expression
2477  ReplaceParamName(fFormula, oldName, name);
2478 
2479 }
2480 
2481 void TFormula::ReplaceParamName(TString & formula, const TString & oldName, const TString & name){
2482  // replace in Formula expression the parameter name
2483  if (!formula.IsNull() ) {
2484  bool found = false;
2485  for(list<TFormulaFunction>::iterator it = fFuncs.begin(); it != fFuncs.end(); ++it)
2486  {
2487  if(oldName == it->GetName())
2488  {
2489  found = true;
2490  it->fName = name;
2491  break;
2492  }
2493  }
2494  if(!found)
2495  {
2496  Error("SetParName","Parameter %s is not defined.",oldName.Data());
2497  return;
2498  }
2499  // change whitespace to \\s avoid problems in parsing
2500  TString newName = name;
2501  newName.ReplaceAll(" ","\\s");
2502  TString pattern = TString::Format("[%s]",oldName.Data());
2503  TString replacement = TString::Format("[%s]",newName.Data());
2504  formula.ReplaceAll(pattern,replacement);
2505  }
2506 
2507 }
2508 
2509 Double_t TFormula::EvalPar(const Double_t *x,const Double_t *params) const
2510 {
2511 
2512  return DoEval(x, params);
2513 }
2515 {
2516  //*-*
2517  //*-* Sets first 4 variables (e.g. x, y, z, t) and evaluate formula.
2518  //*-*
2519  double xxx[4] = {x,y,z,t};
2520  return DoEval(xxx);
2521 }
2523 {
2524  //*-*
2525  //*-* Sets first 3 variables (e.g. x, y, z) and evaluate formula.
2526  //*-*
2527  double xxx[3] = {x,y,z};
2528  return DoEval(xxx);
2529 }
2531 {
2532  //*-*
2533  //*-* Sets first 2 variables (e.g. x and y) and evaluate formula.
2534  //*-*
2535  double xxx[2] = {x,y};
2536  return DoEval(xxx);
2537 }
2539 {
2540  //*-*
2541  //*-* Sets first variable (e.g. x) and evaluate formula.
2542  //*-*
2543  //double xxx[1] = {x};
2544  double * xxx = &x;
2545  return DoEval(xxx);
2546 }
2547 Double_t TFormula::DoEval(const double * x, const double * params) const
2548 {
2549  //*-*
2550  //*-* Evaluate formula.
2551  //*-* If formula is not ready to execute(missing parameters/variables),
2552  //*-* print these which are not known.
2553  //*-* If parameter has default value, and has not been setted, appropriate warning is shown.
2554  //*-*
2555 
2556  if(!fReadyToExecute)
2557  {
2558  Error("Eval","Formula is invalid and not ready to execute ");
2559  for(auto it = fFuncs.begin(); it != fFuncs.end(); ++it)
2560  {
2561  TFormulaFunction fun = *it;
2562  if(!fun.fFound)
2563  {
2564  printf("%s is unknown.\n",fun.GetName());
2565  }
2566  }
2567  return TMath::QuietNaN();
2568  }
2569  if (fLambdaPtr && TestBit(TFormula::kLambda)) {// case of lambda functions
2570  std::function<double(double *, double *)> & fptr = * ( (std::function<double(double *, double *)> *) fLambdaPtr);
2571  assert(x);
2572  //double * v = (x) ? const_cast<double*>(x) : const_cast<double*>(fClingVariables.data());
2573  double * v = const_cast<double*>(x);
2574  double * p = (params) ? const_cast<double*>(params) : const_cast<double*>(fClingParameters.data());
2575  return fptr(v, p);
2576  }
2577  // this is needed when reading from a file
2578  if (!fClingInitialized) {
2579  Error("Eval","Formula is invalid or not properly initialized - try calling TFormula::Compile");
2580  return TMath::QuietNaN();
2581 #ifdef EVAL_IS_NOT_CONST
2582  // need to replace in cling the name of the pointer of this object
2583  TString oldClingName = fClingName;
2585  fClingInput.ReplaceAll(oldClingName, fClingName);
2587 #endif
2588  }
2589 
2590  Double_t result = 0;
2591  void* args[2];
2592  double * vars = (x) ? const_cast<double*>(x) : const_cast<double*>(fClingVariables.data());
2593  args[0] = &vars;
2594  if (fNpar <= 0)
2595  (*fFuncPtr)(0, 1, args, &result);
2596  else {
2597  double * pars = (params) ? const_cast<double*>(params) : const_cast<double*>(fClingParameters.data());
2598  args[1] = &pars;
2599  (*fFuncPtr)(0, 2, args, &result);
2600  }
2601  return result;
2602 }
2603 
2604 ////////////////////////////////////////////////////////////////////////////////
2605 /// return the expression formula
2606 /// If option = "P" replace the parameter names with their values
2607 /// If option = "CLING" return the actual expression used to build the function passed to cling
2608 /// If option = "CLINGP" replace in the CLING expression the parameter with their values
2609 
2611 {
2612  TString opt(option);
2613  if (opt.IsNull() || TestBit(TFormula::kLambda) ) return fFormula;
2614  opt.ToUpper();
2615 
2616  // if (opt.Contains("N") ) {
2617  // TString formula = fFormula;
2618  // ReplaceParName(formula, ....)
2619  // }
2620 
2621  if (opt.Contains("CLING") ) {
2622  std::string clingFunc = fClingInput.Data();
2623  std::size_t found = clingFunc.find("return");
2624  std::size_t found2 = clingFunc.rfind(";");
2625  if (found == std::string::npos || found2 == std::string::npos) {
2626  Error("GetExpFormula","Invalid Cling expression - return default formula expression");
2627  return fFormula;
2628  }
2629  TString clingFormula = fClingInput(found+7,found2-found-7);
2630  // to be implemented
2631  if (!opt.Contains("P")) return clingFormula;
2632  // replace all "p[" with "[parname"
2633  int i = 0;
2634  while (i < clingFormula.Length()-2 ) {
2635  // look for p[number
2636  if (clingFormula[i] == 'p' && clingFormula[i+1] == '[' && isdigit(clingFormula[i+2]) ) {
2637  int j = i+3;
2638  while ( isdigit(clingFormula[j]) ) { j++;}
2639  if (clingFormula[j] != ']') {
2640  Error("GetExpFormula","Parameters not found - invalid expression - return default cling formula");
2641  return clingFormula;
2642  }
2643  TString parNumbName = clingFormula(i+2,j-i-2);
2644  int parNumber = parNumbName.Atoi();
2645  assert(parNumber < fNpar);
2646  TString replacement = TString::Format("%f",GetParameter(parNumber));
2647  clingFormula.Replace(i,j-i+1, replacement );
2648  i += replacement.Length();
2649  }
2650  i++;
2651  }
2652  return clingFormula;
2653  }
2654  if (opt.Contains("P") ) {
2655  // replace parameter names with their values
2656  TString expFormula = fFormula;
2657  int i = 0;
2658  while (i < expFormula.Length()-2 ) {
2659  // look for [parName]
2660  if (expFormula[i] == '[') {
2661  int j = i+1;
2662  while ( expFormula[j] != ']' ) { j++;}
2663  if (expFormula[j] != ']') {
2664  Error("GetExpFormula","Parameter names not found - invalid expression - return default formula");
2665  return expFormula;
2666  }
2667  TString parName = expFormula(i+1,j-i-1);
2668  TString replacement = TString::Format("%g",GetParameter(parName));
2669  expFormula.Replace(i,j-i+1, replacement );
2670  i += replacement.Length();
2671  }
2672  i++;
2673  }
2674  return expFormula;
2675  }
2676  Warning("GetExpFormula","Invalid option - return defult formula expression");
2677  return fFormula;
2678 }
2679 ////////////////////////////////////////////////////////////////////////////////
2680 /// print the formula and its attributes
2681 
2682 void TFormula::Print(Option_t *option) const
2683 {
2684  printf(" %20s : %s Ndim= %d, Npar= %d, Number= %d \n",GetName(),GetTitle(), fNdim,fNpar,fNumber);
2685  printf(" Formula expression: \n");
2686  printf("\t%s \n",fFormula.Data() );
2687  TString opt(option);
2688  opt.ToUpper();
2689  // do an evaluation as a cross-check
2690  //if (fReadyToExecute) Eval();
2691 
2692  if (opt.Contains("V") ) {
2693  if (fNdim > 0 && !TestBit(TFormula::kLambda)) {
2694  printf("List of Variables: \n");
2695  assert(int(fClingVariables.size()) >= fNdim);
2696  for ( int ivar = 0; ivar < fNdim ; ++ivar) {
2697  printf("Var%4d %20s = %10f \n",ivar,GetVarName(ivar).Data(), fClingVariables[ivar]);
2698  }
2699  }
2700  if (fNpar > 0) {
2701  printf("List of Parameters: \n");
2702  if ( int(fClingParameters.size()) < fNpar)
2703  Error("Print","Number of stored parameters in vector %lu in map %lu is different than fNpar %d",fClingParameters.size(), fParams.size(), fNpar);
2704  assert(int(fClingParameters.size()) >= fNpar);
2705  // print with order passed to Cling function
2706  for ( int ipar = 0; ipar < fNpar ; ++ipar) {
2707  printf("Par%4d %20s = %10f \n",ipar,GetParName(ipar), fClingParameters[ipar] );
2708  }
2709  }
2710  printf("Expression passed to Cling:\n");
2711  printf("\t%s\n",fClingInput.Data() );
2712  }
2713  if(!fReadyToExecute)
2714  {
2715  Warning("Print","Formula is not ready to execute. Missing parameters/variables");
2716  for(list<TFormulaFunction>::const_iterator it = fFuncs.begin(); it != fFuncs.end(); ++it)
2717  {
2718  TFormulaFunction fun = *it;
2719  if(!fun.fFound)
2720  {
2721  printf("%s is unknown.\n",fun.GetName());
2722  }
2723  }
2724  }
2726  {
2727  // we can skip this
2728  // Info("Print","Not all parameters are setted.");
2729  // for(map<TString,TFormulaVariable>::const_iterator it = fParams.begin(); it != fParams.end(); ++it)
2730  // {
2731  // pair<TString,TFormulaVariable> param = *it;
2732  // if(!param.second.fFound)
2733  // {
2734  // printf("%s has default value %lf\n",param.first.Data(),param.second.GetInitialValue());
2735  // }
2736  // }
2737 
2738  }
2739 
2740 
2741 }
2742 ////////////////////////////////////////////////////////////////////////////////
2743 /// Stream a class object.
2744 
2745 void TFormula::Streamer(TBuffer &b)
2746 {
2747  if (b.IsReading() ) {
2748  UInt_t R__s, R__c;
2749  Version_t v = b.ReadVersion(&R__s, &R__c);
2750  //std::cout << "version " << v << std::endl;
2751  if (v <= 8 && v > 3 && v != 6) {
2752  // old TFormula class
2753  ROOT::v5::TFormula * fold = new ROOT::v5::TFormula();
2754  // read old TFormula class
2755  fold->Streamer(b, v, R__s, R__c, TFormula::Class());
2756  //std::cout << "read old tformula class " << std::endl;
2757  TFormula fnew(fold->GetName(), fold->GetExpFormula() );
2758 
2759  *this = fnew;
2760 
2761 // printf("copying content in a new TFormula \n");
2762  SetParameters(fold->GetParameters() );
2763  if (!fReadyToExecute ) {
2764  Error("Streamer","Old formula read from file is NOT valid");
2765  Print("v");
2766  }
2767  delete fold;
2768  return;
2769  }
2770  else if (v > 8) {
2771  // new TFormula class
2772  b.ReadClassBuffer(TFormula::Class(), this, v, R__s, R__c);
2773 
2774  //std::cout << "reading npar = " << GetNpar() << std::endl;
2775 
2776  // initialize the formula
2777  // need to set size of fClingVariables which is transient
2778  //fClingVariables.resize(fNdim);
2779 
2780  // case of formula contains only parameters
2781  if (fFormula.IsNull() ) return;
2782 
2783 
2784  // store parameter values, names and order
2785  std::vector<double> parValues = fClingParameters;
2786  auto paramMap = fParams;
2787  fNpar = fParams.size();
2788 
2789  if (!TestBit(TFormula::kLambda) ) {
2790 
2791  //std::cout << "Streamer::Reading preprocess the formula " << fFormula << " ndim = " << fNdim << " npar = " << fNpar << std::endl;
2792  // for ( auto &p : fParams)
2793  // std::cout << "parameter " << p.first << " index " << p.second << std::endl;
2794 
2795  fClingParameters.clear(); // need to be reset before re-initializing it
2796 
2797  FillDefaults();
2798 
2799 
2801 
2802  //std::cout << "Streamer::after pre-process the formula " << fFormula << " ndim = " << fNdim << " npar = " << fNpar << std::endl;
2803 
2805 
2806  //std::cout << "Streamer::after prepared " << fClingInput << " ndim = " << fNdim << " npar = " << fNpar << std::endl;
2807 
2808 
2809  // restore parameter values
2810  if (fNpar != (int) parValues.size() ) {
2811  Error("Streamer","number of parameters computed (%d) is not same as the stored parameters (%d)",fNpar,int(parValues.size()) );
2812  Print("v");
2813  }
2814  }
2815  else {
2816  // case of lamda expressions
2817  bool ret = InitLambdaExpression(fFormula);
2818  if (ret) {
2819  fReadyToExecute = true;
2820  fClingInitialized = true;
2821  }
2822  }
2823  assert(fNpar == (int) parValues.size() );
2824  std::copy( parValues.begin(), parValues.end(), fClingParameters.begin() );
2825  // restore parameter names and order
2826  if (fParams.size() != paramMap.size() ) {
2827  Warning("Streamer","number of parameters list found (%lu) is not same as the stored one (%lu) - use re-created list",fParams.size(),paramMap.size()) ;
2828  //Print("v");
2829  }
2830  else
2831  //assert(fParams.size() == paramMap.size() );
2832  fParams = paramMap;
2833 
2834  // input formula into Cling
2835  // need to replace in cling the name of the pointer of this object
2836  // TString oldClingName = fClingName;
2837  // fClingName.Replace(fClingName.Index("_0x")+1,fClingName.Length(), TString::Format("%p",this) );
2838  // fClingInput.ReplaceAll(oldClingName, fClingName);
2839  // InputFormulaIntoCling();
2840 
2841 
2842 
2843  if (!TestBit(kNotGlobal)) {
2845  gROOT->GetListOfFunctions()->Add(this);
2846  }
2847  if (!fReadyToExecute ) {
2848  Error("Streamer","Formula read from file is NOT ready to execute");
2849  Print("v");
2850  }
2851  //std::cout << "reading 2 npar = " << GetNpar() << std::endl;
2852 
2853  return;
2854  }
2855  else {
2856  Error("Streamer","Reading version %d is not supported",v);
2857  return;
2858  }
2859  }
2860  else {
2861  // case of writing
2863  //std::cout << "writing npar = " << GetNpar() << std::endl;
2864  }
2865 }
Bool_t IsFuncCall() const
Definition: TFormula.h:43
virtual void Clear(Option_t *option="")
Set name and title to empty strings ("").
Definition: TFormula.cxx:591
static Bool_t IsBracket(const char c)
Definition: TFormula.cxx:162
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:52
Ssiz_t Last(char c) const
Find last occurrence of a character c.
Definition: TString.cxx:864
void HandleExponentiation(TString &formula)
Definition: TFormula.cxx:1145
std::map< TString, Double_t > fConsts
Definition: TFormula.h:121
virtual TString GetExpFormula(Option_t *option="") const
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
Double_t G()
Definition: TMath.h:68
virtual CallFuncIFacePtr_t CallFunc_IFacePtr(CallFunc_t *) const
Definition: TInterpreter.h:275
Int_t fNpar
Definition: TFormula.h:125
Double_t Floor(Double_t x)
Definition: TMath.h:473
virtual TFormula * GetFormula()
Definition: TF1.h:339
Double_t H()
Definition: TMath.h:81
std::vector< TObject * > fLinearParts
Definition: TFormula.h:127
static double p3(double t, double a, double b, double c, double d)
Double_t Eval(Double_t x) const
Definition: TFormula.cxx:2538
Bool_t IsReading() const
Definition: TBuffer.h:81
short Version_t
Definition: RtypesCore.h:61
Int_t GetVarNumber(const char *name) const
Definition: TFormula.cxx:2124
ClassImp(TSeqCollection) Int_t TSeqCollection TIter next(this)
Return index of object in collection.
Ssiz_t Length() const
Definition: TString.h:390
void variables(TString fin="TMVA.root", TString dirName="InputVariables_Id", TString title="TMVA Input Variables", Bool_t isRegression=kFALSE, Bool_t useTMVAStyle=kTRUE)
Definition: variables.cxx:10
void AddVariables(const TString *vars, const Int_t size)
Definition: TFormula.cxx:2017
TMethodCall * fMethod
Definition: TFormula.h:101
const char Option_t
Definition: RtypesCore.h:62
void DoSetParameters(const Double_t *p, Int_t size)
Definition: TFormula.cxx:2384
void HandleLinear(TString &formula)
Definition: TFormula.cxx:1245
void HandlePolN(TString &formula)
Definition: TFormula.cxx:788
Double_t QuietNaN()
Definition: TMath.h:635
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:635
#define assert(cond)
Definition: unittest.h:542
TString fName
Definition: TFormula.h:35
Int_t GetNargs() const
Number of function arguments.
Definition: TFunction.cxx:164
void SetParameters(const Double_t *params)
Definition: TFormula.cxx:2401
virtual void SetName(const char *name)
Change (i.e.
Definition: TNamed.cxx:128
Double_t Sqrt2()
Definition: TMath.h:51
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition: TObject.cxx:892
void SetParNames(const char *name0="p0", const char *name1="p1", const char *name2="p2", const char *name3="p3", const char *name4="p4", const char *name5="p5", const char *name6="p6", const char *name7="p7", const char *name8="p8", const char *name9="p9", const char *name10="p10")
Definition: TFormula.cxx:2438
void SetVariables(const std::pair< TString, Double_t > *vars, const Int_t size)
Definition: TFormula.cxx:2087
void ToUpper()
Change string to upper case.
Definition: TString.cxx:1101
Buffer base class used for serializing objects.
Definition: TBuffer.h:40
CallFunc_t * GetCallFunc() const
Definition: TMethodCall.h:97
Helper class for TFormula.
Definition: TFormula.h:32
#define gROOT
Definition: TROOT.h:340
Basic string class.
Definition: TString.h:137
void AddVariable(const TString &name, Double_t value=0)
Definition: TFormula.cxx:1982
void InputFormulaIntoCling()
pointer to the lambda function
Definition: TFormula.cxx:679
void DoAddParameter(const TString &name, Double_t value, bool processFormula)
Definition: TFormula.cxx:2169
int Int_t
Definition: RtypesCore.h:41
virtual void Copy(TObject &f1) const
Copy this to obj.
Definition: TFormula.cxx:518
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
R__EXTERN TVirtualMutex * gROOTMutex
Definition: TROOT.h:63
const TList * GetListOfAllPublicMethods(Bool_t load=kTRUE)
Returns a list of all public methods of this class and its base classes.
Definition: TClass.cxx:3621
Bool_t fFound
Definition: TFormula.h:38
Bool_t fAllParametersSetted
transient to force re-initialization
Definition: TFormula.h:100
#define gInterpreter
Definition: TInterpreter.h:502
Double_t EulerGamma()
Definition: TMath.h:122
Int_t GetParNumber(const char *name) const
Definition: TFormula.cxx:2251
STL namespace.
const TObject * GetLinearPart(Int_t i) const
Definition: TFormula.cxx:1968
size_t
Definition: TBuffer.cxx:28
TString & Insert(Ssiz_t pos, const char *s)
Definition: TString.h:592
void SetPredefinedParamNames()
Definition: TFormula.cxx:1899
void Print(Option_t *option="") const
print the formula and its attributes
Definition: TFormula.cxx:2682
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:732
Bool_t PrepareFormula(TString &formula)
Definition: TFormula.cxx:1315
TString & Replace(Ssiz_t pos, Ssiz_t n, const char *s)
Definition: TString.h:625
const char * Data() const
Definition: TString.h:349
void ReplaceParamName(TString &formula, const TString &oldname, const TString &name)
Definition: TFormula.cxx:2481
virtual ~TFormula()
Definition: TFormula.cxx:274
Double_t x[n]
Definition: legend1.C:17
Double_t EvalPar(const Double_t *x, const Double_t *params=0) const
Definition: TFormula.cxx:2509
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString...
Definition: TString.cxx:2334
Double_t C()
Definition: TMath.h:63
Double_t R()
Definition: TMath.h:109
void Class()
Definition: Class.C:29
Int_t fArrayPos
Definition: TFormula.h:70
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:33
Double_t Ln10()
Definition: TMath.h:57
UChar_t mod R__LOCKGUARD2(gSrvAuthenticateMutex)
virtual void Copy(TObject &named) const
Copy this to obj.
Definition: TNamed.cxx:83
void SetParameter(const char *name, Double_t value)
Definition: TFormula.cxx:2321
TString GetVarName(Int_t ivar) const
Definition: TFormula.cxx:2138
Double_t Log10(Double_t x)
Definition: TMath.h:529
virtual Bool_t Declare(const char *code)=0
static double p2(double t, double a, double b, double c)
void(* Generic_t)(void *, int, void **, void *)
Definition: TInterpreter.h:77
TInterpreter::CallFuncIFacePtr_t::Generic_t fFuncPtr
unique name passed to Cling to define the function ( double clingName(double*x, double*p) ) ...
Definition: TFormula.h:104
static Bool_t IsScientificNotation(const TString &formula, int ipos)
Definition: TFormula.cxx:183
TFormula & operator=(const TFormula &rhs)
Definition: TFormula.cxx:430
static bool IsReservedName(const char *name)
Definition: TFormula.cxx:265
TString & Append(const char *cs)
Definition: TString.h:492
std::vector< std::vector< double > > Data
static const std::string pattern("pattern")
TString GetExpFormula(Option_t *option="") const
return the expression formula If option = "P" replace the parameter names with their values If option...
Definition: TFormula.cxx:2610
Int_t Atoi() const
Return integer value of string.
Definition: TString.cxx:1964
Double_t fValue
Definition: TFormula.h:69
void * fLambdaPtr
function pointer
Definition: TFormula.h:105
Method or function calling interface.
Definition: TMethodCall.h:41
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:918
Int_t GetNdim() const
Definition: TFormula.h:174
void SetParName(Int_t ipar, const char *name)
Definition: TFormula.cxx:2454
static std::unordered_map< std::string, void * > gClingFunctions
Definition: TFormula.cxx:149
Double_t Infinity()
Definition: TMath.h:648
A doubly linked list.
Definition: TList.h:47
Double_t Sigma()
Definition: TMath.h:100
Int_t fNdim
Definition: TFormula.h:124
Int_t GetNargsOpt() const
Number of function optional (default) arguments.
Definition: TFunction.cxx:174
virtual void Delete(Option_t *option="")
Delete this object.
Definition: TObject.cxx:228
static Bool_t IsHexadecimal(const TString &formula, int ipos)
Definition: TFormula.cxx:194
std::vector< Double_t > fClingParameters
cached variables
Definition: TFormula.h:97
Int_t GetNpar() const
Definition: TFormula.h:175
void FillDefaults()
Definition: TFormula.cxx:690
void HandleParametrizedFunctions(TString &formula)
Definition: TFormula.cxx:881
static Bool_t IsOperator(const char c)
Definition: TFormula.cxx:151
void SetName(const char *name)
Change (i.e.
Definition: TFormula.cxx:2065
The F O R M U L A class.
Definition: TFormula.h:89
unsigned int UInt_t
Definition: RtypesCore.h:42
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:173
TMarker * m
Definition: textangle.C:8
static Bool_t IsDefaultVariableName(const TString &name)
Definition: TFormula.cxx:177
Double_t K()
Definition: TMath.h:95
Double_t E()
Definition: TMath.h:54
bool operator()(const TString &a, const TString &b) const
Definition: TFormula.cxx:215
static Bool_t IsFunctionNameChar(const char c)
Definition: TFormula.cxx:172
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
The ROOT global object gROOT contains a list of all defined classes.
Definition: TClass.h:81
const std::string sname
Definition: testIO.cxx:45
static double p1(double t, double a, double b)
Int_t CountChar(Int_t c) const
Return number of times character c occurs in the string.
Definition: TString.cxx:443
Double_t DoEval(const Double_t *x, const Double_t *p=nullptr) const
Definition: TFormula.cxx:2547
Bool_t IsNull() const
Definition: TString.h:387
Double_t GetVariable(const char *name) const
Definition: TFormula.cxx:2111
Bool_t PrepareEvalMethod()
Definition: TFormula.cxx:623
TString fFormula
Definition: TFormula.h:123
void InitWithPrototype(TClass *cl, const char *method, const char *proto, Bool_t objectIsConst=kFALSE, ROOT::EFunctionMatchMode mode=ROOT::kConversionMatch)
Initialize the method invocation environment.
TString fClingName
pointer to methocall
Definition: TFormula.h:102
Double_t Pi()
Definition: TMath.h:44
std::vector< Double_t > fClingVariables
input function passed to Cling
Definition: TFormula.h:96
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
int Ssiz_t
Definition: RtypesCore.h:63
virtual Double_t * GetParameters() const
Definition: TFormula.h:249
std::map< TString, TString > fFunctionsShortcuts
Definition: TFormula.h:122
#define ClassImp(name)
Definition: Rtypes.h:279
double f(double x)
Double_t * GetParameters() const
Definition: TFormula.cxx:2303
double Double_t
Definition: RtypesCore.h:55
Int_t GetNumber() const
Definition: TFormula.h:176
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
Double_t y[n]
Definition: legend1.C:17
Int_t fNumber
Definition: TFormula.h:126
Int_t Compile(const char *expression="")
Definition: TFormula.cxx:479
static TClass * GetClass(const char *name, Bool_t load=kTRUE, Bool_t silent=kFALSE)
Static method returning pointer to TClass of the specified class name.
Definition: TClass.cxx:2881
Int_t GetNargs() const
Definition: TFormula.h:42
#define name(a, b)
Definition: linkTestLib0.cpp:5
void SetVariable(const TString &name, Double_t value)
Definition: TFormula.cxx:2155
const char * GetParName(Int_t ipar) const
Definition: TFormula.cxx:2288
Mother of all ROOT objects.
Definition: TObject.h:58
Bool_t InitLambdaExpression(const char *formula)
Definition: TFormula.cxx:442
TString fName
Definition: TFormula.h:68
Global functions class (global functions are obtained from CINT).
Definition: TFunction.h:30
Bool_t IsDigit() const
Returns true if all characters in string are digits (0-9) or white spaces, i.e.
Definition: TString.cxx:1806
Bool_t IsValid() const
Definition: TFormula.h:186
void Streamer(TBuffer &b, const TClass *onfile_class)
const Ssiz_t kNPOS
Definition: Rtypes.h:115
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
1-Dim function class
Definition: TF1.h:149
void ProcessFormula(TString &formula)
Definition: TFormula.cxx:1549
Each ROOT class (see TClass) has a linked list of methods.
Definition: TMethod.h:40
TF1 * f1
Definition: legend1.C:11
R__EXTERN Int_t gDebug
Definition: Rtypes.h:128
double result[121]
Bool_t fReadyToExecute
Definition: TFormula.h:98
R__EXTERN TInterpreter * gCling
Definition: TInterpreter.h:504
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition: TString.h:582
void ExtractFunctors(TString &formula)
Definition: TFormula.cxx:1346
TString fClingInput
Definition: TFormula.h:95
TObject * obj
virtual TObject * FindObject(const char *name) const
Must be redefined in derived classes.
Definition: TObject.cxx:379
float value
Definition: math.cpp:443
std::map< TString, TFormulaVariable > fVars
Definition: TFormula.h:119
const Int_t n
Definition: legend1.C:16
Bool_t fClingInitialized
trasient to force initialization
Definition: TFormula.h:99
std::map< TString, Int_t, TFormulaParamOrder > fParams
list of variable names
Definition: TFormula.h:120
Double_t GetParameter(const char *name) const
Definition: TFormula.cxx:2263
Bool_t IsValid() const
Return true if the method call has been properly initialized and is usable.
std::list< TFormulaFunction > fFuncs
Definition: TFormula.h:118
Another helper class for TFormula.
Definition: TFormula.h:65
Double_t LogE()
Definition: TMath.h:60
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0
void PreProcessFormula(TString &formula)
Definition: TFormula.cxx:1293
const char * GetName() const
Definition: TFormula.h:40
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:904