Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
df029_SQlitePlatformDistribution.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_dataframe
3/// \notebook -js
4/// Use RDataFrame to display data about ROOT downloads.
5///
6/// In order to display the Platform Distribution of ROOT, we choose to create two TH1F
7/// histograms: one that includes all types of platforms, other filtering and classifying them.
8/// This procedure is using a lambda expression taking as parameter the values
9/// stored in the "Platform" column from the database. At the end, the histograms are filled
10/// with their specific demand regarding the platform's type.
11///
12/// \macro_code
13/// \macro_image
14///
15/// \date August 2018
16/// \author Alexandra-Maria Dobrescu
17
18void df029_SQlitePlatformDistribution() {
19
20 auto rdf = ROOT::RDF::FromSqlite("http://root.cern/files/root_download_stats.sqlite", "SELECT * FROM accesslog;");
21
22 TH1F hRootPlatform("hrootPlatform", "Platform Distribution", 7, 0, -1);
23 TH1F hShortRootPlatform("hShortRootPlatform", "Short Platform Distribution", 7, 0, -1);
24
25 auto fillRootPlatform = [&hRootPlatform, &hShortRootPlatform] ( const std::string &platform ) {
26 TString Platform = platform;
27 TString Platform_0(Platform(0,5));
28 TString Platform_1(Platform(0,6));
29 TString Platform_2(Platform(0,8));
30
31 if ( Platform.Contains("win32") ){
32 hShortRootPlatform.Fill(Platform_0,1);
33 } else if ( Platform.Contains("Linux") ){
34 hShortRootPlatform.Fill(Platform_0,1);
35 } else if ( Platform.Contains("source") ){
36 hShortRootPlatform.Fill(Platform_1,1);
37 } else if ( Platform.Contains("macosx64") ){
38 hShortRootPlatform.Fill(Platform_2,1);
39 } else if ( Platform.Contains("IRIX64") ){
40 hShortRootPlatform.Fill(Platform_1,1);
41 }
42
43 hRootPlatform.Fill(Platform,1);
44 };
45
46 rdf.Foreach( fillRootPlatform, { "Platform" } );
47
48 auto PlatformDistributionHistogram = new TCanvas();
49
50 hRootPlatform.GetXaxis()->LabelsOption("a");
51 hRootPlatform.LabelsDeflate("X");
52 hRootPlatform.DrawClone();
53
54 auto shortPlatformDistributionHistogram = new TCanvas();
55
56 hShortRootPlatform.GetXaxis()->LabelsOption("a");
57 hShortRootPlatform.LabelsDeflate("X");
58 hShortRootPlatform.DrawClone();
59}
The Canvas class.
Definition TCanvas.h:23
1-D histogram with a float per channel (see TH1 documentation)}
Definition TH1.h:577
Basic string class.
Definition TString.h:139
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:636
RDataFrame FromSqlite(std::string_view fileName, std::string_view query)
Factory method to create a SQlite RDataFrame.