TGraphSmooth.cxx

Go to the documentation of this file.
00001 // @(#)root/hist:$Id: TGraphSmooth.cxx 35246 2010-09-13 14:04:31Z brun $
00002 // Author: Christian Stratowa 30/09/2001
00003 
00004 /*************************************************************************
00005  * Copyright (C) 2006, Rene Brun and Fons Rademakers.                    *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 /******************************************************************************
00013 * Copyright(c) 2001-    , Dr. Christian Stratowa, Vienna, Austria.            *
00014 * Author: Christian Stratowa with help from Rene Brun.                        *
00015 *                                                                             *
00016 * Algorithms for smooth regression adapted from:                              *
00017 * R: A Computer Language for Statistical Data Analysis                        *
00018 * Copyright (C) 1996 Robert Gentleman and Ross Ihaka                          *
00019 * Copyright (C) 1999-2001 Robert Gentleman, Ross Ihaka and the                *
00020 * R Development Core Team                                                     *
00021 * R is free software, for licensing see the GNU General Public License        *
00022 * http://www.ci.tuwien.ac.at/R-project/welcome.html                           *
00023 *                                                                             *
00024 ******************************************************************************/
00025 
00026 
00027 #include "Riostream.h"
00028 #include "TMath.h"
00029 #include "TGraphSmooth.h"
00030 #include "TGraphErrors.h"
00031 
00032 ClassImp(TGraphSmooth);
00033 
00034 //______________________________________________________________________
00035 // TGraphSmooth
00036 //
00037 // A helper class to smooth TGraph
00038 // see examples in $ROOTSYS/tutorials/graphs/motorcycle.C and approx.C
00039 //
00040 //______________________________________________________________________
00041 TGraphSmooth::TGraphSmooth(): TNamed()
00042 {
00043 //*-*-*-*-*-*-*-*-*Default GraphSmooth constructor *-*-*-*-*-*-*-*-*-*-*
00044 //                 ===============================
00045 
00046    fNin  = 0;
00047    fNout = 0;
00048    fGin  = 0;
00049    fGout = 0;
00050    fMinX = 0;
00051    fMaxX = 0;
00052 }
00053 
00054 //______________________________________________________________________
00055 TGraphSmooth::TGraphSmooth(const char *name): TNamed(name,"")
00056 {
00057 //*-*-*-*-*-*-*-*-*GraphSmooth constructor *-*-*-*-*-*-*-*-*-*-*-*-*-*-*
00058 //                 =======================
00059 
00060    fNin  = 0;
00061    fNout = 0;
00062    fGin  = 0;
00063    fGout = 0;
00064    fMinX = 0;
00065    fMaxX = 0;
00066 }
00067 
00068 //______________________________________________________________________
00069 TGraphSmooth::~TGraphSmooth()
00070 {
00071 //*-*-*-*-*-*-*-*-*GraphSmooth destructor*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
00072 //                 ======================
00073 
00074    if (fGout) delete fGout;
00075    fGin  = 0;
00076    fGout = 0;
00077 }
00078 
00079 //______________________________________________________________________
00080 void TGraphSmooth::Smoothin(TGraph *grin)
00081 {
00082 //*-*-*-*-*-*-*-*-*Sort input data points*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
00083 //                 ======================
00084 
00085    if (fGout) {delete fGout; fGout = 0;}
00086    fGin = grin;
00087 
00088    fNin = fGin->GetN();
00089    Double_t *xin = new Double_t[fNin];
00090    Double_t *yin = new Double_t[fNin];
00091    Int_t i;
00092    for (i=0;i<fNin;i++) {
00093       xin[i] = fGin->GetX()[i];
00094       yin[i] = fGin->GetY()[i];
00095    }
00096 
00097 // sort input x, y
00098    Int_t *index = new Int_t[fNin];
00099    TMath::Sort(fNin, xin, index, kFALSE);
00100    for (i=0;i<fNin;i++) {
00101       fGin->SetPoint(i, xin[index[i]], yin[index[i]]);
00102    }
00103 
00104    fMinX = fGin->GetX()[0];  //already sorted!
00105    fMaxX = fGin->GetX()[fNin-1];
00106 
00107    delete [] index;
00108    delete [] xin;
00109    delete [] yin;
00110 }
00111 
00112 //______________________________________________________________________
00113 TGraph *TGraphSmooth::SmoothKern(TGraph *grin, Option_t *option,
00114                       Double_t bandwidth, Int_t nout, Double_t *xout)
00115 {
00116 //*-*-*-*-*-*-*-*-*Smooth data with Kernel smoother*-*-*-*-*-*-*-*-*-*-*
00117 //                 ================================
00118 //
00119 // Smooth grin with the Nadaraya-Watson kernel regression estimate.
00120 //
00121 // Arguments:
00122 // grin:      input graph
00123 //
00124 // option:    the kernel to be used: "box", "normal"
00125 // bandwidth: the bandwidth. The kernels are scaled so that their quartiles
00126 //            (viewed as probability densities) are at +/- 0.25*bandwidth.
00127 // nout:      If xout is not specified, interpolation takes place at equally
00128 //            spaced points spanning the interval [min(x), max(x)], where
00129 //            nout = max(nout, number of input data).
00130 // xout:      an optional set of values at which to evaluate the fit
00131 //
00132 
00133    TString opt = option;
00134    opt.ToLower();
00135    Int_t kernel = 1;
00136    if (opt.Contains("normal")) kernel = 2;
00137 
00138    Smoothin(grin);
00139 
00140    Double_t delta = 0;
00141    Int_t *index = 0;
00142    if (xout == 0) {
00143       fNout = TMath::Max(nout, fNin);
00144       delta = (fMaxX - fMinX)/(fNout - 1);
00145    } else {
00146       fNout = nout;
00147       index = new Int_t[nout];
00148       TMath::Sort(nout, xout, index, kFALSE);
00149    }
00150 
00151    fGout = new TGraph(fNout);
00152    for (Int_t i=0;i<fNout;i++) {
00153       if (xout == 0) fGout->SetPoint(i,fMinX + i*delta, 0);
00154       else           fGout->SetPoint(i,xout[index[i]], 0);
00155    }
00156 
00157    BDRksmooth(fGin->GetX(), fGin->GetY(), fNin, fGout->GetX(),
00158                  fGout->GetY(), fNout, kernel, bandwidth);
00159 
00160    if (index) {delete [] index; index = 0;}
00161 
00162    return fGout;
00163 }
00164 
00165 //______________________________________________________________________
00166 void TGraphSmooth::BDRksmooth(Double_t *x, Double_t *y, Int_t n, Double_t *xp,
00167                    Double_t *yp, Int_t np, Int_t kernel, Double_t bw)
00168 {
00169 //*-*-*-*-*-*-*-*-*Smooth data with specified kernel*-*-*-*-*-*-*-*-*-*-*
00170 //*-*              =================================
00171 //
00172 //   Based on R function ksmooth: Translated to C++ by C. Stratowa
00173 //   (R source file: ksmooth.c by B.D.Ripley Copyright (C) 1998)
00174 //
00175 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
00176 
00177    Int_t    imin = 0;
00178    Double_t cutoff = 0.0;
00179 
00180 // bandwidth is in units of half inter-quartile range
00181    if (kernel == 1) {
00182       bw *= 0.5;
00183       cutoff = bw;
00184    }
00185    if (kernel == 2) {
00186       bw *= 0.3706506;
00187       cutoff = 4*bw;
00188    }
00189 
00190    while ((x[imin] < xp[0] - cutoff) && (imin < n)) imin++;
00191 
00192    for (Int_t j=0;j<np;j++) {
00193       Double_t xx, w;
00194       Double_t num = 0.0;
00195       Double_t den = 0.0;
00196       Double_t x0 = xp[j];
00197       for (Int_t i=imin;i<n;i++) {
00198          if (x[i] < x0 - cutoff) imin = i;
00199          if (x[i] > x0 + cutoff) break;
00200          xx = TMath::Abs(x[i] - x0)/bw;
00201          if (kernel == 1) w = 1;
00202          else             w = TMath::Exp(-0.5*xx*xx);
00203          num += w*y[i];
00204          den += w;
00205       }
00206       if (den > 0) {
00207          yp[j] = num/den;
00208       } else {
00209          yp[j] = 0; //should be NA_REAL (see R.h) or nan("NAN")
00210       }
00211    }
00212 }
00213 
00214 
00215 //______________________________________________________________________
00216 TGraph *TGraphSmooth::SmoothLowess(TGraph *grin, Option_t *option ,
00217                       Double_t span, Int_t iter, Double_t delta)
00218 {
00219 //*-*-*-*-*-*-*-*-*Smooth data with Lowess smoother*-*-*-*-*-*-*-*-*-*-*
00220 //                 ================================
00221 //
00222 // This function performs the computations for the LOWESS smoother
00223 // (see the reference below). Lowess returns the output points
00224 // x and y which give the coordinates of the smooth.
00225 //
00226 // Arguments:
00227 // grin:  Input graph
00228 //
00229 // span:  the smoother span. This gives the proportion of points in the plot
00230 //        which influence the smooth at each value.
00231 //        Larger values give more smoothness.
00232 // iter:  the number of robustifying iterations which should be performed.
00233 //        Using smaller values of iter will make lowess run faster.
00234 // delta: values of x which lie within delta of each other replaced by a
00235 //        single value in the output from lowess.
00236 //        For delta = 0, delta will be calculated.
00237 //
00238 // References:
00239 // Cleveland, W. S. (1979) Robust locally weighted regression and smoothing
00240 //        scatterplots. J. Amer. Statist. Assoc. 74, 829-836.
00241 // Cleveland, W. S. (1981) LOWESS: A program for smoothing scatterplots
00242 //        by robust locally weighted regression.
00243 //        The American Statistician, 35, 54.
00244 //                 ==================
00245 
00246    TString opt = option;
00247    opt.ToLower();
00248 
00249    Smoothin(grin);
00250 
00251    if (delta == 0) {delta = 0.01*(TMath::Abs(fMaxX - fMinX));}
00252 
00253 // output X, Y
00254    fNout = fNin;
00255    fGout = new TGraphErrors(fNout);
00256 
00257    for (Int_t i=0;i<fNout;i++) {
00258       fGout->SetPoint(i,fGin->GetX()[i], 0);
00259    }
00260 
00261    Lowess(fGin->GetX(), fGin->GetY(), fNin, fGout->GetY(), span, iter, delta);
00262 
00263    return fGout;
00264 }
00265 
00266 //______________________________________________________________________
00267 void TGraphSmooth::Lowess(Double_t *x, Double_t *y, Int_t n, Double_t *ys,
00268                    Double_t span, Int_t iter, Double_t delta)
00269 {
00270 //*-*-*-*-*-*-*-*-*Lowess regression smoother*-*-*-*-*-*-*-*-*-*-*-*-*-*
00271 //                 ==========================
00272 //   Based on R function clowess: Translated to C++ by C. Stratowa
00273 //   (R source file: lowess.c by R Development Core Team (C) 1999-2001)
00274 //
00275 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
00276 
00277    Int_t    i, iiter, j, last, m1, m2, nleft, nright, ns;
00278    Double_t alpha, c1, c9, cmad, cut, d1, d2, denom, r;
00279    Bool_t   ok;
00280 
00281    if (n < 2) {
00282       ys[0] = y[0];
00283       return;
00284    }
00285 
00286 // nleft, nright, last, etc. must all be shifted to get rid of these:
00287    x--;
00288    y--;
00289    ys--;
00290 
00291    Double_t *rw  = ((TGraphErrors*)fGout)->GetEX();
00292    Double_t *res = ((TGraphErrors*)fGout)->GetEY();
00293 
00294 // at least two, at most n poInt_ts
00295    ns = TMath::Max(2, TMath::Min(n, (Int_t)(span*n + 1e-7)));
00296 
00297 // robustness iterations
00298    iiter = 1;
00299    while (iiter <= iter+1) {
00300       nleft = 1;
00301       nright = ns;
00302       last = 0;   // index of prev estimated poInt_t
00303       i = 1;      // index of current poInt_t
00304 
00305       for(;;) {
00306          if (nright < n) {
00307          // move nleft,  nright to right if radius decreases
00308             d1 = x[i] - x[nleft];
00309             d2 = x[nright+1] - x[i];
00310 
00311          // if d1 <= d2 with x[nright+1] == x[nright], lowest fixes
00312             if (d1 > d2) {
00313             // radius will not decrease by move right
00314                nleft++;
00315                nright++;
00316                continue;
00317             }
00318          }
00319 
00320       // fitted value at x[i]
00321          Bool_t iterg1 = iiter>1;
00322          Lowest(&x[1], &y[1], n, x[i], ys[i], nleft, nright,
00323                       res, iterg1, rw, ok);
00324          if (!ok) ys[i] = y[i];
00325 
00326       // all weights zero copy over value (all rw==0)
00327          if (last < i-1) {
00328             denom = x[i]-x[last];
00329 
00330          // skipped poInt_ts -- Int_terpolate non-zero - proof?
00331             for(j = last+1; j < i; j++) {
00332                alpha = (x[j]-x[last])/denom;
00333                ys[j] = alpha*ys[i] + (1.-alpha)*ys[last];
00334             }
00335       }
00336 
00337       // last poInt_t actually estimated
00338          last = i;
00339 
00340       // x coord of close poInt_ts
00341          cut = x[last] + delta;
00342          for (i = last+1; i <= n; i++) {
00343             if (x[i] > cut)
00344                break;
00345             if (x[i] == x[last]) {
00346                ys[i] = ys[last];
00347                last = i;
00348             }
00349          }
00350          i = TMath::Max(last+1, i-1);
00351          if (last >= n)
00352             break;
00353       }
00354 
00355    // residuals
00356       for(i=0; i < n; i++)
00357          res[i] = y[i+1] - ys[i+1];
00358 
00359    // compute robustness weights except last time
00360       if (iiter > iter)
00361          break;
00362       for(i=0 ; i<n ; i++)
00363          rw[i] = TMath::Abs(res[i]);
00364 
00365    // compute cmad := 6 * median(rw[], n)
00366       m1 = n/2;
00367    // partial sort, for m1 & m2
00368       Psort(rw, n, m1);
00369       if(n % 2 == 0) {
00370          m2 = n-m1-1;
00371          Psort(rw, n, m2);
00372          cmad = 3.*(rw[m1]+rw[m2]);
00373       } else { /* n odd */
00374          cmad = 6.*rw[m1];
00375       }
00376 
00377       c9 = 0.999*cmad;
00378       c1 = 0.001*cmad;
00379       for(i=0 ; i<n ; i++) {
00380          r = TMath::Abs(res[i]);
00381          if (r <= c1)
00382             rw[i] = 1.;
00383          else if (r <= c9)
00384             rw[i] = (1.-(r/cmad)*(r/cmad))*(1.-(r/cmad)*(r/cmad));
00385          else
00386             rw[i] = 0.;
00387       }
00388       iiter++;
00389    }
00390 }
00391 
00392 //______________________________________________________________________
00393 void TGraphSmooth::Lowest(Double_t *x, Double_t *y, Int_t n, Double_t &xs,
00394                    Double_t &ys, Int_t nleft, Int_t nright, Double_t *w,
00395                    Bool_t userw, Double_t *rw, Bool_t &ok)
00396 {
00397 //*-*-*-*-*-*-*-*-*Fit value at x[i] *-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
00398 //                 =================
00399 //  Based on R function lowest: Translated to C++ by C. Stratowa
00400 //  (R source file: lowess.c by R Development Core Team (C) 1999-2001)
00401 //
00402 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
00403 
00404    Int_t    nrt, j;
00405    Double_t a, b, c, d, h, h1, h9, r, range;
00406 
00407    x--;
00408    y--;
00409    w--;
00410    rw--;
00411 
00412    range = x[n]-x[1];
00413    h = TMath::Max(xs-x[nleft], x[nright]-xs);
00414    h9 = 0.999*h;
00415    h1 = 0.001*h;
00416 
00417 // sum of weights
00418    a = 0.;
00419    j = nleft;
00420    while (j <= n) {
00421    // compute weights (pick up all ties on right)
00422       w[j] = 0.;
00423       r = TMath::Abs(x[j] - xs);
00424       if (r <= h9) {
00425          if (r <= h1) {
00426             w[j] = 1.;
00427          } else {
00428             d = (r/h)*(r/h)*(r/h);
00429             w[j] = (1.- d)*(1.- d)*(1.- d);
00430          }
00431          if (userw)
00432             w[j] *= rw[j];
00433          a += w[j];
00434       } else if (x[j] > xs)
00435          break;
00436       j = j+1;
00437    }
00438 
00439 // rightmost pt (may be greater than nright because of ties)
00440    nrt = j-1;
00441    if (a <= 0.)
00442       ok = kFALSE;
00443    else {
00444       ok = kTRUE;
00445    // weighted least squares: make sum of w[j] == 1
00446       for(j=nleft ; j<=nrt ; j++)
00447          w[j] /= a;
00448       if (h > 0.) {
00449          a = 0.;
00450       // use linear fit weighted center of x values
00451          for(j=nleft ; j<=nrt ; j++)
00452             a += w[j] * x[j];
00453          b = xs - a;
00454          c = 0.;
00455          for(j=nleft ; j<=nrt ; j++)
00456             c += w[j]*(x[j]-a)*(x[j]-a);
00457          if (TMath::Sqrt(c) > 0.001*range) {
00458             b /= c;
00459          // poInt_ts are spread out enough to compute slope
00460             for(j=nleft; j <= nrt; j++)
00461                w[j] *= (b*(x[j]-a) + 1.);
00462          }
00463       }
00464       ys = 0.;
00465       for(j=nleft; j <= nrt; j++)
00466          ys += w[j] * y[j];
00467    }
00468 }
00469 
00470 //______________________________________________________________________
00471 TGraph *TGraphSmooth::SmoothSuper(TGraph *grin, Option_t *,
00472         Double_t bass, Double_t span, Bool_t isPeriodic, Double_t *w)
00473 {
00474 //*-*-*-*-*-*-*-*-*Smooth data with Super smoother*-*-*-*-*-*-*-*-*-*-*-*
00475 //                 ===============================
00476 //
00477 // Smooth the (x, y) values by Friedman's ``super smoother''.
00478 //
00479 // Arguments:
00480 // grin: graph for smoothing
00481 //
00482 // span: the fraction of the observations in the span of the running lines
00483 //        smoother, or 0 to choose this by leave-one-out cross-validation.
00484 // bass: controls the smoothness of the fitted curve.
00485 //        Values of up to 10 indicate increasing smoothness.
00486 // isPeriodic: if TRUE, the x values are assumed to be in [0, 1]
00487 //        and of period 1.
00488 // w:     case weights
00489 //
00490 // Details:
00491 // supsmu is a running lines smoother which chooses between three spans for
00492 // the lines. The running lines smoothers are symmetric, with k/2 data points
00493 // each side of the predicted point, and values of k as 0.5 * n, 0.2 * n and
00494 // 0.05 * n, where n is the number of data points. If span is specified,
00495 // a single smoother with span span * n is used.
00496 //
00497 // The best of the three smoothers is chosen by cross-validation for each
00498 // prediction. The best spans are then smoothed by a running lines smoother
00499 // and the final prediction chosen by linear interpolation.
00500 //
00501 // The FORTRAN code says: ``For small samples (n < 40) or if there are
00502 // substantial serial correlations between observations close in x - value,
00503 // then a prespecified fixed span smoother (span > 0) should be used.
00504 // Reasonable span values are 0.2 to 0.4.''
00505 //
00506 // References:
00507 // Friedman, J. H. (1984) SMART User's Guide.
00508 //           Laboratory for Computational Statistics,
00509 //           Stanford University Technical Report No. 1.
00510 //
00511 // Friedman, J. H. (1984) A variable span scatterplot smoother.
00512 //           Laboratory for Computational Statistics,
00513 //           Stanford University Technical Report No. 5.
00514 //                 ==================
00515 
00516    if (span < 0 || span > 1) {
00517       cout << "Error: Span must be between 0 and 1" << endl;
00518       return 0;
00519    }
00520 
00521    Smoothin(grin);
00522 
00523    Int_t iper = 1;
00524    if (isPeriodic) {
00525       iper = 2;
00526       if (fMinX < 0 || fMaxX > 1) {
00527          cout << "Error: x must be between 0 and 1 for periodic smooth" << endl;
00528          return 0;
00529       }
00530    }
00531 
00532 // output X, Y
00533    fNout = fNin;
00534    fGout = new TGraph(fNout);
00535    Int_t i;
00536    for (i=0; i<fNout; i++) {
00537       fGout->SetPoint(i,fGin->GetX()[i], 0);
00538    }
00539 
00540 // weights
00541    Double_t *weight = new Double_t[fNin];
00542    for (i=0; i<fNin; i++) {
00543       if (w == 0) weight[i] = 1;
00544       else        weight[i] = w[i];
00545    }
00546 
00547 // temporary storage array
00548    Int_t nTmp = (fNin+1)*8;
00549    Double_t *tmp = new Double_t[nTmp];
00550    for (i=0; i<nTmp; i++) {
00551       tmp[i] = 0;
00552    }
00553 
00554    BDRsupsmu(fNin, fGin->GetX(), fGin->GetY(), weight, iper, span, bass, fGout->GetY(), tmp);
00555 
00556    delete [] tmp;
00557    delete [] weight;
00558 
00559    return fGout;
00560 }
00561 
00562 //______________________________________________________________________
00563 void TGraphSmooth::BDRsupsmu(Int_t n, Double_t *x, Double_t *y, Double_t *w,
00564      Int_t iper, Double_t span, Double_t alpha, Double_t *smo, Double_t *sc)
00565 {
00566 //*-*-*-*-*-*-*-*-*Friedmann´s super smoother *-*-*-*-*-*-*-*-*-*-*-*-*-*
00567 //                 ==========================
00568 //
00569 //  super smoother (Friedman, 1984).
00570 //
00571 //  version 10/10/84
00572 //
00573 //  coded  and copywrite (c) 1984 by:
00574 //
00575 //                         Jerome H. Friedman
00576 //                      department of statistics
00577 //                                and
00578 //                 stanford linear accelerator center
00579 //                         stanford university
00580 //
00581 //  all rights reserved.
00582 //
00583 //
00584 //  input:
00585 //     n : number of observations (x,y - pairs).
00586 //     x(n) : ordered abscissa values.
00587 //     y(n) : corresponding ordinate (response) values.
00588 //     w(n) : weight for each (x,y) observation.
00589 //     iper : periodic variable flag.
00590 //        iper=1 => x is ordered interval variable.
00591 //        iper=2 => x is a periodic variable with values
00592 //                  in the range (0.0,1.0) and period 1.0.
00593 //     span : smoother span (fraction of observations in window).
00594 //            span=0.0 => automatic (variable) span selection.
00595 //     alpha : controls high frequency (small span) penality
00596 //             used with automatic span selection (bass tone control).
00597 //             (alpha.le.0.0 or alpha.gt.10.0 => no effect.)
00598 //  output:
00599 //    smo(n) : smoothed ordinate (response) values.
00600 //  scratch:
00601 //    sc(n,7) : internal working storage.
00602 //
00603 //  note:
00604 //     for small samples (n < 40) or if there are substantial serial
00605 //     correlations between observations close in x - value, then
00606 //     a prespecified fixed span smoother (span > 0) should be
00607 //     used. reasonable span values are 0.2 to 0.4.
00608 //
00609 // current implementation:
00610 //   Based on R function supsmu: Translated to C++ by C. Stratowa
00611 //   (R source file: ppr.f by B.D.Ripley Copyright (C) 1994-97)
00612 //
00613 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
00614 
00615 // Local variables
00616    Int_t sc_offset;
00617    Int_t i, j, jper;
00618    Double_t a, f;
00619    Double_t sw, sy, resmin, vsmlsq;
00620    Double_t scale;
00621    Double_t d1, d2;
00622 
00623    Double_t spans[3] = { 0.05, 0.2, 0.5 };
00624    Double_t big = 1e20;
00625    Double_t sml = 1e-7;
00626    Double_t eps = 0.001;
00627 
00628 // Parameter adjustments
00629    sc_offset = n + 1;
00630    sc -= sc_offset;
00631    --smo;
00632    --w;
00633    --y;
00634    --x;
00635 
00636 // Function Body
00637    if (x[n] <= x[1]) {
00638       sy = 0.0;
00639       sw = sy;
00640       for (j=1;j<=n;++j) {
00641          sy += w[j] * y[j];
00642          sw += w[j];
00643       }
00644 
00645       a = 0.0;
00646       if (sw > 0.0) a = sy / sw;
00647       for (j=1;j<=n;++j) smo[j] = a;
00648       return;
00649    }
00650 
00651    i = (Int_t)(n / 4);
00652    j = i * 3;
00653    scale = x[j] - x[i];
00654    while (scale <= 0.0) {
00655       if (j < n) ++j;
00656       if (i > 1) --i;
00657       scale = x[j] - x[i];
00658    }
00659 
00660 // Computing 2nd power
00661    d1 = eps * scale;
00662    vsmlsq = d1 * d1;
00663    jper = iper;
00664    if (iper == 2 && (x[1] < 0.0 || x[n] > 1.0)) {
00665       jper = 1;
00666    }
00667    if (jper < 1 || jper > 2) {
00668       jper = 1;
00669    }
00670    if (span > 0.0) {
00671       BDRsmooth(n, &x[1], &y[1], &w[1], span, jper, vsmlsq,
00672                       &smo[1], &sc[sc_offset]);
00673       return;
00674    }
00675 
00676    Double_t *h = new Double_t[n+1];
00677    for (i = 1; i <= 3; ++i) {
00678       BDRsmooth(n, &x[1], &y[1], &w[1], spans[i - 1], jper, vsmlsq,
00679                       &sc[((i<<1)-1)*n + 1], &sc[n*7 + 1]);
00680       BDRsmooth(n, &x[1], &sc[n*7 + 1], &w[1], spans[1], -jper, vsmlsq,
00681                       &sc[(i<<1)*n + 1], &h[1]);
00682    }
00683 
00684    for (j=1; j<=n; ++j) {
00685       resmin = big;
00686       for (i=1; i<=3; ++i) {
00687          if (sc[j + (i<<1)*n] < resmin) {
00688             resmin = sc[j + (i<<1)*n];
00689             sc[j + n*7] = spans[i-1];
00690          }
00691       }
00692 
00693       if (alpha>0.0 && alpha<=10.0 && resmin<sc[j + n*6]  && resmin>0.0) {
00694       // Computing MAX
00695          d1 = TMath::Max(sml,(resmin/sc[j + n*6]));
00696          d2 = 10. - alpha;
00697          sc[j + n*7] += (spans[2] - sc[j + n*7]) * TMath::Power(d1, d2);
00698       }
00699    }
00700 
00701    BDRsmooth(n, &x[1], &sc[n*7 + 1], &w[1], spans[1], -jper, vsmlsq,
00702                    &sc[(n<<1) + 1], &h[1]);
00703 
00704    for (j=1; j<=n; ++j) {
00705       if (sc[j + (n<<1)] <= spans[0]) {
00706          sc[j + (n<<1)] = spans[0];
00707       }
00708       if (sc[j + (n<<1)] >= spans[2]) {
00709          sc[j + (n<<1)] = spans[2];
00710       }
00711       f = sc[j + (n<<1)] - spans[1];
00712       if (f < 0.0) {
00713          f = -f / (spans[1] - spans[0]);
00714          sc[j + (n<<2)] = (1.0 - f) * sc[j + n*3] + f * sc[j + n];
00715       } else {
00716          f /= spans[2] - spans[1];
00717          sc[j + (n<<2)] = (1.0 - f) * sc[j + n*3] + f * sc[j + n*5];
00718       }
00719    }
00720 
00721    BDRsmooth(n, &x[1], &sc[(n<<2) + 1], &w[1], spans[0], -jper, vsmlsq,
00722                    &smo[1], &h[1]);
00723 
00724    delete [] h;
00725    return;
00726 }
00727 
00728 //______________________________________________________________________
00729 void TGraphSmooth::BDRsmooth(Int_t n, Double_t *x, Double_t *y, Double_t *w,
00730      Double_t span, Int_t iper, Double_t vsmlsq, Double_t *smo, Double_t *acvr)
00731 {
00732 //*-*-*-*-*-*-*-*-* Function for super smoother *-*-*-*-*-*-*-*-*-*-*-*
00733 //                 ============================
00734 //
00735 //   Based on R function supsmu: Translated to C++ by C. Stratowa
00736 //   (R source file: ppr.f by B.D.Ripley Copyright (C) 1994-97)
00737 //
00738 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
00739 
00740 // Local variables
00741    Int_t i, j, j0, in, out, it, jper, ibw;
00742    Double_t a, h1, d1;
00743    Double_t xm, ym, wt, sy, fbo, fbw;
00744    Double_t cvar, var, tmp, xti, xto;
00745 
00746 // Parameter adjustments
00747    --acvr;
00748    --smo;
00749    --w;
00750    --y;
00751    --x;
00752 
00753 // Function Body
00754    xm = 0.;
00755    ym = xm;
00756    var = ym;
00757    cvar = var;
00758    fbw = cvar;
00759    jper = TMath::Abs(iper);
00760 
00761    ibw = (Int_t)(span * 0.5 * n + 0.5);
00762    if (ibw < 2) {
00763       ibw = 2;
00764    }
00765 
00766    it = 2*ibw + 1;
00767    for (i=1; i<=it; ++i) {
00768       j = i;
00769       if (jper == 2) {
00770          j = i - ibw - 1;
00771       }
00772       xti = x[j];
00773       if (j < 1) {
00774          j = n + j;
00775          xti = x[j] - 1.0;
00776       }
00777       wt = w[j];
00778       fbo = fbw;
00779       fbw += wt;
00780       if (fbw > 0.0) {
00781          xm = (fbo * xm + wt * xti) / fbw;
00782          ym = (fbo * ym + wt * y[j]) / fbw;
00783       }
00784       tmp = 0.0;
00785       if (fbo > 0.0) {
00786          tmp = fbw * wt * (xti - xm) / fbo;
00787       }
00788       var += tmp * (xti - xm);
00789       cvar += tmp * (y[j] - ym);
00790    }
00791 
00792    for (j=1; j<=n; ++j) {
00793       out = j - ibw - 1;
00794       in = j + ibw;
00795       if (!(jper != 2 && (out < 1 || in > n))) {
00796          if (out < 1) {
00797             out = n + out;
00798             xto = x[out] - 1.0;
00799             xti = x[in];
00800          } else if (in > n) {
00801             in -= n;
00802             xti = x[in] + 1.0;
00803             xto = x[out];
00804          } else {
00805             xto = x[out];
00806             xti = x[in];
00807          }
00808 
00809          wt = w[out];
00810          fbo = fbw;
00811          fbw -= wt;
00812          tmp = 0.0;
00813          if (fbw > 0.0) {
00814             tmp = fbo * wt * (xto - xm) / fbw;
00815          }
00816          var -= tmp * (xto - xm);
00817          cvar -= tmp * (y[out] - ym);
00818          if (fbw > 0.0) {
00819             xm = (fbo * xm - wt * xto) / fbw;
00820             ym = (fbo * ym - wt * y[out]) / fbw;
00821          }
00822          wt = w[in];
00823          fbo = fbw;
00824          fbw += wt;
00825          if (fbw > 0.0) {
00826             xm = (fbo * xm + wt * xti) / fbw;
00827             ym = (fbo * ym + wt * y[in]) / fbw;
00828          }
00829          tmp = 0.0;
00830          if (fbo > 0.0) {
00831             tmp = fbw * wt * (xti - xm) / fbo;
00832          }
00833          var += tmp * (xti - xm);
00834          cvar += tmp * (y[in] - ym);
00835       }
00836 
00837       a = 0.0;
00838       if (var > vsmlsq) {
00839          a = cvar / var;
00840       }
00841       smo[j] = a * (x[j] - xm) + ym;
00842 
00843       if (iper <= 0) {
00844          continue;
00845       }
00846 
00847       h1 = 0.0;
00848       if (fbw > 0.0) {
00849          h1 = 1.0 / fbw;
00850       }
00851       if (var > vsmlsq) {
00852       // Computing 2nd power
00853          d1 = x[j] - xm;
00854          h1 += d1 * d1 / var;
00855       }
00856 
00857       acvr[j] = 0.0;
00858       a = 1.0 - w[j] * h1;
00859       if (a > 0.0) {
00860          acvr[j] = TMath::Abs(y[j] - smo[j]) / a;
00861          continue;
00862       }
00863       if (j > 1) {
00864          acvr[j] = acvr[j-1];
00865       }
00866    }
00867 
00868    j = 1;
00869    do {
00870       j0 = j;
00871       sy = smo[j] * w[j];
00872       fbw = w[j];
00873       if (j < n) {
00874          do {
00875             if (x[j + 1] > x[j]) {
00876                break;
00877             }
00878             ++j;
00879             sy += w[j] * smo[j];
00880             fbw += w[j];
00881          } while (j < n);
00882       }
00883 
00884       if (j > j0) {
00885          a = 0.0;
00886          if (fbw > 0.0) {
00887             a = sy / fbw;
00888          }
00889          for (i=j0; i<=j; ++i) {
00890             smo[i] = a;
00891          }
00892       }
00893       ++j;
00894    } while (j <= n);
00895 
00896    return;
00897 }
00898 
00899 //______________________________________________________________________
00900 void TGraphSmooth::Approxin(TGraph *grin, Int_t /*iKind*/, Double_t &ylow,
00901      Double_t &yhigh, Int_t rule, Int_t iTies)
00902 {
00903 //*-*-*-*-*-*-*-*-*Sort data points and eliminate double x values*-*-*-*
00904 //                 ==============================================
00905 
00906    if (fGout) {delete fGout; fGout = 0;}
00907    fGin = grin;
00908 
00909    fNin = fGin->GetN();
00910    Double_t *xin = new Double_t[fNin];
00911    Double_t *yin = new Double_t[fNin];
00912    Int_t i;
00913    for (i=0;i<fNin;i++) {
00914       xin[i] = fGin->GetX()[i];
00915       yin[i] = fGin->GetY()[i];
00916    }
00917 
00918 // sort/rank input x, y
00919    Int_t *index = new Int_t[fNin];
00920    Int_t *rank  = new Int_t[fNin];
00921    Rank(fNin, xin, index, rank, kFALSE);
00922 
00923 // input X, Y
00924    Int_t vNDup = 0;
00925    Int_t k = 0;
00926    Int_t *dup = new Int_t[fNin];
00927    Double_t *x = new Double_t[fNin];
00928    Double_t *y = new Double_t[fNin];
00929    Double_t vMean, vMin, vMax;
00930    for (i=1;i<fNin+1;i++) {
00931       Int_t ndup = 1;
00932       vMin = vMean = vMax = yin[index[i-1]];
00933       while ((i < fNin) && (rank[index[i]] == rank[index[i-1]])) {
00934          vMean += yin[index[i]];
00935          vMax = (vMax < yin[index[i]]) ? yin[index[i]] : vMax;
00936          vMin = (vMin > yin[index[i]]) ? yin[index[i]] : vMin;
00937          dup[vNDup] = i;
00938          i++;
00939          ndup++;
00940          vNDup++;
00941       }
00942       x[k] = xin[index[i-1]];
00943       if (ndup == 1) {y[k++] = yin[index[i-1]];}
00944       else switch(iTies) {
00945          case 1:
00946             y[k++] = vMean/ndup;
00947             break;
00948          case 2:
00949             y[k++] = vMin;
00950             break;
00951          case 3:
00952             y[k++] = vMax;
00953             break;
00954          default:
00955             y[k++] = vMean/ndup;
00956             break;
00957       }
00958    }
00959    fNin = k;
00960 
00961 // set unique sorted input data x,y as final graph points
00962    fGin->Set(fNin);
00963    for (i=0;i<fNin;i++) {
00964       fGin->SetPoint(i, x[i], y[i]);
00965    }
00966 
00967    fMinX = fGin->GetX()[0];  //already sorted!
00968    fMaxX = fGin->GetX()[fNin-1];
00969 
00970 // interpolate outside interval [min(x),max(x)]
00971    switch(rule) {
00972       case 1:
00973          ylow  = 0;   // = nan("NAN") ??
00974          yhigh = 0;   // = nan("NAN") ??
00975          break;
00976       case 2:
00977          ylow  = fGin->GetY()[0];
00978          yhigh = fGin->GetY()[fNin-1];
00979          break;
00980       default:
00981          break;
00982    }
00983 
00984 // cleanup
00985    delete [] x;
00986    delete [] y;
00987    delete [] dup;
00988    delete [] rank;
00989    delete [] index;
00990    delete [] xin;
00991    delete [] yin;
00992 }
00993 
00994 //______________________________________________________________________
00995 TGraph *TGraphSmooth::Approx(TGraph *grin, Option_t *option, Int_t nout, Double_t *xout,
00996         Double_t yleft, Double_t yright, Int_t rule, Double_t f, Option_t *ties)
00997 {
00998 //*-*-*-*-*-*-*-*-*Approximate data points*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
00999 //                 =======================
01000 //
01001 // Arguments:
01002 // grin:  graph giving the coordinates of the points to be interpolated.
01003 //        Alternatively a single plotting structure can be specified:
01004 //
01005 // option: specifies the interpolation method to be used.
01006 //        Choices are "linear" (iKind = 1) or "constant" (iKind = 2).
01007 // nout:  If xout is not specified, interpolation takes place at n equally
01008 //        spaced points spanning the interval [min(x), max(x)], where
01009 //        nout = max(nout, number of input data).
01010 // xout:  an optional set of values specifying where interpolation is to
01011 //        take place.
01012 // yleft: the value to be returned when input x values less than min(x).
01013 //        The default is defined by the value of rule given below.
01014 // yright: the value to be returned when input x values greater than max(x).
01015 //        The default is defined by the value of rule given below.
01016 // rule:  an integer describing how interpolation is to take place outside
01017 //        the interval [min(x), max(x)]. If rule is 0 then the given yleft
01018 //        and yright values are returned, if it is 1 then 0 is returned
01019 //        for such points and if it is 2, the value at the closest data
01020 //        extreme is used.
01021 // f:     For method="constant" a number between 0 and 1 inclusive,
01022 //        indicating a compromise between left- and right-continuous step
01023 //        functions. If y0 and y1 are the values to the left and right of
01024 //        the point then the value is y0*f+y1*(1-f) so that f=0 is
01025 //        right-continuous and f=1 is left-continuous
01026 // ties:  Handling of tied x values. An integer describing a function with
01027 //        a single vector argument returning a single number result:
01028 //        ties = "ordered" (iTies = 0): input x are "ordered"
01029 //        ties = "mean"    (iTies = 1): function "mean"
01030 //        ties = "min"     (iTies = 2): function "min"
01031 //        ties = "max"     (iTies = 3): function "max"
01032 //
01033 // Details:
01034 // At least two complete (x, y) pairs are required.
01035 // If there are duplicated (tied) x values and ties is a function it is
01036 // applied to the y values for each distinct x value. Useful functions in
01037 // this context include mean, min, and max.
01038 // If ties="ordered" the x values are assumed to be already ordered. The
01039 // first y value will be used for interpolation to the left and the last
01040 // one for interpolation to the right.
01041 //
01042 // Value:
01043 // approx returns a graph with components x and y, containing n coordinates
01044 // which interpolate the given data points according to the method (and rule)
01045 // desired.
01046 
01047    TString opt = option;
01048    opt.ToLower();
01049    Int_t iKind = 0;
01050    if (opt.Contains("linear")) iKind = 1;
01051    else if (opt.Contains("constant")) iKind = 2;
01052 
01053    if (f < 0 || f > 1) {
01054       cout << "Error: Invalid f value" << endl;
01055       return 0;
01056    }
01057 
01058    opt = ties;
01059    opt.ToLower();
01060    Int_t iTies = 0;
01061    if (opt.Contains("ordered")) {
01062       iTies = 0;
01063    } else if (opt.Contains("mean")) {
01064       iTies = 1;
01065    } else if (opt.Contains("min")) {
01066       iTies = 2;
01067    } else if (opt.Contains("max")) {
01068       iTies = 3;
01069    } else {
01070       cout << "Error: Method not known: " << ties << endl;
01071       return 0;
01072    }
01073 
01074 // input X, Y
01075    Double_t ylow  = yleft;
01076    Double_t yhigh = yright;
01077    Approxin(grin, iKind, ylow, yhigh, rule, iTies);
01078 
01079 // output X, Y
01080    Double_t delta = 0;
01081    fNout = nout;
01082    if (xout == 0) {
01083       fNout = TMath::Max(nout, fNin);
01084       delta = (fMaxX - fMinX)/(fNout - 1);
01085    }
01086 
01087    fGout = new TGraph(fNout);
01088 
01089    Double_t x;
01090    for (Int_t i=0;i<fNout;i++) {
01091       if (xout == 0) x = fMinX + i*delta;
01092       else           x = xout[i];
01093       Double_t yout = Approx1(x, f, fGin->GetX(), fGin->GetY(), fNin, iKind, ylow, yhigh);
01094       fGout->SetPoint(i,x, yout);
01095    }
01096 
01097    return fGout;
01098 }
01099 
01100 //______________________________________________________________________
01101 Double_t TGraphSmooth::Approx1(Double_t v, Double_t f, Double_t *x, Double_t *y,
01102          Int_t n, Int_t iKind, Double_t ylow, Double_t yhigh)
01103 {
01104 //*-*-*-*-*-*-*-*-*Approximate one data point*-*-*-*-*-*-*-*-*-*-*-*-*-*
01105 //*-*              ==========================
01106 //
01107 //   Approximate  y(v),  given (x,y)[i], i = 0,..,n-1
01108 //   Based on R function approx1: Translated to C++ by Christian Stratowa
01109 //   (R source file: approx.c by R Development Core Team (C) 1999-2001)
01110 //
01111 
01112    Int_t i = 0;
01113    Int_t j = n - 1;
01114 
01115 // handle out-of-domain points
01116    if(v < x[i]) return ylow;
01117    if(v > x[j]) return yhigh;
01118 
01119 // find the correct interval by bisection
01120    while(i < j - 1) {
01121       Int_t ij = (i + j)/2;
01122       if(v < x[ij]) j = ij;
01123       else i = ij;
01124    }
01125 
01126 // interpolation
01127    if(v == x[j]) return y[j];
01128    if(v == x[i]) return y[i];
01129 
01130    if(iKind == 1) { // linear
01131       return y[i] + (y[j] - y[i]) * ((v - x[i])/(x[j] - x[i]));
01132    } else { // 2 : constant
01133       return y[i] * (1-f) + y[j] * f;
01134    }
01135 }
01136 
01137 // helper functions
01138 //______________________________________________________________________
01139 Int_t TGraphSmooth::Rcmp(Double_t x, Double_t y)
01140 {
01141 //   static function
01142 //   if (ISNAN(x))   return 1;
01143 //   if (ISNAN(y))   return -1;
01144    if (x < y)      return -1;
01145    if (x > y)      return 1;
01146    return 0;
01147 }
01148 
01149 //______________________________________________________________________
01150 void TGraphSmooth::Psort(Double_t *x, Int_t n, Int_t k)
01151 {
01152 //   static function
01153 //   based on R function rPsort: adapted to C++ by Christian Stratowa
01154 //   (R source file: R_sort.c by R Development Core Team (C) 1999-2001)
01155 //
01156 
01157    Double_t v, w;
01158    Int_t pL, pR, i, j;
01159 
01160    for (pL = 0, pR = n - 1; pL < pR; ) {
01161       v = x[k];
01162       for(i = pL, j = pR; i <= j;) {
01163          while (TGraphSmooth::Rcmp(x[i], v) < 0) i++;
01164          while (TGraphSmooth::Rcmp(v, x[j]) < 0) j--;
01165          if (i <= j) { w = x[i]; x[i++] = x[j]; x[j--] = w; }
01166       }
01167       if (j < k) pL = i;
01168       if (k < i) pR = j;
01169    }
01170 }
01171 
01172 //______________________________________________________________________
01173 void TGraphSmooth::Rank(Int_t n, Double_t *a, Int_t *index, Int_t *rank, Bool_t down)
01174 {
01175 //   static function
01176    if (n <= 0) return;
01177    if (n == 1) {
01178       index[0] = 0;
01179       rank[0] = 0;
01180       return;
01181    }
01182 
01183    TMath::Sort(n,a,index,down);
01184 
01185    Int_t k = 0;
01186    for (Int_t i=0;i<n;i++) {
01187       if ((i > 0) && (a[index[i]] == a[index[i-1]])) {
01188          rank[index[i]] = i-1;
01189          k++;
01190       }
01191       rank[index[i]] = i-k;
01192    }
01193 }

Generated on Tue Jul 5 14:24:18 2011 for ROOT_528-00b_version by  doxygen 1.5.1