Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
float16.C File Reference

Detailed Description

View in nbviewer Open in SWAN
Tutorial illustrating use and precision of the Float16_t data type.

See the double32.C tutorial for all the details.

#include "ROOT/TSeq.hxx"
#include "TCanvas.h"
#include "TFile.h"
#include "TGraph.h"
#include "TH1.h"
#include "TLegend.h"
#include "TMath.h"
#include "TRandom3.h"
#include "TTree.h"
class DemoFloat16 {
private:
float fF32; // reference member with full single precision
Float16_t fF16; // saved as a 16 bit floating point number
Float16_t fI16; //[-pi,pi] saved as a 16 bit unsigned int
Float16_t fI14; //[-pi,pi,14] saved as a 14 bit unsigned int
Float16_t fI10; //[-pi,pi,10] saved as a 10 bit unsigned int
Float16_t fI8; //[-pi,pi,8] saved as a 8 bit unsigned int
Float16_t fI6; //[-pi,pi,6] saved as a 6 bit unsigned int
Float16_t fI4; //[-pi,pi,4] saved as a 4 bit unsigned int
Float16_t fR8; //[0, 0, 8] saved as a 16 bit float with a 8 bits mantissa
Float16_t fR6; //[0, 0, 6] saved as a 16 bit float with a 6 bits mantissa
Float16_t fR4; //[0, 0, 4] saved as a 16 bit float with a 4 bits mantissa
Float16_t fR2; //[0, 0, 2] saved as a 16 bit float with a 2 bits mantissa
public:
DemoFloat16() = default;
void Set(float ref) { fF32 = fF16 = fI16 = fI14 = fI10 = fI8 = fI6 = fI4 = fR8 = fR6 = fR4 = fR2 = ref; }
};
void float16()
{
const auto nEntries = 200000;
const auto xmax = TMath::Pi();
const auto xmin = -xmax;
// create a Tree with nEntries objects DemoFloat16
TFile::Open("DemoFloat16.root", "recreate");
TTree tree("tree", "DemoFloat16");
DemoFloat16 demoInstance;
auto demoInstanceBranch = tree.Branch("d", "DemoFloat16", &demoInstance, 4000);
for (auto i : ROOT::TSeqI(nEntries)) {
demoInstance.Set(r.Uniform(xmin, xmax));
tree.Fill();
}
tree.Write();
// Now we can proceed with the analysis of the sizes on disk of all branches
// Create the frame histogram and the graphs
auto branches = demoInstanceBranch->GetListOfBranches();
const auto nb = branches->GetEntries();
auto br = static_cast<TBranch *>(branches->At(0));
const Long64_t zip64 = br->GetZipBytes();
auto h = new TH1F("h", "Float16_t compression and precision", nb, 0, nb);
h->SetMaximum(18);
h->SetStats(false);
auto gcx = new TGraph();
gcx->SetName("gcx");
gcx->SetMarkerStyle(kFullSquare);
gcx->SetMarkerColor(kBlue);
auto gdrange = new TGraph();
gdrange->SetName("gdrange");
gdrange->SetMarkerStyle(kFullCircle);
gdrange->SetMarkerColor(kRed);
auto gdval = new TGraph();
gdval->SetName("gdval");
gdval->SetMarkerStyle(kFullTriangleUp);
gdval->SetMarkerColor(kBlack);
// loop on branches to get the precision and compression factors
for (auto i : ROOT::TSeqI(nb)) {
br = static_cast<TBranch *>(branches->At(i));
const auto brName = br->GetName();
h->GetXaxis()->SetBinLabel(i + 1, brName);
auto const cx = double(zip64) / br->GetZipBytes();
gcx->SetPoint(i, i + 0.5, cx);
if (i == 0) continue;
tree.Draw(Form("(fF32-%s)/(%g)", brName, xmax - xmin), "", "goff");
const auto rmsDrange = TMath::RMS(nEntries, tree.GetV1());
const auto drange = TMath::Max(0., -TMath::Log10(rmsDrange));
gdrange->SetPoint(i - 1, i + 0.5, drange);
tree.Draw(Form("(fF32-%s)/fF32 >> hdval_%s", brName, brName), "", "goff");
const auto rmsDval = TMath::RMS(nEntries, tree.GetV1());
const auto dval = TMath::Max(0., -TMath::Log10(rmsDval));
gdval->SetPoint(i - 1, i + 0.5, dval);
tree.Draw(Form("(fF32-%s) >> hdvalabs_%s", brName, brName), "", "goff");
auto hdval = gDirectory->Get<TH1F>(Form("hdvalabs_%s", brName));
hdval->GetXaxis()->SetTitle("Difference wrt reference value");
auto c = new TCanvas(brName, brName, 800, 600);
c->SetGrid();
c->SetLogy();
hdval->DrawClone();
}
auto c1 = new TCanvas("c1", "c1", 800, 600);
c1->SetGrid();
h->Draw();
h->GetXaxis()->LabelsOption("v");
gcx->Draw("lp");
gdrange->Draw("lp");
gdval->Draw("lp");
// Finally build a legend
auto legend = new TLegend(0.3, 0.6, 0.9, 0.9);
legend->SetHeader(Form("%d entries within the [-#pi, #pi] range", nEntries));
legend->AddEntry(gcx, "Compression factor", "lp");
legend->AddEntry(gdrange, "Log of precision wrt range: p = -Log_{10}( RMS( #frac{Ref - x}{range} ) ) ", "lp");
legend->AddEntry(gdval, "Log of precision wrt value: p = -Log_{10}( RMS( #frac{Ref - x}{Ref} ) ) ", "lp");
legend->Draw();
}
#define c(i)
Definition RSha256.hxx:101
#define h(i)
Definition RSha256.hxx:106
float Float16_t
Definition RtypesCore.h:58
long long Long64_t
Definition RtypesCore.h:80
@ kRed
Definition Rtypes.h:66
@ kBlack
Definition Rtypes.h:65
@ kBlue
Definition Rtypes.h:66
@ kFullSquare
Definition TAttMarker.h:55
@ kFullTriangleUp
Definition TAttMarker.h:55
@ kFullCircle
Definition TAttMarker.h:55
#define gDirectory
Definition TDirectory.h:384
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
float xmin
float xmax
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2467
A TTree is a list of TBranches.
Definition TBranch.h:93
Long64_t GetZipBytes(Option_t *option="") const
Return total number of zip bytes in the branch if option ="*" includes all sub-branches of this branc...
Definition TBranch.cxx:2238
The Canvas class.
Definition TCanvas.h:23
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
Definition TFile.cxx:4075
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
1-D histogram with a float per channel (see TH1 documentation)}
Definition TH1.h:577
TAxis * GetXaxis()
Definition TH1.h:322
This class displays a legend box (TPaveText) containing several legend entries.
Definition TLegend.h:23
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:164
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
Random number generator class based on M.
Definition TRandom3.h:27
A TTree represents a columnar dataset.
Definition TTree.h:79
return c1
Definition legend1.C:41
This file contains a specialised ROOT message handler to test for diagnostic in unit tests.
TSeq< int > TSeqI
Definition TSeq.hxx:203
Double_t RMS(Long64_t n, const T *a, const Double_t *w=nullptr)
Returns the Standard Deviation of an array a with length n.
Definition TMath.h:1188
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:250
constexpr Double_t Pi()
Definition TMath.h:37
Double_t Log10(Double_t x)
Returns the common (base-10) logarithm of x.
Definition TMath.h:762
Definition tree.py:1
Author
Danilo Piparo

Definition in file float16.C.