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

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