Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
df027_SQliteDependencyOverVersion.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_dataframe
3/// \notebook -js
4/// Plot the ROOT downloads based on the version reading a remote sqlite3 file.
5///
6/// This tutorial uses the Reduce method which allows to extract the minimum time
7/// stored in the SQlite3 database.
8/// The next step is to create a TH1F Histogram, which will be filled with the values stored in
9/// two different columns from the database. This procedure is simplified with a lambda
10/// expression that takes as parameters the values stored in the "Time" and "Version" columns.
11///
12/// \macro_code
13/// \macro_image
14///
15/// \date August 2018
16/// \authors Alexandra-Maria Dobrescu, Sergey Linev
17
18
19void df027_SQliteDependencyOverVersion ()
20{
21 auto rdfb = ROOT::RDF::FromSqlite("http://root.cern/files/root_download_stats.sqlite", "SELECT * FROM accesslog;");
22
23 auto minTimeStr = *rdfb.Reduce([](std::string a, std::string b) {return std::min(a, b);}, "Time", std::string("Z"));
24
25 std::cout << "Minimum time is '" << minTimeStr << "'" << std::endl;
26
27 double minTime = TDatime(minTimeStr.c_str()).Convert();
28 double maxTime = minTime + 3600.*24*365.25*4; // cover approx 4 years from minimal time
29
30 auto rdf = rdfb.Define("datime", [](const std::string &time){return TDatime(time.c_str()).Convert();}, {"Time"});
31
32 auto h614 = rdf.Filter([](const std::string &v){ return 0 == v.find("6.14");}, {"Version"})
33 .Histo1D({"h614", "Download time for version 6.14", 64, minTime, maxTime}, {"datime"});
34
35 auto h616 = rdf.Filter([](const std::string &v){ return 0 == v.find("6.16");}, {"Version"})
36 .Histo1D({"h616", "Download time for version 6.16", 64, minTime, maxTime}, {"datime"});
37
38 auto h618 = rdf.Filter([](const std::string &v){ return 0 == v.find("6.18");}, {"Version"})
39 .Histo1D({"h618", "Download time for version 6.18", 64, minTime, maxTime}, {"datime"});
40
41 auto customize_histo = [](TH1D &histo) {
42 auto *xaxis = histo.GetXaxis();
43 xaxis->SetTimeDisplay(1);
44 xaxis->SetLabelSize(0.02);
45 xaxis->SetNdivisions(512, kFALSE);
46 xaxis->SetTimeFormat("%Y-%m-%d%F1970-00-00 00:00:00");
47 histo.SetStats(kFALSE);
48 };
49
50 customize_histo(*h614);
51 customize_histo(*h616);
52 customize_histo(*h618);
53
54 std::vector drawables{h614->Clone(), h616->Clone(), h618->Clone()};
55
56 auto c1 = new TCanvas("c1","Download time", 800, 1500);
57 c1->Divide(1, drawables.size());
58 for (unsigned n = 0; n < drawables.size(); ++n)
59 c1->GetPad(n + 1)->Add(drawables[n]);
60}
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
The Canvas class.
Definition TCanvas.h:23
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:182
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:670
return c1
Definition legend1.C:41
const Int_t n
Definition legend1.C:16
RDataFrame FromSqlite(std::string_view fileName, std::string_view query)
Factory method to create a SQlite RDataFrame.