Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
double32.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 Double32_t data type
5 /// You should run this tutorial with ACLIC: a dictionary will be automatically
6 /// created.
7 /// ~~~{.bash}
8 /// root > .x double32.C+
9 /// ~~~
10 /// The following cases are supported for streaming a Double32_t type
11 /// depending on the range declaration in the comment field of the data member:
12 ///
13 /// Case | Declaration
14 /// -----|------------
15 /// A | Double32_t fNormal;
16 /// B | Double32_t fTemperature; //[0,100]
17 /// C | Double32_t fCharge; //[-1,1,2]
18 /// D | Double32_t fVertex[3]; //[-30,30,10]
19 /// E | Double32_t fChi2; //[0,0,6]
20 /// F | Int_t fNsp;<br>Double32_t* fPointValue; //[fNsp][0,3]
21 ///
22 /// * Case A fNormal is converted from a Double_t to a Float_t
23 /// * Case B fTemperature is converted to a 32 bit unsigned integer
24 /// * Case C fCharge is converted to a 2 bits unsigned integer
25 /// * Case D the array elements of fVertex are converted to an unsigned 10 bits integer
26 /// * Case E fChi2 is converted to a Float_t with truncated precision at 6 bits
27 /// * Case F the fNsp elements of array fPointvalue are converted to an unsigned 32 bit integer. Note that the range specifier must follow the dimension specifier.
28 ///
29 /// Case B has more precision than case A: 9 to 10 significative digits and 6 to 7 digits respectively.
30 /// The range specifier has the general format: [xmin,xmax] or [xmin,xmax,nbits]. Examples
31 /// * [0,1]
32 /// * [-10,100];
33 /// * [-pi,pi], [-pi/2,pi/4],[-2pi,2*pi]
34 /// * [-10,100,16]
35 /// * [0,0,8]
36 /// Note that:
37 /// * If nbits is not specified, or nbits <2 or nbits>32 it is set to 32
38 /// * If (xmin==0 and xmax==0 and nbits <=14) the double word will be converted to a float and its mantissa truncated to nbits significative bits.
39 ///
40 /// ## IMPORTANT NOTE
41 /// Lets assume an original variable double x.
42 /// When using the format [0,0,8] (i.e. range not specified) you get the best
43 /// relative precision when storing and reading back the truncated x, say xt.
44 /// The variance of (x-xt)/x will be better than when specifying a range
45 /// for the same number of bits. However the precision relative to the
46 /// range (x-xt)/(xmax-xmin) will be worse, and vice-versa.
47 /// The format [0,0,8] is also interesting when the range of x is infinite
48 /// or unknown.
49 ///
50 /// \macro_image
51 /// \macro_code
52 ///
53 /// \author Rene Brun
54 
55 #include "ROOT/TSeq.hxx"
56 #include "TCanvas.h"
57 #include "TFile.h"
58 #include "TGraph.h"
59 #include "TH1.h"
60 #include "TLegend.h"
61 #include "TMath.h"
62 #include "TRandom3.h"
63 #include "TTree.h"
64 
65 class DemoDouble32 {
66 private:
67  Double_t fD64; // reference member with full double precision
68  Double32_t fF32; // saved as a 32 bit Float_t
69  Double32_t fI32; //[-pi,pi] saved as a 32 bit unsigned int
70  Double32_t fI30; //[-pi,pi,30] saved as a 30 bit unsigned int
71  Double32_t fI28; //[-pi,pi,28] saved as a 28 bit unsigned int
72  Double32_t fI26; //[-pi,pi,26] saved as a 26 bit unsigned int
73  Double32_t fI24; //[-pi,pi,24] saved as a 24 bit unsigned int
74  Double32_t fI22; //[-pi,pi,22] saved as a 22 bit unsigned int
75  Double32_t fI20; //[-pi,pi,20] saved as a 20 bit unsigned int
76  Double32_t fI18; //[-pi,pi,18] saved as a 18 bit unsigned int
77  Double32_t fI16; //[-pi,pi,16] saved as a 16 bit unsigned int
78  Double32_t fI14; //[-pi,pi,14] saved as a 14 bit unsigned int
79  Double32_t fI12; //[-pi,pi,12] saved as a 12 bit unsigned int
80  Double32_t fI10; //[-pi,pi,10] saved as a 10 bit unsigned int
81  Double32_t fI8; //[-pi,pi, 8] saved as a 8 bit unsigned int
82  Double32_t fI6; //[-pi,pi, 6] saved as a 6 bit unsigned int
83  Double32_t fI4; //[-pi,pi, 4] saved as a 4 bit unsigned int
84  Double32_t fI2; //[-pi,pi, 2] saved as a 2 bit unsigned int
85  Double32_t fR14; //[0, 0, 14] saved as a 32 bit float with a 14 bits mantissa
86  Double32_t fR12; //[0, 0, 12] saved as a 32 bit float with a 12 bits mantissa
87  Double32_t fR10; //[0, 0, 10] saved as a 32 bit float with a 10 bits mantissa
88  Double32_t fR8; //[0, 0, 8] saved as a 32 bit float with a 8 bits mantissa
89  Double32_t fR6; //[0, 0, 6] saved as a 32 bit float with a 6 bits mantissa
90  Double32_t fR4; //[0, 0, 4] saved as a 32 bit float with a 4 bits mantissa
91  Double32_t fR2; //[0, 0, 2] saved as a 32 bit float with a 2 bits mantissa
92 
93 public:
94  DemoDouble32() = default;
95  void Set(Double_t ref)
96  {
97  fD64 = fF32 = fI32 = fI30 = fI28 = fI26 = fI24 = fI22 = fI20 = fI18 = fI16 = fI14 = fI12 = fI10 = fI8 = fI6 =
98  fI4 = fI2 = fR14 = fR12 = fR10 = fR8 = fR6 = fR4 = fR2 = ref;
99  }
100 };
101 
102 void double32()
103 {
104  const auto nEntries = 40000;
105  const auto xmax = TMath::Pi();
106  const auto xmin = -xmax;
107 
108  // create a Tree with nEntries objects DemoDouble32
109  TFile::Open("DemoDouble32.root", "recreate");
110  TTree tree("tree", "DemoDouble32");
111  DemoDouble32 demoInstance;
112  auto demoInstanceBranch = tree.Branch("d", "DemoDouble32", &demoInstance, 4000);
113  TRandom3 r;
114  for (auto i : ROOT::TSeqI(nEntries)) {
115  demoInstance.Set(r.Uniform(xmin, xmax));
116  tree.Fill();
117  }
118  tree.Write();
119 
120  // Now we can proceed with the analysis of the sizes on disk of all branches
121 
122  // Create the frame histogram and the graphs
123  auto branches = demoInstanceBranch->GetListOfBranches();
124  const auto nb = branches->GetEntries();
125  auto br = static_cast<TBranch *>(branches->At(0));
126  const Long64_t zip64 = br->GetZipBytes();
127 
128  auto h = new TH1F("h", "Double32_t compression and precision", nb, 0, nb);
129  h->SetMaximum(18);
130  h->SetStats(0);
131 
132  auto gcx = new TGraph();
133  gcx->SetName("gcx");
134  gcx->SetMarkerStyle(kFullSquare);
135  gcx->SetMarkerColor(kBlue);
136 
137  auto gdrange = new TGraph();
138  gdrange->SetName("gdrange");
139  gdrange->SetMarkerStyle(kFullCircle);
140  gdrange->SetMarkerColor(kRed);
141 
142  auto gdval = new TGraph();
143  gdval->SetName("gdval");
144  gdval->SetMarkerStyle(kFullTriangleUp);
145  gdval->SetMarkerColor(kBlack);
146 
147  // loop on branches to get the precision and compression factors
148  for (auto i : ROOT::TSeqI(nb)) {
149  auto br = static_cast<TBranch *>(branches->At(i));
150  const auto brName = br->GetName();
151 
152  h->GetXaxis()->SetBinLabel(i + 1, brName);
153  const auto cx = double(zip64) / br->GetZipBytes();
154  gcx->SetPoint(i, i + 0.5, cx);
155  if (i == 0 ) continue;
156 
157  tree.Draw(Form("(fD64-%s)/(%g)", brName, xmax - xmin), "", "goff");
158  const auto rmsDrange = TMath::RMS(nEntries, tree.GetV1());
159  const auto drange = TMath::Max(0., -TMath::Log10(rmsDrange));
160  gdrange->SetPoint(i-1, i + 0.5, drange);
161 
162  tree.Draw(Form("(fD64-%s)/fD64", brName), "", "goff");
163  const auto rmsDVal = TMath::RMS(nEntries, tree.GetV1());
164  const auto dval = TMath::Max(0., -TMath::Log10(rmsDVal));
165  gdval->SetPoint(i-1, i + 0.5, dval);
166 
167  tree.Draw(Form("(fD64-%s) >> hdvalabs_%s", brName, brName), "", "goff");
168  auto hdval = gDirectory->Get<TH1F>(Form("hdvalabs_%s", brName));
169  hdval->GetXaxis()->SetTitle("Difference wrt reference value");
170  auto c = new TCanvas(brName, brName, 800, 600);
171  c->SetGrid();
172  c->SetLogy();
173  hdval->DrawClone();
174  }
175 
176  auto c1 = new TCanvas("c1", "c1", 800, 600);
177  c1->SetGrid();
178 
179  h->Draw();
180  h->GetXaxis()->LabelsOption("v");
181  gcx->Draw("lp");
182  gdrange->Draw("lp");
183  gdval->Draw("lp");
184 
185  // Finally build a legend
186  auto legend = new TLegend(0.3, 0.6, 0.9, 0.9);
187  legend->SetHeader(Form("%d entries within the [-#pi, #pi] range", nEntries));
188  legend->AddEntry(gcx, "Compression factor", "lp");
189  legend->AddEntry(gdrange, "Log of precision wrt range: p = -Log_{10}( RMS( #frac{Ref - x}{range} ) ) ", "lp");
190  legend->AddEntry(gdval, "Log of precision wrt value: p = -Log_{10}( RMS( #frac{Ref - x}{Ref} ) ) ", "lp");
191  legend->Draw();
192 }