Logo ROOT  
Reference Guide
df027_SQliteDependencyOverVersion.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_dataframe
3 /// \notebook -jss
4 /// \brief 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 /// \author Alexandra-Maria Dobrescu
17 
18 void df027_SQliteDependencyOverVersion () {
19 
20  auto rdfb = ROOT::RDF::MakeSqliteDataFrame("http://root.cern/files/root_download_stats.sqlite", "SELECT * FROM accesslog;" );
21 
22  auto minTimeStr = *rdfb.Reduce([](std::string a, std::string b) {return std::min(a, b);}, "Time", std::string("Z"));
23 
24  std::cout << "Minimum time is '" << minTimeStr << "'" << std::endl;
25 
26  double minTime = TDatime(minTimeStr.c_str()).Convert();
27  double now = TDatime().Convert();
28 
29  auto rdf = rdfb.Define("datime", [](const std::string &time){return TDatime(time.c_str()).Convert();}, {"Time"});
30 
31  auto h614 = rdf.Filter([](const std::string &v){ return 0 == v.find("6.14");}, {"Version"})
32  .Histo1D({"h614", "Download time for version 6.14", 16, minTime, now}, {"datime"});
33 
34  auto h616 = rdf.Filter([](const std::string &v){ return 0 == v.find("6.16");}, {"Version"})
35  .Histo1D({"h616", "Download time for version 6.16", 16, minTime, now}, {"datime"});
36 
37  auto h618 = rdf.Filter([](const std::string &v){ return 0 == v.find("6.18");}, {"Version"})
38  .Histo1D({"h618", "Download time for version 6.18", 16, minTime, now}, {"datime"});
39 
40  // Add here a newer version!
41 
42  auto histoList = {h614, h616, h618};
43  auto canvases = new std::vector<TCanvas*>(histoList.size());
44 
46  gStyle->SetOptStat(0);
47  auto histoIdx = 0U;
48  for (auto histo : histoList) {
49  canvases->at(histoIdx) = new TCanvas();
50  histo->GetXaxis()->LabelsOption("v");
51  histo->GetXaxis()->SetTimeDisplay(1);
52  histo->GetXaxis()->SetLabelSize(0.02);
53  histo->GetXaxis()->SetNdivisions(512, kFALSE);
54  histo->GetXaxis()->SetTimeFormat("%Y-%m-%d");
55  histo->DrawClone();
56  histoIdx++;
57  }
58 
59 }
ROOT::RDF::MakeSqliteDataFrame
RDataFrame MakeSqliteDataFrame(std::string_view fileName, std::string_view query)
Factory method to create a SQlite RDataFrame.
Definition: RSqliteDS.cxx:545
v
@ v
Definition: rootcling_impl.cxx:3635
b
#define b(i)
Definition: RSha256.hxx:118
gStyle
R__EXTERN TStyle * gStyle
Definition: TStyle.h:412
a
auto * a
Definition: textangle.C:12
kFALSE
const Bool_t kFALSE
Definition: RtypesCore.h:92
TStyle::SetOptStat
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:1592
TStyle::SetTimeOffset
void SetTimeOffset(Double_t toffset)
Change the time offset for time plotting.
Definition: TStyle.cxx:1803
TCanvas
Definition: TCanvas.h:23
TDatime
Definition: TDatime.h:37
TDatime::Convert
UInt_t Convert(Bool_t toGMT=kFALSE) const
Convert fDatime from TDatime format to the standard time_t format.
Definition: TDatime.cxx:182