]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/qaRec/AliTRDcalibration.cxx
Move to THnSparse instead of TH2F to reduce memory size
[u/mrichter/AliRoot.git] / TRD / qaRec / AliTRDcalibration.cxx
1
2 /**************************************************************************
3 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 *                                                                        *
5 * Author: The ALICE Off-line Project.                                    *
6 * Contributors are mentioned in the code where appropriate.              *
7 *                                                                        *
8 * Permission to use, copy, modify and distribute this software and its   *
9 * documentation strictly for non-commercial purposes is hereby granted   *
10 * without fee, provided that the above copyright notice appears in all   *
11 * copies and that both the copyright notice and this permission notice   *
12 * appear in the supporting documentation. The authors make no claims     *
13 * about the suitability of this software for any purpose. It is          *
14 * provided "as is" without express or implied warranty.                  *
15 **************************************************************************/
16
17 /////////////////////////////////////////////////////////////////////////////////
18 //                                                                             
19 // AliTRDcalibration                                                            
20 //                                                                             
21 // Task to run the calibration offline.
22 // Author:
23 //   R. Bailhache (rbailhache@ikf.uni-frankfurt.de, R.Bailhache@gsi.de)
24 //           
25 //////////////////////////////////////////////////////////////////////////////////
26
27
28 #include "Riostream.h"
29 #include "TChain.h"
30 #include "TTree.h"
31 #include "TProfile2D.h"
32 #include "TH2I.h"
33 #include "TH1F.h"
34 #include "TList.h"
35 #include "TMath.h"
36 #include "TCanvas.h"
37 #include "TObject.h"
38 #include "TFile.h"
39 #include "TObjArray.h"
40 #include "TGraph.h"
41 #include "TStyle.h"
42 #include "TLegend.h"
43 #include "TGraphErrors.h"
44 #include "THnSparse.h"
45
46 #include "AliTRDrecoTask.h"
47 #include "AliAnalysisManager.h"
48
49 #include "AliESDInputHandler.h"
50 #include "AliTRDtrackV1.h"
51 #include "AliTRDseedV1.h"
52 #include "AliTRDcluster.h"
53 #include "info/AliTRDtrackInfo.h"
54 #include "AliTRDcalibDB.h"
55
56 #include "AliTRDCalibraFillHisto.h"
57 #include "AliTRDCalibraFit.h"
58 #include "AliTRDCalibraVdriftLinearFit.h"
59 #include "AliTRDCalibraMode.h"
60 #include "AliTRDCalibraVector.h"
61 #include "./Cal/AliTRDCalPad.h"
62 #include "./Cal/AliTRDCalDet.h"
63 #include "AliCDBMetaData.h"
64 #include "AliCDBManager.h"
65 #include "AliCDBStorage.h"
66
67 #include "AliLog.h"
68
69 #include "AliTRDcalibration.h"
70
71
72 ClassImp(AliTRDcalibration)
73
74 //________________________________________________________________________
75 AliTRDcalibration::AliTRDcalibration() 
76   :AliTRDrecoTask("Calibration", "Calibration on tracks")
77   ,fTrackInfo(0)
78   ,ftrdTrack(0)
79   ,fcl(0)
80   ,fTRDCalibraFillHisto(0)
81   ,fNbTRDTrack(0)
82   ,fNbTRDTrackOffline(0)
83   ,fNbTRDTrackStandalone(0)
84   ,fNbTRDTracklet(0)
85   ,fNbTRDTrackletOffline(0)
86   ,fNbTRDTrackletStandalone(0)
87   ,fNbTimeBin(0x0)
88   ,fNbTimeBinOffline(0x0)
89   ,fNbTimeBinStandalone(0x0)
90   ,fNbClusters(0)
91   ,fNbClustersOffline(0)
92   ,fNbClustersStandalone(0)
93   ,fPHSum(0)
94   ,fCHSum(0)
95   ,fDetSum(0)
96   ,fDetSumVector(0)
97   ,fHisto2d(kTRUE)
98   ,fVector2d(kFALSE)
99   ,fVdriftLinear(kTRUE)
100   ,flow(0)
101   ,fhigh(30)
102   ,fNbTimeBins(30)
103   ,ffillZero(kFALSE)
104   ,fnormalizeNbOfCluster(kFALSE)
105   ,fmaxCluster(0)
106   ,fOfflineTracks(kFALSE)
107   ,fStandaloneTracks(kFALSE)
108   ,fCompressPerDetector(kFALSE)
109   ,fRunNumber(0)
110   ,fkNameDirectory("local://.")
111   ,fGraph(0x0)
112   ,fPostProcess(kFALSE)
113 {
114   // Constructor
115   
116   fNRefFigures = 17;
117
118   for(Int_t k = 0; k < 3; k++)
119     {
120       fNz[k]=0;
121       fNrphi[k]=0;
122     }
123
124 }  
125 //________________________________________________________________________
126 AliTRDcalibration::~AliTRDcalibration() 
127 {
128   // Default destructor
129
130   if(fNbTRDTrack) delete fNbTRDTrack;
131   if(fNbTRDTrackOffline) delete fNbTRDTrackOffline;
132   if(fNbTRDTrackStandalone) delete fNbTRDTrackStandalone;
133   if(fNbTRDTracklet) delete fNbTRDTracklet;
134   if(fNbTRDTrackletOffline) delete fNbTRDTrackletOffline;
135   if(fNbTRDTrackletStandalone) delete fNbTRDTrackletStandalone;
136   if(fNbTimeBin) delete fNbTimeBin;
137   if(fNbTimeBinOffline) delete fNbTimeBinOffline;
138   if(fNbTimeBinStandalone) delete fNbTimeBinStandalone;
139   if(fNbClusters) delete fNbClusters;
140   if(fNbClustersOffline) delete fNbClustersOffline;
141   if(fNbClustersStandalone) delete fNbClustersStandalone;
142   if(fPHSum) delete fPHSum;
143   if(fCHSum) delete fCHSum;
144   if(fDetSum) delete fDetSum;
145   if(fDetSumVector) delete fDetSumVector;
146   if(fkNameDirectory) delete fkNameDirectory;
147   if(fGraph){fGraph->Delete(); delete fGraph;}
148
149 }
150 //________________________________________________________________________
151 void AliTRDcalibration::CreateOutputObjects() 
152 {
153   // Create output objects
154
155   OpenFile(0, "RECREATE");
156   
157   // Number of time bins
158   AliTRDcalibDB *cal = AliTRDcalibDB::Instance();
159   fNbTimeBins = cal->GetNumberOfTimeBins();
160   
161   // instance calibration: what to calibrate
162   fTRDCalibraFillHisto = AliTRDCalibraFillHisto::Instance();
163   fTRDCalibraFillHisto->SetHisto2d(fHisto2d); // choose to use histograms
164   fTRDCalibraFillHisto->SetVector2d(fVector2d); // choose to use vectors
165   fTRDCalibraFillHisto->SetCH2dOn();  // choose to calibrate the gain
166   fTRDCalibraFillHisto->SetPH2dOn();  // choose to calibrate the drift velocity
167   fTRDCalibraFillHisto->SetPRF2dOn(); // choose to look at the PRF
168   fTRDCalibraFillHisto->SetLinearFitterOn(fVdriftLinear); // Other possibility vdrift VDRIFT
169   fTRDCalibraFillHisto->SetLinearFitterDebugOn(fVdriftLinear); // Other possibility vdrift
170
171   // segmentation (should be per default the max and add at the end)
172   for(Int_t k = 0; k < 3; k++){
173     if(((fNz[k] != 10) && (fNrphi[k] != 10)) && ((fNz[k] != 100) && (fNrphi[k] != 100))) {
174       fTRDCalibraFillHisto->SetNz(k,fNz[k]);                                    // Mode calibration
175       fTRDCalibraFillHisto->SetNrphi(k,fNrphi[k]);                             // Mode calibration
176     }
177     else {
178       if((fNz[k] == 100) && (fNrphi[k] == 100))  {
179         if(fVector2d) AliInfo("The mode all together is not supported by the vector method");
180         fTRDCalibraFillHisto->SetAllTogether(k);
181       }
182       if((fNz[k] == 10) && (fNrphi[k] == 10))  {
183         if(fVector2d) AliInfo("The mode per supermodule is not supported by the vector method");
184         fTRDCalibraFillHisto->SetPerSuperModule(k);
185       }
186     }
187   }
188
189   // Debug level
190   fTRDCalibraFillHisto->SetDebugLevel(fDebugLevel); //debug stuff
191
192   // Init the stuff
193   fTRDCalibraFillHisto->Init2Dhistos(); // initialise the histos
194
195   // cuts
196   fTRDCalibraFillHisto->SetNumberClusters(flow); // At least flow clusters
197   fTRDCalibraFillHisto->SetNumberClustersf(fhigh); // The more fhigh clusters
198   fTRDCalibraFillHisto->SetFillWithZero(ffillZero); // Fill zeros
199   fTRDCalibraFillHisto->SetNormalizeNbOfCluster(fnormalizeNbOfCluster); // For iterations
200
201   // Add them to the container
202   fContainer = new TObjArray();
203   if(fHisto2d) {
204     fContainer->Add(fTRDCalibraFillHisto->GetCH2d()); //TH2I
205     fContainer->Add(fTRDCalibraFillHisto->GetPH2d()); //TProfile2D
206     fContainer->Add(fTRDCalibraFillHisto->GetPRF2d()); //TProfile2D
207   }
208   if(fVdriftLinear) fContainer->Add(fTRDCalibraFillHisto->GetVdriftLinearFit()); // Other drift velocity 
209   if(fVector2d) fContainer->Add(fTRDCalibraFillHisto->GetCalibraVector()); //calibra vector
210       
211   if(fDebugLevel) {
212     
213     // Init the debug histos
214     fNbTRDTrack = new TH1F("TRDTrack","TRDTrack",500,0,500);
215     fNbTRDTrack->Sumw2();
216     fNbTRDTrackOffline = new TH1F("TRDTrackOffline","TRDTrackOffline",500,0,500);
217     fNbTRDTrackOffline->Sumw2();
218     fNbTRDTrackStandalone = new TH1F("TRDTrackStandalone","TRDTrackStandalone",500,0,500);
219     fNbTRDTrackStandalone->Sumw2();
220     //
221     fNbTRDTracklet = new TH1F("TRDTracklet","TRDTracklet",540,0.,540.);
222     fNbTRDTracklet->Sumw2();
223     fNbTRDTrackletOffline = new TH1F("TRDTrackletOffline","TRDTrackletOffline",540,0.,540.);
224     fNbTRDTrackletOffline->Sumw2();
225     fNbTRDTrackletStandalone = new TH1F("TRDTrackletStandalone","TRDTrackletStandalone",540,0.,540.);
226     fNbTRDTrackletStandalone->Sumw2();
227     //
228     fNbTimeBin = new TH1F("TimeBin","TimeBin",35,0,35);
229     fNbTimeBin->Sumw2();
230     fNbTimeBinOffline = new TH1F("TimeBinOffline","TimeBinOffline",35,0,35);
231     fNbTimeBinOffline->Sumw2();
232     fNbTimeBinStandalone = new TH1F("TimeBinStandalone","TimeBinStandalone",35,0,35);
233     fNbTimeBinStandalone->Sumw2();
234     //
235     fNbClusters = new TH1F("NbClusters","",35,0,35);
236     fNbClusters->Sumw2();
237     fNbClustersOffline = new TH1F("NbClustersOffline","",35,0,35);
238     fNbClustersOffline->Sumw2();
239     fNbClustersStandalone = new TH1F("NbClustersStandalone","",35,0,35);
240     fNbClustersStandalone->Sumw2();
241     //
242     fPHSum = new TProfile2D("PH2dSum","Nz0Nrphi0"
243                             ,fNbTimeBins,-0.05,(Double_t)(fNbTimeBins/10.0)-0.05
244                             ,540,0,540);
245     fPHSum->SetYTitle("Det/pad groups");
246     fPHSum->SetXTitle("time [#mus]");
247     fPHSum->SetZTitle("<PH> [a.u.]");
248     fPHSum->SetStats(0);
249     //
250     fCHSum = new TH2I("CH2dSum","Nz0Nrphi0",100,0,300,540,0,540);
251     fCHSum->SetYTitle("Det/pad groups");
252     fCHSum->SetXTitle("charge deposit [a.u]");
253     fCHSum->SetZTitle("counts");
254     fCHSum->SetStats(0);
255     fCHSum->Sumw2();
256     
257     // Add them
258     fContainer->Add(fNbTRDTrack);
259     fContainer->Add(fNbTRDTrackOffline);
260     fContainer->Add(fNbTRDTrackStandalone);
261     fContainer->Add(fNbTRDTracklet);
262     fContainer->Add(fNbTRDTrackletOffline);
263     fContainer->Add(fNbTRDTrackletStandalone);
264     fContainer->Add(fNbTimeBin);
265     fContainer->Add(fNbTimeBinOffline);
266     fContainer->Add(fNbTimeBinStandalone);
267     fContainer->Add(fNbClusters);
268     fContainer->Add(fNbClustersOffline);
269     fContainer->Add(fNbClustersStandalone);
270     fContainer->Add(fPHSum);
271     fContainer->Add(fCHSum);
272
273   }
274
275 }
276
277 //________________________________________________________________________
278 void AliTRDcalibration::Exec(Option_t *) 
279 {
280   //
281   // Execute function where the reference data are filled
282   //
283
284   if(!fTracks) return;
285   
286   // In total
287   Int_t nbTrdTracks = 0;
288   // standalone
289   Int_t nbTrdTracksStandalone = 0;
290   // offline
291   Int_t nbTrdTracksOffline = 0;
292   
293
294   //
295   // Loop on track in the event
296   //
297   //printf("Total of %d\n",fTracks->GetEntriesFast());
298   for(Int_t itrk=0; itrk < fTracks->GetEntriesFast(); itrk++){
299     
300     //printf("itrk %d\n",itrk);
301
302     fTrackInfo = (AliTRDtrackInfo*)fTracks->UncheckedAt(itrk);
303     ftrdTrack = fTrackInfo->GetTrack();
304     if(!ftrdTrack) continue;
305
306     nbTrdTracks++;
307   
308     fTRDCalibraFillHisto->UpdateHistogramsV1(ftrdTrack);
309
310     if(fDebugLevel) {
311       
312       Bool_t standalonetracklet = kFALSE;  
313       const AliTRDseedV1 *tracklet = 0x0;
314       //
315       // Loop on tracklet in the event
316       //
317       for(Int_t itr = 0; itr < 6; itr++){
318         //printf("itr %d\n",itr);
319         if(!(tracklet = ftrdTrack->GetTracklet(itr))) continue;
320         if(!tracklet->IsOK()) continue;
321         // standalone
322         if(tracklet->IsStandAlone()) standalonetracklet = kTRUE;
323         Int_t nbclusters = 0;
324         // For PH
325         Double_t phtb[AliTRDseedV1::kNtb];
326         memset(phtb, 0, AliTRDseedV1::kNtb*sizeof(Double_t));
327         // For CH
328         Double_t sum = 0.0;
329         // normalisation
330         Float_t normalisation = 6.67;
331         Int_t detector = 0;
332         Int_t crossrow = 0;
333         for(int ic=0; ic<AliTRDseedV1::kNtb; ic++){
334           // Check no shared clusters
335           for(int icc=AliTRDseedV1::kNtb; icc<AliTRDseedV1::kNclusters; icc++){
336             if((fcl = tracklet->GetClusters(icc)))  crossrow = 1;
337           }
338           if(!(fcl = tracklet->GetClusters(ic))) continue;
339           nbclusters++;
340           Int_t time = fcl->GetPadTime();
341           Float_t ch =  tracklet->GetdQdl(ic);
342           Float_t qcl = TMath::Abs(fcl->GetQ());
343           detector = fcl->GetDetector();          
344           if((time>-1) && (time<fNbTimeBins)) phtb[time]=qcl;
345           sum += ch/normalisation;
346           fNbTimeBin->Fill(time);
347           if(tracklet->IsStandAlone()) fNbTimeBinStandalone->Fill(time);
348           else fNbTimeBinOffline->Fill(time);
349         }
350
351         fNbTRDTracklet->Fill(detector);
352         if(tracklet->IsStandAlone()) fNbTRDTrackletStandalone->Fill(detector);
353         else fNbTRDTrackletOffline->Fill(detector);
354         
355         fNbClusters->Fill(nbclusters);
356         if(tracklet->IsStandAlone())  fNbClustersStandalone->Fill(nbclusters);
357         else  fNbClustersOffline->Fill(nbclusters);
358         
359         if((nbclusters > flow) && (nbclusters < fhigh)){
360           fCHSum->Fill(sum/20.0,0.0);
361           for(int ic=0; ic<fNbTimeBins; ic++){
362             //printf("ic %d and time %f and cluster %f \n",ic,(Double_t)(ic/10.0),(Double_t)phtb[ic]);
363             if(ffillZero) fPHSum->Fill((Double_t)(ic/10.0),0.0,(Double_t)phtb[ic]);
364             else {
365               if(phtb[ic] > 0.0) fPHSum->Fill((Double_t)(ic/10.0),0.0,(Double_t)phtb[ic]);
366             }
367           }
368         }
369       }
370     
371     if(standalonetracklet) nbTrdTracksStandalone++;
372     else nbTrdTracksOffline++;
373     
374     }
375     
376   }
377   
378   if(fDebugLevel) {
379     
380     //Fill Histos
381     fNbTRDTrack->Fill(nbTrdTracks);
382     fNbTRDTrackStandalone->Fill(nbTrdTracksStandalone);
383     fNbTRDTrackOffline->Fill(nbTrdTracksOffline);
384     
385   }
386
387   //printf("Nbof tracks %d\n",nbTrdTracks);
388   
389   TFile *file = TFile::Open("test.root","RECREATE");
390   fContainer->Write();
391   file->Close();
392
393   // Post output data
394   PostData(0, fContainer);
395
396   //printf("post container\n");
397   
398 }      
399 //________________________________________________________________________
400 void AliTRDcalibration::Terminate(Option_t *) 
401 {
402   // Draw result to the screen
403   // Called once at the end of the query
404
405   //printf("terminate\n");
406
407   if(fTRDCalibraFillHisto) fTRDCalibraFillHisto->DestroyDebugStreamer();
408  
409 }
410 //________________________________________________________
411 Bool_t AliTRDcalibration::GetRefFigure(Int_t ifig)
412 {
413   //
414   // Draw filled histos
415   //
416   
417   gStyle->SetPalette(1);
418   gStyle->SetOptStat(1111);
419   gStyle->SetPadBorderMode(0);
420   gStyle->SetCanvasColor(10);
421   gStyle->SetPadLeftMargin(0.13);
422   gStyle->SetPadRightMargin(0.13);
423
424   if(!fContainer) return kFALSE;
425   
426   switch(ifig){
427   case kNbTrack:{
428     TCanvas *c0 = new TCanvas("c0","c0",10,10,510,510);
429     TLegend *legNbTrack = new TLegend(.7, .7, .98, .98);
430     legNbTrack->SetBorderSize(1);
431     TH1F *h  = 0x0;
432     TH1F *ho = 0x0;
433     TH1F *hs = 0x0;
434     if(!(h = (TH1F *)fContainer->FindObject("TRDTrack"))) break;
435     if(!(ho = (TH1F *)fContainer->FindObject("TRDTrackOffline"))) break;
436     if(!(hs = (TH1F *)fContainer->FindObject("TRDTrackStandalone"))) break;
437     c0->cd();
438     //gPad->SetLogy();
439     gPad->SetGridy();
440     gPad->SetGridx();
441     h->Draw();
442     ho->Draw("same");
443     hs->Draw("same");
444     legNbTrack->AddEntry(h, "all", "p");
445     legNbTrack->AddEntry(ho, "offline", "p");
446     legNbTrack->AddEntry(hs, "standalone", "p");
447     legNbTrack->Draw("same");
448     return kTRUE;
449   }
450   case kNbTracklet:{
451     TLegend *legNbTracklet = new TLegend(.7, .7, .98, .98);
452     legNbTracklet->SetBorderSize(1);
453     TH1F *h = 0x0;
454     TH1F *ho = 0x0;
455     TH1F *hs = 0x0;
456     if(!(h = (TH1F *)fContainer->FindObject("TRDTracklet"))) break;
457     if(!(ho = (TH1F *)fContainer->FindObject("TRDTrackletOffline"))) break;
458     if(!(hs = (TH1F *)fContainer->FindObject("TRDTrackletStandalone"))) break;
459     h->Draw();
460     ho->Draw("same");
461     hs->Draw("same");
462     legNbTracklet->AddEntry(h, "all", "p");
463     legNbTracklet->AddEntry(ho, "offline", "p");
464     legNbTracklet->AddEntry(hs, "standalone", "p");
465     legNbTracklet->Draw("same");
466     gPad->SetLogy();
467     gPad->SetGridy();
468     gPad->SetGridx();
469     return kTRUE;
470   }
471   case kNbTimeBin:{
472     TLegend *legNbTimeBin = new TLegend(.7, .7, .98, .98);
473     legNbTimeBin->SetBorderSize(1);
474     TH1F *h = 0x0;
475     TH1F *ho = 0x0;
476     TH1F *hs = 0x0;
477     if(!(h = (TH1F *)fContainer->FindObject("TimeBin"))) break;
478     if(!(ho = (TH1F *)fContainer->FindObject("TimeBinOffline"))) break;
479     if(!(hs = (TH1F *)fContainer->FindObject("TimeBinStandalone"))) break;
480     h->Draw();
481     ho->Draw("same");
482     hs->Draw("same");
483     legNbTimeBin->AddEntry(h, "all", "p");
484     legNbTimeBin->AddEntry(ho, "offline", "p");
485     legNbTimeBin->AddEntry(hs, "standalone", "p");
486     legNbTimeBin->Draw("same");
487     //gPad->SetLogy();
488     gPad->SetGridy();
489     gPad->SetGridx();
490     return kTRUE;
491   }
492   case kNbClusters:{
493     TLegend *legNbClusters = new TLegend(.7, .7, .98, .98);
494     legNbClusters->SetBorderSize(1);
495     TH1F *h = 0x0;
496     TH1F *ho = 0x0;
497     TH1F *hs = 0x0;
498     if(!(h = (TH1F *)fContainer->FindObject("NbClusters"))) break;
499     if(!(ho = (TH1F *)fContainer->FindObject("NbClustersOffline"))) break;
500     if(!(hs = (TH1F *)fContainer->FindObject("NbClustersStandalone"))) break;
501     h->Draw();
502     ho->Draw("same");
503     hs->Draw("same");
504     legNbClusters->AddEntry(h, "all", "p");
505     legNbClusters->AddEntry(ho, "offline", "p");
506     legNbClusters->AddEntry(hs, "standalone", "p");
507     legNbClusters->Draw("same");
508     gPad->SetLogy();
509     gPad->SetGridy();
510     gPad->SetGridx();
511     return kTRUE;
512   }
513   case kPHSum:{
514     TProfile2D *h = 0x0;
515     if(!(h = (TProfile2D *)fContainer->FindObject("PH2dSum"))) break;
516     TH1D *projh = h->ProjectionX("projh",1,1,"e");
517     projh->Draw();
518     gPad->SetGridy();
519     gPad->SetGridx();
520     return kTRUE;
521   }
522   case kCHSum:{
523     TH2I *h = 0x0;
524     if(!(h = (TH2I *)fContainer->FindObject("CH2dSum"))) break;
525     TH1D *projh = h->ProjectionX("projhh",1,1,"e");
526     projh->Draw();
527     gPad->SetGridy();
528     gPad->SetGridx();
529     return kTRUE;
530   }
531   case kPH2D:{
532     if(!fHisto2d) {
533       AliInfo("Histo was not filled!");
534       break;
535     }
536     TProfile2D *h = 0x0;
537     if(!(h = (TProfile2D *)fContainer->FindObject("PH2d"))) break;
538     h->Draw("lego");
539     return kTRUE;
540   }
541   case kCH2D:{
542     if(!fHisto2d) {
543       AliInfo("Histo was not filled!");
544       break;
545     }
546     TH2I *h = 0x0;
547     if(!(h = (TH2I *)fContainer->FindObject("CH2d"))) break;
548     h->Draw("lego");
549     return kTRUE;
550   }
551   case kPRF2D:{
552     if(!fHisto2d) {
553       AliInfo("Histo was not filled!");
554       break;
555     }
556     TProfile2D *h = 0x0;
557     if(!(h = (TProfile2D *)fContainer->FindObject("PRF2d"))) break;
558     h->Draw("lego");
559     return kTRUE;
560   }
561   case kPH2DVector:{
562     if(!fVector2d) {
563       AliInfo("vector was not filled!");
564       break;
565     }
566     AliTRDCalibraVector *v = 0x0;
567     TGraphErrors *vdet = 0x0; 
568     if(!(v = (AliTRDCalibraVector *)fContainer->FindObject("AliTRDCalibraVector"))) break;
569     Int_t detectormax = -1;
570     Int_t groupmax    = -1;
571     if(!v->FindTheMaxEntries(1,detectormax,groupmax)) break;
572     if(!(vdet = v->ConvertVectorPHTGraphErrors((Int_t)detectormax,groupmax,"plotPH2dVector"))) break;
573     Int_t nbeentries = 0;
574     TH1F *ko = v->CorrectTheError(vdet,nbeentries);
575     ko->Draw();
576     AliInfo(Form("There are %d entries in the detector %d and group %d",nbeentries,detectormax,groupmax));
577     return kTRUE;
578   }
579 case kCH2DVector:{
580     if(!fVector2d) {
581       AliInfo("vector was not filled!");
582       break;
583     }
584     AliTRDCalibraVector *v = 0x0;
585     TH1F *vdet = 0x0; 
586     if(!(v = (AliTRDCalibraVector *)fContainer->FindObject("AliTRDCalibraVector"))) break;
587     Int_t detectormax = -1;
588     Int_t groupmax    = -1;
589     if(!v->FindTheMaxEntries(0,detectormax,groupmax)) break;
590     if(!(vdet = v->ConvertVectorCHHisto((Int_t)detectormax,groupmax,"plotCH2dVector"))) break;
591     vdet->Draw();
592     AliInfo(Form("The detectormax and groupmax are %d and %d",detectormax,groupmax));
593     return kTRUE;
594   }
595   case kPRF2DVector:{
596     if(!fVector2d) {
597       AliInfo("vector was not filled!");
598       break;
599     }
600     AliTRDCalibraVector *v = 0x0;
601     TGraphErrors *vdet = 0x0; 
602     if(!(v = (AliTRDCalibraVector *)fContainer->FindObject("AliTRDCalibraVector"))) break;
603     Int_t detectormax  = -1;
604     Int_t groupmax     = -1;
605     Int_t nbeentries   = 0;
606     if(!v->FindTheMaxEntries(2,detectormax,groupmax)) break;
607     if(!(vdet = v->ConvertVectorPRFTGraphErrors((Int_t)detectormax,groupmax,"plotPRF2dVector"))) break;
608     TH1F *ko = v->CorrectTheError(vdet,nbeentries);
609     ko->Draw();
610     AliInfo(Form("The detectormax and groupmax are %d and %d",detectormax,groupmax));
611     return kTRUE;
612   }
613   case kLinearFitter:{
614     if(!fVdriftLinear) {
615       AliInfo("vdrift linear was not filled!");
616       break;
617     }
618     AliTRDCalibraVdriftLinearFit *h = 0x0;
619     THnSparseS *hdetector = 0x0; 
620     if(!(h = (AliTRDCalibraVdriftLinearFit *)fContainer->FindObject("AliTRDCalibraVdriftLinearFit"))) break;
621     Double_t entries[540];
622     for(Int_t k = 0; k < 540; k++){
623       entries[k] = 0.0;
624       hdetector = 0x0;
625       if(!(hdetector = (THnSparseS *)h->GetLinearFitterHisto(k,kFALSE))) continue;
626       entries[k] = hdetector->GetEntries();
627     }
628     Double_t max = -10.0;
629     Double_t detectormax = -1;
630     for(Int_t k = 0; k < 540; k++){
631       if(entries[k] > max) {
632         max = entries[k];
633         detectormax = k;
634       }
635     }
636     hdetector = 0x0;
637     if((max == 0.0) || (detectormax <0.0) || (detectormax >=540.0)) break;
638     if(!(hdetector = (THnSparseS *)h->GetLinearFitterHisto((Int_t)detectormax,kFALSE))) break;
639     AliInfo(Form("The detector with the maximum of entries is %d",detectormax));
640     hdetector->Draw();
641     return kTRUE;
642   }
643   case kGainFactor:{
644     if(!fPostProcess){
645       if(!PostProcess()) break;
646     }
647     TGraph *fgain = (TGraph *) fGraph->At(0);
648     if(!fgain) break;
649     fgain->Draw("ALP");
650     return kTRUE;
651   }
652   case kVdriftT0Factor:{
653     if(!fPostProcess){
654       if(!PostProcess()) break;
655     }
656     TCanvas *c = new TCanvas("c","c",10,10,510,510);
657     c->Divide(2,1);
658     TGraph *fvd = (TGraph *) fGraph->At(1);
659     if(fvd){
660       c->cd(1);
661       fvd->Draw("ALP");
662     } 
663     TGraph *ft0 = (TGraph *) fGraph->At(2);
664     if(ft0){
665       c->cd(2);
666       ft0->Draw("ALP");
667     } 
668     return kTRUE;
669   }
670   case kVdriftLorentzAngleFactor:{
671     if(!fVdriftLinear) {
672       AliInfo("vdrift linear was not filled!");
673       break;
674     }
675     if(!fPostProcess){
676       if(!PostProcess()) break;
677     }
678     TCanvas *c = new TCanvas("c","c",10,10,510,510);
679     c->Divide(2,1);
680     TGraph *fvdl = (TGraph *) fGraph->At(3);
681     if(fvdl){
682       c->cd(1);
683       fvdl->Draw("ALP");
684     } 
685     TGraph *flal = (TGraph *) fGraph->At(4);
686     if(flal){
687       c->cd(2);
688       flal->Draw("ALP");
689     } 
690     return kTRUE;
691   }
692   case kPRFFactor:{
693     if(!fPostProcess){
694       if(!PostProcess()) break;
695     }
696     TGraph *fprf = (TGraph *) fGraph->At(5);
697     if(!fprf) break;
698     fprf->Draw("ALP");
699     return kTRUE;
700   }
701   }
702   
703   return kFALSE;
704   
705 }
706 //________________________________________________________________________
707 Bool_t AliTRDcalibration::PostProcess()
708 {
709   // 
710   // Fit the filled histos
711   //
712
713   if(!fGraph){
714     fGraph = new TObjArray(6);
715     fGraph->SetOwner();
716   }
717   else {
718     delete fGraph;
719     PostProcess();
720   }
721
722   Bool_t storage[3] = {kFALSE,kFALSE,kFALSE};
723
724   // storage element
725   AliCDBManager *man = AliCDBManager::Instance();
726   man->SetDefaultStorage("local://$ALICE_ROOT");
727   AliCDBStorage* storLoc = man->GetStorage(fkNameDirectory);
728   if (!storLoc)
729     return kFALSE;
730   man->SetRun(fRunNumber);
731
732   // MetaData
733   AliCDBMetaData mdDet; 
734   mdDet.SetObjectClassName("AliTRDCalDet");
735   AliCDBMetaData mdPad; 
736   mdPad.SetObjectClassName("AliTRDCalPad");
737   
738   // Objects for fitting
739   AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
740   calibra->SetDebugLevel(2); // 0 rien, 1 fitvoir, 2 debug files, 3 one detector  
741   
742   // Take the stuff
743   if (!fContainer) {
744     Printf("ERROR: list not available");
745     return kFALSE;
746   }
747
748   if(fHisto2d && fVector2d) AliInfo("We will only look at histos. Set fHisto2d off if you don't want");
749   AliTRDCalibraVector *calibraVector = 0x0;
750   if(fVector2d) calibraVector = (AliTRDCalibraVector *) fContainer->FindObject("CalibraVector");
751   //
752   // GAIN TH2I
753   //
754   Bool_t pass = kFALSE; 
755   AliTRDCalibraVector *vvect = 0x0;
756   if(fHisto2d) {
757     TH2I *histogain = (TH2I *) fContainer->FindObject("CH2d");  
758     if(histogain) {
759       histogain->SetDirectory(0);
760       calibra->SetMinEntries(20); 
761       if(fCompressPerDetector){
762         if(AddStatsPerDetector(histogain)) pass = calibra->AnalyseCH(fCHSum);
763       }
764       else pass = calibra->AnalyseCH(histogain);
765     }
766   }
767   else {
768     if(fVector2d && calibraVector) {
769       calibra->SetMinEntries(20); 
770       if(fCompressPerDetector){
771         if(!(vvect = calibraVector->AddStatsPerDetectorCH())) return kFALSE;
772         pass = calibra->AnalyseCH(vvect);
773       }
774       else pass = calibra->AnalyseCH(calibraVector);
775     }
776   }
777   
778   if(pass)
779     {
780       Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(0))
781         + 6*  18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(0));
782       Int_t nbfit = calibra->GetNumberFit();  //Number of fits
783       Int_t nbE   = calibra->GetNumberEnt();  //Number of detector mit entries
784       // enough statistics
785       if ((nbtg >                  0) && 
786           (nbfit        >= 0.001*nbE))
787         {
788           // create the cal objects
789           calibra->PutMeanValueOtherVectorFit(1,kTRUE); 
790           TObjArray object           = calibra->GetVectorFit();
791           AliTRDCalDet *objgaindet   = calibra->CreateDetObjectGain(&object);
792           TObject *objgainpad        = calibra->CreatePadObjectGain();
793           // store
794           AliCDBId id1("TRD/Calib/ChamberGainFactor",fRunNumber, AliCDBRunRange::Infinity()); 
795           storLoc->Put((TObject *)objgaindet, id1, &mdDet); 
796           AliCDBId id2("TRD/Calib/LocalGainFactor",fRunNumber, AliCDBRunRange::Infinity()); 
797           storLoc->Put((TObject *)objgainpad, id2, &mdPad);
798           storage[0] = kTRUE;
799           // Make graph
800           TGraph *graph = 0x0;
801           if(FillGraphIndex(&object,graph)){ 
802             fGraph->AddAt(graph,0);
803           }
804         }//if(enough statistics?)
805       calibra->ResetVectorFit();
806     }
807   else return kFALSE;
808   
809   //
810   // VDRIFT average pulse height
811   //
812   pass = kFALSE; 
813   if(fHisto2d) {
814     TProfile2D *histodriftvelocity = (TProfile2D *) fContainer->FindObject("PH2d");  
815     if(histodriftvelocity) {
816       histodriftvelocity->SetDirectory(0);  
817       calibra->SetMinEntries(20*20);  
818       if(fCompressPerDetector){
819         if(AddStatsPerDetector(histodriftvelocity,1)) {
820           pass = calibra->AnalysePH(fDetSumVector);
821         }
822       }
823       else pass = calibra->AnalysePH(histodriftvelocity);
824     }
825   }
826   else {
827     if(fVector2d && calibraVector) {
828       calibra->SetMinEntries(20*20);  
829       if(fCompressPerDetector){
830         if(!(vvect = calibraVector->AddStatsPerDetectorPH())) return kFALSE;
831         pass = calibra->AnalysePH(vvect);
832       }
833       else pass = calibra->AnalysePH(calibraVector);  
834     }
835   }
836
837   if(pass) {
838     Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(1))
839       + 6*  18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(1));
840     Int_t nbfit  = calibra->GetNumberFit();
841     Int_t nbE    = calibra->GetNumberEnt();
842     // enough statistics
843     if ((nbtg >                  0) && 
844         (nbfit        >= 0.001*nbE))
845       {
846         // create the cal objects
847         calibra->PutMeanValueOtherVectorFit(1,kTRUE);
848         calibra->PutMeanValueOtherVectorFit2(1,kTRUE); 
849         TObjArray object  = calibra->GetVectorFit();
850         AliTRDCalDet *objdriftvelocitydet = calibra->CreateDetObjectVdrift(&object,kTRUE);
851         TObject *objdriftvelocitypad      = calibra->CreatePadObjectVdrift();
852         TObjArray objectt          = calibra->GetVectorFit2();
853         AliTRDCalDet *objtime0det  = calibra->CreateDetObjectT0(&object,kTRUE);
854         TObject *objtime0pad       = calibra->CreatePadObjectT0();
855         // store
856         AliCDBId id1("TRD/Calib/ChamberVdrift",fRunNumber, AliCDBRunRange::Infinity()); 
857         storLoc->Put((TObject *)objdriftvelocitydet, id1, &mdDet); 
858         AliCDBId id2("TRD/Calib/LocalVdrift",fRunNumber, AliCDBRunRange::Infinity()); 
859         storLoc->Put((TObject *)objdriftvelocitypad, id2, &mdPad); 
860         //
861         AliCDBId idd1("TRD/Calib/ChamberT0",fRunNumber, AliCDBRunRange::Infinity()); 
862         storLoc->Put((TObject *)objtime0det, idd1, &mdDet); 
863         AliCDBId idd2("TRD/Calib/LocalT0",fRunNumber, AliCDBRunRange::Infinity()); 
864         storLoc->Put((TObject *)objtime0pad, idd2, &mdPad); 
865         // Make graph
866         TGraph *graph = 0x0;
867         if(FillGraphIndex(&object,graph)){ 
868           fGraph->AddAt(graph,1);
869         }
870         TGraph *graphh = 0x0;
871         if(FillGraphIndex(&objectt,graphh)){ 
872           fGraph->AddAt(graphh,2);
873         }
874       }//if(enough statistics)
875     calibra->ResetVectorFit();
876   }
877   else return kFALSE;
878   
879   //
880   // PRF
881   //
882   pass = kFALSE; 
883   if(fHisto2d) {
884     TProfile2D *histoprf = (TProfile2D *) fContainer->FindObject("PRF2d");
885     if (histoprf) {
886       histoprf->SetDirectory(0);  
887       calibra->SetMinEntries(600); 
888       if(fCompressPerDetector){
889         if(AddStatsPerDetector(histoprf,2)) pass = calibra->AnalysePRFMarianFit(fDetSumVector);
890       }
891       else pass = calibra->AnalysePRFMarianFit(histoprf);
892     }
893   }
894   else {
895     if(fVector2d && calibraVector) {
896       calibra->SetMinEntries(600);  
897       if(fCompressPerDetector){
898         if(!(vvect =calibraVector->AddStatsPerDetectorPRF())) return kFALSE;
899         pass = calibra->AnalysePRFMarianFit(vvect);
900       }
901       else pass = calibra->AnalysePRFMarianFit(calibraVector);  
902     }
903   }
904   
905   if(pass){
906     Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(2))
907       + 6*  18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(2));
908     Int_t nbfit        = calibra->GetNumberFit();
909     Int_t nbE          = calibra->GetNumberEnt();
910     // enough statistics
911     if ((nbtg >                  0) && 
912         (nbfit        >= 0.001*nbE)) {
913       TObjArray object            = calibra->GetVectorFit();
914       TObject *objPRFpad          = calibra->CreatePadObjectPRF(&object);
915       // store
916       AliCDBId id2("TRD/Calib/PRFWidth",fRunNumber, AliCDBRunRange::Infinity()); 
917       storLoc->Put((TObject *)objPRFpad, id2, &mdPad); 
918       // Make graph
919       TGraph *graph = 0x0;
920       if(FillGraphIndex(&object,graph)){ 
921         fGraph->AddAt(graph,5);
922       }
923     }
924     calibra->ResetVectorFit();
925   }
926   else return kFALSE;
927   
928   //
929   // VDRIFT linear fit 
930   //
931   AliTRDCalibraVdriftLinearFit *vlinearfit = (AliTRDCalibraVdriftLinearFit *) fContainer->FindObject("LinearVdriftFit"); 
932   if (vlinearfit) {
933     calibra->SetMinEntries(20*20);     
934     if(calibra->AnalyseLinearFitters(vlinearfit)) {
935       
936       Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(2))
937         + 6*  18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(2));
938       Int_t nbfit  = calibra->GetNumberFit();
939       Int_t nbE    = calibra->GetNumberEnt();
940       // enough statistics
941       if ((nbtg >                  0) && 
942           (nbfit        >= 0.001*nbE))
943         {
944           // create the cal objects
945           calibra->PutMeanValueOtherVectorFit(1,kTRUE);
946           calibra->PutMeanValueOtherVectorFit2(1,kTRUE); 
947           TObjArray object  = calibra->GetVectorFit();
948           AliTRDCalDet *objdriftvelocitydet = calibra->CreateDetObjectVdrift(&object,kTRUE);
949           TObject *objdriftvelocitypad      = calibra->CreatePadObjectVdrift();
950           TObjArray objectt          = calibra->GetVectorFit2();
951           AliTRDCalDet *objtime0det  = calibra->CreateDetObjectT0(&object,kTRUE);
952           TObject *objtime0pad       = calibra->CreatePadObjectT0();
953           // store dummy
954           AliCDBId id1("TRD/Calib/ChamberVdriftLinear",fRunNumber, AliCDBRunRange::Infinity()); 
955           storLoc->Put((TObject *)objdriftvelocitydet, id1, &mdDet); 
956           AliCDBId id2("TRD/Calib/LocalVdriftLinear",fRunNumber, AliCDBRunRange::Infinity()); 
957           storLoc->Put((TObject *)objdriftvelocitypad, id2, &mdPad); 
958           //
959           AliCDBId idd1("TRD/Calib/ChamberLorentzAngle",fRunNumber, AliCDBRunRange::Infinity()); 
960           storLoc->Put((TObject *)objtime0det, idd1, &mdDet); 
961           AliCDBId idd2("TRD/Calib/LocalLorentzAngle",fRunNumber, AliCDBRunRange::Infinity()); 
962           storLoc->Put((TObject *)objtime0pad, idd2, &mdPad); 
963           // Make graph
964           TGraph *graph = 0x0;
965           if(FillGraphIndex(&object,graph)){ 
966             fGraph->AddAt(graph,3);
967           }
968           TGraph *graphh = 0x0;
969           if(FillGraphIndex(&objectt,graphh)){ 
970             fGraph->AddAt(graphh,4);
971           }
972         }//if(enough statistics)
973     }// if fit
974     calibra->ResetVectorFit();
975   }
976   else return kFALSE;
977   
978   fPostProcess = kTRUE;
979   
980   return kTRUE;
981   
982 }
983
984 //________________________________________________________________________
985 Bool_t AliTRDcalibration::FillGraphIndex(const TObjArray *vectora,TGraph *graph) const
986 {
987   //
988   // Fill one value (the first one) per detector
989   //
990
991   Int_t loop = (Int_t) vectora->GetEntriesFast();
992   if(loop != 540) {
993     AliInfo("The Vector Fit is not complete!");
994     return kFALSE;
995   }
996   
997   Double_t x[540];
998   Double_t y[540];
999   for (Int_t k = 0; k < loop; k++) {
1000     if(!vectora->At(k)){
1001       x[k] = -1.0;
1002       y[k] = -1.0;
1003       continue;
1004     }
1005     x[k]  = ((AliTRDCalibraFit::AliTRDFitInfo *) vectora->At(k))->GetDetector();
1006     y[k]  = ((Float_t *)((AliTRDCalibraFit::AliTRDFitInfo *) vectora->At(k))->GetCoef())[0];
1007   }
1008
1009   if(!graph) graph = new TGraph(540,&x[0],&y[0]);
1010   else{ 
1011     graph->~TGraph();
1012     new(graph) TGraph(540,&x[0],&y[0]);
1013   }
1014
1015   return kTRUE;
1016
1017 }
1018 //________________________________________________________________________
1019 Bool_t AliTRDcalibration::AddStatsPerDetector(const TH2I *ch) 
1020 {
1021   //
1022   // Add statistic per detector
1023   //
1024   
1025   AliTRDCalibraMode calibMode = AliTRDCalibraMode();
1026   const char *name = ch->GetTitle();
1027   if(!SetNzFromTObject(name,0,&calibMode)) return 0x0;
1028   if(!SetNrphiFromTObject(name,0,&calibMode)) return 0x0;
1029   if(((calibMode.GetNz(0) == 100) && (calibMode.GetNrphi(0) == 100)) || ((calibMode.GetNz(0) == 10) && (calibMode.GetNrphi(0) == 10))) return kFALSE;
1030
1031   Int_t    nybins  = ch->GetNbinsY();// groups number
1032   Int_t    nxbins  = ch->GetNbinsX();// number of bins X
1033   TAxis   *xaxis   = ch->GetXaxis();
1034   Double_t lowedge  = xaxis->GetBinLowEdge(1);
1035   Double_t upedge   = xaxis->GetBinUpEdge(nxbins);
1036
1037   // number per chamber 2
1038   calibMode.ModePadCalibration(2,0);
1039   calibMode.ModePadFragmentation(0,2,0,0);
1040   calibMode.SetDetChamb2(0);
1041   Int_t perChamber2 = (Int_t) calibMode.GetDetChamb2(0);
1042
1043   // number per other chamber
1044   calibMode.ModePadCalibration(0,0);
1045   calibMode.ModePadFragmentation(0,0,0,0);
1046   calibMode.SetDetChamb0(0);
1047   Int_t perChamber0 = (Int_t) calibMode.GetDetChamb0(0);
1048
1049   if(nybins != (6*18*perChamber2+6*4*18*perChamber0)) return kFALSE;
1050
1051   // Create Histo
1052   TString nname((const char *)ch->GetName());
1053   nname  += "PerDetector";
1054   TString title("Nz");
1055   title += 0;
1056   title += "Nrphi";
1057   title += 0;
1058   if(!fCHSum) fCHSum = new TH2I((const char *)nname,(const char *)title
1059                                 ,nxbins,lowedge,upedge,540,0,540);
1060   else{ 
1061     fCHSum->~TH2I();
1062     new(fCHSum) TH2I((const Char_t *) nname,(const char *)title
1063                      ,nxbins,lowedge,upedge,540,0,540);
1064   }
1065   fCHSum->SetYTitle("Detector number");
1066   fCHSum->SetXTitle("charge deposit [a.u]");
1067   fCHSum->SetZTitle("counts");
1068   fCHSum->SetStats(0);
1069   fCHSum->Sumw2();
1070
1071   Int_t counter = 0;
1072   
1073   for(Int_t det = 0; det < 540; det++){
1074
1075     Int_t numberofgroup = 0;
1076     if(AliTRDgeometry::GetStack(det) == 2) numberofgroup = perChamber2;
1077     else numberofgroup = perChamber0;
1078     TH1I *projch = (TH1I *) ch->ProjectionX("projch",counter+1,counter+numberofgroup,(Option_t *)"e");
1079     projch->SetDirectory(0);
1080        
1081     for(Int_t nx = 0; nx <= nxbins; nx++) {
1082       fCHSum->SetBinContent(nx,det+1,projch->GetBinContent(nx));
1083       fCHSum->SetBinError(nx,det+1,projch->GetBinError(nx));
1084     }
1085
1086     counter += numberofgroup;
1087     
1088     delete projch;
1089
1090   }
1091
1092   return kTRUE;
1093
1094 }
1095 //_____________________________________________________________________________________________________________________
1096 Bool_t AliTRDcalibration::AddStatsPerDetector(const TProfile2D *ph,Int_t i)
1097 {
1098   //
1099   // Add statistic per detector
1100   //
1101
1102   AliTRDCalibraMode calibMode = AliTRDCalibraMode();
1103   const char *name = ph->GetTitle();
1104   //printf("name %s\n",name);
1105   if(!SetNzFromTObject(name,0,&calibMode)) return kFALSE;
1106   if(!SetNrphiFromTObject(name,0,&calibMode)) return kFALSE;
1107   if(((calibMode.GetNz(0) == 100) && (calibMode.GetNrphi(0) == 100)) || ((calibMode.GetNz(0) == 10) && (calibMode.GetNrphi(0) == 10))) return kFALSE;
1108   //printf("Found mode Mz %d, Nrphi %d\n",calibMode.GetNz(0),calibMode.GetNrphi(0));  
1109
1110
1111   Int_t    nybins  = ph->GetNbinsY();// groups number
1112   Int_t    nxbins  = ph->GetNbinsX();// number of bins X
1113   TAxis   *xaxis = ph->GetXaxis();
1114   Double_t lowedge  = xaxis->GetBinLowEdge(1);
1115   Double_t upedge   = xaxis->GetBinUpEdge(nxbins);
1116
1117   // number per chamber 2
1118   calibMode.ModePadCalibration(2,0);
1119   calibMode.ModePadFragmentation(0,2,0,0);
1120   calibMode.SetDetChamb2(0);
1121   Int_t perChamber2 = (Int_t) calibMode.GetDetChamb2(0);
1122
1123   // number per other chamber
1124   calibMode.ModePadCalibration(0,0);
1125   calibMode.ModePadFragmentation(0,0,0,0);
1126   calibMode.SetDetChamb0(0);
1127   Int_t perChamber0 = (Int_t) calibMode.GetDetChamb0(0);
1128
1129   if(nybins != (6*18*perChamber2+6*4*18*perChamber0)) return kFALSE;
1130   
1131   // Create calvector 
1132   TString nbname((const char *)ph->GetName());
1133   nbname  += "PerDetectorVector";
1134   if(!fDetSumVector) fDetSumVector = new AliTRDCalibraVector();
1135   else{ 
1136     fDetSumVector->~AliTRDCalibraVector();
1137     new(fDetSumVector) AliTRDCalibraVector();
1138   }
1139   if(i==1){
1140     fDetSumVector->SetTimeMax(nxbins);
1141   }
1142   if(i==2){
1143     fDetSumVector->SetNumberBinPRF(nxbins);
1144     fDetSumVector->SetPRFRange(TMath::Abs(lowedge));
1145   }
1146   fDetSumVector->SetDetCha0(i,1);
1147   fDetSumVector->SetDetCha2(i,1);
1148   fDetSumVector->SetNzNrphi(i,0,0);
1149   if(i==2) {
1150     Int_t nbg = GetNumberOfGroupsPRF((const char *)name);
1151     fDetSumVector->SetNbGroupPRF(nbg);
1152   }
1153
1154   // Create Histo
1155   TString nname((const char *)ph->GetName());
1156   nname  += "PerDetector";
1157   TString title("Nz");
1158   title += 0;
1159   title += "Nrphi";
1160   title += 0;
1161   if(!fDetSum) fDetSum = new TH2D((const char *)nname,(const Char_t *) title
1162                                 ,nxbins,lowedge,upedge,540,0,540);
1163   else{ 
1164     fDetSum->~TH2D();
1165     new(fDetSum) TH2D((const Char_t *) nname,(const Char_t *) title
1166                      ,nxbins,lowedge,upedge,540,0,540);
1167   }
1168   fDetSum->SetYTitle("Detector number");
1169   fDetSum->SetXTitle(xaxis->GetTitle());
1170   fDetSum->SetStats(0);
1171   
1172   Int_t counter = 0;
1173
1174   for(Int_t det = 0; det < 540; det++){
1175
1176     Int_t numberofgroup = 0;
1177     if(AliTRDgeometry::GetStack(det) == 2) numberofgroup = perChamber2;
1178     else numberofgroup = perChamber0;
1179     
1180     for(Int_t nx = 1; nx <= nxbins; nx++) {
1181       
1182       Double_t entries = 0.0;
1183       Double_t sumw2 = 0.0;
1184       Double_t sumw = 0.0;
1185
1186       for(Int_t k = counter+1; k <= (counter+numberofgroup); k++){
1187         Int_t  binnumber = ph->GetBin(nx,k);
1188         entries += ph->GetBinEntries(binnumber);
1189         sumw2 += (ph->GetBinError(binnumber)*ph->GetBinError(binnumber)+ph->GetBinContent(binnumber)*ph->GetBinContent(binnumber))*ph->GetBinEntries(binnumber);
1190         sumw += ph->GetBinContent(binnumber)*ph->GetBinEntries(binnumber);
1191       }
1192
1193       Double_t mean = 0.0;
1194       if(entries > 0.0) mean = sumw/entries;
1195       Double_t squaremean = 0.0;
1196       if(entries > 0.0) squaremean = sumw2/entries;
1197       Double_t errorf = squaremean - mean*mean;
1198       Double_t error = 0.0;
1199       if(entries > 0.0) error = TMath::Sqrt(TMath::Abs(errorf)/entries);
1200       
1201       fDetSum->SetBinContent(nx,det+1,mean);
1202       fDetSum->SetBinError(nx,det+1,error);
1203
1204       if(i==1) fDetSumVector->FillVectorPH(det,0,nx-1,(Int_t)entries,(Float_t)mean,(Float_t)squaremean);
1205       if(i==2) fDetSumVector->FillVectorPRF(det,0,nx-1,(Int_t)entries,(Float_t)mean,(Float_t)squaremean);
1206       
1207     }
1208     
1209     counter += numberofgroup;
1210
1211   }
1212
1213   return kTRUE;
1214
1215   
1216 }
1217 //_____________________________________________________________________________
1218 Bool_t AliTRDcalibration::SetNrphiFromTObject(const char *name, Int_t i, AliTRDCalibraMode *calibMode) const
1219 {
1220   //
1221   // Set the granularity from object
1222   //  
1223   
1224   const Char_t *patternrphi0 = "Nrphi0";
1225   const Char_t *patternrphi1 = "Nrphi1";
1226   const Char_t *patternrphi2 = "Nrphi2";
1227   const Char_t *patternrphi3 = "Nrphi3";
1228   const Char_t *patternrphi4 = "Nrphi4";
1229   const Char_t *patternrphi5 = "Nrphi5";
1230   const Char_t *patternrphi6 = "Nrphi6";
1231
1232   
1233   const Char_t *patternrphi10 = "Nrphi10";
1234   const Char_t *patternrphi100 = "Nrphi100";
1235   const Char_t *patternz10 = "Nz10";
1236   const Char_t *patternz100 = "Nz100";
1237
1238   // Nrphi mode
1239   if ((strstr(name,patternrphi100)) && (strstr(name,patternz100))) {
1240     calibMode->SetAllTogether(i);
1241     return kTRUE;
1242   }
1243   if ((strstr(name,patternrphi10)) && (strstr(name,patternz10))) {
1244     calibMode->SetPerSuperModule(i);
1245     return kTRUE;
1246   }
1247   
1248   if (strstr(name,patternrphi0)) {
1249     calibMode->SetNrphi(i ,0);
1250     return kTRUE;
1251   }
1252   if (strstr(name,patternrphi1)) {
1253     calibMode->SetNrphi(i, 1);
1254     return kTRUE;
1255   }
1256   if (strstr(name,patternrphi2)) {
1257     calibMode->SetNrphi(i, 2);
1258     return kTRUE;
1259   }
1260   if (strstr(name,patternrphi3)) {
1261     calibMode->SetNrphi(i, 3);
1262     return kTRUE;
1263   }
1264   if (strstr(name,patternrphi4)) {
1265     calibMode->SetNrphi(i, 4);
1266     return kTRUE;
1267   }
1268   if (strstr(name,patternrphi5)) {
1269     calibMode->SetNrphi(i, 5);
1270     return kTRUE;
1271   }
1272   if (strstr(name,patternrphi6)) {
1273     calibMode->SetNrphi(i, 6);
1274     return kTRUE;
1275   }
1276   
1277   calibMode->SetNrphi(i ,0);
1278   return kFALSE;
1279   
1280 }
1281 //_____________________________________________________________________________
1282 Bool_t AliTRDcalibration::SetNzFromTObject(const char *name, Int_t i, AliTRDCalibraMode *calibMode) const
1283 {
1284   //
1285   // Set fNz[i] of the AliTRDCalibraFit::Instance()
1286   // corresponding to the given TObject
1287   //
1288
1289   // Some patterns
1290   const Char_t *patternz0    = "Nz0";
1291   const Char_t *patternz1    = "Nz1";
1292   const Char_t *patternz2    = "Nz2";
1293   const Char_t *patternz3    = "Nz3";
1294   const Char_t *patternz4    = "Nz4";
1295
1296   const Char_t *patternrphi10 = "Nrphi10";
1297   const Char_t *patternrphi100 = "Nrphi100";
1298   const Char_t *patternz10 = "Nz10";
1299   const Char_t *patternz100 = "Nz100";
1300
1301   if ((strstr(name,patternrphi100)) && (strstr(name,patternz100))) {
1302     calibMode->SetAllTogether(i);
1303     return kTRUE;
1304   }
1305   if ((strstr(name,patternrphi10)) && (strstr(name,patternz10))) {
1306     calibMode->SetPerSuperModule(i);
1307     return kTRUE;
1308   }
1309   if (strstr(name,patternz0)) {
1310     calibMode->SetNz(i, 0);
1311     return kTRUE;
1312   }
1313   if (strstr(name,patternz1)) {
1314     calibMode->SetNz(i ,1);
1315     return kTRUE;
1316   }
1317   if (strstr(name,patternz2)) {
1318     calibMode->SetNz(i ,2);
1319     return kTRUE;
1320   }
1321   if (strstr(name,patternz3)) {
1322     calibMode->SetNz(i ,3);
1323     return kTRUE;  
1324   }
1325   if (strstr(name,patternz4)) {
1326     calibMode->SetNz(i ,4);
1327     return kTRUE;
1328   }
1329  
1330   calibMode->SetNz(i ,0);
1331   return kFALSE;
1332 }
1333 //____________________________________________________________________________________________________
1334 Int_t AliTRDcalibration::GetNumberOfGroupsPRF(const char* nametitle) const
1335 {
1336   //
1337   // Get numberofgroupsprf
1338   //
1339   
1340   // Some patterns
1341   const Char_t *pattern0 = "Ngp0";
1342   const Char_t *pattern1 = "Ngp1";
1343   const Char_t *pattern2 = "Ngp2";
1344   const Char_t *pattern3 = "Ngp3";
1345   const Char_t *pattern4 = "Ngp4";
1346   const Char_t *pattern5 = "Ngp5";
1347   const Char_t *pattern6 = "Ngp6";
1348
1349   // Nrphi mode
1350   if (strstr(nametitle,pattern0)) {
1351     return 0;
1352   }
1353   if (strstr(nametitle,pattern1)) {
1354     return 1;
1355   }
1356   if (strstr(nametitle,pattern2)) {
1357     return 2;
1358   }
1359   if (strstr(nametitle,pattern3)) {
1360     return 3;
1361   }
1362   if (strstr(nametitle,pattern4)) {
1363     return 4;
1364   }
1365   if (strstr(nametitle,pattern5)) {
1366     return 5;
1367   }
1368   if (strstr(nametitle,pattern6)){
1369     return 6;
1370   }
1371   else return -1;
1372  
1373
1374 }
1375