TCutG.cxx

Go to the documentation of this file.
00001 // @(#)root/graf:$Id: TCutG.cxx 37309 2010-12-06 00:10:29Z pcanal $
00002 // Author: Rene Brun   16/05/97
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2000, 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 //______________________________________________________________________________
00014 /* Begin_Html
00015 <center><h2>Graphical cut class</h2></center>
00016 A Graphical cut.
00017 <p>
00018 A TCutG object is a closed polygon defining a closed region in a x,y plot.
00019 It can be created via the graphics editor option "CutG" or directly by
00020 invoking its constructor. The first and last points should be the same.
00021 <p>
00022 To create a TCutG via the graphics editor, use the left button to select the
00023 points building the contour of the cut. Click on the right button to close the
00024 TCutG. When it is created via the graphics editor, the TCutG object is named
00025 "CUTG". It is recommended to immediatly change the name by using the context
00026 menu item "SetName". When the graphics editor is used, the names of the
00027 variables X,Y are automatically taken from the current pad title.
00028 <p>
00029 Example:
00030 <p>
00031 Assume a TTree object T and:
00032 <pre>
00033      Root > T.Draw("abs(fMomemtum)%fEtot")
00034 </pre>
00035 the TCutG members fVarX, fVary will be set to:
00036 <pre>
00037      fVarx = fEtot
00038      fVary = abs(fMomemtum)
00039 </pre>
00040 
00041 A graphical cut can be used in a TTree selection expression:
00042 <pre>
00043     Root > T.Draw("fEtot","cutg1")
00044 </pre>
00045 where "cutg1" is the name of an existing graphical cut.
00046 <p>
00047 Note that, as shown in the example above, a graphical cut may be used in a
00048 selection expression when drawing TTrees expressions of 1-d, 2-d or
00049 3-dimensions. The expressions used in TTree::Draw can reference the variables in
00050 the fVarX, fVarY of the graphical cut plus other variables.
00051 <p>
00052 When the TCutG object is created by TTree::Draw, it is added to the list of special objects in
00053 the main TROOT object pointed by gROOT. To retrieve a pointer to this object
00054 from the code or command line, do:
00055 <pre>
00056     TCutG *mycutg;
00057     mycutg = (TCutG*)gROOT->GetListOfSpecials()->FindObject("CUTG")
00058     mycutg->SetName("mycutg");
00059 </pre>
00060 <p>
00061 When the TCutG is not created via TTree::Draw, one must set the variable names
00062 corresponding to x,y if one wants to use the cut as input to TTree::Draw,eg
00063 <pre>
00064    TCutG *cutg = new TCutG("mycut",5);
00065    cutg->SetVarX("y");
00066    cutg->SetVarY("x");
00067    cutg->SetPoint(0,-0.3586207,1.509534);
00068    cutg->SetPoint(1,-1.894181,-0.529661);
00069    cutg->SetPoint(2,0.07780173,-1.21822);
00070    cutg->SetPoint(3,-1.0375,-0.07944915);
00071    cutg->SetPoint(4,0.756681,0.1853814);
00072    cutg->SetPoint(5,-0.3586207,1.509534);
00073 </pre>
00074    
00075 Example of use of a TCutG in TTree::Draw:
00076 <pre>
00077        tree.Draw("x:y","mycutg && z>0 %% sqrt(x)>1")
00078 </pre>
00079 <p>
00080 A Graphical cut may be drawn via TGraph::Draw. It can be edited like a normal
00081 TGraph.
00082 <p>
00083 A Graphical cut may be saved to a file via TCutG::Write.
00084 End_Html */
00085 
00086 #include <string.h>
00087 
00088 #include "Riostream.h"
00089 #include "TROOT.h"
00090 #include "TCutG.h"
00091 #include "TVirtualPad.h"
00092 #include "TPaveText.h"
00093 #include "TH2.h"
00094 #include "TClass.h"
00095 #include "TMath.h"
00096 
00097 ClassImp(TCutG)
00098 
00099 
00100 //______________________________________________________________________________
00101 TCutG::TCutG() : TGraph()
00102 {
00103    // TCutG default constructor.
00104 
00105    fObjectX  = 0;
00106    fObjectY  = 0;
00107 }
00108 
00109 
00110 //______________________________________________________________________________
00111 TCutG::TCutG(const TCutG &cutg)
00112       :TGraph(cutg)
00113 {
00114    // TCutG copy constructor
00115 
00116    fVarX    = cutg.fVarX;
00117    fVarY    = cutg.fVarY;
00118    fObjectX = cutg.fObjectX;
00119    fObjectY = cutg.fObjectY;
00120 }
00121 
00122 
00123 //______________________________________________________________________________
00124 TCutG::TCutG(const char *name, Int_t n)
00125       :TGraph(n)
00126 {
00127    // TCutG normal constructor.
00128 
00129    fObjectX  = 0;
00130    fObjectY  = 0;
00131    SetName(name);
00132    delete gROOT->GetListOfSpecials()->FindObject(name);
00133    gROOT->GetListOfSpecials()->Add(this);
00134 
00135    // Take name of cut variables from pad title if title contains ":"
00136    if (gPad) {
00137       TPaveText *ptitle = (TPaveText*)gPad->FindObject("title");
00138       if (!ptitle) return;
00139       TText *ttitle = ptitle->GetLineWith(":");
00140       if (!ttitle) ttitle = ptitle->GetLineWith("{");
00141       if (!ttitle) ttitle = ptitle->GetLine(0);
00142       if (!ttitle) return;
00143       const char *title = ttitle->GetTitle();
00144       Int_t nch = strlen(title);
00145       char *vars = new char[nch+1];
00146       strlcpy(vars,title,nch+1);
00147       char *col = strstr(vars,":");
00148       if (col) {
00149          *col = 0;
00150          col++;
00151          char *brak = strstr(col," {");
00152          if (brak) *brak = 0;
00153          fVarY = vars;
00154          fVarX = col;
00155       } else {
00156          char *brak = strstr(vars," {");
00157          if (brak) *brak = 0;
00158          fVarX = vars;
00159       }
00160       delete [] vars;
00161    }
00162 }
00163 
00164 
00165 //______________________________________________________________________________
00166 TCutG::TCutG(const char *name, Int_t n, const Float_t *x, const Float_t *y)
00167       :TGraph(n,x,y)
00168 {
00169    // TCutG normal constructor.
00170 
00171    fObjectX  = 0;
00172    fObjectY  = 0;
00173    SetName(name);
00174    delete gROOT->GetListOfSpecials()->FindObject(name);
00175    gROOT->GetListOfSpecials()->Add(this);
00176 
00177    // Take name of cut variables from pad title if title contains ":"
00178    if (gPad) {
00179       TPaveText *ptitle = (TPaveText*)gPad->FindObject("title");
00180       if (!ptitle) return;
00181       TText *ttitle = ptitle->GetLineWith(":");
00182       if (!ttitle) ttitle = ptitle->GetLineWith("{");
00183       if (!ttitle) ttitle = ptitle->GetLine(0);
00184       if (!ttitle) return;
00185       const char *title = ttitle->GetTitle();
00186       Int_t nch = strlen(title);
00187       char *vars = new char[nch+1];
00188       strlcpy(vars,title,nch+1);
00189       char *col = strstr(vars,":");
00190       if (col) {
00191          *col = 0;
00192          col++;
00193          char *brak = strstr(col," {");
00194          if (brak) *brak = 0;
00195          fVarY = vars;
00196          fVarX = col;
00197       } else {
00198          char *brak = strstr(vars," {");
00199          if (brak) *brak = 0;
00200          fVarX = vars;
00201       }
00202       delete [] vars;
00203    }
00204 }
00205 
00206 
00207 //______________________________________________________________________________
00208 TCutG::TCutG(const char *name, Int_t n, const Double_t *x, const Double_t *y)
00209       :TGraph(n,x,y)
00210 {
00211    // TCutG normal constructor.
00212 
00213    fObjectX  = 0;
00214    fObjectY  = 0;
00215    SetName(name);
00216    delete gROOT->GetListOfSpecials()->FindObject(name);
00217    gROOT->GetListOfSpecials()->Add(this);
00218 
00219    // Take name of cut variables from pad title if title contains ":"
00220    if (gPad) {
00221       TPaveText *ptitle = (TPaveText*)gPad->FindObject("title");
00222       if (!ptitle) return;
00223       TText *ttitle = ptitle->GetLineWith(":");
00224       if (!ttitle) ttitle = ptitle->GetLineWith("{");
00225       if (!ttitle) ttitle = ptitle->GetLine(0);
00226       if (!ttitle) return;
00227       const char *title = ttitle->GetTitle();
00228       Int_t nch = strlen(title);
00229       char *vars = new char[nch+1];
00230       strlcpy(vars,title,nch+1);
00231       char *col = strstr(vars,":");
00232       if (col) {
00233          *col = 0;
00234          col++;
00235          char *brak = strstr(col," {");
00236          if (brak) *brak = 0;
00237          fVarY = vars;
00238          fVarX = col;
00239       } else {
00240          char *brak = strstr(vars," {");
00241          if (brak) *brak = 0;
00242          fVarX = vars;
00243       }
00244       delete [] vars;
00245    }
00246 }
00247 
00248 
00249 //______________________________________________________________________________
00250 TCutG::~TCutG()
00251 {
00252    // TCutG destructor.
00253 
00254    delete fObjectX;
00255    delete fObjectY;
00256    gROOT->GetListOfSpecials()->Remove(this);
00257 }
00258 
00259 //______________________________________________________________________________
00260 Double_t TCutG::Area() const
00261 {
00262    // Compute the area inside this TCutG
00263    // The algorithm uses Stoke's theorem over the border of the closed polygon.
00264    // Just as a reminder: Stoke's theorem reduces a surface integral
00265    // to a line integral over the border of the surface integral.
00266    Double_t a = 0;
00267    Int_t n = GetN();
00268    for (Int_t i=0;i<n-1;i++) {
00269       a += (fX[i]-fX[i+1])*(fY[i]+fY[i+1]);
00270    }
00271    a *= 0.5;
00272    return a;   
00273 }
00274 
00275 //______________________________________________________________________________
00276 void TCutG::Center(Double_t &cx, Double_t &cy) const
00277 {
00278    // Compute the center x,y of this TCutG
00279    // The algorithm uses Stoke's theorem over the border of the closed polygon.
00280    // Just as a reminder: Stoke's theorem reduces a surface integral
00281    // to a line integral over the border of the surface integral.
00282    Int_t n = GetN();
00283    Double_t a  = 0;
00284    cx = cy = 0;
00285    Double_t t;
00286    for (Int_t i=0;i<n-1;i++) {
00287       t   = 2*fX[i]*fY[i] + fY[i]*fX[i+1] + fX[i]*fY[i+1] + 2*fX[i+1]*fY[i+1];
00288       cx += (fX[i]-fX[i+1])*t;
00289       cy += (-fY[i]+fY[i+1])*t;
00290       a  += (fX[i]-fX[i+1])*(fY[i]+fY[i+1]);
00291    }
00292    a  *= 0.5;
00293    cx *= 1./(6*a);
00294    cy *= 1./(6*a);
00295 }
00296 
00297 //______________________________________________________________________________
00298 Double_t TCutG::IntegralHist(TH2 *h, Option_t *option) const
00299 {
00300    // Compute the integral of 2-d histogram h for all bins inside the cut
00301    // if option "width" is specified, the integral is the sum of
00302    // the bin contents multiplied by the bin width in x and in y.
00303 
00304    if (!h) return 0;
00305    Int_t n = GetN();
00306    Double_t xmin= 1e200;
00307    Double_t xmax = -xmin;
00308    Double_t ymin = xmin;
00309    Double_t ymax = xmax;
00310    for (Int_t i=0;i<n;i++) {
00311       if (fX[i] < xmin) xmin = fX[i];
00312       if (fX[i] > xmax) xmax = fX[i];
00313       if (fY[i] < ymin) ymin = fY[i];
00314       if (fY[i] > ymax) ymax = fY[i];
00315    }
00316    TAxis *xaxis = h->GetXaxis();
00317    TAxis *yaxis = h->GetYaxis();
00318    Int_t binx1 = xaxis->FindBin(xmin);
00319    Int_t binx2 = xaxis->FindBin(xmax);
00320    Int_t biny1 = yaxis->FindBin(ymin);
00321    Int_t biny2 = yaxis->FindBin(ymax);
00322    Int_t nbinsx = h->GetNbinsX();
00323    Stat_t integral = 0;
00324 
00325    // Loop on bins for which the bin center is in the cut
00326    TString opt = option;
00327    opt.ToLower();
00328    Bool_t width = kFALSE;
00329    if (opt.Contains("width")) width = kTRUE;
00330    Int_t bin, binx, biny;
00331    for (biny=biny1;biny<=biny2;biny++) {
00332       Double_t y = yaxis->GetBinCenter(biny);
00333       for (binx=binx1;binx<=binx2;binx++) {
00334          Double_t x = xaxis->GetBinCenter(binx);
00335          if (!IsInside(x,y)) continue;
00336          bin = binx +(nbinsx+2)*biny;
00337          if (width) integral += h->GetBinContent(bin)*xaxis->GetBinWidth(binx)*yaxis->GetBinWidth(biny);
00338          else       integral += h->GetBinContent(bin);
00339       }
00340    }
00341    return integral;
00342 }
00343 
00344 
00345 //______________________________________________________________________________
00346 void TCutG::SavePrimitive(ostream &out, Option_t *option /*= ""*/)
00347 {
00348    // Save primitive as a C++ statement(s) on output stream out.
00349 
00350    char quote = '"';
00351    out<<"   "<<endl;
00352    if (gROOT->ClassSaved(TCutG::Class())) {
00353       out<<"   ";
00354    } else {
00355       out<<"   TCutG *";
00356    }
00357    out<<"cutg = new TCutG("<<quote<<GetName()<<quote<<","<<fNpoints<<");"<<endl;
00358    out<<"   cutg->SetVarX("<<quote<<GetVarX()<<quote<<");"<<endl;
00359    out<<"   cutg->SetVarY("<<quote<<GetVarY()<<quote<<");"<<endl;
00360    out<<"   cutg->SetTitle("<<quote<<GetTitle()<<quote<<");"<<endl;
00361 
00362    SaveFillAttributes(out,"cutg",0,1001);
00363    SaveLineAttributes(out,"cutg",1,1,1);
00364    SaveMarkerAttributes(out,"cutg",1,1,1);
00365 
00366    for (Int_t i=0;i<fNpoints;i++) {
00367       out<<"   cutg->SetPoint("<<i<<","<<fX[i]<<","<<fY[i]<<");"<<endl;
00368    }
00369    out<<"   cutg->Draw("
00370       <<quote<<option<<quote<<");"<<endl;
00371 }
00372 
00373 //______________________________________________________________________________
00374 void TCutG::SetObjectX(TObject *obj) 
00375 {
00376    // Set the X object (and delete the previous one if any).
00377 
00378    delete fObjectX;
00379    fObjectX = obj;
00380 }
00381 
00382 //______________________________________________________________________________
00383 void TCutG::SetObjectY(TObject *obj)
00384 {
00385    // Set the Y object (and delete the previous one if any).
00386    
00387    delete fObjectY;
00388    fObjectY = obj;
00389 }
00390 
00391 //______________________________________________________________________________
00392 void TCutG::SetVarX(const char *varx)
00393 {
00394    // Set X variable.
00395 
00396    fVarX = varx;
00397    delete fObjectX;
00398    fObjectX = 0;
00399 }
00400 
00401 
00402 //______________________________________________________________________________
00403 void TCutG::SetVarY(const char *vary)
00404 {
00405    // Set Y variable.
00406 
00407    fVarY = vary;
00408    delete fObjectY;
00409    fObjectY = 0;
00410 }
00411 
00412 
00413 //______________________________________________________________________________
00414 void TCutG::Streamer(TBuffer &R__b)
00415 {
00416    // Stream an object of class TCutG.
00417 
00418    if (R__b.IsReading()) {
00419       R__b.ReadClassBuffer(TCutG::Class(), this);
00420       gROOT->GetListOfSpecials()->Add(this);
00421    } else {
00422       R__b.WriteClassBuffer(TCutG::Class(), this);
00423    }
00424 }

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