]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/CalibMacros/CalibTimeVdrift.C
Adding macro for ExB global fit
[u/mrichter/AliRoot.git] / TPC / CalibMacros / CalibTimeVdrift.C
1 /* 
2    Small macro to generate and update OCDB entries for a given run:
3    To be used together with pother macros to create OCDB entries
4    // How to use it...
5    //
6    // 1. go to the directory with the "CalibObjectsTrain1.root" file and start your aliroot
7    //
8    // 2. copy and paste the following commands
9    //
10    gSystem->Load("libSTEER");
11    gSystem->Load("libANALYSIS");
12    gSystem->Load("libSTAT");
13    gSystem->Load("libTPCcalib");
14    gSystem->AddIncludePath("-I$ALICE_ROOT/STEER");
15    gSystem->AddIncludePath("-I$ALICE_ROOT/TPC");
16    .L $ALICE_ROOT/TPC/CalibMacros/CalibTimeVdrift.C+
17    
18    //
19    // 3. setup your OCDB output path e.g:
20    //
21    ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
22    //
23    // 4. simple execution
24    //
25    CalibTimeVdrift("CalibObjectsTrain1.root", 110000, AliCDBRunRange::Infinity())
26    //
27    // 5. try to visualize new entry
28    //   
29    TObjArray *picArray = new TObjArray;
30    MakeDefaultPlots(vdriftArray,picArray);
31 */
32
33
34 #if !defined(__CINT__) || defined(__MAKECINT__)
35 #include "TMap.h"
36 #include "TGraphErrors.h"
37 #include "AliExternalTrackParam.h"
38 #include "TFile.h"
39 #include "TGraph.h"
40 #include "TMultiGraph.h"
41 #include "TCanvas.h"
42 #include "THnSparse.h"
43 #include "TLegend.h"
44 #include "TPad.h"
45 #include "TH2D.h"
46
47 #include "AliESDfriend.h"
48
49
50 #include "AliTPCcalibTime.h"
51 #include "AliSplineFit.h"
52 #include "AliCDBMetaData.h"
53 #include "AliCDBId.h"
54 #include "AliCDBManager.h"
55 #include "AliCDBStorage.h"
56 #include "AliTPCcalibBase.h"
57 #include "AliTPCcalibDB.h"
58 #include "AliTPCcalibDButil.h"
59 #include "AliRelAlignerKalman.h"
60 #include "AliTPCParamSR.h"
61 #endif
62
63
64 TString ocdbStorage="dummy";
65 const Int_t    kMinEntries=100;     // minimal number of entries
66 Int_t startRun=0, endRun=0;
67 Int_t startTime=0, endTime=0;
68 TObjArray * vdriftArray = 0;
69 AliTPCcalibTime * timeDrift=0;
70 //
71 // functions:
72 //
73 void CalibTimeVdrift(Char_t* file="CalibObjectsTrain1.root", Int_t ustartRun=0, Int_t uendRun=AliCDBRunRange::Infinity(),TString ocdbStorage="");
74 void UpdateOCDBDrift(Int_t ustartRun, Int_t uendRun);
75 void AddHistoGraphs(  TObjArray * vdriftArray, AliTPCcalibTime *timeDrift, Int_t minEntries);
76 void AddAlignmentGraphs(  TObjArray * vdriftArray, AliTPCcalibTime *timeDrift);
77 void AddLaserGraphs(  TObjArray * vdriftArray, AliTPCcalibTime *timeDrift);
78 void SetDefaultGraphDrift(TGraph *graph, Int_t color, Int_t style);
79 TGraphErrors* FilterGraphMedianAbs(TGraphErrors * graph, Float_t cut,Double_t &medianY);
80 TGraphErrors* FilterGraphDrift(TGraphErrors * graph, Float_t errSigmaCut, Float_t medianCutAbs);
81 TGraphErrors * MakeGraphFilter0(THnSparse *hisN, Int_t itime, Int_t ival, Int_t minEntries, Double_t offset=0);
82
83
84 void GetRunRange(AliTPCcalibTime* timeDrift){
85   //
86   // find the fist and last run
87   //
88   TObjArray *hisArray =timeDrift->GetHistoDrift();
89   {for (Int_t i=0; i<hisArray->GetEntriesFast(); i++){
90     THnSparse* addHist=(THnSparse*)hisArray->UncheckedAt(i);
91     if (addHist->GetEntries()<kMinEntries) continue;
92     if (!addHist) continue;
93     TH1D* histo    =addHist->Projection(3);
94     TH1D* histoTime=addHist->Projection(0);
95     printf("%s\t%f\t%d\t%d\n",histo->GetName(), histo->GetEntries(),histo->FindFirstBinAbove(0),histo->FindLastBinAbove(0));
96
97     if (startRun<=0){ 
98       startRun=histo->FindFirstBinAbove(0);
99       endRun  =histo->FindLastBinAbove(0);
100     }else{
101       startRun=TMath::Min(histo->FindFirstBinAbove(0),startRun);
102       endRun  =TMath::Max(histo->FindLastBinAbove(0),endRun);
103     }
104     if (startTime==0){ 
105       startTime=histoTime->FindFirstBinAbove(0);
106       endTime  =histoTime->FindLastBinAbove(0);
107     }else{
108       startTime=TMath::Min(histoTime->FindFirstBinAbove(0),startTime);
109       endTime  =TMath::Max(histoTime->FindLastBinAbove(0),endTime);
110     }
111     delete histo;
112     delete histoTime;
113   }}
114   if (startRun<0) startRun=0;
115   if (endRun<0) endRun=100000000;
116   printf("Run range  :\t%d-%d\n", startRun, endRun);
117   printf("Time range :\t%d-%d\n", startTime, endTime);
118
119 }
120
121
122
123 void CalibTimeVdrift(Char_t* file, Int_t ustartRun, Int_t uendRun, TString pocdbStorage){
124   //
125   //
126   //
127   const Int_t    kMinEntries=500;     // minimal number of entries
128   if (pocdbStorage.Length()>0) ocdbStorage=pocdbStorage;
129   //
130   // 1. Initialization and run range setting
131   TFile fcalib(file);
132   timeDrift=(AliTPCcalibTime*)fcalib.Get("calibTime");
133   startRun=ustartRun;
134   endRun=ustartRun; 
135   TObjArray *hisArray =timeDrift->GetHistoDrift();  
136   GetRunRange(timeDrift);
137   for (Int_t i=0; i<hisArray->GetEntriesFast(); i++){
138     THnSparse* addHist=(THnSparse*)hisArray->At(i);
139     if (!addHist) continue;
140     if (startTime<endTime) addHist->GetAxis(0)->SetRange(startTime-1,endTime+1);
141     if (startRun<endRun) addHist->GetAxis(3)->SetRange(startRun-1,endRun+1);
142   }
143   //
144   //
145   // 2. extraction of the information
146   //
147   vdriftArray = new TObjArray();
148   AddAlignmentGraphs(vdriftArray,timeDrift);
149   AddHistoGraphs(vdriftArray,timeDrift,kMinEntries);
150   AddLaserGraphs(vdriftArray,timeDrift);
151   //
152   // 3. Append QA plots
153   //
154   MakeDefaultPlots(vdriftArray,vdriftArray);
155   //
156   //
157   // 4. update of OCDB
158   //
159   //
160   UpdateOCDBDrift(ustartRun,uendRun);
161 }
162
163 void UpdateOCDBDrift( Int_t ustartRun, Int_t uendRun){
164   //
165   // Update OCDB 
166   //
167   AliCDBMetaData *metaData= new AliCDBMetaData();
168   metaData->SetObjectClassName("TObjArray");
169   metaData->SetResponsible("Marian Ivanov");
170   metaData->SetBeamPeriod(1);
171   metaData->SetAliRootVersion("05-25-01"); //root version
172   metaData->SetComment("Calibration of the time dependence of the drift velocity");
173   AliCDBId* id1=NULL;
174   id1=new AliCDBId("TPC/Calib/TimeDrift", ustartRun, uendRun);
175   AliCDBStorage* gStorage = AliCDBManager::Instance()->GetStorage(ocdbStorage);
176   gStorage->Put(vdriftArray, (*id1), metaData);
177 }
178
179
180
181 void UpdateDriftParam(AliTPCParam *param, TObjArray *arr, Int_t startRun){
182   //
183   //  update the OCDB entry for the nominal time0
184   //
185   //
186   //  AliTPCParam * param = AliTPCcalibDB::Instance()->GetParameters();
187   AliTPCParam *paramNew = (AliTPCParam *)param->Clone();
188   TGraphErrors *grT =  (TGraphErrors *)arr->FindObject("ALIGN_ITSM_TPC_T0");
189   Double_t deltaTcm = TMath::Median(grT->GetN(),grT->GetY());
190   Double_t deltaT   = deltaTcm/param->GetDriftV();
191   paramNew->SetL1Delay(param->GetL1Delay()-deltaT);
192   paramNew->Update();
193
194   AliCDBMetaData *metaData= new AliCDBMetaData();
195   metaData->SetObjectClassName("TObjArray");
196   metaData->SetResponsible("Marian Ivanov");
197   metaData->SetBeamPeriod(1);
198   metaData->SetAliRootVersion("05-25-02"); //root version
199   metaData->SetComment("Updated calibration of nominal time 0");
200   AliCDBId* id1=NULL;
201   id1=new AliCDBId("TPC/Calib/Parameters", startRun, AliCDBRunRange::Infinity());
202   AliCDBStorage* gStorage = AliCDBManager::Instance()->GetStorage(ocdbStorage);
203   gStorage->Put(param, (*id1), metaData);
204
205 }
206
207
208 void PrintArray(TObjArray *array){
209   //
210   //
211   //
212   Int_t entries = array->GetEntries();
213   for (Int_t i=0; i<entries; i++){
214     printf("%d\t %s\n", i,  array->At(i)->GetName());
215   }
216 }
217
218
219
220 TGraphErrors* FilterGraphDrift(TGraphErrors * graph, Float_t errSigmaCut, Float_t medianCutAbs){
221   // 2 filters:
222   //    1. filter graph - error cut errSigmaCut
223   //    2. filter graph - medianCutAbs around median
224   //
225   // errSigmaCut   - cut on error
226   // medianCutAbs  - cut on value around median
227   Double_t dummy=0;               //   
228   //
229   // 1. filter graph - error cut errSigmaCut
230   //              
231   TGraphErrors *graphF; 
232   graphF = AliTPCcalibDButil::FilterGraphMedianErr(graph,errSigmaCut,dummy);
233   delete graph;
234   if (!graphF) return 0;
235   graph = AliTPCcalibDButil::FilterGraphMedianErr(graphF,errSigmaCut,dummy);
236   delete graphF;
237   if (!graph) return 0;
238   //
239   // filter graph - kMedianCutAbs around median
240   // 
241   graphF=FilterGraphMedianAbs(graph, medianCutAbs,dummy);
242   delete graph;
243   if (!graphF) return 0;
244   graph=FilterGraphMedianAbs(graphF, medianCutAbs,dummy);
245   delete graphF;
246   if (!graph) return 0;
247   return graph;
248 }
249
250
251
252 TGraphErrors* FilterGraphMedianAbs(TGraphErrors * graph, Float_t cut,Double_t &medianY){
253   //
254   // filter outlyer measurement
255   // Only points around median +- cut filtered 
256   //
257   if (!graph) return  0;
258   Int_t kMinPoints=2;
259   Int_t npoints0 = graph->GetN();
260   Int_t npoints=0;
261   Float_t  rmsY=0;
262   Double_t *outx=new Double_t[npoints0];
263   Double_t *outy=new Double_t[npoints0];
264   Double_t *errx=new Double_t[npoints0];
265   Double_t *erry=new Double_t[npoints0];
266   //
267   //
268   if (npoints0<kMinPoints) return 0;
269   for (Int_t iter=0; iter<3; iter++){
270     npoints=0;
271     for (Int_t ipoint=0; ipoint<npoints0; ipoint++){
272       if (graph->GetY()[ipoint]==0) continue;
273       if (iter>0 &&TMath::Abs(graph->GetY()[ipoint]-medianY)>cut) continue;  
274       outx[npoints]  = graph->GetX()[ipoint];
275       outy[npoints]  = graph->GetY()[ipoint];
276       errx[npoints]  = graph->GetErrorX(ipoint);
277       erry[npoints]  = graph->GetErrorY(ipoint);
278       npoints++;
279     }
280     if (npoints<=1) break;
281     medianY  =TMath::Median(npoints,outy);
282     rmsY   =TMath::RMS(npoints,outy);
283   }
284   TGraphErrors *graphOut=0;
285   if (npoints>1) graphOut= new TGraphErrors(npoints,outx,outy,errx,erry); 
286   return graphOut;
287 }
288
289
290 void AddHistoGraphs(  TObjArray * vdriftArray, AliTPCcalibTime *timeDrift, Int_t minEntries){
291   //
292   // Add graphs corresponding to the alignment
293   //
294   const Double_t kErrSigmaCut=5;      // error sigma cut - for filtering
295   const Double_t kMedianCutAbs=0.03;  // error sigma cut - for filtering
296   //
297   TObjArray * array=timeDrift->GetHistoDrift();
298   if (array){
299     THnSparse* hist=NULL;
300     // 2.a) cosmics with different triggers
301     for (Int_t i=0; i<array->GetEntriesFast();i++){
302       hist=(THnSparseF*)array->UncheckedAt(i);
303       if(!hist) continue;
304       if (hist->GetEntries()<minEntries) continue;
305       //hist->Print();
306       TString name=hist->GetName();
307       Int_t dim[4]={0,1,2,3};
308       THnSparse* newHist=hist->Projection(4,dim);
309       newHist->SetName(name);
310       TGraphErrors* graph=AliTPCcalibBase::FitSlices(newHist,2,0,400,100,0.05,0.95, kTRUE);
311       printf("name=%s graph=%i, N=%i\n", name.Data(), graph==0, graph->GetN());
312       Int_t pos=name.Index("_");
313       name=name(pos,name.Capacity()-pos);
314       TString graphName=graph->ClassName();
315       graphName+=name;
316       graphName.ToUpper();
317       //
318       graph = FilterGraphDrift(graph, kErrSigmaCut, kMedianCutAbs);
319       if (!graph) {
320         printf("Graph =%s filtered out\n", name.Data());
321         continue;
322       }
323       //
324       graph->SetMarkerStyle(i%8+20);
325       graph->SetMarkerColor(i%7);
326       graph->GetXaxis()->SetTitle("Time");
327       graph->GetYaxis()->SetTitle("v_{dcor}");
328       graph->SetName(graphName);
329       graph->SetTitle(graphName);
330       printf("Graph %d\t=\t%s\n", i, graphName.Data());
331       vdriftArray->Add(graph);
332     }
333   }
334 }
335
336
337
338
339 void AddAlignmentGraphs(  TObjArray * vdriftArray, AliTPCcalibTime *timeDrift){
340   //
341   // Add graphs corresponding to alignment to the object array
342   //
343   TObjArray *arrayITS=0;
344   TObjArray *arrayTOF=0;
345   TObjArray *arrayTRD=0;
346   TMatrixD *mstatITS=0;
347   TMatrixD *mstatTOF=0;
348   TMatrixD *mstatTRD=0;
349   //
350   arrayITS=timeDrift->GetAlignITSTPC();
351   arrayTRD=timeDrift->GetAlignTRDTPC();
352   arrayTOF=timeDrift->GetAlignTOFTPC();
353
354   if (arrayITS->GetEntries()>0) mstatITS= AliTPCcalibDButil::MakeStatRelKalman(arrayITS,0.9,50,0.025);
355   if (arrayTOF->GetEntries()>0) mstatTOF= AliTPCcalibDButil::MakeStatRelKalman(arrayTOF,0.9,1000,0.025);
356   if (arrayTRD->GetEntries()>0) mstatTRD= AliTPCcalibDButil::MakeStatRelKalman(arrayTRD,0.9,50,0.025);
357   //
358   TObjArray * arrayITSP= AliTPCcalibDButil::SmoothRelKalman(arrayITS,*mstatITS, 0, 5.);
359   TObjArray * arrayITSM= AliTPCcalibDButil::SmoothRelKalman(arrayITS,*mstatITS, 1, 5.);
360   TObjArray * arrayITSB= AliTPCcalibDButil::SmoothRelKalman(arrayITSP,arrayITSM);
361   TObjArray * arrayTOFP= AliTPCcalibDButil::SmoothRelKalman(arrayTOF,*mstatTOF, 0, 5.);
362   TObjArray * arrayTOFM= AliTPCcalibDButil::SmoothRelKalman(arrayTOF,*mstatTOF, 1, 5.);
363   TObjArray * arrayTOFB= AliTPCcalibDButil::SmoothRelKalman(arrayTOFP,arrayTOFM);
364
365   TObjArray * arrayTRDP= 0x0;
366   TObjArray * arrayTRDM= 0x0;
367   TObjArray * arrayTRDB= 0x0;
368   arrayTRDP= AliTPCcalibDButil::SmoothRelKalman(arrayTRD,*mstatTRD, 0, 5.);
369   arrayTRDM= AliTPCcalibDButil::SmoothRelKalman(arrayTRD,*mstatTRD, 1, 5.);
370   arrayTRDB= AliTPCcalibDButil::SmoothRelKalman(arrayTRDP,arrayTRDM);
371   //
372   //
373   Int_t entries=TMath::Max(arrayITS->GetEntriesFast(),arrayTOF->GetEntriesFast());
374   TObjArray *arrays[12]={arrayITS, arrayITSP, arrayITSM, arrayITSB,
375                          arrayTRD, arrayTRDP, arrayTRDM, arrayTRDB,
376                          arrayTOF, arrayTOFP, arrayTOFM, arrayTOFB};
377   TString   grnames[12]={"ALIGN_ITS", "ALIGN_ITSP", "ALIGN_ITSM", "ALIGN_ITSB",
378                          "ALIGN_TRD", "ALIGN_TRDP", "ALIGN_TRDM","ALIGN_TRDB",
379                          "ALIGN_TOF", "ALIGN_TOFP", "ALIGN_TOFM","ALIGN_TOFB"};
380   TString   grpar[9]={"DELTAPSI", "DELTATHETA", "DELTAPHI",
381                       "DELTAX", "DELTAY", "DELTAZ",
382                       "DRIFTVD", "T0", "VDGY"};
383
384   
385   TVectorD vX(entries);
386   TVectorD vY(entries);
387   TVectorD vEx(entries);
388   TVectorD vEy(entries);
389   TObjArray *arr=0;
390   for (Int_t iarray=0; iarray<12; iarray++){
391     arr = arrays[iarray];
392     if (arr==0) continue;
393     for (Int_t ipar=0; ipar<9; ipar++){      
394       Int_t counter=0;
395       for (Int_t itime=0; itime<arr->GetEntriesFast(); itime++){
396         AliRelAlignerKalman * kalman = (AliRelAlignerKalman *) arr->UncheckedAt(itime);
397         if (!kalman) continue;
398         vX[counter]=kalman->GetTimeStamp();
399         vY[counter]=(*(kalman->GetState()))[ipar];
400         if (ipar==6) vY[counter]=1./(*(kalman->GetState()))[ipar]-1;
401         vEx[counter]=0;
402         vEy[counter]=TMath::Sqrt((*(kalman->GetStateCov()))(ipar,ipar));
403         counter++;
404       }
405     
406       TGraphErrors * graph=new TGraphErrors(counter, vX.GetMatrixArray(),
407                                           vY.GetMatrixArray(),
408                                           vEx.GetMatrixArray(),
409                                           vEy.GetMatrixArray());
410       TString grName=grnames[iarray];
411       grName+="_TPC_";
412       grName+=grpar[ipar];
413       graph->SetName(grName.Data());
414       vdriftArray->AddLast(graph);
415     }
416   }  
417 }
418
419
420
421
422 void AddLaserGraphs(  TObjArray * vdriftArray, AliTPCcalibTime *timeDrift){
423   //
424   // add graphs for laser
425   //
426   const Double_t delayL0L1 = 0.071;  //this is hack for 1/2 weeks
427   THnSparse *hisN=0;
428   TGraphErrors *grLaser[6]={0,0,0,0,0,0};
429   hisN = timeDrift->GetHistVdriftLaserA(0);
430   if (timeDrift->GetHistVdriftLaserA(0)){
431     grLaser[0]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(0),0,2,5,delayL0L1);
432     grLaser[0]->SetName("GRAPH_MEAN_DELAY_LASER_ALL_A");
433     vdriftArray->AddLast(grLaser[0]);
434   }    
435   if (timeDrift->GetHistVdriftLaserA(1)){
436     grLaser[1]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(1),0,2,5);
437     grLaser[1]->SetName("GRAPH_MEAN_DRIFT_LASER_ALL_A");
438     vdriftArray->AddLast(grLaser[1]);
439   }    
440   if (timeDrift->GetHistVdriftLaserA(2)){
441     grLaser[2]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserA(2),0,2,5);
442     grLaser[2]->SetName("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
443     vdriftArray->AddLast(grLaser[2]);
444   }    
445   if (timeDrift->GetHistVdriftLaserC(0)){
446     grLaser[3]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(0),0,2,5,delayL0L1);
447     grLaser[3]->SetName("GRAPH_MEAN_DELAY_LASER_ALL_C");
448     vdriftArray->AddLast(grLaser[3]);
449   }    
450   if (timeDrift->GetHistVdriftLaserC(1)){
451     grLaser[4]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(1),0,2,5);
452     grLaser[4]->SetName("GRAPH_MEAN_DRIFT_LASER_ALL_C");
453     vdriftArray->AddLast(grLaser[4]);
454   }    
455   if (timeDrift->GetHistVdriftLaserC(2)){
456     grLaser[5]=MakeGraphFilter0(timeDrift->GetHistVdriftLaserC(2),0,2,5);
457     grLaser[5]->SetName("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");    
458     vdriftArray->AddLast(grLaser[5]);
459   }    
460   for (Int_t i=0; i<6;i++){
461     if (grLaser[i]) {
462       SetDefaultGraphDrift(grLaser[i], 1,(i+20));
463       grLaser[i]->GetYaxis()->SetTitle("Laser Correction");
464     }
465   }
466 }
467  
468  
469 TGraphErrors * MakeGraphFilter0(THnSparse *hisN, Int_t itime, Int_t ival, Int_t minEntries, Double_t offset){
470   //
471   // Make graph with mean values and rms
472   //
473   hisN->GetAxis(itime)->SetRange(0,100000000);
474   hisN->GetAxis(ival)->SetRange(0,100000000);
475   TH1 * hisT      = hisN->Projection(itime);
476   TH1 * hisV      = hisN->Projection(ival);
477   //
478   Int_t firstBinA = hisT->FindFirstBinAbove(2);
479   Int_t lastBinA  = hisT->FindLastBinAbove(2);    
480   Int_t firstBinV = hisV->FindFirstBinAbove(0);
481   Int_t lastBinV  = hisV->FindLastBinAbove(0);    
482   hisN->GetAxis(itime)->SetRange(firstBinA,lastBinA);
483   hisN->GetAxis(ival)->SetRange(firstBinV,lastBinV);
484   Int_t entries=0;
485   for (Int_t ibin=firstBinA; ibin<lastBinA; ibin++){
486     Double_t cont = hisT->GetBinContent(ibin);
487     if (cont<minEntries) continue;
488     entries++;
489   }
490   TVectorD vecTime(entries);
491   TVectorD vecMean0(entries);
492   TVectorD vecRMS0(entries);
493   TVectorD vecMean1(entries);
494   TVectorD vecRMS1(entries);
495   entries=0;
496   {for (Int_t ibin=firstBinA; ibin<lastBinA; ibin++){
497       Double_t cont = hisT->GetBinContent(ibin);
498       if (cont<minEntries) continue;
499       hisN->GetAxis(itime)->SetRange(ibin-1,ibin+1);
500       Double_t time = hisT->GetBinCenter(ibin);
501       TH1 * his = hisN->Projection(ival);
502       Double_t nentries0= his->GetBinContent(his->FindBin(0));
503       if (cont-nentries0<minEntries) continue;
504       //
505       his->SetBinContent(his->FindBin(0),0);
506       vecTime[entries]=time;
507       vecMean0[entries]=his->GetMean()+offset;
508       vecMean1[entries]=his->GetMeanError();
509       vecRMS0[entries] =his->GetRMS();
510       vecRMS1[entries] =his->GetRMSError();
511       delete his;  
512       entries++;
513     }}
514   delete hisT;
515   delete hisV;
516   TGraphErrors * graph =  new TGraphErrors(entries,vecTime.GetMatrixArray(), vecMean0.GetMatrixArray(),                                    0, vecMean1.GetMatrixArray());
517   return graph;
518 }
519
520
521
522
523
524
525
526
527 void SetDefaultGraphDrift(TGraph *graph, Int_t color, Int_t style){
528   //
529   //
530   //
531   graph->GetXaxis()->SetTimeDisplay(kTRUE);
532   graph->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
533   graph->SetMaximum( 0.025);
534   graph->SetMinimum(-0.025);
535   graph->GetXaxis()->SetTitle("Time");
536   graph->GetYaxis()->SetTitle("v_{dcorr}");
537   //
538   graph->GetYaxis()->SetLabelSize(0.03);
539   graph->GetXaxis()->SetLabelSize(0.03);
540   //
541   graph->GetXaxis()->SetNdivisions(10,5,0);
542   graph->GetYaxis()->SetNdivisions(10,5,0);
543   //
544   graph->GetXaxis()->SetLabelOffset(0.02);
545   graph->GetYaxis()->SetLabelOffset(0.005);
546   //
547   graph->GetXaxis()->SetTitleOffset(1.3);
548   graph->GetYaxis()->SetTitleOffset(1.2);
549   //
550   graph->SetMarkerColor(color);
551   graph->SetLineColor(color);
552   graph->SetMarkerStyle(style);
553 }
554
555 void SetPadStyle(TPad *pad, Float_t mx0, Float_t mx1, Float_t my0, Float_t my1){
556   //
557   //
558   pad->SetTicks(1,1);
559   pad->SetMargin(mx0,mx1,my0,my1);
560 }
561
562
563 void MakeDefaultPlots(TObjArray * arr, TObjArray *picArray){
564   //
565   //
566   //
567   // margins
568   Float_t mx0=0.12, mx1=0.1, my0=0.15, my1=0.1;
569   //
570   TGraphErrors* laserA       =(TGraphErrors*)arr->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_A");
571   TGraphErrors* laserC       =(TGraphErrors*)arr->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_C");
572   TGraphErrors* cosmic       =(TGraphErrors*)arr->FindObject("TGRAPHERRORS_MEAN_VDRIFT_COSMICS_ALL");
573   TGraphErrors* cross        =(TGraphErrors*)arr->FindObject("TGRAPHERRORS_VDRIFT_CROSS_ALL");
574   TGraphErrors* itstpcP       =(TGraphErrors*)arr->FindObject("ALIGN_ITSP_TPC_DRIFTVD");
575   TGraphErrors* itstpcM       =(TGraphErrors*)arr->FindObject("ALIGN_ITSM_TPC_DRIFTVD");
576   TGraphErrors* itstpcB       =(TGraphErrors*)arr->FindObject("ALIGN_ITSB_TPC_DRIFTVD");
577   //
578   if (laserA)  SetDefaultGraphDrift(laserA,2,25);
579   if (laserC)  SetDefaultGraphDrift(laserC,4,26);
580   if (cosmic)  SetDefaultGraphDrift(cosmic,3,27);
581   if (cross)   SetDefaultGraphDrift(cross,4,28);
582   if (itstpcP) SetDefaultGraphDrift(itstpcP,2,29);
583   if (itstpcM) SetDefaultGraphDrift(itstpcM,4,30);
584   if (itstpcB) SetDefaultGraphDrift(itstpcB,1,31);
585   //
586   //
587   TPad *pad=0;
588   //
589   // Laser-Laser
590   //
591   if (laserA&&laserC){
592     pad = new TCanvas("TPCLaserVDrift","TPCLaserVDrift");
593     laserA->Draw("alp");
594     SetPadStyle(pad,mx0,mx1,my0,my1);
595     laserA->Draw("apl");
596     laserC->Draw("p");
597     TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
598     legend->AddEntry(laserA,"Laser A side");
599     legend->AddEntry(laserC,"Laser C side");
600     legend->Draw();    
601     picArray->AddLast(pad);
602   }
603
604   if (itstpcP&&itstpcM){
605     pad = new TCanvas("ITSTPC","ITSTPC");
606     itstpcP->Draw("alp");
607     SetPadStyle(pad,mx0,mx1,my0,my1);    
608     itstpcP->Draw("alp");
609     gPad->Clear();
610     itstpcM->Draw("apl");
611     itstpcP->Draw("p");
612     itstpcB->Draw("p");
613     TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
614     legend->AddEntry(itstpcP,"ITS-TPC smooth plus");
615     legend->AddEntry(itstpcM,"ITS-TPC smooth minus");
616     legend->AddEntry(itstpcB,"ITS-TPC smooth ");
617     legend->Draw();    
618     picArray->AddLast(pad);
619   }
620
621   if (itstpcB&&laserA){
622     pad = new TCanvas("ITSTPC_LASER","ITSTPC_LASER");
623     SetPadStyle(pad,mx0,mx1,my0,my1);    
624     laserA->Draw("alp");
625     itstpcP->Draw("p");
626     itstpcM->Draw("p");
627     itstpcB->Draw("p");
628     TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
629     legend->AddEntry(laserA,"TPC laser");
630     legend->AddEntry(itstpcP,"ITS-TPC smooth plus");   
631     legend->AddEntry(itstpcM,"ITS-TPC smooth minus");   
632     legend->AddEntry(itstpcB,"ITS-TPC smooth ");
633     legend->Draw();
634     picArray->AddLast(pad);
635   }
636
637   if (itstpcP&&cross){ 
638     pad = new TCanvas("ITSTPC_CROSS","ITSTPC_CROSS");
639     SetPadStyle(pad,mx0,mx1,my0,my1);    
640     itstpcP->Draw("alp");
641     pad->Clear();
642     cross->Draw("ap");
643     itstpcP->Draw("p");
644     //
645     TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
646
647     legend->AddEntry(cross,"TPC cross tracks");
648     legend->AddEntry(itstpcB,"ITS-TPC smooth");
649     legend->Draw();        
650     picArray->AddLast(pad);
651   }
652   if (itstpcP&&cosmic){ 
653     pad = new TCanvas("ITSTPC_COSMIC","ITSTPC_COSMIC");
654     SetPadStyle(pad,mx0,mx1,my0,my1);    
655     itstpcP->Draw("alp");
656     pad->Clear();
657     cosmic->Draw("ap");
658     itstpcP->Draw("p");
659     //
660     TLegend *legend = new TLegend(mx0+0.01,1-my1-0.2, 0.5, 1-my1-0.01, "Drift velocity correction");
661
662     legend->AddEntry(cosmic,"TPC cross tracks0 up-down");
663     legend->AddEntry(itstpcB,"ITS-TPC smooth");
664     legend->Draw();        
665     picArray->AddLast(pad);
666   }
667 }
668
669
670