Logo ROOT   6.07/09
Reference Guide
stressVdt.cxx
Go to the documentation of this file.
1 /// Simple program to benchmark vdt accuracy and cpu performance.
2 #include <vector>
3 #include <iostream>
4 #include <cmath> //for log2
5 #include <assert.h>
6 #include <limits>
7 #include <iostream>
8 #include <iomanip>
9 #include <map>
10 #include <string>
11 
12 #include "vdt/vdtMath.h"
13 
14 #include "TStopwatch.h"
15 #include "TRandom3.h"
16 #include "TH1F.h"
17 #include "TCanvas.h"
18 #include "TError.h"
19 
20 const double RANGE=3000.;
21 const uint32_t SIZE= 16777216;
22 
23 //------------------------------------------------------------------------------
24 // Not good for floating point, but just to get the bit difference
25 template <class T>
26 uint64_t fp2uint (T /*x*/)
27 {
28  T::not_implemented; // "Static assert" in C++03
29  return 0;
30 }
31 
32 template <>
33 uint64_t fp2uint<double> (double x)
34 {
35  return vdt::details::dp2uint64(x);
36 }
37 
38 template <>
39 uint64_t fp2uint<float> (float x)
40 {
41  return vdt::details::sp2uint32(x);
42 }
43 
44 //------------------------------------------------------------------------------
45 /// Returns most significative different bit
46 template <class T>
47 inline uint32_t diffbit(const T a,const T b )
48 {
49  uint64_t ia = fp2uint<T>(a);
50  uint64_t ib = fp2uint<T>(b);
51  uint64_t c = ia>ib? ia-ib : ib-ia;
52  return log2(c)+1;
53 }
54 
55 //------------------------------------------------------------------------------
56 // This allows to vectorise on very modern compilers (>=gcc 4.7)
57 // Note how the templating mechanism allows the compiler to *inline* the
58 // function. It is much more efficient than a void ptr or std::function!
59 template <typename T, typename F>
60 inline void calculateValues(F mathFunc,
61  const std::vector<T>& inputVector,
62  std::vector<T>& outputVector)
63 {
64 
65  const uint32_t size = inputVector.size();
66 
67  for (unsigned int i=0;i<size;++i){
68  outputVector[i]=mathFunc(inputVector[i]);
69  }
70 
71 }
72 
73 ////////////////////////////////////////////////////////////////////////////////
74 
75 template <typename T>
76 void compareOutputs(const std::vector<T>& inputVector1,
77  const std::vector<T>& inputVector2,
78  std::vector<uint32_t>& outputVector)
79 {
80  assert(inputVector1.size()==inputVector2.size() &&
81  inputVector1.size()==outputVector.size());
82 
83  const uint32_t size = inputVector1.size();
84 
85  for (unsigned int i=0;i<size;++i)
86  outputVector[i]=diffbit(inputVector1[i],inputVector2[i]);
87 }
88 
89 //------------------------------------------------------------------------------
90 
96 
97 template <typename T>
98 void fillRandom(std::vector<T>& randomV,
99  const rangeType theRangeType)
100 {
101  // Yeah, well, maybe it can be done better. But this is not a tutorial about
102  // random generation!
103  const uint32_t size=randomV.size();
104  static TRandom3 rndmGenerator(123);
105  T* arr = &(randomV[0]);
106  rndmGenerator.RndmArray(size,arr);
107  if (kReal == theRangeType ) for (uint32_t i=0;i<size;++i) randomV[i]=(randomV[i]-0.5)*2*RANGE;
108  if (kExp == theRangeType ) for (uint32_t i=0;i<size;++i) randomV[i]=(randomV[i]-0.5)*2*705.;
109  if (kExpf == theRangeType ) for (uint32_t i=0;i<size;++i) randomV[i]=(randomV[i]-0.5)*2*85.;
110  if (kRealPlus == theRangeType ) for (uint32_t i=0;i<size;++i) randomV[i]=randomV[i]*RANGE+0.000001;
111  if (km1p1 == theRangeType ) for (uint32_t i=0;i<size;++i) randomV[i]=(randomV[i]-0.5)*2;
112 
113 }
114 
115 //------------------------------------------------------------------------------
116 
117 template<typename T>
119  const std::vector<T>& VDTVals,
120  const std::vector<T>& SystemVals)
121 {
122  const uint32_t size = VDTVals.size();
123  std::vector<uint32_t> diff(size);
124  compareOutputs(VDTVals,SystemVals,diff);
125  uint32_t theDiff=0;
126  for (uint32_t i =0;i<size;++i){
127  theDiff=diff[i];
128  histo.Fill(theDiff);
129  }
130 }
131 
132 ////////////////////////////////////////////////////////////////////////////////
133 
134 template <typename T, typename F>
135 inline double measureTiming(F mathFunc,
136  const std::vector<T>& inputVector,
137  std::vector<T>& outputVector)
138 {
140  timer.Start();
141  calculateValues<T>(mathFunc,inputVector,outputVector);
142  timer.Stop();
143  return timer.RealTime();
144 }
145 
146 //------------------------------------------------------------------------------
147 
148  /*
149  Reference values on a Intel(R) Core(TM) i7 CPU 950 @ 3.07GHz
150  gcc 4.8.2 -Ofast
151 
152  Expf 0.110623 0.821863 7.4294 3.35636 6
153  Sinf 0.151893 9.83798 64.7692 0.235055 9
154  Cosf 0.11277 9.87508 87.5683 0.234271 8
155  Tanf 0.167273 9.84792 58.8733 0.519784 9
156  Atanf 0.0855272 0.529288 6.18854 0.366963 2
157  Logf 0.114023 0.465541 4.08287 0.270999 2
158  Isqrtf 0.0328619 0.275043 8.36965 4.35237 7
159  Asinf 0.0958891 0.415733 4.33556 0.60167 3
160  Acosf 0.099067 0.470179 4.74607 0.480427 10
161  Exp 0.204585 0.64904 3.17247 0.137142 2
162  Sin 0.327537 1.5099 4.60986 0.253579 2
163  Cos 0.299601 1.5038 5.01933 0.253664 2
164  Tan 0.276369 2.13009 7.7074 0.351065 5
165  Atan 0.355532 0.902413 2.5382 0.326134 2
166  Log 0.244172 1.25513 5.14034 0.385112 2
167  Isqrt 0.167836 0.283752 1.69065 0.453692 2
168  Asin 0.40379 0.869315 2.15289 0.318644 2
169  Acos 0.392566 0.864706 2.2027 0.391922 11
170  */
171 
172  // The reference values: speedup and accuracy. Some contingency is given
173 struct staticInitHelper{
174  std::map<std::string, std::pair<float,uint32_t> > referenceValues;
175 
176  staticInitHelper()
177  {
178 #if defined(R__B64)
179  const int sincosTolerance = 4;
180  const int tanTolerance = 5;
181 #else
182  const int sincosTolerance = 19;
183  const int tanTolerance = 19;
184 #endif
185  referenceValues["Expf"] = std::make_pair(1.f,8);
186  referenceValues["Sinf"] = std::make_pair(1.f,11);
187  referenceValues["Cosf"] = std::make_pair(1.f,10);
188  referenceValues["Tanf"] = std::make_pair(1.f,11);
189  referenceValues["Atanf"] = std::make_pair(1.f,4);
190  referenceValues["Logf"] = std::make_pair(1.f,4);
191  referenceValues["Isqrtf"]= std::make_pair(1.f,9);
192  referenceValues["Asinf"] = std::make_pair(1.f,5);
193  referenceValues["Acosf"] = std::make_pair(1.f,12);
194  referenceValues["Exp"] = std::make_pair(1.f,4);
195  referenceValues["Sin"] = std::make_pair(1.f,sincosTolerance);
196  referenceValues["Cos"] = std::make_pair(1.f,sincosTolerance);
197  referenceValues["Tan"] = std::make_pair(1.f,tanTolerance);
198  referenceValues["Atan"] = std::make_pair(1.f,4);
199  referenceValues["Log"] = std::make_pair(1.f,4);
200  referenceValues["Isqrt"] = std::make_pair(.4f,4); // Fix fluctuation on x86_64-slc5-gcc47
201  referenceValues["Asin"] = std::make_pair(1.f,4);
202  referenceValues["Acos"] = std::make_pair(1.f,13);
203  }
204 } gbl;
205 
206 template <typename T, typename F1, typename F2>
207 inline void compareFunctions(const std::string& label,
208  F1 vdtFunc,
209  F2 systemFunc,
210  const std::vector<T>& inputVector,
211  std::vector<T>& outputVectorVDT,
212  std::vector<T>& outputVectorSystem,
213  float& speedup,
214  uint32_t& maxdiffBit,
215  TH1F& histo)
216 {
217  double timeVdt = measureTiming<T>(vdtFunc,inputVector,outputVectorVDT);
218  double timeSystem = measureTiming<T>(systemFunc,inputVector,outputVectorSystem);
219  std::string name(label);
220  std::string title(label);
221  title+=";Diff. Bit;#";
222  histo.Reset();
223  histo.SetName(label.c_str());
224  histo.SetTitle(title.c_str());
225  treatBinDiffHisto(histo,outputVectorVDT,outputVectorSystem);
226  double meandiffBit = histo.GetMean();
227  maxdiffBit = 0;
228  const uint32_t xmax=histo.GetXaxis()->GetXmax();
229 
230  for (uint32_t i=1;i<=xmax;i++){
231  if ( histo.GetBinContent(i) > 0.f )
232  maxdiffBit=i-1;
233  }
234 
235  speedup = timeSystem/timeVdt ;
236 
237  std::cout << std::setw(8)
238  << label << std::setw(10)
239  << timeVdt << std::setw(10)
240  << timeSystem << std::setw(15)
241  << speedup << std::setw(15)
242  << meandiffBit << std::setw(15)
243  << maxdiffBit << std::endl;
244 
245  // Draw it
246  TCanvas c;
247  c.cd();
248  histo.SetLineColor(kBlue);
249  histo.SetLineWidth(2);
250  histo.Draw("hist");
251  c.SetLogy();
252  name+=".png";
253  Int_t oldErrorIgnoreLevel = gErrorIgnoreLevel; // we know we are printing..
254  gErrorIgnoreLevel=1001;
255  c.Print(name.c_str());
256  gErrorIgnoreLevel=oldErrorIgnoreLevel;
257 }
258 
259 //
260 void checkFunction(const std::string& label,float /*speedup*/, uint32_t maxdiffBit)
261 {
262 // Remove check on the speed as routinely this program is ran on virtual build nodes
263 // and several factors may cause fluctuations in the result.
264 // if (gbl.referenceValues[label].first > speedup)
265 // std::cerr << "Note " << label << " was slower than the system library.\n";
266  if (gbl.referenceValues[label].second < maxdiffBit)
267  std::cerr << "WARNING " << label << " is too inaccurate!\n";
268 }
269 
270 //------------------------------------------------------------------------------
271 // Some helpers
272 inline float isqrtf(float x) {return 1.f/sqrt(x);};
273 inline double isqrt(double x) {return 1./sqrt(x);};
274 
275 //------------------------------------------------------------------------------
276 // Artificially chunck the work to help the compiler manage everything.
277 // If all the content is moved to a single function, the vectorization breaks.
278 // Same holds for the check of speed and accuracy. Technically it could be part
279 // of compareFunctions.
280 void spStep1()
281 {
282  std::vector<float> VDTVals(SIZE);
283  std::vector<float> SystemVals(SIZE);
284  std::vector<float> realNumbers(SIZE);
285  TH1F histo("bitDiffHisto","willbechanged",32,0,32);
286 
287  float speedup;
288  uint32_t maxdiffBit;
289 
290  fillRandom(realNumbers,kExpf);
291  compareFunctions<float>("Expf", vdt::fast_expf, expf, realNumbers, VDTVals, SystemVals, speedup, maxdiffBit, histo);
292  checkFunction("Expf",speedup, maxdiffBit);
293 
294  fillRandom(realNumbers,kReal);
295  compareFunctions<float>("Sinf", vdt::fast_sinf, sinf, realNumbers, VDTVals, SystemVals, speedup, maxdiffBit, histo);
296  checkFunction("Sinf",speedup, maxdiffBit);
297  compareFunctions<float>("Cosf", vdt::fast_cosf, cosf, realNumbers, VDTVals, SystemVals, speedup, maxdiffBit, histo);
298  checkFunction("Cosf",speedup, maxdiffBit);
299  compareFunctions<float>("Tanf", vdt::fast_tanf, tanf, realNumbers, VDTVals, SystemVals, speedup, maxdiffBit, histo);
300  checkFunction("Tanf",speedup, maxdiffBit);
301  compareFunctions<float>("Atanf", vdt::fast_atanf, atanf, realNumbers, VDTVals, SystemVals, speedup, maxdiffBit, histo);
302  checkFunction("Atanf",speedup, maxdiffBit);
303 }
304 
305 void spStep2()
306 {
307  std::vector<float> VDTVals(SIZE);
308  std::vector<float> SystemVals(SIZE);
309  std::vector<float> realNumbers(SIZE);
310  TH1F histo("bitDiffHisto","willbechanged",32,0,32);
311 
312  float speedup;
313  uint32_t maxdiffBit;
314 
315  fillRandom(realNumbers,kRealPlus);
316  compareFunctions<float>("Logf", vdt::fast_logf, logf, realNumbers, VDTVals, SystemVals, speedup, maxdiffBit, histo);
317  checkFunction("Logf",speedup, maxdiffBit);
318  compareFunctions<float>("Isqrtf", vdt::fast_isqrtf, isqrtf, realNumbers, VDTVals, SystemVals, speedup, maxdiffBit, histo);
319  checkFunction("Isqrtf",speedup, maxdiffBit);
320 }
321 
322 void spStep3()
323 {
324  std::vector<float> VDTVals(SIZE);
325  std::vector<float> SystemVals(SIZE);
326  std::vector<float> realNumbers(SIZE);
327  TH1F histo("bitDiffHisto","willbechanged",32,0,32);
328 
329  float speedup;
330  uint32_t maxdiffBit;
331 
332  fillRandom(realNumbers,km1p1);
333  compareFunctions<float>("Asinf", vdt::fast_asinf, asinf, realNumbers, VDTVals, SystemVals, speedup, maxdiffBit, histo);
334  checkFunction("Asinf",speedup, maxdiffBit);
335  compareFunctions<float>("Acosf", vdt::fast_acosf, acosf, realNumbers, VDTVals, SystemVals, speedup, maxdiffBit, histo);
336  checkFunction("Acosf",speedup, maxdiffBit);
337 }
338 
339 void dpStep1()
340 {
341  std::vector<double> VDTVals(SIZE);
342  std::vector<double> SystemVals(SIZE);
343  std::vector<double> realNumbers(SIZE);
344  TH1F histo("bitDiffHisto","willbechanged",64,0,64);
345 
346  float speedup;
347  uint32_t maxdiffBit;
348 
349  fillRandom(realNumbers,kExp);
350  compareFunctions<double>("Exp", vdt::fast_exp, exp, realNumbers, VDTVals, SystemVals, speedup, maxdiffBit, histo);
351  checkFunction("Exp",speedup, maxdiffBit);
352  fillRandom(realNumbers,kReal);
353  compareFunctions<double>("Sin", vdt::fast_sin, sin, realNumbers, VDTVals, SystemVals, speedup, maxdiffBit, histo);
354  checkFunction("Sin",speedup, maxdiffBit);
355  compareFunctions<double>("Cos", vdt::fast_cos, cos, realNumbers, VDTVals, SystemVals, speedup, maxdiffBit, histo);
356  checkFunction("Cos",speedup, maxdiffBit);
357  compareFunctions<double>("Tan", vdt::fast_tan, tan, realNumbers, VDTVals, SystemVals, speedup, maxdiffBit, histo);
358  checkFunction("Tan",speedup, maxdiffBit);
359  compareFunctions<double>("Atan", vdt::fast_atan, atan, realNumbers, VDTVals, SystemVals, speedup, maxdiffBit, histo);
360  checkFunction("Atan",speedup, maxdiffBit);
361 }
362 
363 void dpStep2()
364 {
365  std::vector<double> VDTVals(SIZE);
366  std::vector<double> SystemVals(SIZE);
367  std::vector<double> realNumbers(SIZE);
368  TH1F histo("bitDiffHisto","willbechanged",64,0,64);
369 
370  float speedup;
371  uint32_t maxdiffBit;
372 
373  fillRandom(realNumbers,kRealPlus);
374  compareFunctions<double>("Log", vdt::fast_log, log, realNumbers, VDTVals, SystemVals, speedup, maxdiffBit, histo);
375  checkFunction("Log",speedup, maxdiffBit);
376  compareFunctions<double>("Isqrt", vdt::fast_isqrt, isqrt, realNumbers, VDTVals, SystemVals, speedup, maxdiffBit, histo);
377  checkFunction("Isqrt",speedup, maxdiffBit);
378 }
379 
380 void dpStep3()
381 {
382  std::vector<double> VDTVals(SIZE);
383  std::vector<double> SystemVals(SIZE);
384  std::vector<double> realNumbers(SIZE);
385  TH1F histo("bitDiffHisto","willbechanged",64,0,64);
386 
387  float speedup;
388  uint32_t maxdiffBit;
389 
390  fillRandom(realNumbers,km1p1);
391  compareFunctions<double>("Asin", vdt::fast_asin, asin, realNumbers, VDTVals, SystemVals, speedup, maxdiffBit, histo);
392  checkFunction("Asin",speedup, maxdiffBit);
393  compareFunctions<double>("Acos", vdt::fast_acos, acos, realNumbers, VDTVals, SystemVals, speedup, maxdiffBit, histo);
394  checkFunction("Acos",speedup, maxdiffBit);
395 }
396 //------------------------------------------------------------------------------
397 
398 int main(){
399 
400  std::cout << "Test performed on " << SIZE << " random numbers\n"
401  << std::setw(8)
402  << "Name" << std::setw(10)
403  << "VDT (s)" << std::setw(10)
404  << "Sys (s)" << std::setw(15)
405  << "Speedup" << std::setw(15)
406  << "<diff Bit>" << std::setw(15)
407  << "max(diff Bit)" << std::endl;
408 
409  // Single precision ----
410  spStep1();
411  spStep2();
412  spStep3();
413 
414  // Double precision ----
415  dpStep1();
416  dpStep2();
417  dpStep3();
418 
419  return 0;
420 
421 }
422 
423 
424 
425 
426 
427 
428 
429 
430 
431 
432 
433 
434 
435 
436 
437 
438 
439 
440 
441 
442 
443 
444 
445 
446 
447 
448 
449 
450 
451 
452 
453 
454 
455 
456 
457 
458 
459 
460 
461 
462 
463 
464 
465 
466 
467 
468 
469 
470 
471 
void spStep3()
Definition: stressVdt.cxx:322
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition: TAttLine.h:49
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3127
float fast_expf(float initial_x)
Exponential Function single precision.
Definition: exp.h:118
float fast_asinf(float x)
Single Precision asin.
Definition: asin.h:157
float fast_tanf(float x)
Definition: tan.h:130
Random number generator class based on M.
Definition: TRandom3.h:29
Double_t RealTime()
Stop the stopwatch (if it is running) and return the realtime (in seconds) passed between the start a...
Definition: TStopwatch.cxx:110
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH1.cxx:4640
float fast_sinf(float x)
Definition: sin.h:41
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Definition: TStopwatch.cxx:58
R__EXTERN Int_t gErrorIgnoreLevel
Definition: TError.h:107
double fast_exp(double initial_x)
Exponential Function double precision.
Definition: exp.h:70
return c
double T(double x)
Definition: ChebyshevPol.h:34
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:659
uint32_t sp2uint32(float x)
Converts a float to an int.
const double RANGE
Simple program to benchmark vdt accuracy and cpu performance.
Definition: stressVdt.cxx:20
tomato 1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:575
double measureTiming(F mathFunc, const std::vector< T > &inputVector, std::vector< T > &outputVector)
Definition: stressVdt.cxx:135
struct staticInitHelper gbl
int Int_t
Definition: RtypesCore.h:41
TArc * a
Definition: textangle.C:12
double fast_asin(double x)
Double Precision asin.
Definition: asin.h:120
virtual void Print(const char *filename="") const
Save Pad contents in a file in one of various formats.
Definition: TPad.cxx:4160
uint32_t diffbit(const T a, const T b)
Returns most significative different bit.
Definition: stressVdt.cxx:47
float isqrtf(float x)
Definition: stressVdt.cxx:272
double cos(double)
void treatBinDiffHisto(TH1F &histo, const std::vector< T > &VDTVals, const std::vector< T > &SystemVals)
Definition: stressVdt.cxx:118
virtual void Reset(Option_t *option="")
Reset.
Definition: TH1.cxx:9164
const uint32_t SIZE
Definition: stressVdt.cxx:21
double sqrt(double)
double fast_acos(double x)
Definition: asin.h:207
TStopwatch timer
Definition: pirndm.C:37
void Stop()
Stop the stopwatch.
Definition: TStopwatch.cxx:77
Double_t x[n]
Definition: legend1.C:17
double acos(double)
void dpStep2()
Definition: stressVdt.cxx:363
double fast_atan(double x)
Fast Atan implementation double precision.
Definition: atan.h:86
double fast_cos(double x)
Double precision cosine: just call sincos.
Definition: cos.h:22
void compareFunctions(const std::string &label, F1 vdtFunc, F2 systemFunc, const std::vector< T > &inputVector, std::vector< T > &outputVectorVDT, std::vector< T > &outputVectorSystem, float &speedup, uint32_t &maxdiffBit, TH1F &histo)
Definition: stressVdt.cxx:207
#define F1(x, y, z)
Definition: TMD5.cxx:267
void dpStep3()
Definition: stressVdt.cxx:380
double sin(double)
void checkFunction(const std::string &label, float, uint32_t maxdiffBit)
Definition: stressVdt.cxx:260
uint64_t fp2uint< double >(double x)
Definition: stressVdt.cxx:33
double fast_log(double x)
Definition: log.h:89
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:46
#define F(x, y, z)
virtual Double_t GetMean(Int_t axis=1) const
For axis = 1,2 or 3 returns the mean value of the histogram along X,Y or Z axis.
Definition: TH1.cxx:6717
void spStep2()
Definition: stressVdt.cxx:305
double fast_sin(double x)
Double precision sine: just call sincos.
Definition: sin.h:37
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2853
float fast_isqrtf(float x)
Two iterations.
Definition: sqrt.h:85
uint64_t dp2uint64(double x)
Converts a double to an unsigned long long.
#define F2(x, y, z)
Definition: TMD5.cxx:268
void spStep1()
Definition: stressVdt.cxx:280
float xmax
Definition: THbookFile.cxx:93
double asin(double)
void dpStep1()
Definition: stressVdt.cxx:339
virtual void RndmArray(Int_t n, Float_t *array)
Return an array of n random numbers uniformly distributed in ]0,1].
Definition: TRandom3.cxx:138
virtual void SetName(const char *name)
Change the name of this histogram.
Definition: TH1.cxx:8028
The Canvas class.
Definition: TCanvas.h:41
void fillRandom(std::vector< T > &randomV, const rangeType theRangeType)
Definition: stressVdt.cxx:98
rangeType
Definition: stressVdt.cxx:91
double f(double x)
uint64_t fp2uint(T)
Definition: stressVdt.cxx:26
double isqrt(double x)
Definition: stressVdt.cxx:273
uint64_t fp2uint< float >(float x)
Definition: stressVdt.cxx:39
Double_t GetXmax() const
Definition: TAxis.h:140
double atan(double)
int main()
Definition: stressVdt.cxx:398
double tan(double)
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
float fast_logf(float x)
Definition: log.h:172
Definition: Rtypes.h:61
virtual void SetTitle(const char *title)
Change (i.e.
Definition: TH1.cxx:5984
float fast_cosf(float x)
Definition: cos.h:26
void compareOutputs(const std::vector< T > &inputVector1, const std::vector< T > &inputVector2, std::vector< uint32_t > &outputVector)
Definition: stressVdt.cxx:76
double exp(double)
double fast_tan(double x)
Double precision tangent implementation.
Definition: tan.h:84
void calculateValues(F mathFunc, const std::vector< T > &inputVector, std::vector< T > &outputVector)
Definition: stressVdt.cxx:60
char name[80]
Definition: TGX11.cxx:109
double log(double)
TAxis * GetXaxis()
Definition: TH1.h:324
float fast_atanf(float xx)
Fast Atan implementation single precision.
Definition: atan.h:129
double fast_isqrt(double x)
Four iterations.
Definition: sqrt.h:55
float fast_acosf(float x)
Definition: asin.h:211
Stopwatch class.
Definition: TStopwatch.h:30
virtual void SetLogy(Int_t value=1)
Set Lin/Log scale for Y.
Definition: TPad.cxx:5347