Logo ROOT   6.12/07
Reference Guide
tdf003_profiles.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_tdataframe
3 /// \notebook -nodraw
4 /// This tutorial illustrates how to use TProfiles in combination with the
5 /// TDataFrame. See the documentation of TProfile and TProfile2D to better
6 /// understand the analogy of this code with the example one.
7 ///
8 /// \macro_code
9 ///
10 /// \date February 2017
11 /// \author Danilo Piparo
12 
13 // A simple helper function to fill a test tree: this makes the example
14 // stand-alone.
15 void fill_tree(const char *treeName, const char *fileName)
16 {
18  d.Define("px", []() { return gRandom->Gaus(); })
19  .Define("py", []() { return gRandom->Gaus(); })
20  .Define("pz", [](double px, double py) { return sqrt(px * px + py * py); }, {"px", "py"})
21  .Snapshot(treeName, fileName);
22 }
23 
24 void tdf003_profiles()
25 {
26  // We prepare an input tree to run on
27  auto fileName = "tdf003_profiles.root";
28  auto treeName = "myTree";
29  fill_tree(treeName, fileName);
30 
31  // We read the tree from the file and create a TDataFrame.
32  ROOT::Experimental::TDataFrame d(treeName, fileName, {"px", "py", "pz"});
33 
34  // Create the profiles
35  auto hprof1d = d.Profile1D({"hprof1d", "Profile of pz versus px", 64, -4, 4});
36  auto hprof2d = d.Profile2D({"hprof2d", "Profile of pz versus px and py", 40, -4, 4, 40, -4, 4, 0, 20});
37 
38  // And Draw
39  auto c1 = new TCanvas("c1", "Profile histogram example", 200, 10, 700, 500);
40  hprof1d->DrawClone();
41  auto c2 = new TCanvas("c2", "Profile2D histogram example", 200, 10, 700, 500);
42  c2->cd();
43  hprof2d->DrawClone();
44 }
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
Definition: TRandom.cxx:256
return c1
Definition: legend1.C:41
double sqrt(double)
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
The Canvas class.
Definition: TCanvas.h:31
TResultProxy<::TProfile > Profile1D(const TProfile1DModel &model, std::string_view v1Name="", std::string_view v2Name="")
Fill and return a one-dimensional profile (lazy action)
return c2
Definition: legend2.C:14
ROOT&#39;s TDataFrame offers a high level interface for analyses of data stored in TTrees.
Definition: TDataFrame.hxx:39