Logo ROOT  
Reference Guide
df027_SQliteDependencyOverVersion.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_dataframe
3/// \notebook -jss
4/// Plot the ROOT downloads based on the version reading a remote sqlite3 file with RSqliteDS.
5/// This tutorial uses the Reduce method which allows to extract the minimum time
6/// stored in the SQlite3 database.
7/// The next step is to create a TH1F Histogram, which will be filled with the values stored in
8/// two different columns from the database. This procedure is simplified with a lambda
9/// expression that takes as parameters the values stored in the "Time" and "Version" columns.
10///
11/// \macro_code
12/// \macro_image
13///
14/// \date August 2018
15/// \author Alexandra-Maria Dobrescu
16
17void df027_SQliteDependencyOverVersion () {
18
19 auto rdfb = ROOT::RDF::MakeSqliteDataFrame("http://root.cern/files/root_download_stats.sqlite", "SELECT * FROM accesslog;" );
20
21 auto minTimeStr = *rdfb.Reduce([](std::string a, std::string b) {return std::min(a, b);}, "Time", std::string("Z"));
22
23 std::cout << "Minimum time is '" << minTimeStr << "'" << std::endl;
24
25 double minTime = TDatime(minTimeStr.c_str()).Convert();
26 double now = TDatime().Convert();
27
28 auto rdf = rdfb.Define("datime", [](const std::string &time){return TDatime(time.c_str()).Convert();}, {"Time"});
29
30 auto h614 = rdf.Filter([](const std::string &v){ return 0 == v.find("6.14");}, {"Version"})
31 .Histo1D({"h614", "Download time for version 6.14", 16, minTime, now}, {"datime"});
32
33 auto h616 = rdf.Filter([](const std::string &v){ return 0 == v.find("6.16");}, {"Version"})
34 .Histo1D({"h616", "Download time for version 6.16", 16, minTime, now}, {"datime"});
35
36 auto h618 = rdf.Filter([](const std::string &v){ return 0 == v.find("6.18");}, {"Version"})
37 .Histo1D({"h618", "Download time for version 6.18", 16, minTime, now}, {"datime"});
38
39 // Add here a newer version!
40
41 auto histoList = {h614, h616, h618};
42 auto canvases = new std::vector<TCanvas*>(histoList.size());
43
46 auto histoIdx = 0U;
47 for (auto histo : histoList) {
48 canvases->at(histoIdx) = new TCanvas();
49 histo->GetXaxis()->LabelsOption("v");
50 histo->GetXaxis()->SetTimeDisplay(1);
51 histo->GetXaxis()->SetLabelSize(0.02);
52 histo->GetXaxis()->SetNdivisions(512, kFALSE);
53 histo->GetXaxis()->SetTimeFormat("%Y-%m-%d");
54 histo->DrawClone();
55 histoIdx++;
56 }
57
58}
#define b(i)
Definition: RSha256.hxx:100
const Bool_t kFALSE
Definition: RtypesCore.h:90
R__EXTERN TStyle * gStyle
Definition: TStyle.h:410
The Canvas class.
Definition: TCanvas.h:27
This class stores the date and time with a precision of one second in an unsigned 32 bit word (950130...
Definition: TDatime.h:37
UInt_t Convert(Bool_t toGMT=kFALSE) const
Convert fDatime from TDatime format to the standard time_t format.
Definition: TDatime.cxx:181
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
Definition: TStyle.cxx:1590
void SetTimeOffset(Double_t toffset)
Change the time offset for time plotting.
Definition: TStyle.cxx:1801
RDataFrame MakeSqliteDataFrame(std::string_view fileName, std::string_view query)
Factory method to create a SQlite RDataFrame.
Definition: RSqliteDS.cxx:541
auto * a
Definition: textangle.C:12