]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/scripts/ESDDataSize.C
renaming function to avoid library conflict (discovered in test/generators/TUHKMgen)
[u/mrichter/AliRoot.git] / FMD / scripts / ESDDataSize.C
1 void
2 FillEtas(AliFMDFloatMap& etas)
3 {
4   AliCDBManager*  cdb  = AliCDBManager::Instance();
5   cdb->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
6   cdb->SetRun(0);
7   AliGeomManager::LoadGeometry("geometry.root");
8   AliFMDGeometry* geom = AliFMDGeometry::Instance();
9   geom->Init();
10   geom->InitTransformations();
11   
12   etas.Clear();
13   for (UShort_t d = 1; d <= 3; d++) { 
14     UShort_t nRing = (d == 1 ? 1 : 2);
15     for (UShort_t q = 0; q < nRing; q++) { 
16       Char_t   r    = (q == 0 ? 'I' : 'O');
17       UShort_t nStr = (q == 0 ? 512 : 256);
18       for (UShort_t t = 0; t < nStr; t++) { 
19         Double_t x, y, z;
20         geom->Detector2XYZ(d, r, 0, t, x, y, z);
21         Double_t phi   =  TMath::ATan2(y, x);
22         Double_t rr    =  TMath::Sqrt(y * y + x * x);
23         Double_t theta =  TMath::ATan2(rr, z);
24         Double_t eta   = -TMath::Log(TMath::Tan(theta / 2));
25         etas(d, r, 0, t) =  eta;
26       }
27     }
28   }
29 }
30
31 Double_t
32 MakeMult()
33 {
34   return TMath::Min(20., gRandom->Landau(1, .1));
35 }
36
37 void
38 FillESDFMD(AliESDFMD& fmd, AliFMDFloatMap& etas, Double_t ratio)
39 {
40   fmd.Clear();
41   
42   Int_t nTotal  = 51200;
43   Int_t nFilled = 0;
44
45   for (UShort_t d = 1; d <= 3; d++) { 
46     UShort_t nRing = (d == 1 ? 1 : 2);
47     for (UShort_t q = 0; q < nRing; q++) { 
48       Char_t   r    = (q == 0 ? 'I' : 'O');
49       UShort_t nSec = (q == 0 ?  20 :  40);
50       UShort_t nStr = (q == 0 ? 512 : 256);
51       for (UShort_t s = 0; s < nSec; s++) {
52         for (UShort_t t = 0; t < nStr; t++) { 
53           if (s == 0) { 
54             Double_t eta = etas(d, r, s, t);
55             fmd.SetEta(d, r, s, t, eta);
56           }
57           Double_t now  = Double_t(nFilled) / nTotal;
58           Bool_t   fill = now < ratio;
59           fmd.SetMultiplicity(d, r, s, t, (!fill ? 0 : MakeMult()));
60           if (fill) nFilled++;
61         }
62       }
63     }
64   }
65 }
66
67 Float_t
68 MakeTree(Int_t ratio, Int_t split, AliFMDFloatMap& etas)
69 {
70   // Split is 
71   // 
72   //  0 - No split of members
73   //  1 - Split on data and sub-objects.
74   //  2 - Same as 1
75   // 99 - Default - same as 2
76   
77   AliESDFMD* fmd = new AliESDFMD;
78   TString    fname(Form("test_%03d_%02d.root", ratio, split));
79   TString    tname(Form("T_%03d_%02d", ratio, split));
80   TFile*     file = TFile::Open(fname.Data(),"RECREATE");
81   TTree*     tree = new TTree(tname.Data(), 
82                               Form("Ratio=%3d%%, Split=%2d", ratio, split));
83   TBranch* branch = tree->Branch(split == 0 ? "FMDESD" : "FMDESD.", 
84                                  &fmd, 32000, split);
85   // branch->SetCompressionLevel(9);
86   
87   std::cout << "  Doing " << tree->GetTitle() << " " << std::flush;
88   Int_t nEvents = 10;
89   for (Int_t i = 0; i < nEvents; i++) { 
90     std::cout << "." << std::flush;
91     FillESDFMD(*fmd, etas, Double_t(ratio) / 100);
92     tree->Fill();
93   }
94   
95   
96   // branch->Print();
97   tree->GetCurrentFile()->cd();
98   tree->Write();
99   tree->GetCurrentFile()->Close();
100   
101   file = TFile::Open(fname.Data(), "UPDATE");
102   tree = static_cast<TTree*>(file->Get(tname.Data()));
103   if (split == 0) 
104     branch = tree->FindBranch("FMDESD");
105   else if (split == 1) { 
106     TLeaf* leaf = tree->FindLeaf("FMDESD.fMultiplicity");
107     if (leaf) branch = leaf->GetBranch();
108     else Warning("FillTree", "Failed to find FMDESD.fMultiplicity");
109   }
110   else if (split > 1) { 
111     TLeaf* leaf = tree->FindLeaf("FMDESD.fMultiplicity.fData");
112     if (leaf) branch = leaf->GetBranch();
113     else Warning("FillTree", "Failed to find FMDESD.fMultiplicity.fData");
114   }
115   
116   TH1* h = new TH1F("h", tree->GetTitle(), 100, 0, 20);
117   tree->Draw("FMDESD.fMultiplicity.fData>>h", 
118              "FMDESD.fMultiplicity.fData<20&&FMDESD.fMultiplicity.fData>0");
119   Int_t   nEntries = h->GetEntries();
120   Float_t fratio   = Float_t(nEntries) / 51200 / nEvents;
121   
122   Int_t    comp    = branch->GetCompressionLevel();
123   Long64_t totSize = branch->GetTotBytes();
124   Long64_t zipSize = branch->GetZipBytes();
125   Float_t  compz   = (zipSize ? Float_t(totSize+0.00001)/zipSize : 1);
126   std::cout << " done - " 
127             << "Compression: " << compz << " (" << comp << ", " 
128             << Int_t(fratio*1000+.5)/10. << "%)" << std::endl;
129   // tree->GetCurrentFile()->Write();
130   // tree->GetCurrentFile()->Close();
131
132   // delete tree;
133   // delete fmd;
134   
135   return compz;
136 }
137
138 void
139 ESDDataSize()
140 {
141   gStyle->SetPalette(1);
142   
143   AliFMDFloatMap etas(3, 2, 1, 512);
144   FillEtas(etas);
145   
146   Int_t  ratios[] = { 1, 10, 25, 50, 100, -1 };
147   Int_t  splits[] = { 0, 1, 99, -1 };
148   Int_t* pratio   = ratios;
149
150   Double_t arat[] = { 0.5, 1.5, 9.5, 10.5, 24.5, 25.5, 49.5, 50.5, 99.5, 100.5};
151   Double_t aspl[] = { -0.5, 0.5, 1.5, 98.5, 99.5 }; 
152   TH2* h = new TH2F("size", "ESD size", 9, arat, 4, aspl);
153   h->GetXaxis()->SetTitle("Fill ratio");
154   h->GetYaxis()->SetTitle("Split level");
155   h->SetDirectory(0);
156   h->Draw("COLZ");
157   std::cout << "Loop over ratios" << std::endl;
158   while (*pratio > 0) { 
159     Int_t  ratio  = *pratio++;
160     Int_t* psplit = splits;
161
162     std::cout << " Loop over splits" << std::endl;
163     while (*psplit >= 0) { 
164       Int_t split = *psplit++;
165       
166       Double_t comp = MakeTree(ratio, split, etas);
167       h->Fill(ratio, split, comp);
168       // h->Draw("COLZ");
169     }
170   }
171   h->Draw("LEGO2");
172 }
173
174