Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
float16.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_io
3 /// \notebook -js
4 /// Tutorial illustrating use and precision of the Float16_t data type.
5 /// See the double32.C tutorial for all the details.
6 /// \macro_image
7 /// \macro_code
8 ///
9 /// \author Danilo Piparo
10 
11 #include "ROOT/TSeq.hxx"
12 #include "TCanvas.h"
13 #include "TFile.h"
14 #include "TGraph.h"
15 #include "TH1.h"
16 #include "TLegend.h"
17 #include "TMath.h"
18 #include "TRandom3.h"
19 #include "TTree.h"
20 
21 class DemoFloat16 {
22 private:
23  float fF32; // reference member with full single precision
24  Float16_t fF16; // saved as a 16 bit floating point number
25  Float16_t fI16; //[-pi,pi] saved as a 16 bit unsigned int
26  Float16_t fI14; //[-pi,pi,14] saved as a 14 bit unsigned int
27  Float16_t fI10; //[-pi,pi,10] saved as a 10 bit unsigned int
28  Float16_t fI8; //[-pi,pi,8] saved as a 8 bit unsigned int
29  Float16_t fI6; //[-pi,pi,6] saved as a 6 bit unsigned int
30  Float16_t fI4; //[-pi,pi,4] saved as a 4 bit unsigned int
31  Float16_t fR8; //[0, 0, 8] saved as a 16 bit float with a 8 bits mantissa
32  Float16_t fR6; //[0, 0, 6] saved as a 16 bit float with a 6 bits mantissa
33  Float16_t fR4; //[0, 0, 4] saved as a 16 bit float with a 4 bits mantissa
34  Float16_t fR2; //[0, 0, 2] saved as a 16 bit float with a 2 bits mantissa
35 
36 public:
37  DemoFloat16() = default;
38  void Set(float ref) { fF32 = fF16 = fI16 = fI14 = fI10 = fI8 = fI6 = fI4 = fR8 = fR6 = fR4 = fR2 = ref; }
39 };
40 
41 void float16()
42 {
43  const auto nEntries = 200000;
44  const auto xmax = TMath::Pi();
45  const auto xmin = -xmax;
46 
47  // create a Tree with nEntries objects DemoFloat16
48  TFile::Open("DemoFloat16.root", "recreate");
49  TTree tree("tree", "DemoFloat16");
50  DemoFloat16 demoInstance;
51  auto demoInstanceBranch = tree.Branch("d", "DemoFloat16", &demoInstance, 4000);
52  TRandom3 r;
53  for (auto i : ROOT::TSeqI(nEntries)) {
54  demoInstance.Set(r.Uniform(xmin, xmax));
55  tree.Fill();
56  }
57  tree.Write();
58 
59  // Now we can proceed with the analysis of the sizes on disk of all branches
60 
61  // Create the frame histogram and the graphs
62  auto branches = demoInstanceBranch->GetListOfBranches();
63  const auto nb = branches->GetEntries();
64  auto br = static_cast<TBranch *>(branches->At(0));
65  const Long64_t zip64 = br->GetZipBytes();
66 
67  auto h = new TH1F("h", "Float16_t compression and precision", nb, 0, nb);
68  h->SetMaximum(18);
69  h->SetStats(0);
70 
71  auto gcx = new TGraph();
72  gcx->SetName("gcx");
73  gcx->SetMarkerStyle(kFullSquare);
74  gcx->SetMarkerColor(kBlue);
75 
76  auto gdrange = new TGraph();
77  gdrange->SetName("gdrange");
78  gdrange->SetMarkerStyle(kFullCircle);
79  gdrange->SetMarkerColor(kRed);
80 
81  auto gdval = new TGraph();
82  gdval->SetName("gdval");
83  gdval->SetMarkerStyle(kFullTriangleUp);
84  gdval->SetMarkerColor(kBlack);
85 
86  // loop on branches to get the precision and compression factors
87  for (auto i : ROOT::TSeqI(nb)) {
88  br = static_cast<TBranch *>(branches->At(i));
89  const auto brName = br->GetName();
90  h->GetXaxis()->SetBinLabel(i + 1, brName);
91  auto const cx = double(zip64) / br->GetZipBytes();
92  gcx->SetPoint(i, i + 0.5, cx);
93  if (i == 0) continue;
94 
95  tree.Draw(Form("(fF32-%s)/(%g)", brName, xmax - xmin), "", "goff");
96  const auto rmsDrange = TMath::RMS(nEntries, tree.GetV1());
97  const auto drange = TMath::Max(0., -TMath::Log10(rmsDrange));
98  gdrange->SetPoint(i - 1, i + 0.5, drange);
99 
100  tree.Draw(Form("(fF32-%s)/fF32 >> hdval_%s", brName, brName), "", "goff");
101  const auto rmsDval = TMath::RMS(nEntries, tree.GetV1());
102  const auto dval = TMath::Max(0., -TMath::Log10(rmsDval));
103  gdval->SetPoint(i - 1, i + 0.5, dval);
104 
105  tree.Draw(Form("(fF32-%s) >> hdvalabs_%s", brName, brName), "", "goff");
106  auto hdval = gDirectory->Get<TH1F>(Form("hdvalabs_%s", brName));
107  hdval->GetXaxis()->SetTitle("Difference wrt reference value");
108  auto c = new TCanvas(brName, brName, 800, 600);
109  c->SetGrid();
110  c->SetLogy();
111  hdval->DrawClone();
112  }
113 
114  auto c1 = new TCanvas("c1", "c1", 800, 600);
115  c1->SetGrid();
116 
117  h->Draw();
118  h->GetXaxis()->LabelsOption("v");
119  gcx->Draw("lp");
120  gdrange->Draw("lp");
121  gdval->Draw("lp");
122 
123  // Finally build a legend
124  auto legend = new TLegend(0.3, 0.6, 0.9, 0.9);
125  legend->SetHeader(Form("%d entries within the [-#pi, #pi] range", nEntries));
126  legend->AddEntry(gcx, "Compression factor", "lp");
127  legend->AddEntry(gdrange, "Log of precision wrt range: p = -Log_{10}( RMS( #frac{Ref - x}{range} ) ) ", "lp");
128  legend->AddEntry(gdval, "Log of precision wrt value: p = -Log_{10}( RMS( #frac{Ref - x}{Ref} ) ) ", "lp");
129  legend->Draw();
130 }