]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/TRD/AliTRDpidRefMakerLQ.cxx
add variable TRD segmentation levels. Possible options are
[u/mrichter/AliRoot.git] / PWG1 / TRD / AliTRDpidRefMakerLQ.cxx
1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 *                                                                        *
4 * Author: The ALICE Off-line Project.                                    *
5 * Contributors are mentioned in the code where appropriate.              *
6 *                                                                        *
7 * Permission to use, copy, modify and distribute this software and its   *
8 * documentation strictly for non-commercial purposes is hereby granted   *
9 * without fee, provided that the above copyright notice appears in all   *
10 * copies and that both the copyright notice and this permission notice   *
11 * appear in the supporting documentation. The authors make no claims     *
12 * about the suitability of this software for any purpose. It is          *
13 * provided "as is" without express or implied warranty.                  *
14 **************************************************************************/
15
16 /* $Id: AliTRDpidRefMakerLQ.cxx 34163 2009-08-07 11:28:51Z cblume $ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //
20 //
21 //  TRD calibration class for building reference data for PID
22 //  - 2D reference histograms (responsible A.Bercuci) 
23 //  - 3D reference histograms (not yet implemented) (responsible A.Bercuci)
24 //
25 //   Origin
26 //   Alex Bercuci  (A.Bercuci@gsi.de)
27 //
28 ///////////////////////////////////////////////////////////////////////////////
29
30 #include <TSystem.h>
31 #include <TROOT.h>
32 #include <TFile.h>
33 #include <TTree.h>
34 #include <TMath.h>
35 #include <TEventList.h>
36 #include <TH2D.h>
37 #include <TH2I.h>
38 #include <TCanvas.h>
39
40 #include "AliLog.h"
41 #include "../STAT/TKDPDF.h"
42 #include "../STAT/TKDInterpolator.h"
43 #include "AliTRDpidRefMakerLQ.h"
44 #include "Cal/AliTRDCalPID.h"
45 #include "Cal/AliTRDCalPIDLQ.h"
46 #include "AliTRDseedV1.h"
47 #include "AliTRDcalibDB.h"
48 #include "AliTRDgeometry.h"
49 #include "info/AliTRDpidInfo.h"
50
51 ClassImp(AliTRDpidRefMakerLQ)
52
53 //__________________________________________________________________
54 AliTRDpidRefMakerLQ::AliTRDpidRefMakerLQ() 
55   : AliTRDpidRefMaker()
56   , fPDF(NULL)
57 {
58   //
59   // AliTRDpidRefMakerLQ default constructor
60   //
61   SetNameTitle("refMakerLQ", "PID(LQ) Reference Maker");
62 }
63
64 //__________________________________________________________________
65 AliTRDpidRefMakerLQ::AliTRDpidRefMakerLQ(const char *name)
66   : AliTRDpidRefMaker(name, "PID(LQ) Reference Maker")
67   , fPDF(NULL)
68 {
69   //
70   // AliTRDpidRefMakerLQ default constructor
71   //
72   DefineOutput(3, TObjArray::Class());
73 }
74
75 //__________________________________________________________________
76 AliTRDpidRefMakerLQ::~AliTRDpidRefMakerLQ()
77 {
78   //
79   // AliTRDCalPIDQRef destructor
80   //
81   if(fPDF){
82     //fPDF->Write("PDF_LQ", TObject::kSingleKey);
83     fPDF->Delete();
84     delete fPDF;
85   }
86 }
87
88 // //________________________________________________________________________
89 void AliTRDpidRefMakerLQ::UserCreateOutputObjects()
90 {
91   // Create histograms
92   // Called once
93
94   AliTRDpidRefMaker::UserCreateOutputObjects();
95   fContainer = Histos();
96
97   OpenFile(2, "RECREATE");
98   fPDF = new TObjArray(AliTRDCalPIDLQ::GetNrefs());
99   fPDF->SetOwner();fPDF->SetName("PDF_LQ");
100 }
101
102
103 //__________________________________________________________________
104 TObjArray* AliTRDpidRefMakerLQ::Histos()
105 {
106   // Create histograms
107
108   if(!fContainer){
109     AliError("Monitor histograms missing.");
110     return NULL;
111   }
112   Int_t n(fContainer->GetEntriesFast());
113
114   // save dE/dx references
115   TH1 *h(NULL);
116   for(Int_t ip=AliTRDCalPID::kNMom; ip--;){ 
117     TObjArray *arr = new TObjArray(2*AliPID::kSPECIES);
118     arr->SetName(Form("Pbin%02d", ip)); arr->SetOwner();
119     for(Int_t is=AliPID::kSPECIES; is--;) {
120       h = new TH1D(Form("h1%s%d", AliPID::ParticleShortName(is), ip), Form("1D %s @ Pbin[%d]", AliPID::ParticleName(is), ip), 50, 7., 12.);
121       h->GetXaxis()->SetTitle("log(dE/dx) [au]");
122       h->GetYaxis()->SetTitle("#");
123       h->SetLineColor(AliTRDCalPIDLQ::GetPartColor(is));
124       arr->AddAt(h, is);
125     }
126     for(Int_t is=AliPID::kSPECIES; is--;) {
127       h = new TH2D(Form("h2%s%d", AliPID::ParticleShortName(is), ip), Form("2D %s @ Pbin[%d]", AliPID::ParticleName(is), ip), 50, 7., 12., 50, 6.5, 11.);
128       h->GetXaxis()->SetTitle("log(dE/dx_{am}) [au]");
129       h->GetYaxis()->SetTitle("log(dE/dx_{dr}) [au]");
130       h->GetZaxis()->SetTitle("#");
131       arr->AddAt(h, AliPID::kSPECIES+is);
132     }
133     fContainer->AddAt(arr, n+ip);
134   }
135   fNRefFigures+=AliTRDCalPID::kNMom;
136   return fContainer;
137 }
138
139
140 //__________________________________________________________________
141 TObject* AliTRDpidRefMakerLQ::GetOCDBEntry(Option_t */*opt*/)
142 {
143 // Steer loading of OCDB LQ PID
144
145   if(gSystem->AccessPathName(Form("TRD.Calib%s.root", GetName()), kReadPermission)){
146     AliError(Form("File TRD.Calib%s.root not readable", GetName()));
147     return NULL;
148   }
149   AliTRDCalPIDLQ *cal = new AliTRDCalPIDLQ("pidLQ", "LQ TRD PID object");
150   cal->LoadReferences(Form("TRD.Calib%s.root", GetName()));
151   return cal;
152 }
153
154 //__________________________________________________________________
155 Bool_t AliTRDpidRefMakerLQ::GetRefFigure(Int_t ifig)
156 {
157 // Steering reference picture
158
159   if(ifig<0 || ifig>=fNRefFigures){ 
160     AliError("Ref fig requested outside definition.");
161     return kFALSE;
162   }
163   if(!gPad){
164     AliWarning("Please provide a canvas to draw results.");
165     return kFALSE;
166   }
167
168   TObjArray *arr(NULL);TList *l(NULL);TH1 *h(NULL);
169   if(!(arr = (TObjArray*)fContainer->At(ifig))){
170     AliError(Form("PDF container @ pBin[%d] missing.", ifig));
171     return kFALSE;
172   }
173   gPad->Divide(5, 2, 1.e-5, 1.e-5);l=gPad->GetListOfPrimitives(); 
174   for(Int_t is=0; is<AliPID::kSPECIES; is++){
175     ((TVirtualPad*)l->At(is))->cd();
176     if(!(h=(TH1*)arr->At(is))){
177       AliError(Form("1D for %s @ pBin[%d] missing.", AliPID::ParticleShortName(is), ifig));
178       return kFALSE;
179     }
180     h->GetYaxis()->SetRangeUser(0., 1.2);
181     h->Draw("l");
182
183     ((TVirtualPad*)l->At(AliPID::kSPECIES+is))->cd();
184     if(!(h=(TH1*)arr->At(AliPID::kSPECIES+is))){
185       AliError(Form("2D for %s @ pBin[%d] missing.", AliPID::ParticleShortName(is), ifig));
186       return kFALSE;
187     }
188     h->Draw("cont4z");
189   }
190
191   return kTRUE;
192 }
193
194
195 //________________________________________________________________________
196 Bool_t AliTRDpidRefMakerLQ::Load(const Char_t */*fname*/)
197 {
198   const Char_t *name("PIDrefMaker");
199   if(gSystem->AccessPathName(Form("TRD.Calib%s.root", name), kReadPermission)){
200     AliError(Form("File TRD.Calib%s.root not readable", name));
201     return kFALSE;
202   }
203   if(!TFile::Open(Form("TRD.Calib%s.root", name))){
204     AliError(Form("File TRD.Calib%s.root corrupted", name));
205     return kFALSE;
206   }
207   if (!(fData = dynamic_cast<TTree*>(gFile->Get(name)))) {
208     AliError(Form("Tree %s not available", name));
209     return kFALSE;
210   }
211   LinkPIDdata();
212
213   TObjArray *o(NULL);
214 /*  if(!(o = (TObjArray*)gFile->Get(Form("Moni%s", GetName())))){
215     AliWarning(Form("Monitor container Moni%s not available.", name));
216     return kFALSE;
217   }
218   fContainer = (TObjArray*)o->Clone("monitor");
219   fNRefFigures=AliTRDCalPID::kNMom;
220 */
221   // temporary until new calibration data are being produced
222   Histos();
223
224
225   if(!TFile::Open(Form("TRD.Calib%s.root", GetName()), "UPDATE")){
226     AliError(Form("File TRD.Calib%s.root corrupted", GetName()));
227     return kFALSE;
228   }
229   if(!(o = (TObjArray*)gFile->Get(Form("PDF_LQ")))) {
230     AliInfo("PDF container not available. Create.");
231     fPDF = new TObjArray(AliTRDCalPIDLQ::GetNrefs());
232     fPDF->SetOwner();fPDF->SetName("PDF_LQ");
233   } else fPDF = (TObjArray*)o->Clone("PDF_LQ");
234
235   return kTRUE;
236 }
237
238
239 //________________________________________________________________________
240 void AliTRDpidRefMakerLQ::UserExec(Option_t */*opt*/)
241 {
242 // Load PID data into local data storage
243
244   AliTRDpidInfo *pid(NULL);
245   const AliTRDpidInfo::AliTRDpidData *data(NULL);
246   Char_t s(-1);
247   for(Int_t itrk=fInfo->GetEntriesFast(); itrk--;){
248     if(!(pid=(AliTRDpidInfo*)fInfo->At(itrk))) continue;
249     if((s=pid->GetPID())<0) continue;
250     for(Int_t itrklt=pid->GetNtracklets();itrklt--;){
251       data=pid->GetData(itrklt);
252       Int_t ip(data->Momentum());
253       
254
255       Double_t dedx[] = {0., 0.};
256       if(!AliTRDCalPIDLQ::CookdEdx(data->fdEdx, dedx)) continue;
257       ((TH2*)((TObjArray*)fContainer->At(ip))->At(s+Int_t(AliPID::kSPECIES)))->Fill(dedx[0], dedx[1]);
258       AliTRDCalPIDLQ::CookdEdx(data->fdEdx, dedx, kFALSE);
259       ((TH1*)((TObjArray*)fContainer->At(ip))->At(s))->Fill(dedx[0]);
260     }
261   }
262
263   PostData(1, fContainer);
264   PostData(3, fPDF);
265 }
266
267
268 //________________________________________________________________________
269 Bool_t AliTRDpidRefMakerLQ::PostProcess()
270 {
271 // Analyse merged dedx = f(p) distributions.
272 //   - select momentum - species bins
273 //   - rotate to principal components
274 //   - locally interpolate with TKDPDF
275 //   - save interpolation to monitoring histograms
276 //   - write pdf to file for loading to OCDB
277 // 
278
279   TCanvas *fMonitor(NULL);
280   // allocate working storage
281   const Int_t kWS(AliPID::kSPECIES*AliTRDCalPID::kNMom);
282   Float_t *data[2*kWS], *data1D[kWS];
283   for(Int_t i(0); i<kWS; i++){ 
284     data1D[i]   = new Float_t[kMaxStat];
285     data[i]     = new Float_t[kMaxStat];
286     data[kWS+i] = new Float_t[kMaxStat];
287   }
288   Int_t ndata[kWS]; memset(ndata, 0, kWS*sizeof(Int_t));
289
290   AliDebug(1, Form("Loading data[%d]", fData->GetEntries()));
291   for(Int_t itrk=0; itrk < fData->GetEntries(); itrk++){
292     if(!(fData->GetEntry(itrk))) continue;
293     Int_t sbin(fPIDbin);
294     for(Int_t ily=fPIDdataArray->fNtracklets; ily--;){
295       Int_t pbin(fPIDdataArray->fData[ily].fPLbin & 0xf);
296       
297       Double_t dedx[] = {0., 0.}, 
298                dedx1D[] = {0., 0.};
299       if(!AliTRDCalPIDLQ::CookdEdx(fPIDdataArray->fData[ily].fdEdx, dedx)) continue;
300       AliTRDCalPIDLQ::CookdEdx(fPIDdataArray->fData[ily].fdEdx, dedx1D, kFALSE);
301       Int_t idx=AliTRDCalPIDLQ::GetModelID(pbin,sbin);
302       if(ndata[idx]==kMaxStat) continue;
303
304       // store data
305       data1D[idx][ndata[idx]]   = dedx1D[0];
306       data[idx][ndata[idx]]     = dedx[0];
307       data[idx+kWS][ndata[idx]] = dedx[1];
308       ndata[idx]++;
309     }
310   }
311
312   TKDPDF *pdf(NULL);
313   Int_t in(0); Float_t par[6], *pp(NULL);
314   for(Int_t ip=0;ip<AliTRDCalPID::kNMom;ip++){ 
315     for(Int_t is=AliPID::kSPECIES; is--;) {
316       // estimate bucket statistics
317       Int_t idx(AliTRDCalPIDLQ::GetModelID(ip,is)),
318             nb(kMinBuckets), // number of buckets
319         ns((Int_t)(((Float_t)(ndata[idx]))/nb));    //statistics/bucket
320             
321       AliDebug(2, Form("pBin[%d] sBin[%d] n[%d] ns[%d] nb[%d]", ip, is, ndata[idx], ns, nb));
322       if(ns<Int_t(kMinStat)){
323         AliWarning(Form("Not enough entries [%d] for %s[%d].", ndata[idx], AliPID::ParticleShortName(is), ip));
324         continue;
325       }
326
327       // build helper 1D PDF
328       pdf = new TKDPDF(ndata[idx], 1, ns, &data1D[idx]);
329       pdf->SetCOG(kFALSE);
330       pdf->SetWeights();
331       //pdf->SetStore();
332       pdf->SetAlpha(15.);
333       //pdf.GetStatus();
334       fPDF->AddAt(pdf, idx);
335       in=pdf->GetNTNodes(); pp=&par[0];
336       while(in--){
337         const TKDNodeInfo *nn = pdf->GetNodeInfo(in);
338         nn->GetCOG(pp);
339         Double_t p(par[0]),r,e;
340         pdf->Eval(&p,r,e,1);
341       }
342
343       // build helper 2D PDF
344       Float_t *ldata[2]={data[idx], data[kWS+idx]};
345       pdf = new TKDPDF(ndata[idx], 2, ns, ldata);
346       pdf->SetCOG(kFALSE);
347       pdf->SetWeights();
348       //pdf->SetStore();
349       pdf->SetAlpha(5.);
350       //pdf.GetStatus();
351       fPDF->AddAt(pdf, AliTRDCalPIDLQ::GetModelID(ip,is, 2));
352       in=pdf->GetNTNodes(); pp=&par[0];
353       while(in--){
354         const TKDNodeInfo *nn = pdf->GetNodeInfo(in);
355         nn->GetCOG(pp);
356         Double_t p[] = {par[0], par[1]}, r,e;
357         pdf->Eval(p,r,e,1);
358       }
359 /*
360       Int_t nnodes = pdf.GetNTNodes(),
361             nside = Int_t(0.05*nnodes),
362             nzeros = 4*(nside+1);
363       printf("nnodes[%d] nside[%d] nzeros[%d]\n", nnodes, nside, nzeros);
364     
365
366       // Build interpolator on the pdf skeleton
367       TKDInterpolator interpolator(2, nnodes+nzeros); 
368       for(Int_t in=nnodes; in--;)
369         interpolator.SetNode(in, *pdf.GetNodeInfo(in));
370       TKDNodeInfo *nodes = new TKDNodeInfo[nzeros], *node = &nodes[0];
371       Float_t ax0min, ax0max, ax1min, ax1max;
372       pdf.GetRange(0, ax0min, ax0max); Float_t dx = (ax0max-ax0min)/nside;
373       pdf.GetRange(1, ax1min, ax1max); Float_t dy = (ax1max-ax1min)/nside;
374       printf("x=[%f %f] y[%f %f]\n", ax0min, ax0max, ax1min, ax1max);
375
376       Int_t jn = nnodes; 
377       SetZeroes(&interpolator, node, nside, jn, ax0min, dx, ax1min, -dy, 'x');
378       SetZeroes(&interpolator, node, nside, jn, ax1min, dy, ax0max, dx, 'y');
379       SetZeroes(&interpolator, node, nside, jn, ax0max,-dx, ax1max, dy, 'x');
380       SetZeroes(&interpolator, node, nside, jn ,ax1max, -dy, ax0min, -dx, 'y');
381       delete [] nodes;
382       Int_t in=nnodes; Float_t par[6], *pp=&par[0];
383       while(in--){
384         const TKDNodeInfo *nn = interpolator.GetNodeInfo(in);
385         nn->GetCOG(pp);
386         //printf("evaluate for node[%d] @ [%f %f]\n",in, par[0], par[1]);
387         Double_t p[] = {par[0], par[1]}, r,e;
388         interpolator.Eval(p,r,e,1);
389       }
390 */
391
392       // visual on-line monitoring
393       if(HasOnlineMonitor()){
394         if(!fMonitor) fMonitor = new TCanvas("cc", "PDF 2D LQ", 500, 500);
395         pdf->DrawProjection();
396         fMonitor->Modified(); fMonitor->Update(); 
397         fMonitor->SaveAs(Form("pdf_%s%02d.gif", AliPID::ParticleShortName(is), ip));
398       }
399
400       //fContainer->ls();
401       // save a discretization of the PDF for result monitoring
402       Double_t rxy[]={0.,0.};
403       TH2 *h2s = (TH2D*)((TObjArray*)fContainer->At(ip))->At(AliPID::kSPECIES+is);
404       TAxis *ax = h2s->GetXaxis(), *ay = h2s->GetYaxis();
405       h2s->Clear();
406       for(int ix=1; ix<=ax->GetNbins(); ix++){
407         rxy[0] = ax->GetBinCenter(ix);
408         for(int iy=1; iy<=ay->GetNbins(); iy++){
409           rxy[1] = ay->GetBinCenter(iy);
410       
411           Double_t rr,ee;
412           pdf->Eval(rxy, rr, ee, kFALSE);
413           if(rr<0. || ee/rr>.15) continue; // 15% relative error
414           //printf("x[%2d] x[%2d] r[%f] e[%f]\n", ix, iy, rr, ee);
415           h2s->SetBinContent(ix, iy, rr);
416         }
417       }
418
419       pdf=dynamic_cast<TKDPDF*>(fPDF->At(idx));
420       TH1 *h1 = (TH1D*)((TObjArray*)fContainer->At(ip))->At(is);
421       ax = h1->GetXaxis();
422       h1->Clear();
423       for(int ix=1; ix<=ax->GetNbins(); ix++){
424         rxy[0] = ax->GetBinCenter(ix);
425       
426         Double_t rr,ee;
427         pdf->Eval(rxy, rr, ee, kFALSE);
428         if(rr<0. || ee/rr>.15) continue; // 15% relative error
429         //printf("x[%2d] x[%2d] r[%f] e[%f]\n", ix, iy, rr, ee);
430         h1->SetBinContent(ix, rr);
431       }
432     }
433   }
434   for(Int_t i(0); i<kWS; i++){ 
435     delete [] data1D[i];
436     delete [] data[i];
437     delete [] data[i+kWS];
438   }
439   return kTRUE;
440 }
441
442
443 //__________________________________________________________________
444 void AliTRDpidRefMakerLQ::SetZeroes(TKDInterpolator *interpolator, TKDNodeInfo *node, Int_t n, Int_t& idx, Float_t x, Float_t dx, Float_t y, Float_t dy, const Char_t opt)
445 {
446 // Set extra nodes to ensure boundary conditions
447   
448   printf("SetZeroes(%c)\n", opt);
449   Float_t par[6], val[] = {0., 1.};
450   Int_t a[6];
451   if(opt=='x'){
452     a[0]=0; a[1]=1; a[2]=2; a[3]=3; a[4]=4; a[5]=5;
453   } else if(opt=='y'){
454     a[0]=1; a[1]=0; a[2]=4; a[3]=5; a[4]=2; a[5]=3;
455   } else return;
456   Float_t tmp;
457   par[a[1]] = y;
458   par[a[4]] = y; par[a[5]] = y+dy;
459   if(dy<0.){tmp=par[a[4]]; par[a[4]]=par[a[5]]; par[a[5]]=tmp;}
460   for(Int_t in=n; in--; node++, idx++, x+=dx){
461     par[a[0]] = x+.5*dx;
462     par[a[2]] = x;  par[a[3]] = x+dx;
463     if(dx<0.){tmp=par[a[2]]; par[a[2]]=par[a[3]]; par[a[3]]=tmp;}
464     node->SetNode(2, par, val);
465     printf("\n\tnode[%d]\n", idx); node->Print();
466     interpolator->SetNode(idx, *node);
467   }
468   par[a[0]] = x;
469   par[a[2]] = x;  par[a[3]] = x+dx;
470   if(dx<0.){tmp=par[a[2]]; par[a[2]]=par[a[3]]; par[a[3]]=tmp;}
471   node->SetNode(2, par, val);
472   printf("\n\tnode[%d]\n", idx); node->Print();
473   interpolator->SetNode(idx, *node);node++;idx++;
474 }
475