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