source: trunk/MagicSoft/Mars/mbase/MMath.cc@ 9195

Last change on this file since 9195 was 9195, checked in by tbretz, 16 years ago
*** empty log message ***
File size: 24.5 KB
Line 
1/* ======================================================================== *\
2! $Name: not supported by cvs2svn $:$Id: MMath.cc,v 1.42 2008-12-21 18:09:49 tbretz Exp $
3! --------------------------------------------------------------------------
4!
5! *
6! * This file is part of MARS, the MAGIC Analysis and Reconstruction
7! * Software. It is distributed to you in the hope that it can be a useful
8! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
9! * It is distributed WITHOUT ANY WARRANTY.
10! *
11! * Permission to use, copy, modify and distribute this software and its
12! * documentation for any purpose is hereby granted without fee,
13! * provided that the above copyright notice appear in all copies and
14! * that both that copyright notice and this permission notice appear
15! * in supporting documentation. It is provided "as is" without express
16! * or implied warranty.
17! *
18!
19!
20! Author(s): Thomas Bretz 3/2004 <mailto:tbretz@astro.uni-wuerzburg.de>
21!
22! Copyright: MAGIC Software Development, 2000-2005
23!
24!
25\* ======================================================================== */
26
27/////////////////////////////////////////////////////////////////////////////
28//
29// MMath
30//
31// Mars - Math package (eg Significances, etc)
32//
33/////////////////////////////////////////////////////////////////////////////
34#include "MMath.h"
35
36#include <stdlib.h> // atof (Ubuntu 8.10)
37
38#ifndef ROOT_TVector2
39#include <TVector2.h>
40#endif
41
42#ifndef ROOT_TVector3
43#include <TVector3.h>
44#endif
45
46#ifndef ROOT_TArrayD
47#include <TArrayD.h>
48#endif
49
50#ifndef ROOT_TComplex
51#include <TComplex.h>
52#endif
53
54//NamespaceImp(MMath);
55
56// --------------------------------------------------------------------------
57//
58// Calculate Significance as
59// significance = (s-b)/sqrt(s+k*k*b) mit k=s/b
60//
61// s: total number of events in signal region
62// b: number of background events in signal region
63//
64Double_t MMath::Significance(Double_t s, Double_t b)
65{
66 const Double_t k = b==0 ? 0 : s/b;
67 const Double_t f = s+k*k*b;
68
69 return f==0 ? 0 : (s-b)/TMath::Sqrt(f);
70}
71
72// --------------------------------------------------------------------------
73//
74// Symmetrized significance - this is somehow analog to
75// SignificanceLiMaSigned
76//
77// Returns Significance(s,b) if s>b otherwise -Significance(b, s);
78//
79Double_t MMath::SignificanceSym(Double_t s, Double_t b)
80{
81 return s>b ? Significance(s, b) : -Significance(b, s);
82}
83
84// --------------------------------------------------------------------------
85//
86// calculates the significance according to Li & Ma
87// ApJ 272 (1983) 317, Formula 17
88//
89// s // s: number of on events
90// b // b: number of off events
91// alpha = t_on/t_off; // t: observation time
92//
93// The significance has the same (positive!) value for s>b and b>s.
94//
95// Returns -1 if s<0 or b<0 or alpha<0 or the argument of sqrt<0
96//
97// Here is some eMail written by Daniel Mazin about the meaning of the arguments:
98//
99// > Ok. Here is my understanding:
100// > According to Li&Ma paper (correctly cited in MMath.cc) alpha is the
101// > scaling factor. The mathematics behind the formula 17 (and/or 9) implies
102// > exactly this. If you scale OFF to ON first (using time or using any other
103// > method), then you cannot use formula 17 (9) anymore. You can just try
104// > the formula before scaling (alpha!=1) and after scaling (alpha=1), you
105// > will see the result will be different.
106//
107// > Here are less mathematical arguments:
108//
109// > 1) the better background determination you have (smaller alpha) the more
110// > significant is your excess, thus your analysis is more sensitive. If you
111// > normalize OFF to ON first, you loose this sensitivity.
112//
113// > 2) the normalization OFF to ON has an error, which naturally depends on
114// > the OFF and ON. This error is propagating to the significance of your
115// > excess if you use the Li&Ma formula 17 correctly. But if you normalize
116// > first and use then alpha=1, the error gets lost completely, you loose
117// > somehow the criteria of goodness of the normalization.
118//
119Double_t MMath::SignificanceLiMa(Double_t s, Double_t b, Double_t alpha)
120{
121 const Double_t sum = s+b;
122
123 if (s<0 || b<0 || alpha<=0)
124 return -1;
125
126 const Double_t l = s==0 ? 0 : s*TMath::Log(s/sum*(alpha+1)/alpha);
127 const Double_t m = b==0 ? 0 : b*TMath::Log(b/sum*(alpha+1) );
128
129 return l+m<0 ? -1 : TMath::Sqrt((l+m)*2);
130}
131
132/*
133Double_t MMath::SignificanceLiMaErr(Double_t s, Double_t b, Double_t alpha)
134{
135 Double_t S = SignificanceLiMa(s, b, alpha);
136 if (S<0)
137 return -1;
138
139 const Double_t sum = s+b;
140
141
142 Double_t l = TMath::Log(s/sum*(alpha+1)/alpha)/TMath::Sqrt(2*S);
143 Double_t m = TMath::Log(s/sum*(alpha+1)/alpha)/TMath::Sqrt(2*S);
144
145
146 const Double_t sum = s+b;
147
148 if (s<0 || b<0 || alpha<=0)
149 return -1;
150
151 const Double_t l = s==0 ? 0 : s*TMath::Log(s/sum*(alpha+1)/alpha);
152 const Double_t m = b==0 ? 0 : b*TMath::Log(b/sum*(alpha+1) );
153
154 return l+m<0 ? -1 : TMath::Sqrt((l+m)*2);
155}
156*/
157
158// --------------------------------------------------------------------------
159//
160// Calculates MMath::SignificanceLiMa(s, b, alpha). Returns 0 if the
161// calculation has failed. Otherwise the Li/Ma significance which was
162// calculated. If s<b a negative value is returned.
163//
164Double_t MMath::SignificanceLiMaSigned(Double_t s, Double_t b, Double_t alpha)
165{
166 const Double_t sig = SignificanceLiMa(s, b, alpha);
167 if (sig<=0)
168 return 0;
169
170 return TMath::Sign(sig, s-alpha*b);
171}
172
173// --------------------------------------------------------------------------
174//
175// Return Li/Ma (5) for the error of the excess, under the assumption that
176// the existance of a signal is already known. (basically signal/error
177// calculated by error propagation)
178//
179Double_t MMath::SignificanceExc(Double_t s, Double_t b, Double_t alpha)
180{
181 const Double_t error = ErrorExc(s, b, alpha);
182 if (error==0)
183 return 0;
184
185 const Double_t Ns = s - alpha*b;
186
187 return Ns/error;
188}
189
190// --------------------------------------------------------------------------
191//
192// Calculate the error of s-alpha*b by error propagation
193//
194Double_t MMath::ErrorExc(Double_t s, Double_t b, Double_t alpha)
195{
196 const Double_t sN = s + alpha*alpha*b;
197 return sN<0 ? 0 : TMath::Sqrt(sN);
198}
199
200// --------------------------------------------------------------------------
201//
202// Returns: 2/(sigma*sqrt(2))*integral[0,x](exp(-(x-mu)^2/(2*sigma^2)))
203//
204Double_t MMath::GaussProb(Double_t x, Double_t sigma, Double_t mean)
205{
206 if (x<mean)
207 return 0;
208
209 static const Double_t sqrt2 = TMath::Sqrt(2.);
210
211 const Double_t rc = TMath::Erf((x-mean)/(sigma*sqrt2));
212
213 if (rc<0)
214 return 0;
215 if (rc>1)
216 return 1;
217
218 return rc;
219}
220
221// ------------------------------------------------------------------------
222//
223// Return the "median" (at 68.3%) value of the distribution of
224// abs(a[i]-Median)
225//
226template <class Size, class Element>
227Double_t MMath::MedianDevImp(Size n, const Element *a, Double_t &med)
228{
229 static const Double_t prob = 0.682689477208650697; //MMath::GaussProb(1.0);
230
231 // Sanity check
232 if (n <= 0 || !a)
233 return 0;
234
235 // Get median of distribution
236 med = TMath::Median(n, a);
237
238 // Create the abs(a[i]-med) distribution
239 Double_t arr[n];
240 for (int i=0; i<n; i++)
241 arr[i] = TMath::Abs(a[i]-med);
242
243 //return TMath::Median(n, arr)/0.67449896936; //MMath::GaussProb(x)=0.5
244
245 // Define where to divide (floor because the highest possible is n-1)
246 const Size div = TMath::FloorNint(n*prob);
247
248 // Calculate result
249 Double_t dev = TMath::KOrdStat(n, arr, div);
250 if (n%2 == 0)
251 {
252 dev += TMath::KOrdStat(n, arr, div-1);
253 dev /= 2;
254 }
255
256 return dev;
257}
258
259// ------------------------------------------------------------------------
260//
261// Return the "median" (at 68.3%) value of the distribution of
262// abs(a[i]-Median)
263//
264Double_t MMath::MedianDev(Long64_t n, const Short_t *a, Double_t &med)
265{
266 return MedianDevImp(n, a, med);
267}
268
269// ------------------------------------------------------------------------
270//
271// Return the "median" (at 68.3%) value of the distribution of
272// abs(a[i]-Median)
273//
274Double_t MMath::MedianDev(Long64_t n, const Int_t *a, Double_t &med)
275{
276 return MedianDevImp(n, a, med);
277}
278
279// ------------------------------------------------------------------------
280//
281// Return the "median" (at 68.3%) value of the distribution of
282// abs(a[i]-Median)
283//
284Double_t MMath::MedianDev(Long64_t n, const Float_t *a, Double_t &med)
285{
286 return MedianDevImp(n, a, med);
287}
288
289// ------------------------------------------------------------------------
290//
291// Return the "median" (at 68.3%) value of the distribution of
292// abs(a[i]-Median)
293//
294Double_t MMath::MedianDev(Long64_t n, const Double_t *a, Double_t &med)
295{
296 return MedianDevImp(n, a, med);
297}
298
299// ------------------------------------------------------------------------
300//
301// Return the "median" (at 68.3%) value of the distribution of
302// abs(a[i]-Median)
303//
304Double_t MMath::MedianDev(Long64_t n, const Long_t *a, Double_t &med)
305{
306 return MedianDevImp(n, a, med);
307}
308
309// ------------------------------------------------------------------------
310//
311// Return the "median" (at 68.3%) value of the distribution of
312// abs(a[i]-Median)
313//
314Double_t MMath::MedianDev(Long64_t n, const Long64_t *a, Double_t &med)
315{
316 return MedianDevImp(n, a, med);
317}
318
319Double_t MMath::MedianDev(Long64_t n, const Short_t *a) { Double_t med; return MedianDevImp(n, a, med); }
320Double_t MMath::MedianDev(Long64_t n, const Int_t *a) { Double_t med; return MedianDevImp(n, a, med); }
321Double_t MMath::MedianDev(Long64_t n, const Float_t *a) { Double_t med; return MedianDevImp(n, a, med); }
322Double_t MMath::MedianDev(Long64_t n, const Double_t *a) { Double_t med; return MedianDevImp(n, a, med); }
323Double_t MMath::MedianDev(Long64_t n, const Long_t *a) { Double_t med; return MedianDevImp(n, a, med); }
324Double_t MMath::MedianDev(Long64_t n, const Long64_t *a) { Double_t med; return MedianDevImp(n, a, med); }
325
326// --------------------------------------------------------------------------
327//
328// This function reduces the precision to roughly 0.5% of a Float_t by
329// changing its bit-pattern (Be carefull, in rare cases this function must
330// be adapted to different machines!). This is usefull to enforce better
331// compression by eg. gzip.
332//
333void MMath::ReducePrecision(Float_t &val)
334{
335 UInt_t &f = (UInt_t&)val;
336
337 f += 0x00004000;
338 f &= 0xffff8000;
339}
340
341// -------------------------------------------------------------------------
342//
343// Quadratic interpolation
344//
345// calculate the parameters of a parabula such that
346// y(i) = a + b*x(i) + c*x(i)^2
347//
348// If the determinant==0 an empty TVector3 is returned.
349//
350TVector3 MMath::GetParab(const TVector3 &x, const TVector3 &y)
351{
352 const Double_t x1 = x(0);
353 const Double_t x2 = x(1);
354 const Double_t x3 = x(2);
355
356 const Double_t y1 = y(0);
357 const Double_t y2 = y(1);
358 const Double_t y3 = y(2);
359
360 const double det =
361 + x2*x3*x3 + x1*x2*x2 + x3*x1*x1
362 - x2*x1*x1 - x3*x2*x2 - x1*x3*x3;
363
364
365 if (det==0)
366 return TVector3();
367
368 const double det1 = 1.0/det;
369
370 const double ai11 = x2*x3*x3 - x3*x2*x2;
371 const double ai12 = x3*x1*x1 - x1*x3*x3;
372 const double ai13 = x1*x2*x2 - x2*x1*x1;
373
374 const double ai21 = x2*x2 - x3*x3;
375 const double ai22 = x3*x3 - x1*x1;
376 const double ai23 = x1*x1 - x2*x2;
377
378 const double ai31 = x3 - x2;
379 const double ai32 = x1 - x3;
380 const double ai33 = x2 - x1;
381
382 return TVector3((ai11*y1 + ai12*y2 + ai13*y3) * det1,
383 (ai21*y1 + ai22*y2 + ai23*y3) * det1,
384 (ai31*y1 + ai32*y2 + ai33*y3) * det1);
385}
386
387// --------------------------------------------------------------------------
388//
389// Interpolate the points with x-coordinates vx and y-coordinates vy
390// by a parabola (second order polynomial) and return the value at x.
391//
392Double_t MMath::InterpolParabLin(const TVector3 &vx, const TVector3 &vy, Double_t x)
393{
394 const TVector3 c = GetParab(vx, vy);
395 return c(0) + c(1)*x + c(2)*x*x;
396}
397
398// --------------------------------------------------------------------------
399//
400// Interpolate the points with x-coordinates vx=(-1,0,1) and
401// y-coordinates vy by a parabola (second order polynomial) and return
402// the value at x.
403//
404Double_t MMath::InterpolParabLin(const TVector3 &vy, Double_t x)
405{
406 const TVector3 c(vy(1), (vy(2)-vy(0))/2, vy(0)/2 - vy(1) + vy(2)/2);
407 return c(0) + c(1)*x + c(2)*x*x;
408}
409
410Double_t MMath::InterpolParabLog(const TVector3 &vx, const TVector3 &vy, Double_t x)
411{
412 const Double_t l0 = TMath::Log10(vx(0));
413 const Double_t l1 = TMath::Log10(vx(1));
414 const Double_t l2 = TMath::Log10(vx(2));
415
416 const TVector3 vx0(l0, l1, l2);
417 return InterpolParabLin(vx0, vy, TMath::Log10(x));
418}
419
420Double_t MMath::InterpolParabCos(const TVector3 &vx, const TVector3 &vy, Double_t x)
421{
422 const Double_t l0 = TMath::Cos(vx(0));
423 const Double_t l1 = TMath::Cos(vx(1));
424 const Double_t l2 = TMath::Cos(vx(2));
425
426 const TVector3 vx0(l0, l1, l2);
427 return InterpolParabLin(vx0, vy, TMath::Cos(x));
428}
429
430// --------------------------------------------------------------------------
431//
432// Analytically calculated result of a least square fit of:
433// y = A*e^(B*x)
434// Equal weights
435//
436// It returns TArrayD(2) = { A, B };
437//
438// see: http://mathworld.wolfram.com/LeastSquaresFittingExponential.html
439//
440TArrayD MMath::LeastSqFitExpW1(Int_t n, Double_t *x, Double_t *y)
441{
442 Double_t sumxsqy = 0;
443 Double_t sumylny = 0;
444 Double_t sumxy = 0;
445 Double_t sumy = 0;
446 Double_t sumxylny = 0;
447 for (int i=0; i<n; i++)
448 {
449 sumylny += y[i]*TMath::Log(y[i]);
450 sumxy += x[i]*y[i];
451 sumxsqy += x[i]*x[i]*y[i];
452 sumxylny += x[i]*y[i]*TMath::Log(y[i]);
453 sumy += y[i];
454 }
455
456 const Double_t dev = sumy*sumxsqy - sumxy*sumxy;
457
458 const Double_t a = (sumxsqy*sumylny - sumxy*sumxylny)/dev;
459 const Double_t b = (sumy*sumxylny - sumxy*sumylny)/dev;
460
461 TArrayD rc(2);
462 rc[0] = TMath::Exp(a);
463 rc[1] = b;
464 return rc;
465}
466
467// --------------------------------------------------------------------------
468//
469// Analytically calculated result of a least square fit of:
470// y = A*e^(B*x)
471// Greater weights to smaller values
472//
473// It returns TArrayD(2) = { A, B };
474//
475// see: http://mathworld.wolfram.com/LeastSquaresFittingExponential.html
476//
477TArrayD MMath::LeastSqFitExp(Int_t n, Double_t *x, Double_t *y)
478{
479 // -------- Greater weights to smaller values ---------
480 Double_t sumlny = 0;
481 Double_t sumxlny = 0;
482 Double_t sumxsq = 0;
483 Double_t sumx = 0;
484 for (int i=0; i<n; i++)
485 {
486 sumlny += TMath::Log(y[i]);
487 sumxlny += x[i]*TMath::Log(y[i]);
488
489 sumxsq += x[i]*x[i];
490 sumx += x[i];
491 }
492
493 const Double_t dev = n*sumxsq-sumx*sumx;
494
495 const Double_t a = (sumlny*sumxsq - sumx*sumxlny)/dev;
496 const Double_t b = (n*sumxlny - sumx*sumlny)/dev;
497
498 TArrayD rc(2);
499 rc[0] = TMath::Exp(a);
500 rc[1] = b;
501 return rc;
502}
503
504// --------------------------------------------------------------------------
505//
506// Analytically calculated result of a least square fit of:
507// y = A+B*ln(x)
508//
509// It returns TArrayD(2) = { A, B };
510//
511// see: http://mathworld.wolfram.com/LeastSquaresFittingLogarithmic.html
512//
513TArrayD MMath::LeastSqFitLog(Int_t n, Double_t *x, Double_t *y)
514{
515 Double_t sumylnx = 0;
516 Double_t sumy = 0;
517 Double_t sumlnx = 0;
518 Double_t sumlnxsq = 0;
519 for (int i=0; i<n; i++)
520 {
521 sumylnx += y[i]*TMath::Log(x[i]);
522 sumy += y[i];
523 sumlnx += TMath::Log(x[i]);
524 sumlnxsq += TMath::Log(x[i])*TMath::Log(x[i]);
525 }
526
527 const Double_t b = (n*sumylnx-sumy*sumlnx)/(n*sumlnxsq-sumlnx*sumlnx);
528 const Double_t a = (sumy-b*sumlnx)/n;
529
530 TArrayD rc(2);
531 rc[0] = a;
532 rc[1] = b;
533 return rc;
534}
535
536// --------------------------------------------------------------------------
537//
538// Analytically calculated result of a least square fit of:
539// y = A*x^B
540//
541// It returns TArrayD(2) = { A, B };
542//
543// see: http://mathworld.wolfram.com/LeastSquaresFittingPowerLaw.html
544//
545TArrayD MMath::LeastSqFitPowerLaw(Int_t n, Double_t *x, Double_t *y)
546{
547 Double_t sumlnxlny = 0;
548 Double_t sumlnx = 0;
549 Double_t sumlny = 0;
550 Double_t sumlnxsq = 0;
551 for (int i=0; i<n; i++)
552 {
553 sumlnxlny += TMath::Log(x[i])*TMath::Log(y[i]);
554 sumlnx += TMath::Log(x[i]);
555 sumlny += TMath::Log(y[i]);
556 sumlnxsq += TMath::Log(x[i])*TMath::Log(x[i]);
557 }
558
559 const Double_t b = (n*sumlnxlny-sumlnx*sumlny)/(n*sumlnxsq-sumlnx*sumlnx);
560 const Double_t a = (sumlny-b*sumlnx)/n;
561
562 TArrayD rc(2);
563 rc[0] = TMath::Exp(a);
564 rc[1] = b;
565 return rc;
566}
567
568// --------------------------------------------------------------------------
569//
570// Calculate the intersection of two lines defined by (x1;y1) and (x2;x2)
571// Returns the intersection point.
572//
573// It is assumed that the lines intersect. If there is no intersection
574// TVector2() is returned (which is not destinguishable from
575// TVector2(0,0) if the intersection is at the coordinate source)
576//
577// Formula from: http://mathworld.wolfram.com/Line-LineIntersection.html
578//
579TVector2 MMath::GetIntersectionPoint(const TVector2 &x1, const TVector2 &y1, const TVector2 &x2, const TVector2 &y2)
580{
581 TMatrix d(2,2);
582 d[0][0] = x1.X()-y1.X();
583 d[0][1] = x2.X()-y2.X();
584 d[1][0] = x1.Y()-y1.Y();
585 d[1][1] = x2.Y()-y2.Y();
586
587 const Double_t denom = d.Determinant();
588 if (denom==0)
589 return TVector2();
590
591 TMatrix l1(2,2);
592 TMatrix l2(2,2);
593
594 l1[0][0] = x1.X();
595 l1[0][1] = y1.X();
596 l2[0][0] = x2.X();
597 l2[0][1] = y2.X();
598
599 l1[1][0] = x1.Y();
600 l1[1][1] = y1.Y();
601 l2[1][0] = x2.Y();
602 l2[1][1] = y2.Y();
603
604 TMatrix a(2,2);
605 a[0][0] = l1.Determinant();
606 a[0][1] = l2.Determinant();
607 a[1][0] = x1.X()-y1.X();
608 a[1][1] = x2.X()-y2.X();
609
610 const Double_t X = a.Determinant()/denom;
611
612 a[1][0] = x1.Y()-y1.Y();
613 a[1][1] = x2.Y()-y2.Y();
614
615 const Double_t Y = a.Determinant()/denom;
616
617 return TVector2(X, Y);
618}
619
620// --------------------------------------------------------------------------
621//
622// Solves: x^2 + ax + b = 0;
623// Return number of solutions returned as x1, x2
624//
625Int_t MMath::SolvePol2(Double_t a, Double_t b, Double_t &x1, Double_t &x2)
626{
627 const Double_t r = a*a - 4*b;
628 if (r<0)
629 return 0;
630
631 if (r==0)
632 {
633 x1 = x2 = -a/2;
634 return 1;
635 }
636
637 const Double_t s = TMath::Sqrt(r);
638
639 x1 = (-a+s)/2;
640 x2 = (-a-s)/2;
641
642 return 2;
643}
644
645// --------------------------------------------------------------------------
646//
647// This is a helper function making the execution of SolverPol3 a bit faster
648//
649static inline Double_t ReMul(const TComplex &c1, const TComplex &th)
650{
651 const TComplex c2 = TComplex::Cos(th/3.);
652 return c1.Re() * c2.Re() - c1.Im() * c2.Im();
653}
654
655// --------------------------------------------------------------------------
656//
657// Solves: x^3 + ax^2 + bx + c = 0;
658// Return number of the real solutions, returned as z1, z2, z3
659//
660// Algorithm adapted from http://home.att.net/~srschmitt/cubizen.heml
661// Which is based on the solution given in
662// http://mathworld.wolfram.com/CubicEquation.html
663//
664// -------------------------------------------------------------------------
665//
666// Exact solutions of cubic polynomial equations
667// by Stephen R. Schmitt Algorithm
668//
669// An exact solution of the cubic polynomial equation:
670//
671// x^3 + a*x^2 + b*x + c = 0
672//
673// was first published by Gerolamo Cardano (1501-1576) in his treatise,
674// Ars Magna. He did not discoverer of the solution; a professor of
675// mathematics at the University of Bologna named Scipione del Ferro (ca.
676// 1465-1526) is credited as the first to find an exact solution. In the
677// years since, several improvements to the original solution have been
678// discovered. Zeno source code
679//
680// http://home.att.net/~srschmitt/cubizen.html
681//
682// % compute real or complex roots of cubic polynomial
683// function cubic( var z1, z2, z3 : real, a, b, c : real ) : real
684//
685// var Q, R, D, S, T : real
686// var im, th : real
687//
688// Q := (3*b - a^2)/9
689// R := (9*b*a - 27*c - 2*a^3)/54
690// D := Q^3 + R^2 % polynomial discriminant
691//
692// if (D >= 0) then % complex or duplicate roots
693//
694// S := sgn(R + sqrt(D))*abs(R + sqrt(D))^(1/3)
695// T := sgn(R - sqrt(D))*abs(R - sqrt(D))^(1/3)
696//
697// z1 := -a/3 + (S + T) % real root
698// z2 := -a/3 - (S + T)/2 % real part of complex root
699// z3 := -a/3 - (S + T)/2 % real part of complex root
700// im := abs(sqrt(3)*(S - T)/2) % complex part of root pair
701//
702// else % distinct real roots
703//
704// th := arccos(R/sqrt( -Q^3))
705//
706// z1 := 2*sqrt(-Q)*cos(th/3) - a/3
707// z2 := 2*sqrt(-Q)*cos((th + 2*pi)/3) - a/3
708// z3 := 2*sqrt(-Q)*cos((th + 4*pi)/3) - a/3
709// im := 0
710//
711// end if
712//
713// return im % imaginary part
714//
715// end function
716//
717// see also http://en.wikipedia.org/wiki/Cubic_equation
718//
719Int_t MMath::SolvePol3(Double_t a, Double_t b, Double_t c,
720 Double_t &x1, Double_t &x2, Double_t &x3)
721{
722 // Double_t coeff[4] = { 1, a, b, c };
723 // return TMath::RootsCubic(coeff, x1, x2, x3) ? 1 : 3;
724
725 const Double_t Q = (a*a - 3*b)/9;
726 const Double_t R = (9*b*a - 27*c - 2*a*a*a)/54;
727 const Double_t D = R*R - Q*Q*Q; // polynomial discriminant
728
729 // ----- The single-real / duplicate-roots solution -----
730
731 // D<0: three real roots
732 // D>0: one real root
733 // D==0: maximum two real roots (two identical roots)
734
735 // R==0: only one unique root
736 // R!=0: two roots
737
738 if (D==0)
739 {
740 const Double_t r = MMath::Sqrt3(R);
741
742 x1 = r - a/3.; // real root
743 if (R==0)
744 return 1;
745
746 x2 = 2*r - a/3.; // real root
747 return 2;
748 }
749
750 if (D>0) // complex or duplicate roots
751 {
752 const Double_t sqrtd = TMath::Sqrt(D);
753
754 const Double_t A = TMath::Sign(1., R)*MMath::Sqrt3(TMath::Abs(R)+sqrtd);
755
756 // The case A==0 cannot happen. This would imply D==0
757 // if (A==0)
758 // {
759 // x1 = -a/3;
760 // return 1;
761 // }
762
763 x1 = (A+Q/A)-a/3;
764
765 //const Double_t S = MMath::Sqrt3(R + sqrtd);
766 //const Double_t T = MMath::Sqrt3(R - sqrtd);
767 //x1 = (S+T) - a/3.; // real root
768
769 return 1;
770
771 //z2 = (S + T)/2 - a/3.; // real part of complex root
772 //z3 = (S + T)/2 - a/3.; // real part of complex root
773 //im = fabs(sqrt(3)*(S - T)/2) // complex part of root pair
774 }
775
776 // ----- The general solution with three roots ---
777
778 if (Q==0)
779 return 0;
780
781 if (Q>0) // This is here for speed reasons
782 {
783 const Double_t sqrtq = TMath::Sqrt(Q);
784 const Double_t rq = R/TMath::Abs(Q);
785
786 const Double_t t = TMath::ACos(rq/sqrtq)/3;
787
788 static const Double_t sqrt3 = TMath::Sqrt(3.);
789
790 const Double_t sn = TMath::Sin(t)*sqrt3;
791 const Double_t cs = TMath::Cos(t);
792
793 x1 = 2*sqrtq * cs - a/3;
794 x2 = -sqrtq * (sn + cs) - a/3;
795 x3 = sqrtq * (sn - cs) - a/3;
796
797 /* --- Easier to understand but slower ---
798 const Double_t th1 = TMath::ACos(rq/sqrtq);
799 const Double_t th2 = th1 + TMath::TwoPi();
800 const Double_t th3 = th2 + TMath::TwoPi();
801
802 x1 = 2.*sqrtq * TMath::Cos(th1/3.) - a/3.;
803 x2 = 2.*sqrtq * TMath::Cos(th2/3.) - a/3.;
804 x3 = 2.*sqrtq * TMath::Cos(th3/3.) - a/3.;
805 */
806 return 3;
807 }
808
809 const TComplex sqrtq = TComplex::Sqrt(Q);
810 const Double_t rq = R/TMath::Abs(Q);
811
812 const TComplex th1 = TComplex::ACos(rq/sqrtq);
813 const TComplex th2 = th1 + TMath::TwoPi();
814 const TComplex th3 = th2 + TMath::TwoPi();
815
816 // For ReMul, see bove
817 x1 = ReMul(2.*sqrtq, th1) - a/3.;
818 x2 = ReMul(2.*sqrtq, th2) - a/3.;
819 x3 = ReMul(2.*sqrtq, th3) - a/3.;
820
821 return 3;
822}
823
824// --------------------------------------------------------------------------
825//
826// Format a value and its error corresponding to the rules (note
827// this won't work if the error is more then eight orders smaller than
828// the value)
829//
830void MMath::Format(Double_t &v, Double_t &e)
831{
832 // Valid digits
833 Int_t i = TMath::FloorNint(TMath::Log10(v))-TMath::FloorNint(TMath::Log10(e));
834
835 // Check if error starts with 1 or 2. In this case use one
836 // more valid digit
837 TString error = Form("%.0e", e);
838 if (error[0]=='1' || error[0]=='2')
839 {
840 i++;
841 error = Form("%.1e", e);
842 }
843
844 const char *fmt = Form("%%.%de", i);
845
846 v = atof(Form(fmt, v));
847 e = error.Atof();
848}
Note: See TracBrowser for help on using the repository browser.