]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCPreprocessorOffline.cxx
macro to run on alien
[u/mrichter/AliRoot.git] / TPC / AliTPCPreprocessorOffline.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
17
18 /*
19   Responsible: marian.ivanov@cern.ch 
20   Code to analyze the TPC calibration and to produce OCDB entries  
21
22
23    .x ~/rootlogon.C
24    gSystem->Load("libANALYSIS");
25    gSystem->Load("libTPCcalib");
26
27    AliTPCPreprocessorOffline proces;
28    TString ocdbPath="local:////"
29    ocdbPath+=gSystem->GetFromPipe("pwd");
30
31    proces.CalibTimeGain("CalibObjects.root",run0,run1,ocdbPath);
32    proces.CalibTimeVdrift("CalibObjects.root",run0,run1,ocdbPath);
33   // take the raw calibration data from the file CalibObjects.root 
34   // and make a OCDB entry with run  validity run0-run1
35   // results are stored at the ocdbPath - local or alien ...
36   // default storage ""- data stored at current working directory 
37  
38 */
39 #include "Riostream.h"
40 #include <fstream>
41 #include "TMap.h"
42 #include "TGraphErrors.h"
43 #include "AliExternalTrackParam.h"
44 #include "TFile.h"
45 #include "TGraph.h"
46 #include "TMultiGraph.h"
47 #include "TCanvas.h"
48 #include "THnSparse.h"
49 #include "TLegend.h"
50 #include "TPad.h"
51 #include "TH2D.h"
52 #include "TH3D.h"
53 #include "AliTPCROC.h"
54 #include "AliTPCCalROC.h"
55 #include "AliESDfriend.h"
56 #include "AliTPCcalibTime.h"
57 #include "AliSplineFit.h"
58 #include "AliCDBMetaData.h"
59 #include "AliCDBId.h"
60 #include "AliCDBManager.h"
61 #include "AliCDBStorage.h"
62 #include "AliTPCcalibBase.h"
63 #include "AliTPCcalibDB.h"
64 #include "AliTPCcalibDButil.h"
65 #include "AliRelAlignerKalman.h"
66 #include "AliTPCParamSR.h"
67 #include "AliTPCcalibTimeGain.h"
68 #include "AliSplineFit.h"
69 #include "AliTPCPreprocessorOffline.h"
70
71
72 ClassImp(AliTPCPreprocessorOffline)
73
74 AliTPCPreprocessorOffline::AliTPCPreprocessorOffline():
75   TNamed("TPCPreprocessorOffline","TPCPreprocessorOffline"),
76   fMinEntries(500),                      // minimal number of entries for fit
77   startRun(0),                         // start Run - used to make fast selection in THnSparse
78   endRun(0),                           // end   Run - used to make fast selection in THnSparse
79   startTime(0),                        // startTime - used to make fast selection in THnSparse
80   endTime(0),                          // endTime   - used to make fast selection in THnSparse
81   ocdbStorage(""),                   // path to the OCDB storage
82   fVdriftArray(new TObjArray),
83   fTimeDrift(0),
84   fGraphMIP(0),                // graph time dependence of MIP
85   fGraphCosmic(0),             // graph time dependence at Plateu
86   fGraphAttachmentMIP(0),
87   fFitMIP(0),                  // fit of dependence - MIP
88   fFitCosmic(0),               // fit of dependence - Plateu
89   fGainArray(new TObjArray),               // array to be stored in the OCDB
90   fGainMIP(0),          // calibration component for MIP
91   fGainCosmic(0)       // calibration component for cosmic
92 {
93   //
94   // default constructor
95   //
96 }
97
98 AliTPCPreprocessorOffline::~AliTPCPreprocessorOffline() {
99   //
100   // Destructor
101   //
102 }
103
104
105
106
107 void AliTPCPreprocessorOffline::GetRunRange(AliTPCcalibTime * const  timeDrift){
108   //
109   // find the fist and last run
110   //
111   TObjArray *hisArray =timeDrift->GetHistoDrift();
112   {for (Int_t i=0; i<hisArray->GetEntriesFast(); i++){
113     THnSparse* addHist=(THnSparse*)hisArray->UncheckedAt(i);
114     if (addHist->GetEntries()<fMinEntries) continue;
115     if (!addHist) continue;
116     TH1D* histo    =addHist->Projection(3);
117     TH1D* histoTime=addHist->Projection(0);
118     printf("%s\t%f\t%d\t%d\n",histo->GetName(), histo->GetEntries(),histo->FindFirstBinAbove(0),histo->FindLastBinAbove(0));
119
120     if (startRun<=0){ 
121       startRun=histo->FindFirstBinAbove(0);
122       endRun  =histo->FindLastBinAbove(0);
123     }else{
124       startRun=TMath::Min(histo->FindFirstBinAbove(0),startRun);
125       endRun  =TMath::Max(histo->FindLastBinAbove(0),endRun);
126     }
127     if (startTime==0){ 
128       startTime=histoTime->FindFirstBinAbove(0);
129       endTime  =histoTime->FindLastBinAbove(0);
130     }else{
131       startTime=TMath::Min(histoTime->FindFirstBinAbove(0),startTime);
132       endTime  =TMath::Max(histoTime->FindLastBinAbove(0),endTime);
133     }
134     delete histo;
135     delete histoTime;
136   }}
137   if (startRun<0) startRun=0;
138   if (endRun<0) endRun=100000000;
139   printf("Run range  :\t%d-%d\n", startRun, endRun);
140   printf("Time range :\t%d-%d\n", startTime, endTime);
141
142 }
143
144
145
146 void AliTPCPreprocessorOffline::CalibTimeVdrift(const Char_t* file, Int_t ustartRun, Int_t uendRun, TString pocdbStorage){
147   //
148   // make calibration of the drift velocity
149   // Input parameters:
150   //      file                   - the location of input file
151   //      ustartRun, uendrun     - run validity period 
152   //      pocdbStorage           - path to hte OCDB storage
153   //                             - if empty - local storage 'pwd' uesed
154   if (pocdbStorage.Length()>0) ocdbStorage=pocdbStorage;
155   else
156   ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
157   //
158   // 1. Initialization and run range setting
159   TFile fcalib(file);
160   fTimeDrift=(AliTPCcalibTime*)fcalib.Get("calibTime");
161   startRun=ustartRun;
162   endRun=ustartRun; 
163   TObjArray *hisArray =fTimeDrift->GetHistoDrift();  
164   GetRunRange(fTimeDrift);
165   for (Int_t i=0; i<hisArray->GetEntriesFast(); i++){
166     THnSparse* addHist=(THnSparse*)hisArray->At(i);
167     if (!addHist) continue;
168     if (startTime<endTime) addHist->GetAxis(0)->SetRange(startTime-1,endTime+1);
169     if (startRun<endRun) addHist->GetAxis(3)->SetRange(startRun-1,endRun+1);
170   }
171   //
172   //
173   // 2. extraction of the information
174   //
175   fVdriftArray = new TObjArray();
176   AddAlignmentGraphs(fVdriftArray,fTimeDrift);
177   AddHistoGraphs(fVdriftArray,fTimeDrift,fMinEntries);
178   AddLaserGraphs(fVdriftArray,fTimeDrift);
179   //
180   // 3. Append QA plots
181   //
182   MakeDefaultPlots(fVdriftArray,fVdriftArray);
183   //
184   //
185   // 4. update of OCDB
186   //
187   //
188   
189   UpdateOCDBDrift(ustartRun,uendRun,ocdbStorage);
190 }
191
192 void AliTPCPreprocessorOffline::UpdateOCDBDrift( Int_t ustartRun, Int_t uendRun,  const char* storagePath ){
193   //
194   // Update OCDB 
195   //
196   AliCDBMetaData *metaData= new AliCDBMetaData();
197   metaData->SetObjectClassName("TObjArray");
198   metaData->SetResponsible("Marian Ivanov");
199   metaData->SetBeamPeriod(1);
200   metaData->SetAliRootVersion("05-25-01"); //root version
201   metaData->SetComment("Calibration of the time dependence of the drift velocity");
202   AliCDBId* id1=NULL;
203   id1=new AliCDBId("TPC/Calib/TimeDrift", ustartRun, uendRun);
204   AliCDBStorage* gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
205   gStorage->Put(fVdriftArray, (*id1), metaData);
206 }
207
208
209
210 void AliTPCPreprocessorOffline::UpdateDriftParam(AliTPCParam *param, TObjArray *const arr, Int_t lstartRun){
211   //
212   //  update the OCDB entry for the nominal time0
213   //
214   //
215   //  AliTPCParam * param = AliTPCcalibDB::Instance()->GetParameters();
216   AliTPCParam *paramNew = (AliTPCParam *)param->Clone();
217   TGraphErrors *grT =  (TGraphErrors *)arr->FindObject("ALIGN_ITSM_TPC_T0");
218   Double_t deltaTcm = TMath::Median(grT->GetN(),grT->GetY());
219   Double_t deltaT   = deltaTcm/param->GetDriftV();
220   paramNew->SetL1Delay(param->GetL1Delay()-deltaT);
221   paramNew->Update();
222
223   AliCDBMetaData *metaData= new AliCDBMetaData();
224   metaData->SetObjectClassName("TObjArray");
225   metaData->SetResponsible("Marian Ivanov");
226   metaData->SetBeamPeriod(1);
227   metaData->SetAliRootVersion("05-25-02"); //root version
228   metaData->SetComment("Updated calibration of nominal time 0");
229   AliCDBId* id1=NULL;
230   id1=new AliCDBId("TPC/Calib/Parameters", lstartRun, AliCDBRunRange::Infinity());
231   AliCDBStorage* gStorage = AliCDBManager::Instance()->GetStorage(ocdbStorage);
232   gStorage->Put(param, (*id1), metaData);
233
234 }
235
236
237 void AliTPCPreprocessorOffline::PrintArray(TObjArray *array){
238   //
239   // Print the names of the entries in array
240   //
241   Int_t entries = array->GetEntries();
242   for (Int_t i=0; i<entries; i++){
243     if (!array->At(i)) continue;
244     printf("%d\t %s\n", i,  array->At(i)->GetName());
245   }
246 }
247
248
249
250 TGraphErrors* AliTPCPreprocessorOffline::FilterGraphDrift(TGraphErrors * graph, Float_t errSigmaCut, Float_t medianCutAbs){
251   // 2 filters:
252   //    1. filter graph - error cut errSigmaCut
253   //    2. filter graph - medianCutAbs around median
254   //
255   // errSigmaCut   - cut on error
256   // medianCutAbs  - cut on value around median
257   Double_t dummy=0;               //   
258   //
259   // 1. filter graph - error cut errSigmaCut
260   //              
261   TGraphErrors *graphF; 
262   graphF = AliTPCcalibDButil::FilterGraphMedianErr(graph,errSigmaCut,dummy);
263   delete graph;
264   if (!graphF) return 0;
265   graph = AliTPCcalibDButil::FilterGraphMedianErr(graphF,errSigmaCut,dummy);
266   delete graphF;
267   if (!graph) return 0;
268   //
269   // filter graph - kMedianCutAbs around median
270   // 
271   graphF=FilterGraphMedianAbs(graph, medianCutAbs,dummy);
272   delete graph;
273   if (!graphF) return 0;
274   graph=FilterGraphMedianAbs(graphF, medianCutAbs,dummy);
275   delete graphF;
276   if (!graph) return 0;
277   return graph;
278 }
279
280
281
282 TGraphErrors* AliTPCPreprocessorOffline::FilterGraphMedianAbs(TGraphErrors * graph, Float_t cut,Double_t &medianY){
283   //
284   // filter outlyer measurement
285   // Only points around median +- cut filtered 
286   //
287   if (!graph) return  0;
288   Int_t kMinPoints=2;
289   Int_t npoints0 = graph->GetN();
290   Int_t npoints=0;
291   Float_t  rmsY=0;
292   Double_t *outx=new Double_t[npoints0];
293   Double_t *outy=new Double_t[npoints0];
294   Double_t *errx=new Double_t[npoints0];
295   Double_t *erry=new Double_t[npoints0];
296   //
297   //
298   if (npoints0<kMinPoints) return 0;
299   for (Int_t iter=0; iter<3; iter++){
300     npoints=0;
301     for (Int_t ipoint=0; ipoint<npoints0; ipoint++){
302       if (graph->GetY()[ipoint]==0) continue;
303       if (iter>0 &&TMath::Abs(graph->GetY()[ipoint]-medianY)>cut) continue;  
304       outx[npoints]  = graph->GetX()[ipoint];
305       outy[npoints]  = graph->GetY()[ipoint];
306       errx[npoints]  = graph->GetErrorX(ipoint);
307       erry[npoints]  = graph->GetErrorY(ipoint);
308       npoints++;
309     }
310     if (npoints<=1) break;
311     medianY  =TMath::Median(npoints,outy);
312     rmsY   =TMath::RMS(npoints,outy);
313   }
314   TGraphErrors *graphOut=0;
315   if (npoints>1) graphOut= new TGraphErrors(npoints,outx,outy,errx,erry); 
316   return graphOut;
317 }
318
319
320 void AliTPCPreprocessorOffline::AddHistoGraphs(  TObjArray * vdriftArray, AliTPCcalibTime * const timeDrift, Int_t minEntries){
321   //
322   // Add graphs corresponding to the alignment
323   //
324   const Double_t kErrSigmaCut=5;      // error sigma cut - for filtering
325   const Double_t kMedianCutAbs=0.03;  // error sigma cut - for filtering
326   //
327   TObjArray * array=timeDrift->GetHistoDrift();
328   if (array){
329     THnSparse* hist=NULL;
330     // 2.a) cosmics with different triggers
331     for (Int_t i=0; i<array->GetEntriesFast();i++){
332       hist=(THnSparseF*)array->UncheckedAt(i);
333       if(!hist) continue;
334       if (hist->GetEntries()<minEntries) continue;
335       //hist->Print();
336       TString name=hist->GetName();
337       Int_t dim[4]={0,1,2,3};
338       THnSparse* newHist=hist->Projection(4,dim);
339       newHist->SetName(name);
340       TGraphErrors* graph=AliTPCcalibBase::FitSlices(newHist,2,0,400,100,0.05,0.95, kTRUE);
341       printf("name=%s graph=%i, N=%i\n", name.Data(), graph==0, graph->GetN());
342       Int_t pos=name.Index("_");
343       name=name(pos,name.Capacity()-pos);
344       TString graphName=graph->ClassName();
345       graphName+=name;
346       graphName.ToUpper();
347       //
348       graph = FilterGraphDrift(graph, kErrSigmaCut, kMedianCutAbs);
349       if (!graph) {
350         printf("Graph =%s filtered out\n", name.Data());
351         continue;
352       }
353       //
354       graph->SetMarkerStyle(i%8+20);
355       graph->SetMarkerColor(i%7);
356       graph->GetXaxis()->SetTitle("Time");
357       graph->GetYaxis()->SetTitle("v_{dcor}");
358       graph->SetName(graphName);
359       graph->SetTitle(graphName);
360       printf("Graph %d\t=\t%s\n", i, graphName.Data());
361       vdriftArray->Add(graph);
362     }
363   }
364 }
365
366
367
368
369 void AliTPCPreprocessorOffline::AddAlignmentGraphs(  TObjArray * vdriftArray, AliTPCcalibTime *const timeDrift){
370   //
371   // Add graphs corresponding to alignment to the object array
372   //
373   TObjArray *arrayITS=0;
374   TObjArray *arrayTOF=0;
375   TObjArray *arrayTRD=0;
376   TMatrixD *mstatITS=0;
377   TMatrixD *mstatTOF=0;
378   TMatrixD *mstatTRD=0;
379   //
380   arrayITS=timeDrift->GetAlignITSTPC();
381   arrayTRD=timeDrift->GetAlignTRDTPC();
382   arrayTOF=timeDrift->GetAlignTOFTPC();
383
384   if (arrayITS->GetEntries()>0) mstatITS= AliTPCcalibDButil::MakeStatRelKalman(arrayITS,0.9,50,0.025);
385   if (arrayTOF->GetEntries()>0) mstatTOF= AliTPCcalibDButil::MakeStatRelKalman(arrayTOF,0.9,1000,0.025);
386   if (arrayTRD->GetEntries()>0) mstatTRD= AliTPCcalibDButil::MakeStatRelKalman(arrayTRD,0.9,50,0.025);
387   //
388   TObjArray * arrayITSP= AliTPCcalibDButil::SmoothRelKalman(arrayITS,*mstatITS, 0, 5.);
389   TObjArray * arrayITSM= AliTPCcalibDButil::SmoothRelKalman(arrayITS,*mstatITS, 1, 5.);
390   TObjArray * arrayITSB= AliTPCcalibDButil::SmoothRelKalman(arrayITSP,arrayITSM);
391   TObjArray * arrayTOFP= AliTPCcalibDButil::SmoothRelKalman(arrayTOF,*mstatTOF, 0, 5.);
392   TObjArray * arrayTOFM= AliTPCcalibDButil::SmoothRelKalman(arrayTOF,*mstatTOF, 1, 5.);
393   TObjArray * arrayTOFB= AliTPCcalibDButil::SmoothRelKalman(arrayTOFP,arrayTOFM);
394
395   TObjArray * arrayTRDP= 0x0;
396   TObjArray * arrayTRDM= 0x0;
397   TObjArray * arrayTRDB= 0x0;
398   arrayTRDP= AliTPCcalibDButil::SmoothRelKalman(arrayTRD,*mstatTRD, 0, 5.);
399   arrayTRDM= AliTPCcalibDButil::SmoothRelKalman(arrayTRD,*mstatTRD, 1, 5.);
400   arrayTRDB= AliTPCcalibDButil::SmoothRelKalman(arrayTRDP,arrayTRDM);
401   //
402   //
403   Int_t entries=TMath::Max(arrayITS->GetEntriesFast(),arrayTOF->GetEntriesFast());
404   TObjArray *arrays[12]={arrayITS, arrayITSP, arrayITSM, arrayITSB,
405                          arrayTRD, arrayTRDP, arrayTRDM, arrayTRDB,
406                          arrayTOF, arrayTOFP, arrayTOFM, arrayTOFB};
407   TString   grnames[12]={"ALIGN_ITS", "ALIGN_ITSP", "ALIGN_ITSM", "ALIGN_ITSB",
408                          "ALIGN_TRD", "ALIGN_TRDP", "ALIGN_TRDM","ALIGN_TRDB",
409                          "ALIGN_TOF", "ALIGN_TOFP", "ALIGN_TOFM","ALIGN_TOFB"};
410   TString   grpar[9]={"DELTAPSI", "DELTATHETA", "DELTAPHI",
411                       "DELTAX", "DELTAY", "DELTAZ",
412                       "DRIFTVD", "T0", "VDGY"};
413
414   
415   TVectorD vX(entries);
416   TVectorD vY(entries);
417   TVectorD vEx(entries);
418   TVectorD vEy(entries);
419   TObjArray *arr=0;
420   for (Int_t iarray=0; iarray<12; iarray++){
421     arr = arrays[iarray];
422     if (arr==0) continue;
423     for (Int_t ipar=0; ipar<9; ipar++){      
424       Int_t counter=0;
425       for (Int_t itime=0; itime<arr->GetEntriesFast(); itime++){
426         AliRelAlignerKalman * kalman = (AliRelAlignerKalman *) arr->UncheckedAt(itime);
427         if (!kalman) continue;
428         vX[counter]=kalman->GetTimeStamp();
429         vY[counter]=(*(kalman->GetState()))[ipar];
430         if (ipar==6) vY[counter]=1./(*(kalman->GetState()))[ipar]-1;
431         vEx[counter]=0;
432         vEy[counter]=TMath::Sqrt((*(kalman->GetStateCov()))(ipar,ipar));
433         counter++;
434       }
435     
436       TGraphErrors * graph=new TGraphErrors(counter, vX.GetMatrixArray(),
437                                           vY.GetMatrixArray(),
438                                           vEx.GetMatrixArray(),
439                                           vEy.GetMatrixArray());
440       TString grName=grnames[iarray];
441       grName+="_TPC_";
442       grName+=grpar[ipar];
443       graph->SetName(grName.Data());
444       vdriftArray->AddLast(graph);
445     }
446   }  
447 }
448
449
450
451
452 void AliTPCPreprocessorOffline::AddLaserGraphs(  TObjArray * vdriftArray, AliTPCcalibTime *timeDrift){
453   //
454   // add graphs for laser
455   //
456   const Double_t delayL0L1 = 0.071;  //this is hack for 1/2 weeks
457   THnSparse *hisN=0;
458   TGraphErrors *grLaser[6]={0,0,0,0,0,0};
459   hisN = timeDrift->GetHistVdriftLaserA(0);
460   if (timeDrift->GetHistVdriftLaserA(0)){
461     grLaser[0]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(0),0,2,5,delayL0L1);
462     grLaser[0]->SetName("GRAPH_MEAN_DELAY_LASER_ALL_A");
463     vdriftArray->AddLast(grLaser[0]);
464   }    
465   if (timeDrift->GetHistVdriftLaserA(1)){
466     grLaser[1]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(1),0,2,5);
467     grLaser[1]->SetName("GRAPH_MEAN_DRIFT_LASER_ALL_A");
468     vdriftArray->AddLast(grLaser[1]);
469   }    
470   if (timeDrift->GetHistVdriftLaserA(2)){
471     grLaser[2]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(2),0,2,5);
472     grLaser[2]->SetName("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
473     vdriftArray->AddLast(grLaser[2]);
474   }    
475   if (timeDrift->GetHistVdriftLaserC(0)){
476     grLaser[3]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(0),0,2,5,delayL0L1);
477     grLaser[3]->SetName("GRAPH_MEAN_DELAY_LASER_ALL_C");
478     vdriftArray->AddLast(grLaser[3]);
479   }    
480   if (timeDrift->GetHistVdriftLaserC(1)){
481     grLaser[4]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(1),0,2,5);
482     grLaser[4]->SetName("GRAPH_MEAN_DRIFT_LASER_ALL_C");
483     vdriftArray->AddLast(grLaser[4]);
484   }    
485   if (timeDrift->GetHistVdriftLaserC(2)){
486     grLaser[5]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(2),0,2,5);
487     grLaser[5]->SetName("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");    
488     vdriftArray->AddLast(grLaser[5]);
489   }    
490   for (Int_t i=0; i<6;i++){
491     if (grLaser[i]) {
492       SetDefaultGraphDrift(grLaser[i], 1,(i+20));
493       grLaser[i]->GetYaxis()->SetTitle("Laser Correction");
494     }
495   }
496 }
497  
498  
499 TGraphErrors * AliTPCPreprocessorOffline::MakeGraphFilter0(THnSparse *hisN, Int_t itime, Int_t ival, Int_t minEntries, Double_t offset){
500   //
501   // Make graph with mean values and rms
502   //
503   hisN->GetAxis(itime)->SetRange(0,100000000);
504   hisN->GetAxis(ival)->SetRange(0,100000000);
505   TH1 * hisT      = hisN->Projection(itime);
506   TH1 * hisV      = hisN->Projection(ival);
507   //
508   Int_t firstBinA = hisT->FindFirstBinAbove(2);
509   Int_t lastBinA  = hisT->FindLastBinAbove(2);    
510   Int_t firstBinV = hisV->FindFirstBinAbove(0);
511   Int_t lastBinV  = hisV->FindLastBinAbove(0);    
512   hisN->GetAxis(itime)->SetRange(firstBinA,lastBinA);
513   hisN->GetAxis(ival)->SetRange(firstBinV,lastBinV);
514   Int_t entries=0;
515   for (Int_t ibin=firstBinA; ibin<lastBinA; ibin++){
516     Double_t cont = hisT->GetBinContent(ibin);
517     if (cont<minEntries) continue;
518     entries++;
519   }
520   TVectorD vecTime(entries);
521   TVectorD vecMean0(entries);
522   TVectorD vecRMS0(entries);
523   TVectorD vecMean1(entries);
524   TVectorD vecRMS1(entries);
525   entries=0;
526   {for (Int_t ibin=firstBinA; ibin<=lastBinA; ibin++){
527       Double_t cont = hisT->GetBinContent(ibin);
528       if (cont<minEntries) continue;
529       hisN->GetAxis(itime)->SetRange(ibin-1,ibin+1);
530       Double_t time = hisT->GetBinCenter(ibin);
531       TH1 * his = hisN->Projection(ival);
532       Double_t nentries0= his->GetBinContent(his->FindBin(0));
533       if (cont-nentries0<minEntries) continue;
534       //
535       his->SetBinContent(his->FindBin(0),0);
536       vecTime[entries]=time;
537       vecMean0[entries]=his->GetMean()+offset;
538       vecMean1[entries]=his->GetMeanError();
539       vecRMS0[entries] =his->GetRMS();
540       vecRMS1[entries] =his->GetRMSError();
541       delete his;  
542       entries++;
543     }}
544   delete hisT;
545   delete hisV;
546   TGraphErrors * graph =  new TGraphErrors(entries,vecTime.GetMatrixArray(), vecMean0.GetMatrixArray(),                                    0, vecMean1.GetMatrixArray());
547   return graph;
548 }
549
550
551
552
553
554
555
556
557 void AliTPCPreprocessorOffline::SetDefaultGraphDrift(TGraph *graph, Int_t color, Int_t style){
558   //
559   // Set default style for QA views
560   //
561   graph->GetXaxis()->SetTimeDisplay(kTRUE);
562   graph->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
563   graph->SetMaximum( 0.025);
564   graph->SetMinimum(-0.025);
565   graph->GetXaxis()->SetTitle("Time");
566   graph->GetYaxis()->SetTitle("v_{dcorr}");
567   //
568   graph->GetYaxis()->SetLabelSize(0.03);
569   graph->GetXaxis()->SetLabelSize(0.03);
570   //
571   graph->GetXaxis()->SetNdivisions(10,5,0);
572   graph->GetYaxis()->SetNdivisions(10,5,0);
573   //
574   graph->GetXaxis()->SetLabelOffset(0.02);
575   graph->GetYaxis()->SetLabelOffset(0.005);
576   //
577   graph->GetXaxis()->SetTitleOffset(1.3);
578   graph->GetYaxis()->SetTitleOffset(1.2);
579   //
580   graph->SetMarkerColor(color);
581   graph->SetLineColor(color);
582   graph->SetMarkerStyle(style);
583 }
584
585 void AliTPCPreprocessorOffline::SetPadStyle(TPad *pad, Float_t mx0, Float_t mx1, Float_t my0, Float_t my1){
586   // 
587   // Set default pad style for QA
588   // 
589   pad->SetTicks(1,1);
590   pad->SetMargin(mx0,mx1,my0,my1);
591 }
592
593
594 void AliTPCPreprocessorOffline::MakeDefaultPlots(TObjArray * const arr, TObjArray *picArray){
595   //
596   // 0. make a default QA plots
597   // 1. Store them in the array
598   //
599   //
600   Float_t mx0=0.12, mx1=0.1, my0=0.15, my1=0.1;
601   //
602   TGraphErrors* laserA       =(TGraphErrors*)arr->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_A");
603   TGraphErrors* laserC       =(TGraphErrors*)arr->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_C");
604   TGraphErrors* cosmic       =(TGraphErrors*)arr->FindObject("TGRAPHERRORS_MEAN_VDRIFT_COSMICS_ALL");
605   TGraphErrors* cross        =(TGraphErrors*)arr->FindObject("TGRAPHERRORS_VDRIFT_CROSS_ALL");
606   TGraphErrors* itstpcP       =(TGraphErrors*)arr->FindObject("ALIGN_ITSP_TPC_DRIFTVD");
607   TGraphErrors* itstpcM       =(TGraphErrors*)arr->FindObject("ALIGN_ITSM_TPC_DRIFTVD");
608   TGraphErrors* itstpcB       =(TGraphErrors*)arr->FindObject("ALIGN_ITSB_TPC_DRIFTVD");
609   //
610   if (laserA)  SetDefaultGraphDrift(laserA,2,25);
611   if (laserC)  SetDefaultGraphDrift(laserC,4,26);
612   if (cosmic)  SetDefaultGraphDrift(cosmic,3,27);
613   if (cross)   SetDefaultGraphDrift(cross,4,28);
614   if (itstpcP) SetDefaultGraphDrift(itstpcP,2,29);
615   if (itstpcM) SetDefaultGraphDrift(itstpcM,4,30);
616   if (itstpcB) SetDefaultGraphDrift(itstpcB,1,31);
617   //
618   //
619   TPad *pad=0;
620   //
621   // Laser-Laser
622   //
623   if (laserA&&laserC){
624     pad = new TCanvas("TPCLaserVDrift","TPCLaserVDrift");
625     laserA->Draw("alp");
626     SetPadStyle(pad,mx0,mx1,my0,my1);
627     laserA->Draw("apl");
628     laserC->Draw("p");
629     TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
630     legend->AddEntry(laserA,"Laser A side");
631     legend->AddEntry(laserC,"Laser C side");
632     legend->Draw();    
633     picArray->AddLast(pad);
634   }
635
636   if (itstpcP&&itstpcM){
637     pad = new TCanvas("ITSTPC","ITSTPC");
638     itstpcP->Draw("alp");
639     SetPadStyle(pad,mx0,mx1,my0,my1);    
640     itstpcP->Draw("alp");
641     gPad->Clear();
642     itstpcM->Draw("apl");
643     itstpcP->Draw("p");
644     itstpcB->Draw("p");
645     TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
646     legend->AddEntry(itstpcP,"ITS-TPC smooth plus");
647     legend->AddEntry(itstpcM,"ITS-TPC smooth minus");
648     legend->AddEntry(itstpcB,"ITS-TPC smooth ");
649     legend->Draw();    
650     picArray->AddLast(pad);
651   }
652
653   if (itstpcB&&laserA){
654     pad = new TCanvas("ITSTPC_LASER","ITSTPC_LASER");
655     SetPadStyle(pad,mx0,mx1,my0,my1);    
656     laserA->Draw("alp");
657     itstpcP->Draw("p");
658     itstpcM->Draw("p");
659     itstpcB->Draw("p");
660     TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
661     legend->AddEntry(laserA,"TPC laser");
662     legend->AddEntry(itstpcP,"ITS-TPC smooth plus");   
663     legend->AddEntry(itstpcM,"ITS-TPC smooth minus");   
664     legend->AddEntry(itstpcB,"ITS-TPC smooth ");
665     legend->Draw();
666     picArray->AddLast(pad);
667   }
668
669   if (itstpcP&&cross){ 
670     pad = new TCanvas("ITSTPC_CROSS","ITSTPC_CROSS");
671     SetPadStyle(pad,mx0,mx1,my0,my1);    
672     itstpcP->Draw("alp");
673     pad->Clear();
674     cross->Draw("ap");
675     itstpcP->Draw("p");
676     //
677     TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
678
679     legend->AddEntry(cross,"TPC cross tracks");
680     legend->AddEntry(itstpcB,"ITS-TPC smooth");
681     legend->Draw();        
682     picArray->AddLast(pad);
683   }
684   if (itstpcP&&cosmic){ 
685     pad = new TCanvas("ITSTPC_COSMIC","ITSTPC_COSMIC");
686     SetPadStyle(pad,mx0,mx1,my0,my1);    
687     itstpcP->Draw("alp");
688     pad->Clear();
689     cosmic->Draw("ap");
690     itstpcP->Draw("p");
691     //
692     TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
693
694     legend->AddEntry(cosmic,"TPC cross tracks0 up-down");
695     legend->AddEntry(itstpcB,"ITS-TPC smooth");
696     legend->Draw();        
697     picArray->AddLast(pad);
698   }
699 }
700
701
702
703
704 void AliTPCPreprocessorOffline::CalibTimeGain(const Char_t* fileName, Int_t startRunNumber, Int_t endRunNumber,  TString  pocdbStorage){
705   //
706   // Update OCDB gain
707   //
708   ReadGainGlobal(fileName);
709   AnalyzeGain(startRunNumber,endRunNumber, 1000,1.43);
710   AnalyzeAttachment(startRunNumber,endRunNumber);
711   MakeQAPlot(1.43);  
712   if (pocdbStorage.Length()==0) pocdbStorage+="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
713   UpdateOCDBGain( startRunNumber, endRunNumber, pocdbStorage.Data());
714 }
715
716
717
718
719 void AliTPCPreprocessorOffline::ReadGainGlobal(const Char_t* fileName){
720   //
721   // read calibration entries from file
722   // 
723   TFile fcalib(fileName);
724   TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib");
725   if (array){
726     fGainMIP    = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGain");
727     fGainCosmic = ( AliTPCcalibTimeGain *)array->FindObject("calibTimeGainCosmic");
728   }else{
729     fGainMIP    = ( AliTPCcalibTimeGain *)fcalib.Get("calibTimeGain");
730     fGainCosmic = ( AliTPCcalibTimeGain *)fcalib.Get("calibTimeGainCosmic");
731   }
732   TH1 * hisT=0;
733   Int_t firstBinA =0, lastBinA=0;
734
735   if (fGainCosmic){ 
736     hisT= fGainCosmic->GetHistGainTime()->Projection(1);
737     firstBinA = hisT->FindFirstBinAbove(2);
738     lastBinA  = hisT->FindLastBinAbove(2);    
739     fGainCosmic->GetHistGainTime()->GetAxis(1)->SetRange(firstBinA,lastBinA);
740     delete hisT;
741   }
742
743   if (fGainMIP){ 
744     hisT= fGainMIP->GetHistGainTime()->Projection(1);
745     firstBinA = hisT->FindFirstBinAbove(2);
746     lastBinA  = hisT->FindLastBinAbove(2);    
747     fGainMIP->GetHistGainTime()->GetAxis(1)->SetRange(firstBinA,lastBinA);
748     delete hisT;
749   }
750
751 }
752
753
754
755 Bool_t AliTPCPreprocessorOffline::AnalyzeGain(Int_t startRunNumber, Int_t endRunNumber, Int_t minEntriesGaussFit,  Float_t FPtoMIPratio){
756   //
757   // Analyze gain - produce the calibration graphs
758   //
759   fGainMIP->GetHistGainTime()->GetAxis(5)->SetRangeUser(startRunNumber, endRunNumber);
760   // 1.) try to create MIP spline
761   fGainMIP->GetHistGainTime()->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data
762   fGainMIP->GetHistGainTime()->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions
763   //
764   fGraphMIP = AliTPCcalibBase::FitSlices(fGainMIP->GetHistGainTime(),0,1,minEntriesGaussFit,10,0.1,0.7);
765   if (fGraphMIP->GetN()==0) fGraphMIP = 0x0;
766   if (fGraphMIP) fFitMIP = AliTPCcalibTimeGain::MakeSplineFit(fGraphMIP);
767   if (fGraphMIP) fGraphMIP->SetName("TGRAPHERRORS_MEAN_GAIN_BEAM_ALL");// set proper names according to naming convention
768   fGainArray->AddAt(fFitMIP,0);
769   
770
771   // 2.) try to create Cosmic spline
772   if (fGainCosmic){
773     fGainCosmic->GetHistGainTime()->GetAxis(2)->SetRangeUser(0.51,1.49); // only cosmics
774     fGainCosmic->GetHistGainTime()->GetAxis(4)->SetRangeUser(20,100);    // only Fermi-Plateau muons
775     //
776     fGraphCosmic = AliTPCcalibBase::FitSlices(fGainCosmic->GetHistGainTime(),0,1,minEntriesGaussFit,10);
777     if (fGraphCosmic->GetN()==0) fGraphCosmic = 0x0;
778     //
779     if (fGraphCosmic) {
780       for(Int_t i=0; i < fGraphCosmic->GetN(); i++) {
781         fGraphCosmic->GetY()[i] /= FPtoMIPratio;
782         fGraphCosmic->GetEY()[i] /= FPtoMIPratio;
783       }
784     }
785     //
786     if (fGraphCosmic) fFitCosmic = AliTPCcalibTimeGain::MakeSplineFit(fGraphCosmic);
787     if (fGraphCosmic) fGraphCosmic->SetName("TGRAPHERRORS_MEAN_GAIN_COSMIC_ALL"); // set proper names according to naming convention
788     fGainArray->AddAt(fFitCosmic,1);
789   }
790   // with naming convention and backward compatibility
791   fGainArray->AddAt(fGraphMIP,2);
792   fGainArray->AddAt(fGraphCosmic,3);
793   cout << "fGraphCosmic: " << fGraphCosmic << " fGraphMIP " << fGraphMIP << endl;
794   return kTRUE;
795
796 }
797
798
799 Bool_t AliTPCPreprocessorOffline::AnalyzeAttachment(Int_t startRunNumber, Int_t endRunNumber, Int_t minEntriesFit) {
800   //
801   // determine slope as a function of mean driftlength
802   //
803   fGainMIP->GetHistGainTime()->GetAxis(5)->SetRangeUser(startRunNumber, endRunNumber);
804   //
805   fGainMIP->GetHistGainTime()->GetAxis(2)->SetRangeUser(1.51,2.49); // only beam data
806   fGainMIP->GetHistGainTime()->GetAxis(4)->SetRangeUser(0.39,0.51); // only MIP pions
807   //
808   TH3D * hist = fGainMIP->GetHistGainTime()->Projection(1, 0, 3);
809   //
810   Double_t *xvec = new Double_t[hist->GetNbinsX()];
811   Double_t *yvec = new Double_t[hist->GetNbinsX()];
812   Double_t *xerr = new Double_t[hist->GetNbinsX()];
813   Double_t *yerr = new Double_t[hist->GetNbinsX()];
814   Int_t counter  = 0;
815   //
816   for(Int_t i=1; i < hist->GetNbinsX(); i++) {
817     Int_t nsum=0;
818     Int_t imin   =  i;
819     Int_t imax   =  i;    
820     for (Int_t idelta=0; idelta<10; idelta++){
821       //
822       imin   =  TMath::Max(i-idelta,1);
823       imax   =  TMath::Min(i+idelta,hist->GetNbinsX());
824       nsum = TMath::Nint(hist->Integral(imin,imax,1,hist->GetNbinsY()-1,1,hist->GetNbinsZ()-1));
825       //if (nsum==0) break;
826       if (nsum>minEntriesFit) break;
827     }
828     if (nsum<minEntriesFit) continue;
829     //
830     fGainMIP->GetHistGainTime()->GetAxis(1)->SetRangeUser(hist->GetXaxis()->GetBinCenter(imin),hist->GetXaxis()->GetBinCenter(imax)); // define time range
831     TGraphErrors * driftDep = AliTPCcalibBase::FitSlices(fGainMIP->GetHistGainTime(),0,3,100,10,0.1,0.7);
832     if (driftDep->GetN() < 4) {
833       delete driftDep;
834          continue;
835     }
836     //
837     TObjArray arr;
838     //
839     TF1 pol1("polynom1","pol1",10,240);
840     //driftDep->Fit(&pol1,"QNRROB=0.8");
841     driftDep->Fit(&pol1,"QNR");
842     xvec[counter] = 0.5*(hist->GetXaxis()->GetBinCenter(imin)+hist->GetXaxis()->GetBinCenter(imax));
843     yvec[counter] = pol1.GetParameter(1)/pol1.GetParameter(0);
844     xerr[counter] = 0;
845     yerr[counter] = pol1.GetParError(1)/pol1.GetParameter(0);
846     counter++;
847     //
848     delete driftDep;
849   }
850   //
851   fGraphAttachmentMIP = new TGraphErrors(counter, xvec, yvec, xerr, yerr);
852   if (fGraphAttachmentMIP) fGraphAttachmentMIP->SetName("TGRAPHERRORS_MEAN_ATTACHMENT_BEAM_ALL");// set proper names according to naming convention
853   fGainArray->AddLast(fGraphAttachmentMIP);
854   //
855   delete [] xvec;
856   delete [] yvec;
857   delete [] xerr;
858   delete [] yerr;
859   delete hist;
860   //
861   if (counter < 1) return kFALSE;
862   return kTRUE;
863 }
864
865
866
867
868
869 void AliTPCPreprocessorOffline::UpdateOCDBGain(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
870   //
871   // Update OCDB entry
872   //
873   AliCDBMetaData *metaData= new AliCDBMetaData();
874   metaData->SetObjectClassName("TObjArray");
875   metaData->SetResponsible("Alexander Kalweit");
876   metaData->SetBeamPeriod(1);
877   metaData->SetAliRootVersion("05-24-00"); //root version
878   metaData->SetComment("Calibration of the time dependence of the gain due to pressure and temperature changes.");
879   AliCDBId id1("TPC/Calib/TimeGain", startRunNumber, endRunNumber);
880   AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
881   gStorage->Put(fGainArray, id1, metaData);    
882 }
883
884 void AliTPCPreprocessorOffline::MakeQAPlot(Float_t  FPtoMIPratio) {
885   //
886   // Make QA plot to visualize results
887   //
888   //
889   //
890   if (fGraphCosmic) {
891     TCanvas * canvasCosmic = new TCanvas("gain Cosmic", "time dependent gain QA histogram cosmic");
892     canvasCosmic->cd();
893     TH2D * gainHistoCosmic = fGainCosmic->GetHistGainTime()->Projection(0,1);
894     gainHistoCosmic->SetDirectory(0);
895     gainHistoCosmic->SetName("GainHistoCosmic");
896     gainHistoCosmic->GetXaxis()->SetTimeDisplay(kTRUE);
897     gainHistoCosmic->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
898     gainHistoCosmic->Draw("colz");
899     fGraphCosmic->SetMarkerStyle(25);
900     fGraphCosmic->Draw("lp");
901     fGraphCosmic->SetMarkerStyle(25);
902     TGraph * grfFitCosmic = fFitCosmic->MakeGraph(fGraphCosmic->GetX()[0],fGraphCosmic->GetX()[fGraphCosmic->GetN()-1],50000,0);
903     if (grfFitCosmic) {
904       for(Int_t i=0; i < grfFitCosmic->GetN(); i++) {
905         grfFitCosmic->GetY()[i] *= FPtoMIPratio;        
906       }
907       for(Int_t i=0; i < fGraphCosmic->GetN(); i++) {
908         fGraphCosmic->GetY()[i] *= FPtoMIPratio;        
909       }
910     }
911     fGraphCosmic->Draw("lp");
912     grfFitCosmic->SetLineColor(2);
913     grfFitCosmic->Draw("lu");
914     fGainArray->AddLast(gainHistoCosmic);
915     fGainArray->AddLast(canvasCosmic->Clone());
916     delete canvasCosmic;    
917   }
918   if (fFitMIP) {
919     TCanvas * canvasMIP = new TCanvas("gain MIP", "time dependent gain QA histogram MIP");
920     canvasMIP->cd();
921     TH2D * gainHistoMIP    = fGainMIP->GetHistGainTime()->Projection(0,1);
922     gainHistoMIP->SetName("GainHistoCosmic");
923     gainHistoMIP->SetDirectory(0);
924     gainHistoMIP->GetXaxis()->SetTimeDisplay(kTRUE);
925     gainHistoMIP->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
926     gainHistoMIP->Draw("colz");
927     fGraphMIP->SetMarkerStyle(25);
928     fGraphMIP->Draw("lp");
929     TGraph * grfFitMIP = fFitMIP->MakeGraph(fGraphMIP->GetX()[0],fGraphMIP->GetX()[fGraphMIP->GetN()-1],50000,0);
930     grfFitMIP->SetLineColor(2);
931     grfFitMIP->Draw("lu");    
932     fGainArray->AddLast(gainHistoMIP);
933     fGainArray->AddLast(canvasMIP->Clone());
934     delete canvasMIP;    
935   }  
936 }
937
938
939