Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
df040_RResultPtr_lifetimeManagement.C File Reference

Detailed Description

View in nbviewer Open in SWAN
Usage of RResultPtr: Lifetime management.

This tutorial illustrates how to manage the lifetime of RDataFrame results. When RDataFrame results are declared in functions (or scopes in general), they are destroyed at the end of the scope. To prevent this, one needs to copy the RResultPtr or obtain a copy of its underlying shared_ptr.

#include <TCanvas.h>
#include <THStack.h>
#include <random>
#include <vector>
// Canvas that should survive the running of this macro:
{
// Create a simple dataframe that fills random numbers into histograms.
std::mt19937 generator{1};
std::normal_distribution gaus{5., 1.};
auto rdf = bare_rdf.Define("x", [&]() -> double {
return gaus(generator);
}, {});
ROOT::RDF::TH1DModel histoModel{"Histo", "Histo;x", 10, 0, 10};
// Keeping the results alive is vital when they are passed to other entities or when they are drawn.
// Compare the following situations:
// 1. The wrong way (the result ht is destroyed at the end of the loop body):
THStack histStack1("histStack1", "Stacking result histograms (wrong way)");
for(int i=0; i<2; i++) {
auto ht = rdf.Histo1D(histoModel, {"x"});
ht->SetFillColor(kBlue+i);
histStack1.Add(ht.GetPtr()); // Wrong, this histogram will not survive
}
c1 = new TCanvas("c1", "THStack without obtaining a shared_ptr (wrong)");
histStack1.DrawClone();
c1->Draw();
// 2. The right way: Results survive because we copy the shared_ptr:
THStack histStack2("histStack2", "THStack with shared_ptr (correct way)");
std::vector<std::shared_ptr<TH1D>> results;
for(int i=0; i<2; i++) {
auto ht = rdf.Histo1D(histoModel, {"x"});
ht->SetFillColor(kBlue+2*i);
histStack2.Add(ht.GetPtr());
results.push_back(ht.GetSharedPtr()); // Makes the histogram survive
}
c2 = new TCanvas("c2", "Drawing with obtaining a shared_ptr (right)");
histStack2.DrawClone();
c2->Draw();
}
@ 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.
Date
2025
Author
Stephan Hageboeck (CERN)

Definition in file df040_RResultPtr_lifetimeManagement.C.