Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
df040_RResultPtr_lifetimeManagement.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_dataframe
3/// \notebook -nodraw
4/// Usage of RResultPtr: Lifetime management.
5///
6/// This tutorial illustrates how to manage the lifetime of RDataFrame results.
7/// When RDataFrame results are declared in functions (or scopes in general),
8/// they are destroyed at the end of the scope.
9/// To prevent this, one needs to copy the RResultPtr or obtain a copy of its
10/// underlying shared_ptr.
11///
12/// \macro_code
13/// \macro_output
14///
15/// \date 2025
16/// \author Stephan Hageboeck (CERN)
17
18#include <ROOT/RDataFrame.hxx>
19#include <ROOT/RDFHelpers.hxx>
20#include <TCanvas.h>
21#include <THStack.h>
22#include <random>
23#include <vector>
24
25// Canvas that should survive the running of this macro:
26TCanvas *c1, *c2;
27
29{
30 // Create a simple dataframe that fills random numbers into histograms.
32 std::mt19937 generator{1};
33 std::normal_distribution gaus{5., 1.};
34 auto rdf = bare_rdf.Define("x", [&]() -> double {
35 return gaus(generator);
36 }, {});
37
38 ROOT::RDF::TH1DModel histoModel{"Histo", "Histo;x", 10, 0, 10};
39
40 // Keeping the results alive is vital when they are passed to other entities or when they are drawn.
41 // Compare the following situations:
42
43 // 1. The wrong way (the result ht is destroyed at the end of the loop body):
44 THStack histStack1("histStack1", "Stacking result histograms (wrong way)");
45 for(int i=0; i<2; i++) {
46 auto ht = rdf.Histo1D(histoModel, {"x"});
47 ht->SetFillColor(kBlue+i);
48 histStack1.Add(ht.GetPtr()); // Wrong, this histogram will not survive
49 }
50 c1 = new TCanvas("c1", "THStack without obtaining a shared_ptr (wrong)");
51 histStack1.DrawClone();
52 c1->Draw();
53
54 // 2. The right way: Results survive because we copy the shared_ptr:
55 THStack histStack2("histStack2", "THStack with shared_ptr (correct way)");
56 std::vector<std::shared_ptr<TH1D>> results;
57 for(int i=0; i<2; i++) {
58 auto ht = rdf.Histo1D(histoModel, {"x"});
59 ht->SetFillColor(kBlue+2*i);
60 histStack2.Add(ht.GetPtr());
61 results.push_back(ht.GetSharedPtr()); // Makes the histogram survive
62 }
63 c2 = new TCanvas("c2", "Drawing with obtaining a shared_ptr (right)");
64 histStack2.DrawClone();
65 c2->Draw();
66}
@ kBlue
Definition Rtypes.h:66
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree ,...
The Canvas class.
Definition TCanvas.h:23
The Histogram stack class.
Definition THStack.h:40
return c1
Definition legend1.C:41
return c2
Definition legend2.C:14
A struct which stores the parameters of a TH1D.