Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
motorcycle.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_graphs
3 /// \notebook
4 /// Macro to test scatterplot smoothers: ksmooth, lowess, supsmu
5 /// as described in:
6 ///
7 /// Modern Applied Statistics with S-Plus, 3rd Edition
8 /// W.N. Venables and B.D. Ripley
9 /// Chapter 9: Smooth Regression, Figure 9.1
10 ///
11 /// Example is a set of data on 133 observations of acceleration against time
12 /// for a simulated motorcycle accident, taken from Silverman (1985).
13 ///
14 /// \macro_image
15 /// \macro_code
16 ///
17 /// \author Christian Stratowa, Vienna, Austria
18 
19 #include "TString.h"
20 #include "TInterpreter.h"
21 #include <fstream>
22 #include "TH1.h"
23 #include "TGraphSmooth.h"
24 #include "TCanvas.h"
25 #include "TSystem.h"
26 
27 
28 TCanvas *vC1;
29 TGraph *grin, *grout;
30 
31 void DrawSmooth(Int_t pad, const char *title, const char *xt, const char *yt)
32 {
33  vC1->cd(pad);
34  TH1F *vFrame = gPad->DrawFrame(0,-130,60,70);
35  vFrame->SetTitle(title);
36  vFrame->SetTitleSize(0.2);
37  vFrame->SetXTitle(xt);
38  vFrame->SetYTitle(yt);
39  grin->Draw("P");
40  grout->DrawClone("LPX");
41 }
42 
43 void motorcycle()
44 {
45 // data taken from R library MASS: mcycle.txt
46  TString dir = gROOT->GetTutorialDir();
47  dir.Append("/graphs/");
48  dir.ReplaceAll("/./","/");
49 
50 // read file and add to fit object
51  Double_t *x = new Double_t[133];
52  Double_t *y = new Double_t[133];
53  Double_t vX, vY;
54  Int_t vNData = 0;
55  ifstream vInput;
56  vInput.open(Form("%smotorcycle.dat",dir.Data()));
57  while (1) {
58  vInput >> vX >> vY;
59  if (!vInput.good()) break;
60  x[vNData] = vX;
61  y[vNData] = vY;
62  vNData++;
63  }//while
64  vInput.close();
65  grin = new TGraph(vNData,x,y);
66 
67 // draw graph
68  vC1 = new TCanvas("vC1","Smooth Regression",200,10,900,700);
69  vC1->Divide(2,3);
70 
71 // Kernel Smoother
72 // create new kernel smoother and smooth data with bandwidth = 2.0
73  TGraphSmooth *gs = new TGraphSmooth("normal");
74  grout = gs->SmoothKern(grin,"normal",2.0);
75  DrawSmooth(1,"Kernel Smoother: bandwidth = 2.0","times","accel");
76 
77 // redraw ksmooth with bandwidth = 5.0
78  grout = gs->SmoothKern(grin,"normal",5.0);
79  DrawSmooth(2,"Kernel Smoother: bandwidth = 5.0","","");
80 
81 // Lowess Smoother
82 // create new lowess smoother and smooth data with fraction f = 2/3
83  grout = gs->SmoothLowess(grin,"",0.67);
84  DrawSmooth(3,"Lowess: f = 2/3","","");
85 
86 // redraw lowess with fraction f = 0.2
87  grout = gs->SmoothLowess(grin,"",0.2);
88  DrawSmooth(4,"Lowess: f = 0.2","","");
89 
90 // Super Smoother
91 // create new super smoother and smooth data with default bass = 0 and span = 0
92  grout = gs->SmoothSuper(grin,"",0,0);
93  DrawSmooth(5,"Super Smoother: bass = 0","","");
94 
95 // redraw supsmu with bass = 3 (smoother curve)
96  grout = gs->SmoothSuper(grin,"",3);
97  DrawSmooth(6,"Super Smoother: bass = 3","","");
98 
99 // cleanup
100  delete [] x;
101  delete [] y;
102  delete gs;
103 }