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