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