00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052 #include "TFile.h"
00053 #include "TCanvas.h"
00054 #include "TTree.h"
00055 #include "TH1.h"
00056 #include "TMath.h"
00057 #include "TRandom3.h"
00058 #include "TGraph.h"
00059 #include "TLegend.h"
00060 #include "TFrame.h"
00061 #include "TPaveLabel.h"
00062
00063 class DemoDouble32 {
00064 private:
00065 Double_t fD64;
00066 Double32_t fF32;
00067 Double32_t fI32;
00068 Double32_t fI30;
00069 Double32_t fI28;
00070 Double32_t fI26;
00071 Double32_t fI24;
00072 Double32_t fI22;
00073 Double32_t fI20;
00074 Double32_t fI18;
00075 Double32_t fI16;
00076 Double32_t fI14;
00077 Double32_t fI12;
00078 Double32_t fI10;
00079 Double32_t fI8;
00080 Double32_t fI6;
00081 Double32_t fI4;
00082 Double32_t fI2;
00083 Double32_t fR14;
00084 Double32_t fR12;
00085 Double32_t fR10;
00086 Double32_t fR8;
00087 Double32_t fR6;
00088 Double32_t fR4;
00089 Double32_t fR2;
00090
00091 public:
00092 DemoDouble32() {;}
00093 void Set(Double_t ref);
00094 };
00095
00096 void DemoDouble32::Set(Double_t ref) {
00097 fD64 = fF32 = fI32 = fI30 = fI28 = fI26 = fI24 = fI22 = fI20 = ref;
00098 fI18 = fI16 = fI14 = fI12 = fI10 = fI8 = fI6 = fI4 = fI2 = ref;
00099 fR14 = fR12 = fR10 = fR8 = fR6 = fR4 = fR2 = ref;
00100 }
00101
00102 void double32() {
00103
00104
00105 DemoDouble32 *d = new DemoDouble32();
00106
00107
00108 TFile::Open("DemoDouble32.root","recreate");
00109 TTree *T = new TTree("T","DemoDouble32");
00110 TBranch *bd = T->Branch("d","DemoDouble32",&d,4000);
00111 TRandom3 r;
00112 Double_t xmax = TMath::Pi();
00113 Double_t xmin = -xmax;
00114 Int_t i, n = 40000;
00115 for (i=0;i<n;i++) {
00116 d->Set(r.Uniform(xmin,xmax));
00117 T->Fill();
00118 }
00119 T->Write();
00120
00121
00122 TObjArray *branches = bd->GetListOfBranches();
00123 Int_t nb = branches->GetEntries();
00124 TBranch *br = (TBranch*)branches->At(0);
00125 Long64_t zip64 = br->GetZipBytes();
00126 Double_t cx = 1;
00127 Double_t drange = 15;
00128 Double_t dval = 15;
00129 TCanvas *c1 = new TCanvas("c1","c1",800,600);
00130 c1->SetGrid();
00131 c1->SetHighLightColor(0);
00132 c1->SetFillColor(17);
00133 c1->SetFrameFillColor(20);
00134 c1->SetFrameBorderSize(10);
00135 TH1F *h = new TH1F("h","",nb,0,nb);
00136 h->SetMaximum(16);
00137 h->SetStats(0);
00138 h->Draw();
00139 c1->GetFrame()->SetFillColor(21);
00140 c1->GetFrame()->SetBorderSize(12);
00141 TGraph *gcx = new TGraph(nb); gcx->SetName("gcx");
00142 gcx->SetMarkerStyle(21);
00143 gcx->SetMarkerColor(kBlue);
00144 TGraph *gdrange = new TGraph(nb); gdrange->SetName("gdrange");
00145 gdrange->SetMarkerStyle(20);
00146 gdrange->SetMarkerColor(kRed);
00147 TGraph *gdval = new TGraph(nb); gdval->SetName("gdval");
00148 gdval->SetMarkerStyle(20);
00149 gdval->SetMarkerColor(kBlack);
00150 TPaveLabel *title = new TPaveLabel(.15,.92,.85,.97,"Double32_t compression and precision","brNDC");
00151 title->Draw();
00152
00153
00154 for (i=0;i<nb;i++) {
00155 br = (TBranch*)branches->At(i);
00156 h->GetXaxis()->SetBinLabel(i+1,br->GetName());
00157 cx = Double_t(zip64)/Double_t(br->GetZipBytes());
00158 gcx->SetPoint(i,i+0.5,cx);
00159 if (i > 0) {
00160 T->Draw(Form("(fD64-%s)/(%g)",br->GetName(),xmax-xmin),"","goff");
00161 Double_t rms = TMath::RMS(n,T->GetV1());
00162 drange = TMath::Max(0.,-TMath::Log10(rms));
00163 }
00164 gdrange->SetPoint(i,i+0.5,drange);
00165 if (i > 0) {
00166 T->Draw(Form("(fD64-%s)/(fD64+0.01)",br->GetName()),"","goff");
00167 Double_t rms = TMath::RMS(n,T->GetV1());
00168 dval = TMath::Max(0.,-TMath::Log10(rms));
00169 }
00170 gdval->SetPoint(i,i+0.5,dval);
00171 }
00172 gcx->Draw("lp");
00173 gdrange->Draw("lp");
00174 gdval->Draw("lp");
00175 TLegend *legend = new TLegend(0.2,0.7,0.7,0.85);
00176 legend->SetTextFont(72);
00177 legend->SetTextSize(0.04);
00178 legend->AddEntry(gcx,"Compression factor","lp");
00179 legend->AddEntry(gdrange,"Log of precision wrt range","lp");
00180 legend->AddEntry(gdval,"Log of precision wrt value","lp");
00181 legend->Draw();
00182 TPaveLabel *rang = new TPaveLabel(.75,.75,.88,.80,"[-pi,pi]","brNDC");
00183 rang->Draw();
00184 }
00185
00186
00187
00188