]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibTime.cxx
Making drift velocity calibration more robust
[u/mrichter/AliRoot.git] / TPC / AliTPCcalibTime.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 Comments to be written here:
18
19 1. What do we calibrate.
20
21   Time dependence of gain and drift velocity in order to account for changes in: temperature, pressure, gas composition.
22
23   AliTPCcalibTime *calibTime = new AliTPCcalibTime("cosmicTime","cosmicTime",0, 1213.9e+06, 1213.96e+06, 0.04e+04, 0.04e+04);
24
25 */
26
27 #include "Riostream.h"
28 #include "TDatabasePDG.h"
29 #include "TGraphErrors.h"
30 #include "TH1F.h"
31 #include "THnSparse.h"
32 #include "TList.h"
33 #include "TMath.h"
34 #include "TTimeStamp.h"
35 #include "TTree.h"
36 #include "TVectorD.h"
37 //#include "TChain.h"
38 //#include "TFile.h"
39
40 #include "AliDCSSensor.h"
41 #include "AliDCSSensorArray.h"
42 #include "AliESDEvent.h"
43 #include "AliESDInputHandler.h"
44 #include "AliESDVertex.h"
45 #include "AliESDfriend.h"
46 #include "AliLog.h"
47 #include "AliRelAlignerKalman.h"
48 #include "AliTPCCalROC.h"
49 #include "AliTPCParam.h"
50 #include "AliTPCTracklet.h"
51 #include "AliTPCcalibDB.h"
52 #include "AliTPCcalibLaser.h"
53 #include "AliTPCcalibTime.h"
54 #include "AliTPCclusterMI.h"
55 #include "AliTPCseed.h"
56 #include "AliTrackPointArray.h"
57 #include "AliTracker.h"
58 #include "AliKFVertex.h"
59 #include <AliLog.h>
60
61 ClassImp(AliTPCcalibTime)
62
63
64 AliTPCcalibTime::AliTPCcalibTime() 
65   :AliTPCcalibBase(),  
66    fMemoryMode(1), // 0 -do not fill THnSparse with residuals  1- fill only important QA THn 2 - Fill all THnsparse for calibration
67    fLaser(0),       // pointer to laser calibration
68    fDz(0),          // current delta z
69    fCutMaxD(3),        // maximal distance in rfi ditection
70    fCutMaxDz(25),      // maximal distance in rfi ditection
71    fCutTheta(0.03),    // maximal distan theta
72    fCutMinDir(-0.99),  // direction vector products
73    fCutTracks(2500),
74    fArrayLaserA(0),      //laser  fit parameters C
75    fArrayLaserC(0),      //laser  fit parameters A
76    fArrayDz(0),          //NEW! Tmap of V drifts for different triggers
77    fAlignITSTPC(0),      //alignemnt array ITS TPC match
78    fAlignTRDTPC(0),      //alignemnt array TRD TPC match 
79    fAlignTOFTPC(0),      //alignemnt array TOF TPC match
80    fTimeKalmanBin(60*15), //time bin width for kalman - 15 minutes default
81    fTimeBins(0),
82    fTimeStart(0),
83    fTimeEnd(0),
84    fPtBins(0),
85    fPtStart(0),
86    fPtEnd(0),
87    fVdriftBins(0),
88    fVdriftStart(0),
89    fVdriftEnd(0),
90    fRunBins(0),
91    fRunStart(0),
92    fRunEnd(0)
93 {  
94   //
95   // default constructor
96   //
97   AliInfo("Default Constructor");  
98   for (Int_t i=0;i<3;i++) {
99     fHistVdriftLaserA[i]=0;
100     fHistVdriftLaserC[i]=0;
101   }
102   for (Int_t i=0;i<10;i++) {
103     fCosmiMatchingHisto[i]=0;
104   }
105   //
106   for (Int_t i=0;i<5;i++) {
107     fResHistoTPCCE[i]=0;
108     fResHistoTPCITS[i]=0;
109     fResHistoTPCTRD[i]=0;
110     fResHistoTPCTOF[i]=0;
111     fResHistoTPCvertex[i]=0;
112     fTPCVertex[i]=0;
113   }
114   for (Int_t i=0;i<12;i++) {
115     fTPCVertex[i]=0;
116   }
117   for (Int_t i=0;i<5;i++) {
118     fTPCVertexCorrelation[i]=0;
119   }
120   static Int_t counter=0;
121   if (1) {
122     TTimeStamp s;
123     Int_t time=s;
124     AliInfo(Form("Counter Constructor\t%d\t%d",counter,time));
125     counter++;
126   }
127
128 }
129
130 AliTPCcalibTime::AliTPCcalibTime(const Text_t *name, const Text_t *title, UInt_t StartTime, UInt_t EndTime, Int_t deltaIntegrationTimeVdrift, Int_t memoryMode)
131   :AliTPCcalibBase(),
132    fMemoryMode(memoryMode), // 0 -do not fill THnSparse with residuals  1- fill only important QA THn 2 - Fill all THnsparse for calibration
133    fLaser(0),            // pointer to laser calibration
134    fDz(0),               // current delta z
135    fCutMaxD(5*0.5356),   // maximal distance in rfi ditection
136    fCutMaxDz(40),   // maximal distance in rfi ditection
137    fCutTheta(5*0.004644),// maximal distan theta
138    fCutMinDir(-0.99),    // direction vector products
139    fCutTracks(2500),
140    fArrayLaserA(new TObjArray(1000)),      //laser  fit parameters C
141    fArrayLaserC(new TObjArray(1000)),      //laser  fit parameters A
142    fArrayDz(0),            //Tmap of V drifts for different triggers
143    fAlignITSTPC(0),      //alignemnt array ITS TPC match
144    fAlignTRDTPC(0),      //alignemnt array TRD TPC match 
145    fAlignTOFTPC(0),      //alignemnt array TOF TPC match
146    fTimeKalmanBin(60*15), //time bin width for kalman - 15 minutes default
147    fTimeBins(0),
148    fTimeStart(0),
149    fTimeEnd(0),
150    fPtBins(0),
151    fPtStart(0),
152    fPtEnd(0),
153    fVdriftBins(0),
154    fVdriftStart(0),
155    fVdriftEnd(0),
156    fRunBins(0),
157    fRunStart(0),
158    fRunEnd(0)
159 {
160   //
161   // Non deafaul constructor - to be used in the Calibration setups 
162   //
163
164   SetName(name);
165   SetTitle(title);
166   for (Int_t i=0;i<3;i++) {
167     fHistVdriftLaserA[i]=0;
168     fHistVdriftLaserC[i]=0;
169   }
170
171   for (Int_t i=0;i<5;i++) {
172     fResHistoTPCCE[i]=0;
173     fResHistoTPCITS[i]=0;
174     fResHistoTPCTRD[i]=0;
175     fResHistoTPCTOF[i]=0;
176     fResHistoTPCvertex[i]=0;
177   }
178
179
180   AliInfo("Non Default Constructor");
181   fTimeBins   =(EndTime-StartTime)/deltaIntegrationTimeVdrift;
182   fTimeStart  =StartTime; //(((TObjString*)(mapGRP->GetValue("fAliceStartTime")))->GetString()).Atoi();
183   fTimeEnd    =EndTime;   //(((TObjString*)(mapGRP->GetValue("fAliceStopTime")))->GetString()).Atoi();
184   fPtBins     = 400;
185   fPtStart    = -0.04;
186   fPtEnd      =  0.04;
187   fVdriftBins = 500;
188   fVdriftStart= -0.1;
189   fVdriftEnd  =  0.1;
190   fRunBins    = 1000001;
191   fRunStart   = -1.5;
192   fRunEnd     = 999999.5;
193
194   Int_t    binsVdriftLaser[4] = {fTimeBins , fPtBins , fVdriftBins*20, fRunBins };
195   Double_t xminVdriftLaser[4] = {fTimeStart, fPtStart, fVdriftStart  , fRunStart};
196   Double_t xmaxVdriftLaser[4] = {fTimeEnd  , fPtEnd  , fVdriftEnd    , fRunEnd  };
197   TString axisTitle[4]={
198     "T",
199     "#delta_{P/T}",
200     "value",
201     "run"
202   };
203   TString histoName[3]={
204     "Loffset",
205     "Lcorr",
206     "Lgy"
207   };
208
209   
210   for (Int_t i=0;i<3;i++) {
211     fHistVdriftLaserA[i] = new THnSparseF("HistVdriftLaser","HistVdriftLaser;time;p/T ratio;Vdrift;run",4,binsVdriftLaser,xminVdriftLaser,xmaxVdriftLaser);
212     fHistVdriftLaserC[i] = new THnSparseF("HistVdriftLaser","HistVdriftLaser;time;p/T ratio;Vdrift;run",4,binsVdriftLaser,xminVdriftLaser,xmaxVdriftLaser);
213     fHistVdriftLaserA[i]->SetName(histoName[i]);
214     fHistVdriftLaserC[i]->SetName(histoName[i]);
215     for (Int_t iaxis=0; iaxis<4;iaxis++){
216       fHistVdriftLaserA[i]->GetAxis(iaxis)->SetName(axisTitle[iaxis]);
217       fHistVdriftLaserC[i]->GetAxis(iaxis)->SetName(axisTitle[iaxis]);
218     }
219   }
220   fBinsVdrift[0] = fTimeBins;
221   fBinsVdrift[1] = fPtBins;
222   fBinsVdrift[2] = fVdriftBins;
223   fBinsVdrift[3] = fRunBins;
224   fXminVdrift[0] = fTimeStart;
225   fXminVdrift[1] = fPtStart;
226   fXminVdrift[2] = fVdriftStart;
227   fXminVdrift[3] = fRunStart;
228   fXmaxVdrift[0] = fTimeEnd;
229   fXmaxVdrift[1] = fPtEnd;
230   fXmaxVdrift[2] = fVdriftEnd;
231   fXmaxVdrift[3] = fRunEnd;
232
233   fArrayDz=new TObjArray();
234   fAlignITSTPC = new TObjArray;      //alignemnt array ITS TPC match
235   fAlignTRDTPC = new TObjArray;      //alignemnt array ITS TPC match
236   fAlignTOFTPC = new TObjArray;      //alignemnt array ITS TPC match
237   fAlignITSTPC->SetOwner(kTRUE);
238   fAlignTRDTPC->SetOwner(kTRUE);
239   fAlignTOFTPC->SetOwner(kTRUE);
240   
241
242   fCosmiMatchingHisto[0]=new TH1F("Cosmics matching","p0-all"   ,100,-10*0.5356  ,10*0.5356  );
243   fCosmiMatchingHisto[1]=new TH1F("Cosmics matching","p1-all"   ,100,-10*4.541   ,10*4.541   );
244   fCosmiMatchingHisto[2]=new TH1F("Cosmics matching","p2-all"   ,100,-10*0.01134 ,10*0.01134 );
245   fCosmiMatchingHisto[3]=new TH1F("Cosmics matching","p3-all"   ,100,-10*0.004644,10*0.004644);
246   fCosmiMatchingHisto[4]=new TH1F("Cosmics matching","p4-all"   ,100,-10*0.03773 ,10*0.03773 );
247   fCosmiMatchingHisto[5]=new TH1F("Cosmics matching","p0-isPair",100,-10*0.5356  ,10*0.5356  );
248   fCosmiMatchingHisto[6]=new TH1F("Cosmics matching","p1-isPair",100,-10*4.541   ,10*4.541   );
249   fCosmiMatchingHisto[7]=new TH1F("Cosmics matching","p2-isPair",100,-10*0.01134 ,10*0.01134 );
250   fCosmiMatchingHisto[8]=new TH1F("Cosmics matching","p3-isPair",100,-10*0.004644,10*0.004644);
251   fCosmiMatchingHisto[9]=new TH1F("Cosmics matching","p4-isPair",100,-10*0.03773 ,10*0.03773 );
252   for (Int_t i=0;i<12;i++) {
253     fTPCVertex[i]=0;
254   }
255   for (Int_t i=0;i<5;i++) {
256     fTPCVertexCorrelation[i]=0;
257   }
258   BookDistortionMaps();
259   
260 }
261
262 AliTPCcalibTime::~AliTPCcalibTime(){
263   //
264   // Virtual Destructor
265   //
266   static Int_t counter=0;
267   if (1) {
268     TTimeStamp s;
269     Int_t time=s;
270     AliInfo(Form("Counter Destructor\t%s\t%d\t%d",GetName(),counter,time));
271     counter++;
272   }
273   for(Int_t i=0;i<3;i++){
274     if(fHistVdriftLaserA[i]){
275       delete fHistVdriftLaserA[i];
276       fHistVdriftLaserA[i]=NULL;
277     }
278     if(fHistVdriftLaserC[i]){
279       delete fHistVdriftLaserC[i];
280       fHistVdriftLaserC[i]=NULL;
281     }
282   }
283   if(fArrayDz){
284     fArrayDz->SetOwner();
285     fArrayDz->Delete();
286     delete fArrayDz;
287     fArrayDz=NULL;
288   }
289   for(Int_t i=0;i<5;i++){
290     if(fCosmiMatchingHisto[i]){
291       delete fCosmiMatchingHisto[i];
292       fCosmiMatchingHisto[i]=NULL;
293     }
294   }
295
296   for (Int_t i=0;i<5;i++) {
297     delete fResHistoTPCCE[i];
298     delete fResHistoTPCITS[i];
299     delete fResHistoTPCTRD[i];
300     delete fResHistoTPCTOF[i];
301     delete fResHistoTPCvertex[i];
302     fResHistoTPCCE[i]=0;
303     fResHistoTPCITS[i]=0;
304     fResHistoTPCTRD[i]=0;
305     fResHistoTPCTOF[i]=0;
306     fResHistoTPCvertex[i]=0;
307   }
308
309   if (fTPCVertex[0]) {
310     for (Int_t i=0;i<12;i++)  delete fTPCVertex[i];
311   }
312   if (fTPCVertexCorrelation[0]) {
313     for (Int_t i=0;i<5;i++)  delete fTPCVertexCorrelation[i];
314   }
315
316   if (fAlignITSTPC){
317     fAlignITSTPC->SetOwner(kTRUE);
318     fAlignTRDTPC->SetOwner(kTRUE);
319     fAlignTOFTPC->SetOwner(kTRUE);
320     
321     fAlignITSTPC->Delete();
322     fAlignTRDTPC->Delete();
323     fAlignTOFTPC->Delete();
324     delete fAlignITSTPC;
325     delete fAlignTRDTPC;
326     delete fAlignTOFTPC;
327   }
328 }
329
330 // Bool_t AliTPCcalibTime::IsLaser(const AliESDEvent *const /*event*/) const{
331 //   //
332 //   // Indicator is laser event not yet implemented  - to be done using trigger info or event specie
333 //   //
334 //   return kTRUE; //More accurate creteria to be added
335 // }
336 // Bool_t AliTPCcalibTime::IsCosmics(const AliESDEvent *const /*event*/){
337 //   //
338 //   // Indicator is cosmic event not yet implemented - to be done using trigger info or event specie
339 //   //
340
341 //   return kTRUE; //More accurate creteria to be added
342 // }
343 // Bool_t AliTPCcalibTime::IsBeam(const AliESDEvent *const /*event*/) const{
344 //   //
345 //   // Indicator is physic event not yet implemented - to be done using trigger info or event specie
346 //   //
347
348 //   return kTRUE; //More accurate creteria to be added
349 // }
350 void AliTPCcalibTime::ResetCurrent(){
351   //
352   //ResetCurrent
353   //
354   fDz=0; //Reset current dz
355 }
356
357
358
359 void AliTPCcalibTime::Process(AliESDEvent *event){
360   //
361   // main function to make calibration
362   //
363   if(!event) return;
364   if (event->GetNumberOfTracks()<2) return; 
365   AliESDfriend *ESDfriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
366   if (!ESDfriend) {
367     return;
368   }
369   if (ESDfriend->TestSkipBit()) return;
370   
371   ResetCurrent();
372   //if(IsLaser  (event)) 
373   ProcessLaser (event);
374   //if(IsCosmics(event)) 
375   ProcessCosmic(event);
376   //if(IsBeam   (event)) 
377   ProcessBeam  (event);
378 }
379
380 void AliTPCcalibTime::ProcessLaser(AliESDEvent *event){
381   //
382   // Fit drift velocity using laser 
383   // 
384   // 0. cuts
385   const Int_t    kMinTracks     = 40;    // minimal number of laser tracks
386   const Int_t    kMinTracksSide = 20;    // minimal number of tracks per side
387   const Float_t  kMaxDeltaZ     = 30.;   // maximal trigger delay
388   const Float_t  kMaxDeltaV     = 0.05;  // maximal deltaV 
389   const Float_t  kMaxRMS        = 0.1;   // maximal RMS of tracks
390   //
391   /*
392     TCut cutRMS("sqrt(laserA.fElements[4])<0.1&&sqrt(laserC.fElements[4])<0.1");
393     TCut cutZ("abs(laserA.fElements[0]-laserC.fElements[0])<3");
394     TCut cutV("abs(laserA.fElements[1]-laserC.fElements[1])<0.01");
395     TCut cutY("abs(laserA.fElements[2]-laserC.fElements[2])<2");
396     TCut cutAll = cutRMS+cutZ+cutV+cutY;
397   */
398   if (event->GetNumberOfTracks()<kMinTracks) return;
399   //
400   if(!fLaser) fLaser = new AliTPCcalibLaser("laserTPC","laserTPC",kFALSE);
401   fLaser->Process(event);
402   if (fLaser->GetNtracks()<kMinTracks) return;   // small amount of tracks cut
403   if (fLaser->fFitAside->GetNrows()==0  && fLaser->fFitCside->GetNrows()==0) return;  // no fit neither a or C side
404   //
405   // debug streamer  - activate stream level
406   // Use it for tuning of the cuts
407   //
408   // cuts to be applied
409   //
410   Int_t isReject[2]={0,0};
411   //
412   // not enough tracks 
413   if (TMath::Abs((*fLaser->fFitAside)[3]) < kMinTracksSide) isReject[0]|=1; 
414   if (TMath::Abs((*fLaser->fFitCside)[3]) < kMinTracksSide) isReject[1]|=1; 
415   // unreasonable z offset
416   if (TMath::Abs((*fLaser->fFitAside)[0])>kMaxDeltaZ)  isReject[0]|=2;
417   if (TMath::Abs((*fLaser->fFitCside)[0])>kMaxDeltaZ)  isReject[1]|=2;
418   // unreasonable drift velocity
419   if (TMath::Abs((*fLaser->fFitAside)[1]-1)>kMaxDeltaV)  isReject[0]|=4;
420   if (TMath::Abs((*fLaser->fFitCside)[1]-1)>kMaxDeltaV)  isReject[1]|=4;
421   // big chi2
422   if (TMath::Sqrt(TMath::Abs((*fLaser->fFitAside)[4]))>kMaxRMS ) isReject[0]|=8;
423   if (TMath::Sqrt(TMath::Abs((*fLaser->fFitCside)[4]))>kMaxRMS ) isReject[1]|=8;
424
425
426
427
428   if (fStreamLevel>0){
429     printf("Trigger: %s\n",event->GetFiredTriggerClasses().Data());
430
431     TTreeSRedirector *cstream = GetDebugStreamer();
432     if (cstream){
433       TTimeStamp tstamp(fTime);
434       (*cstream)<<"laserInfo"<<
435         "run="<<fRun<<              //  run number
436         "event="<<fEvent<<          //  event number
437         "time="<<fTime<<            //  time stamp of event
438         "trigger="<<fTrigger<<      //  trigger
439         "mag="<<fMagF<<             //  magnetic field
440         //laser
441         "rejectA="<<isReject[0]<<
442         "rejectC="<<isReject[1]<<
443         "laserA.="<<fLaser->fFitAside<<
444         "laserC.="<<fLaser->fFitCside<<
445         "laserAC.="<<fLaser->fFitACside<<
446         "trigger="<<event->GetFiredTriggerClasses()<<
447         "\n";
448     }
449   }
450   //
451   // fill histos
452   //
453   TVectorD vdriftA(5), vdriftC(5),vdriftAC(6);
454   vdriftA=*(fLaser->fFitAside);
455   vdriftC=*(fLaser->fFitCside);
456   vdriftAC=*(fLaser->fFitACside);
457   Int_t npointsA=0, npointsC=0;
458   Float_t chi2A=0, chi2C=0;
459   npointsA= TMath::Nint(vdriftA[3]);
460   chi2A= vdriftA[4];
461   npointsC= TMath::Nint(vdriftC[3]);
462   chi2C= vdriftC[4];
463
464   if (npointsA>kMinTracksSide || npointsC>kMinTracksSide){
465     TVectorD *fitA = new TVectorD(6);
466     TVectorD *fitC = new TVectorD(6);
467     for (Int_t ipar=0; ipar<5; ipar++){
468       (*fitA)[ipar]=vdriftA[ipar];
469       (*fitC)[ipar]=vdriftC[ipar];
470     }
471     (*fitA)[5]=fTime;
472     (*fitC)[5]=fTime;
473     fArrayLaserA->AddLast(fitA);
474     fArrayLaserC->AddLast(fitC);
475   }
476   //
477   
478   TTimeStamp tstamp(fTime);
479   Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
480   Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
481   Double_t driftA=0, driftC=0;
482   if (vdriftA[1]>1.-kMaxDeltaV) driftA = 1./vdriftA[1]-1.;
483   if (vdriftC[1]>1.-kMaxDeltaV) driftC = 1./vdriftC[1]-1.;
484   //
485   Double_t vecDriftLaserA[4]={fTime,(ptrelative0+ptrelative1)/2.0,driftA,event->GetRunNumber()};
486   Double_t vecDriftLaserC[4]={fTime,(ptrelative0+ptrelative1)/2.0,driftC,event->GetRunNumber()};
487   //  Double_t vecDrift[4]      ={fTime,(ptrelative0+ptrelative1)/2.0,1./((*(fLaser->fFitACside))[1])-1,event->GetRunNumber()};
488
489   for (Int_t icalib=0;icalib<3;icalib++){
490     if (icalib==0){ //z0 shift
491       vecDriftLaserA[2]=vdriftA[0]/250.;
492       vecDriftLaserC[2]=vdriftC[0]/250.;
493     }
494     if (icalib==1){ //vdrel shift
495       vecDriftLaserA[2]=driftA;
496       vecDriftLaserC[2]=driftC;
497     }
498     if (icalib==2){ //gy shift - full gy - full drift
499       vecDriftLaserA[2]=vdriftA[2]/250.;
500       vecDriftLaserC[2]=vdriftC[2]/250.;
501     }
502     //if (isReject[0]==0) fHistVdriftLaserA[icalib]->Fill(vecDriftLaserA);
503     //if (isReject[1]==0) fHistVdriftLaserC[icalib]->Fill(vecDriftLaserC);
504     fHistVdriftLaserA[icalib]->Fill(vecDriftLaserA);
505     fHistVdriftLaserC[icalib]->Fill(vecDriftLaserC);
506   }
507 }
508
509 void AliTPCcalibTime::ProcessCosmic(const AliESDEvent *const event){
510   //
511   // process Cosmic event - track matching A side C side
512   //
513   if (!event) {
514     Printf("ERROR: ESD not available");
515     return;
516   }  
517   if (event->GetTimeStamp() == 0 ) {
518     Printf("no time stamp!");
519     return;
520   }
521   
522   //fd
523   // Find cosmic pairs
524   // 
525   // Track0 is choosen in upper TPC part
526   // Track1 is choosen in lower TPC part
527   //
528   const Int_t kMinClustersCross =30;
529   const Int_t kMinClusters      =80;
530   Int_t ntracks=event->GetNumberOfTracks();
531   if (ntracks==0) return;
532   if (ntracks > fCutTracks) return;
533   
534   if (GetDebugLevel()>20) printf("Hallo world: Im here\n");
535   AliESDfriend *esdFriend=(AliESDfriend*)(((AliESDEvent*)event)->FindListObject("AliESDfriend"));
536   
537   TObjArray  tpcSeeds(ntracks);
538   Double_t vtxx[3]={0,0,0};
539   Double_t svtxx[3]={0.000001,0.000001,100.};
540   AliESDVertex vtx(vtxx,svtxx);
541   //
542   // track loop
543   //
544   TArrayI clusterSideA(ntracks);
545   TArrayI clusterSideC(ntracks);
546   for (Int_t i=0;i<ntracks;++i) {
547     clusterSideA[i]=0;
548     clusterSideC[i]=0;
549     AliESDtrack *track = event->GetTrack(i);
550     
551     const AliExternalTrackParam * trackIn = track->GetInnerParam();
552     const AliExternalTrackParam * trackOut = track->GetOuterParam();
553     if (!trackIn) continue;
554     if (!trackOut) continue;
555     
556     AliESDfriendTrack *friendTrack = esdFriend->GetTrack(i);
557     if (!friendTrack) continue;
558     if (friendTrack) ProcessSame(track,friendTrack,event);
559     if (friendTrack) ProcessAlignITS(track,friendTrack,event,esdFriend);
560     if (friendTrack) ProcessAlignTRD(track,friendTrack);
561     if (friendTrack) ProcessAlignTOF(track,friendTrack);
562     TObject *calibObject;
563     AliTPCseed *seed = 0;
564     for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
565     if (seed) {
566       tpcSeeds.AddAt(seed,i);
567       Int_t nA=0, nC=0;
568       for (Int_t irow=159;irow>0;irow--) {
569         AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
570         if (!cl) continue;
571         if ((cl->GetDetector()%36)<18) nA++;
572         if ((cl->GetDetector()%36)>=18) nC++;
573       }
574       clusterSideA[i]=nA;
575       clusterSideC[i]=nC;
576     }
577   }
578   if (ntracks<2) return;
579   //
580   // Find pairs
581   //
582
583   for (Int_t i=0;i<ntracks;++i) {
584     AliESDtrack *track0 = event->GetTrack(i);
585     // track0 - choosen upper part
586     if (!track0) continue;
587     if (!track0->GetOuterParam()) continue;
588     if (track0->GetOuterParam()->GetAlpha()<0) continue;
589     Double_t d1[3];
590     track0->GetDirection(d1);    
591     for (Int_t j=0;j<ntracks;++j) {
592       if (i==j) continue;
593       AliESDtrack *track1 = event->GetTrack(j);   
594       //track 1 lower part
595       if (!track1) continue;
596       if (!track1->GetOuterParam()) continue;
597       if (track0->GetTPCNcls()+ track1->GetTPCNcls()< kMinClusters) continue;
598       Int_t nAC = TMath::Max( TMath::Min(clusterSideA[i], clusterSideC[j]), 
599                               TMath::Min(clusterSideC[i], clusterSideA[j]));
600       if (nAC<kMinClustersCross) continue; 
601       Int_t nA0=clusterSideA[i];
602       Int_t nC0=clusterSideC[i];
603       Int_t nA1=clusterSideA[j];
604       Int_t nC1=clusterSideC[j];
605       //      if (track1->GetOuterParam()->GetAlpha()>0) continue;
606       //
607       Double_t d2[3];
608       track1->GetDirection(d2);
609       
610       AliTPCseed * seed0 = (AliTPCseed*) tpcSeeds.At(i);
611       AliTPCseed * seed1 = (AliTPCseed*) tpcSeeds.At(j);
612       if (! seed0) continue;
613       if (! seed1) continue;
614       Float_t dir = (d1[0]*d2[0] + d1[1]*d2[1] + d1[2]*d2[2]);
615       Float_t dist0  = track0->GetLinearD(0,0);
616       Float_t dist1  = track1->GetLinearD(0,0);
617       //
618       // conservative cuts - convergence to be guarantied
619       // applying before track propagation
620       if (TMath::Abs(TMath::Abs(dist0)-TMath::Abs(dist1))>fCutMaxD) continue;   // distance to the 0,0
621       if (TMath::Abs(dir)<TMath::Abs(fCutMinDir)) continue;               // direction vector product
622       Float_t bz = AliTracker::GetBz();
623       Float_t dvertex0[2];   //distance to 0,0
624       Float_t dvertex1[2];   //distance to 0,0 
625       track0->GetDZ(0,0,0,bz,dvertex0);
626       track1->GetDZ(0,0,0,bz,dvertex1);
627       if (TMath::Abs(dvertex0[1])>250) continue;
628       if (TMath::Abs(dvertex1[1])>250) continue;
629       //
630       //
631       //
632       Float_t dmax = TMath::Max(TMath::Abs(dist0),TMath::Abs(dist1));
633       AliExternalTrackParam param0(*track0);
634       AliExternalTrackParam param1(*track1);
635       //
636       // Propagate using Magnetic field and correct fo material budget
637       //
638       AliTracker::PropagateTrackTo(&param0,dmax+1,TDatabasePDG::Instance()->GetParticle("e-")->Mass(),3,kTRUE);
639       AliTracker::PropagateTrackTo(&param1,dmax+1,TDatabasePDG::Instance()->GetParticle("e-")->Mass(),3,kTRUE);
640       //
641       // Propagate rest to the 0,0 DCA - z should be ignored
642       //
643       //Bool_t b0 = ;
644       param0.PropagateToDCA(&vtx,bz,1000);
645       //Bool_t b1 = 
646       param1.PropagateToDCA(&vtx,bz,1000);
647       param0.GetDZ(0,0,0,bz,dvertex0);
648       param1.GetDZ(0,0,0,bz,dvertex1);
649       Double_t xyz0[3];
650       Double_t xyz1[3];
651       param0.GetXYZ(xyz0);
652       param1.GetXYZ(xyz1);
653       Bool_t isPair = IsPair(&param0,&param1);
654       Bool_t isCross = IsCross(track0, track1);
655       Bool_t isSame = IsSame(track0, track1);
656
657       THnSparse* hist=new THnSparseF("","HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
658       TString shortName=hist->ClassName();
659       shortName+="_MEAN_VDRIFT_COSMICS_";
660       delete hist;
661       hist=NULL;
662
663       if((isSame) || (isCross && isPair)){
664         if (track0->GetTPCNcls()+ track1->GetTPCNcls()> 80) {
665           fDz = param0.GetZ() - param1.GetZ();
666           Double_t sign=(nA0>nA1)? 1:-1; 
667           fDz*=sign;
668           TTimeStamp tstamp(fTime);
669           Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
670           Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
671           Double_t vecDrift[4]={fTime,(ptrelative0+ptrelative1)/2.0,fDz/500.0,event->GetRunNumber()};
672           THnSparse* curHist=NULL;
673           TString name="";
674
675           name=shortName;
676           name+=event->GetFiredTriggerClasses();
677           name.ToUpper();
678           curHist=(THnSparseF*)fArrayDz->FindObject(name);
679           if(!curHist){
680             curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
681             fArrayDz->AddLast(curHist);
682           }
683 //        curHist=(THnSparseF*)(fMapDz->GetValue(event->GetFiredTriggerClasses()));
684 //        if(!curHist){
685 //          curHist=new THnSparseF(event->GetFiredTriggerClasses(),"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
686 //          fMapDz->Add(new TObjString(event->GetFiredTriggerClasses()),curHist);
687 //        }
688           curHist->Fill(vecDrift);
689           
690           name=shortName;
691           name+="ALL";
692           name.ToUpper();
693           curHist=(THnSparseF*)fArrayDz->FindObject(name);
694           if(!curHist){
695             curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
696             fArrayDz->AddLast(curHist);
697           }
698 //        curHist=(THnSparseF*)(fMapDz->GetValue("all"));
699 //        if(!curHist){
700 //          curHist=new THnSparseF("all","HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
701 //          fMapDz->Add(new TObjString("all"),curHist);
702 //        }
703           curHist->Fill(vecDrift);
704         }
705       }
706       TTreeSRedirector *cstream = GetDebugStreamer();
707       if (fStreamLevel>0){
708         if (cstream){
709         (*cstream)<<"trackInfo"<<
710           "tr0.="<<track0<<
711           "tr1.="<<track1<<
712           "p0.="<<&param0<<
713           "p1.="<<&param1<<
714           "nAC="<<nAC<<
715           "nA0="<<nA0<<
716           "nA1="<<nA1<<
717           "nC0="<<nC0<<
718           "nC1="<<nC1<<
719           "isPair="<<isPair<<
720           "isCross="<<isCross<<
721           "isSame="<<isSame<<
722           "fDz="<<fDz<<
723           "fRun="<<fRun<<
724           "fTime="<<fTime<<
725           "\n";
726         }
727       }
728     } // end 2nd order loop        
729   } // end 1st order loop
730   
731   if (fStreamLevel>0){
732     TTreeSRedirector *cstream = GetDebugStreamer();
733     if (cstream){
734       (*cstream)<<"timeInfo"<<
735         "run="<<fRun<<              //  run number
736         "event="<<fEvent<<          //  event number
737         "time="<<fTime<<            //  time stamp of event
738         "trigger="<<fTrigger<<      //  trigger
739         "mag="<<fMagF<<             //  magnetic field
740         // Environment values
741         //
742         // accumulated values
743         //
744         "fDz="<<fDz<<          //! current delta z
745         "trigger="<<event->GetFiredTriggerClasses()<<
746         "\n";
747     }
748   }
749   if (GetDebugLevel()>20) printf("Trigger: %s\n",event->GetFiredTriggerClasses().Data());
750 }
751
752 void AliTPCcalibTime::ProcessBeam(const AliESDEvent *const event){
753   //
754   // Process beam data - calculates vartex
755   //                     from A side and C side
756   // Histogram the differences                 
757   // 
758   const Int_t kMinClusters  =80;
759   const Int_t kMinTracks    =2;      // minimal number of tracks to define the vertex
760   const Int_t kMinTracksVertex=30;   // minimal number of tracks to define the cumulative vertex
761   const Double_t kMaxTgl    =1.2;    // maximal Tgl (z angle)
762   const Double_t kMinPt     =0.2;    // minimal pt
763   const Double_t kMaxD0     =5.;     // cut on distance to the primary vertex first guess
764   const Double_t kMaxZ0     =20; 
765   const Double_t kMaxD      =2.5;    // cut on distance to the primary vertex 
766   const Double_t kMaxZ      =4;      // maximal z distance between tracks form the same side
767   const Double_t kMaxChi2   =15;     // maximal chi2 of the TPCvertex 
768   const Double_t kCumulCovarXY=0.003; //increase the error of cumul vertex 30 microns profile
769   const Double_t kCumulCovarZ=250.;  //increase the error of cumul vertex
770   const Double_t kMaxDvertex = 1.0;  // cut to accept the vertex; 
771   //
772   Int_t  flags=0;
773   const  Int_t  kBuffSize=100;
774   static Double_t deltaZ[kBuffSize]={0};
775   static Int_t counterZ=0;
776   static AliKFVertex cumulVertexA, cumulVertexC, cumulVertexAC; // cumulative vertex 
777   AliKFVertex vertexA, vertexC;
778
779   Float_t  dca0[2]={0,0};
780   Double_t dcaVertex[2]={0,0};
781   Int_t ntracks=event->GetNumberOfTracks();
782   if (ntracks==0) return;
783   if (ntracks > fCutTracks) return;
784   //
785   AliESDfriend *esdFriend=(AliESDfriend*)(((AliESDEvent*)event)->FindListObject("AliESDfriend"));
786   //
787   // Divide tracks to A and C side tracks - using the cluster indexes
788   TObjArray tracksA(ntracks);  
789   TObjArray tracksC(ntracks);  
790   //
791   AliESDVertex *vertexSPD =  (AliESDVertex *)event->GetPrimaryVertexSPD();
792   AliESDVertex *vertex    =  (AliESDVertex *)event->GetPrimaryVertex();
793   AliESDVertex *vertexTracks =  (AliESDVertex *)event->GetPrimaryVertexTracks();
794   Double_t vertexZA[10000], vertexZC[10000];
795   //
796   Int_t ntracksA= 0;
797   Int_t ntracksC= 0;
798   //
799   for (Int_t itrack=0;itrack<ntracks;itrack++) {
800     AliESDtrack *track = event->GetTrack(itrack);
801     AliESDfriendTrack *friendTrack = esdFriend->GetTrack(itrack);
802     if (!friendTrack) continue;
803     if (TMath::Abs(track->GetTgl())>kMaxTgl) continue;
804     if (TMath::Abs(track->Pt())<kMinPt) continue;
805     const AliExternalTrackParam * trackIn  = track->GetInnerParam();
806     TObject *calibObject=0;
807     AliTPCseed *seed = 0;
808     Int_t nA=0, nC=0;
809     for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
810     if (seed) {
811       for (Int_t irow=159;irow>0;irow--) {
812         AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
813         if (!cl) continue;
814         if ((cl->GetDetector()%36)<18) nA++;
815         if ((cl->GetDetector()%36)>=18) nC++;
816       }
817       if ((nA>kMinClusters || nC>kMinClusters) && (nA*nC==0) ){
818         track->GetImpactParameters(dca0[0],dca0[1]);
819         if (TMath::Abs(dca0[0])>kMaxD0) continue;
820         if (TMath::Abs(dca0[1])>kMaxZ0) continue;
821         AliExternalTrackParam pTPCvertex(*trackIn);
822         if (!AliTracker::PropagateTrackToBxByBz(&pTPCvertex,4.+4.*TMath::Abs(dca0[0]),0.1,2,kTRUE)) continue;
823         pTPCvertex.PropagateToDCA(vertex,AliTracker::GetBz(), kMaxD, dcaVertex,0);
824         if (TMath::Abs(dcaVertex[0])>kMaxD) continue;
825         if (nA>kMinClusters &&nC==0) { tracksA.AddLast(pTPCvertex.Clone()); vertexZA[ntracksA++] = pTPCvertex.GetZ();}
826         if (nC>kMinClusters &&nA==0) {tracksC.AddLast(pTPCvertex.Clone());  vertexZC[ntracksC++] = pTPCvertex.GetZ();}
827       }
828     }
829   }
830   Double_t medianZA=TMath::Median(ntracksA, vertexZA);  // tracks median
831   Double_t medianZC=TMath::Median(ntracksC, vertexZC);  // tracks median
832   //
833   ntracksA= tracksA.GetEntriesFast();
834   ntracksC= tracksC.GetEntriesFast();
835   if (ntracksA>kMinTracks && ntracksC>kMinTracks){
836     deltaZ[counterZ%kBuffSize]=medianZA-medianZC;
837     counterZ+=1;
838     Double_t medianDelta=(counterZ>=kBuffSize)? TMath::Median(kBuffSize, deltaZ): TMath::Median(counterZ, deltaZ);
839     if (TMath::Abs(medianDelta-(medianZA-medianZC))>kMaxZ) flags+=16;
840     // increse the error of cumulative vertex at the beginning of event
841     cumulVertexA.Covariance(0,0)+=kCumulCovarXY*kCumulCovarXY;
842     cumulVertexA.Covariance(1,1)+=kCumulCovarXY*kCumulCovarXY;
843     cumulVertexA.Covariance(2,2)+=kCumulCovarZ*kCumulCovarZ;
844     cumulVertexC.Covariance(0,0)+=kCumulCovarXY*kCumulCovarXY;
845     cumulVertexC.Covariance(1,1)+=kCumulCovarXY*kCumulCovarXY;
846     cumulVertexC.Covariance(2,2)+=kCumulCovarZ*kCumulCovarZ;
847     cumulVertexAC.Covariance(0,0)+=kCumulCovarXY*kCumulCovarXY;
848     cumulVertexAC.Covariance(1,1)+=kCumulCovarXY*kCumulCovarXY;
849     cumulVertexAC.Covariance(2,2)+=kCumulCovarZ*kCumulCovarZ;
850     //
851     for (Int_t iA=0; iA<ntracksA; iA++){
852       if (flags!=0) continue;
853       AliExternalTrackParam *aliTrack =  (AliExternalTrackParam *)tracksA.At(iA);
854       if (TMath::Abs(aliTrack->GetZ()-medianZA)>kMaxZ) continue;
855       AliKFParticle part(*aliTrack,211);
856       vertexA+=part;
857     }   
858     for (Int_t iC=0; iC<ntracksC; iC++){
859       if (flags!=0) continue;
860       AliExternalTrackParam *aliTrack =  (AliExternalTrackParam *)tracksC.At(iC);
861       if (TMath::Abs(aliTrack->GetZ()-medianZC)>kMaxZ) continue;
862       AliKFParticle part(*aliTrack,211);
863       vertexC+=part;
864     }    
865     //
866     if (vertexA.GetNDF()<kMinTracks) flags+=32;
867     if (vertexC.GetNDF()<kMinTracks) flags+=32;
868     if (TMath::Abs(vertexA.Z()-medianZA)>kMaxZ) flags+=1;   //apply cuts
869     if (TMath::Abs(vertexC.Z()-medianZC)>kMaxZ) flags+=2;
870     if (TMath::Abs(vertexA.GetChi2()/vertexA.GetNDF()+vertexC.GetChi2()/vertexC.GetNDF())> kMaxChi2) flags+=4;
871     //
872     if (flags==0){
873       for (Int_t iA=0; iA<ntracksA; iA++){
874         if (flags!=0) continue;
875         AliExternalTrackParam *aliTrack =  (AliExternalTrackParam *)tracksA.At(iA);
876         if (TMath::Abs(aliTrack->GetZ()-medianZA)>kMaxZ) continue;
877         AliKFParticle part(*aliTrack,211);
878         cumulVertexA+=part;
879         cumulVertexAC+=part;
880       } 
881       for (Int_t iC=0; iC<ntracksC; iC++){
882         if (flags!=0) continue;
883         AliExternalTrackParam *aliTrack =  (AliExternalTrackParam *)tracksC.At(iC);
884         if (TMath::Abs(aliTrack->GetZ()-medianZC)>kMaxZ) continue;
885         AliKFParticle part(*aliTrack,211);
886         cumulVertexC+=part;
887         cumulVertexAC+=part;
888       }      
889       //
890       if (TMath::Abs(cumulVertexA.X()-vertexA.X())>kMaxDvertex) flags+=64;
891       if (TMath::Abs(cumulVertexA.Y()-vertexA.Y())>kMaxDvertex) flags+=64;
892       if (TMath::Abs(cumulVertexA.Z()-vertexA.Z())>kMaxDvertex) flags+=64;
893       //
894       if (TMath::Abs(cumulVertexC.X()-vertexC.X())>kMaxDvertex) flags+=64;
895       if (TMath::Abs(cumulVertexC.Y()-vertexC.Y())>kMaxDvertex) flags+=64;
896       if (TMath::Abs(cumulVertexC.Z()-vertexC.Z())>kMaxDvertex) flags+=64;
897       
898       
899       if ( flags==0 && cumulVertexC.GetNDF()>kMinTracksVertex&&cumulVertexA.GetNDF()>kMinTracksVertex){
900         Double_t cont[2]={0,fTime};
901         //
902         cont[0]= cumulVertexA.X();
903         fTPCVertex[0]->Fill(cont);
904         cont[0]= cumulVertexC.X();
905         fTPCVertex[1]->Fill(cont);
906         cont[0]= 0.5*(cumulVertexA.X()-cumulVertexC.X());
907         fTPCVertex[2]->Fill(cont);
908         cont[0]= 0.5*(cumulVertexA.X()+cumulVertexC.X())-vertexSPD->GetX();
909         fTPCVertex[3]->Fill(cont);
910         //
911         cont[0]= cumulVertexA.Y();
912         fTPCVertex[4]->Fill(cont);
913         cont[0]= cumulVertexC.Y();
914         fTPCVertex[5]->Fill(cont);
915         cont[0]= 0.5*(cumulVertexA.Y()-cumulVertexC.Y());
916         fTPCVertex[6]->Fill(cont);
917         cont[0]= 0.5*(cumulVertexA.Y()+cumulVertexC.Y())-vertexSPD->GetY();
918         fTPCVertex[7]->Fill(cont);
919         //
920         //
921         cont[0]= 0.5*(cumulVertexA.Z()+cumulVertexC.Z());
922         fTPCVertex[8]->Fill(cont);
923         cont[0]= 0.5*(cumulVertexA.Z()-cumulVertexC.Z());
924         fTPCVertex[9]->Fill(cont);
925         cont[0]= 0.5*(cumulVertexA.Z()-cumulVertexC.Z());
926         fTPCVertex[10]->Fill(cont);      
927         cont[0]= 0.5*(cumulVertexA.Z()+cumulVertexC.Z())-vertexSPD->GetZ();
928         fTPCVertex[11]->Fill(cont);
929         //
930         Double_t correl[2]={0,0};
931         //
932         correl[0]=cumulVertexC.Z();
933         correl[1]=cumulVertexA.Z();
934         fTPCVertexCorrelation[0]->Fill(correl);   // fill A side :TPC
935         correl[0]=cumulVertexA.Z();
936         correl[1]=cumulVertexC.Z(); 
937         fTPCVertexCorrelation[1]->Fill(correl);   // fill C side :TPC
938         //
939         correl[0]=vertexSPD->GetZ();
940         correl[1]=cumulVertexA.Z()-correl[0];
941         fTPCVertexCorrelation[2]->Fill(correl);   // fill A side :ITS
942         correl[1]=cumulVertexC.Z()-correl[0]; 
943         fTPCVertexCorrelation[3]->Fill(correl);   // fill C side :ITS
944         correl[1]=0.5*(cumulVertexA.Z()+cumulVertexC.Z())-correl[0]; 
945         fTPCVertexCorrelation[4]->Fill(correl);   // fill C side :ITS
946       }
947     }        
948     TTreeSRedirector *cstream = GetDebugStreamer();
949     if (cstream){
950       /*
951         TCut cutChi2= "sqrt(vA.fChi2/vA.fNDF+vC.fChi2/vC.fNDF)<10";  // chi2 Cut e.g 10         
952         TCut cutXY= "sqrt((vA.fP[0]-vC.fP[0])^2+(vA.fP[0]-vC.fP[1])^2)<5";   // vertex Cut      
953         TCut cutZ= "abs(vA.fP[2]-mZA)<3&&abs(vC.fP[2]-mZC)<5";           // vertex Cut  
954         tree->Draw("sqrt(vA.fChi2/vA.fNDF)","sqrt(vA.fChi2/vA.fNDF)<100","")
955         
956       */
957       //vertexA.Print();
958       //vertexC.Print();      
959       (*cstream)<<"vertexTPC"<<
960         "flags="<<flags<<        // rejection flags
961         "vSPD.="<<vertexSPD<<    // SPD vertex
962         "vT.="<<vertexTracks<<   // track vertex
963         "v.="<<vertex<<          // esd vertex
964         "mZA="<<medianZA<<       // median Z position at vertex A side
965         "mZC="<<medianZC<<       // median Z position at vertex C side
966         "mDelta="<<medianDelta<< // median delta A side -C side
967         "counter="<<counterZ<<    // counter Z
968         //
969         "vA.="<<&vertexA<<       // vertex A side
970         "vC.="<<&vertexC<<       // vertex C side
971         "cvA.="<<&cumulVertexA<<       // cumulative vertex A side
972         "cvC.="<<&cumulVertexC<<       // cumulative vertex C side
973         "cvAC.="<<&cumulVertexAC<<       // cumulative vertex A+C side
974         "nA="<<ntracksA<<        // contributors
975         "nC="<<ntracksC<<        // contributors
976         "\n";
977     }      
978   }
979   tracksA.Delete();
980   tracksC.Delete();
981 }
982
983 void AliTPCcalibTime::Analyze(){
984   //
985   // Special macro to analyze result of calibration and extract calibration entries
986   // Not yet ported to the Analyze function yet
987   //
988 }
989
990 THnSparse* AliTPCcalibTime::GetHistoDrift(const char* name) const
991 {
992   //
993   // Get histogram for given trigger mask
994   //
995   TIterator* iterator = fArrayDz->MakeIterator();
996   iterator->Reset();
997   TString newName=name;
998   newName.ToUpper();
999   THnSparse* newHist=new THnSparseF(newName,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
1000   THnSparse* addHist=NULL;
1001   while((addHist=(THnSparseF*)iterator->Next())){
1002     //  if(!addHist) continue;
1003     TString histName=addHist->GetName();
1004     if(!histName.Contains(newName)) continue;
1005     addHist->Print();
1006     newHist->Add(addHist);
1007   }
1008   return newHist;
1009 }
1010
1011 TObjArray* AliTPCcalibTime::GetHistoDrift() const
1012 {
1013   //
1014   // return array of histograms
1015   //
1016   return fArrayDz;
1017 }
1018
1019 TGraphErrors* AliTPCcalibTime::GetGraphDrift(const char* name){
1020   //
1021   // Make a drift velocity (delta Z) graph
1022   //
1023   THnSparse* histoDrift=GetHistoDrift(name);
1024   TGraphErrors* graphDrift=NULL;
1025   if(histoDrift){
1026     graphDrift=FitSlices(histoDrift,2,0,400,100,0.05,0.95, kTRUE);
1027     TString end=histoDrift->GetName();
1028     Int_t pos=end.Index("_");
1029     end=end(pos,end.Capacity()-pos);
1030     TString graphName=graphDrift->ClassName();
1031     graphName+=end;
1032     graphName.ToUpper();
1033     graphDrift->SetName(graphName);
1034   }
1035   return graphDrift;
1036 }
1037
1038 TObjArray* AliTPCcalibTime::GetGraphDrift(){
1039   //
1040   // make a array of drift graphs
1041   //
1042   TObjArray* arrayGraphDrift=new TObjArray();
1043   TIterator* iterator=fArrayDz->MakeIterator();
1044   iterator->Reset();
1045   THnSparse* addHist=NULL;
1046   while((addHist=(THnSparseF*)iterator->Next())) arrayGraphDrift->AddLast(GetGraphDrift(addHist->GetName()));
1047   return arrayGraphDrift;
1048 }
1049
1050 AliSplineFit* AliTPCcalibTime::GetFitDrift(const char* name){
1051   //
1052   // Make a fit AliSplinefit  of drift velocity
1053   //
1054   TGraph* graphDrift=GetGraphDrift(name);
1055   AliSplineFit* fitDrift=NULL;
1056   if(graphDrift && graphDrift->GetN()){
1057     fitDrift=new AliSplineFit();
1058     fitDrift->SetGraph(graphDrift);
1059     fitDrift->SetMinPoints(graphDrift->GetN()+1);
1060     fitDrift->InitKnots(graphDrift,2,0,0.001);
1061     fitDrift->SplineFit(0);
1062     TString end=graphDrift->GetName();
1063     Int_t pos=end.Index("_");
1064     end=end(pos,end.Capacity()-pos);
1065     TString fitName=fitDrift->ClassName();
1066     fitName+=end;
1067     fitName.ToUpper();
1068     //fitDrift->SetName(fitName);
1069     delete graphDrift;
1070     graphDrift=NULL;
1071   }
1072   return fitDrift;
1073 }
1074
1075
1076 Long64_t AliTPCcalibTime::Merge(TCollection *const li) {
1077   //
1078   // Object specific merging procedure
1079   //
1080   TIterator* iter = li->MakeIterator();
1081   AliTPCcalibTime* cal = 0;
1082   //
1083   while ((cal = (AliTPCcalibTime*)iter->Next())) {
1084     if (!cal->InheritsFrom(AliTPCcalibTime::Class())) {
1085       Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName());
1086       return -1;
1087     }
1088     for (Int_t imeas=0; imeas<3; imeas++){
1089       if (cal->GetHistVdriftLaserA(imeas) && cal->GetHistVdriftLaserA(imeas)){
1090         fHistVdriftLaserA[imeas]->Add(cal->GetHistVdriftLaserA(imeas));
1091         fHistVdriftLaserC[imeas]->Add(cal->GetHistVdriftLaserC(imeas));
1092       }
1093     }
1094     //
1095     if (fTPCVertexCorrelation[0] && cal->fTPCVertexCorrelation[0]){
1096       for (Int_t imeas=0; imeas<5; imeas++){
1097         if (fTPCVertexCorrelation[imeas] && cal->fTPCVertexCorrelation[imeas]) fTPCVertexCorrelation[imeas]->Add(cal->fTPCVertexCorrelation[imeas]);
1098       }
1099     }
1100       
1101     if (fTPCVertex[0] && cal->fTPCVertex[0]) 
1102       for (Int_t imeas=0; imeas<12; imeas++){
1103         if (fTPCVertex[imeas] && cal->fTPCVertex[imeas]) fTPCVertex[imeas]->Add(cal->fTPCVertex[imeas]);
1104       }
1105     
1106     if (fMemoryMode>0) for (Int_t imeas=0; imeas<5; imeas++){
1107       if (fMemoryMode>1){
1108         if ( cal->GetResHistoTPCCE(imeas) && cal->GetResHistoTPCCE(imeas)){
1109           fResHistoTPCCE[imeas]->Add(cal->fResHistoTPCCE[imeas]);
1110         }else{
1111           fResHistoTPCCE[imeas]=(THnSparse*)cal->fResHistoTPCCE[imeas]->Clone();
1112         }
1113       }
1114       //
1115       if ((fMemoryMode>0) &&cal->GetResHistoTPCITS(imeas) && cal->GetResHistoTPCITS(imeas)){
1116         if (fMemoryMode>1 || (imeas%2)==1) fResHistoTPCITS[imeas]->Add(cal->fResHistoTPCITS[imeas]);
1117         if (fMemoryMode>1) fResHistoTPCvertex[imeas]->Add(cal->fResHistoTPCvertex[imeas]);
1118       }
1119       //
1120       if ((fMemoryMode>1) && cal->fResHistoTPCTRD[imeas]){
1121         if (fResHistoTPCTRD[imeas])
1122           fResHistoTPCTRD[imeas]->Add(cal->fResHistoTPCTRD[imeas]);
1123         else
1124           fResHistoTPCTRD[imeas]=(THnSparse*)cal->fResHistoTPCTRD[imeas]->Clone();
1125       }
1126       //
1127       if  ((fMemoryMode>1) && cal->fResHistoTPCTOF[imeas]){
1128         if (fResHistoTPCTOF[imeas])
1129           fResHistoTPCTOF[imeas]->Add(cal->fResHistoTPCTOF[imeas]);
1130         else
1131           fResHistoTPCTOF[imeas]=(THnSparse*)cal->fResHistoTPCTOF[imeas]->Clone();      
1132       }
1133       //
1134       if (cal->fArrayLaserA){
1135         fArrayLaserA->Expand(fArrayLaserA->GetEntriesFast()+cal->fArrayLaserA->GetEntriesFast());
1136         fArrayLaserC->Expand(fArrayLaserC->GetEntriesFast()+cal->fArrayLaserC->GetEntriesFast());
1137         for (Int_t ical=0; ical<cal->fArrayLaserA->GetEntriesFast(); ical++){
1138           if (cal->fArrayLaserA->UncheckedAt(ical)) fArrayLaserA->AddLast(cal->fArrayLaserA->UncheckedAt(ical)->Clone());
1139           if (cal->fArrayLaserC->UncheckedAt(ical)) fArrayLaserC->AddLast(cal->fArrayLaserC->UncheckedAt(ical)->Clone());
1140         }
1141       }
1142
1143     }
1144     TObjArray* addArray=cal->GetHistoDrift();
1145     if(!addArray) return 0;
1146     TIterator* iterator = addArray->MakeIterator();
1147     iterator->Reset();
1148     THnSparse* addHist=NULL;
1149     if ((fMemoryMode>1)) while((addHist=(THnSparseF*)iterator->Next())){
1150       //      if(!addHist) continue;
1151       addHist->Print();
1152       THnSparse* localHist=(THnSparseF*)fArrayDz->FindObject(addHist->GetName());
1153       if(!localHist){
1154         localHist=new THnSparseF(addHist->GetName(),"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
1155         fArrayDz->AddLast(localHist);
1156       }
1157       localHist->Add(addHist);
1158     }
1159
1160     for(Int_t i=0;i<10;i++) if (cal->GetCosmiMatchingHisto(i)) fCosmiMatchingHisto[i]->Add(cal->GetCosmiMatchingHisto(i));
1161     //
1162     // Merge alignment
1163     //
1164     for (Int_t itype=0; itype<3; itype++){
1165       //
1166       //
1167       TObjArray *arr0= 0;
1168       TObjArray *arr1= 0;
1169       if (itype==0) {arr0=fAlignITSTPC; arr1=cal->fAlignITSTPC;}
1170       if (itype==1) {arr0=fAlignTRDTPC; arr1=cal->fAlignTRDTPC;}
1171       if (itype==2) {arr0=fAlignTOFTPC; arr1=cal->fAlignTOFTPC;}
1172       if (!arr1) continue;
1173       if (!arr0) arr0=new TObjArray(arr1->GetEntriesFast());
1174       if (arr1->GetEntriesFast()>arr0->GetEntriesFast()){
1175         arr0->Expand(arr1->GetEntriesFast());
1176       }
1177       for (Int_t i=0;i<arr1->GetEntriesFast(); i++){
1178         AliRelAlignerKalman *kalman1 = (AliRelAlignerKalman *)arr1->UncheckedAt(i);
1179         AliRelAlignerKalman *kalman0 = (AliRelAlignerKalman *)arr0->UncheckedAt(i);
1180         if (!kalman1)  continue;
1181         if (!kalman0) {arr0->AddAt(new AliRelAlignerKalman(*kalman1),i); continue;}
1182         kalman0->SetRejectOutliers(kFALSE);
1183         kalman0->Merge(kalman1);
1184       }
1185     }
1186
1187   }
1188   return 0;
1189 }
1190
1191 Bool_t  AliTPCcalibTime::IsPair(const AliExternalTrackParam *tr0, const AliExternalTrackParam *tr1){
1192   /*
1193   // 0. Same direction - OPOSITE  - cutDir +cutT    
1194   TCut cutDir("cutDir","dir<-0.99")
1195   // 1. 
1196   TCut cutT("cutT","abs(Tr1.fP[3]+Tr0.fP[3])<0.03")
1197   //
1198   // 2. The same rphi 
1199   TCut cutD("cutD","abs(Tr0.fP[0]+Tr1.fP[0])<5")
1200   //
1201   TCut cutPt("cutPt","abs(Tr1.fP[4]+Tr0.fP[4])<1&&abs(Tr0.fP[4])+abs(Tr1.fP[4])<10");  
1202   // 1/Pt diff cut
1203   */
1204   const Double_t *p0 = tr0->GetParameter();
1205   const Double_t *p1 = tr1->GetParameter();
1206   fCosmiMatchingHisto[0]->Fill(p0[0]+p1[0]);
1207   fCosmiMatchingHisto[1]->Fill(p0[1]-p1[1]);
1208   fCosmiMatchingHisto[2]->Fill(tr1->GetAlpha()-tr0->GetAlpha()+TMath::Pi());
1209   fCosmiMatchingHisto[3]->Fill(p0[3]+p1[3]);
1210   fCosmiMatchingHisto[4]->Fill(p0[4]+p1[4]);
1211   
1212   if (TMath::Abs(p0[3]+p1[3])>fCutTheta) return kFALSE;
1213   if (TMath::Abs(p0[0]+p1[0])>fCutMaxD)  return kFALSE;
1214   if (TMath::Abs(p0[1]-p1[1])>fCutMaxDz)  return kFALSE;
1215   Double_t d0[3], d1[3];
1216   tr0->GetDirection(d0);    
1217   tr1->GetDirection(d1);       
1218   if (d0[0]*d1[0] + d0[1]*d1[1] + d0[2]*d1[2] >fCutMinDir) return kFALSE;
1219
1220   fCosmiMatchingHisto[5]->Fill(p0[0]+p1[0]);
1221   fCosmiMatchingHisto[6]->Fill(p0[1]-p1[1]);
1222   fCosmiMatchingHisto[7]->Fill(tr1->GetAlpha()-tr0->GetAlpha()+TMath::Pi());
1223   fCosmiMatchingHisto[8]->Fill(p0[3]+p1[3]);
1224   fCosmiMatchingHisto[9]->Fill(p0[4]+p1[4]);
1225
1226   return kTRUE;  
1227 }
1228 Bool_t AliTPCcalibTime::IsCross(const AliESDtrack *const tr0, const AliESDtrack *const tr1){
1229   //
1230   // check if the cosmic pair of tracks crossed A/C side
1231   // 
1232   Bool_t result= tr0->GetOuterParam()->GetZ()*tr1->GetOuterParam()->GetZ()<0;
1233   if (result==kFALSE) return result;
1234   result=kTRUE;
1235   return result;
1236 }
1237
1238 Bool_t AliTPCcalibTime::IsSame(const AliESDtrack *const tr0, const AliESDtrack *const tr1){
1239   // 
1240   // track crossing the CE
1241   // 0. minimal number of clusters 
1242   // 1. Same sector +-1
1243   // 2. Inner and outer track param on opposite side
1244   // 3. Outer and inner track parameter close each to other
1245   // 3. 
1246   Bool_t result=kTRUE;
1247   //
1248   // inner and outer on opposite sides in z
1249   //
1250   const Int_t knclCut0  = 30;
1251   const Double_t kalphaCut = 0.4;
1252   //
1253   // 0. minimal number of clusters
1254   //
1255   if (tr0->GetTPCNcls()<knclCut0) return kFALSE;
1256   if (tr1->GetTPCNcls()<knclCut0) return kFALSE;
1257   //
1258   // 1. alpha cut - sector+-1
1259   //
1260   if (TMath::Abs(tr0->GetOuterParam()->GetAlpha()-tr1->GetOuterParam()->GetAlpha())>kalphaCut) return kFALSE;
1261   //
1262   // 2. Z crossing
1263   //
1264   if (tr0->GetOuterParam()->GetZ()*tr0->GetInnerParam()->GetZ()>0) result&=kFALSE;
1265   if (tr1->GetOuterParam()->GetZ()*tr1->GetInnerParam()->GetZ()>0) result&=kFALSE;
1266   if (result==kFALSE){
1267     return result;
1268   }
1269   //
1270   //
1271   const Double_t *p0I = tr0->GetInnerParam()->GetParameter();
1272   const Double_t *p1I = tr1->GetInnerParam()->GetParameter();
1273   const Double_t *p0O = tr0->GetOuterParam()->GetParameter();
1274   const Double_t *p1O = tr1->GetOuterParam()->GetParameter();
1275   //
1276   if (TMath::Abs(p0I[0]-p1I[0])>fCutMaxD)  result&=kFALSE;
1277   if (TMath::Abs(p0I[1]-p1I[1])>fCutMaxDz) result&=kFALSE;
1278   if (TMath::Abs(p0I[2]-p1I[2])>fCutTheta) result&=kFALSE;
1279   if (TMath::Abs(p0I[3]-p1I[3])>fCutTheta) result&=kFALSE;
1280   if (TMath::Abs(p0O[0]-p1O[0])>fCutMaxD)  result&=kFALSE;
1281   if (TMath::Abs(p0O[1]-p1O[1])>fCutMaxDz) result&=kFALSE;
1282   if (TMath::Abs(p0O[2]-p1O[2])>fCutTheta) result&=kFALSE;
1283   if (TMath::Abs(p0O[3]-p1O[3])>fCutTheta) result&=kFALSE;
1284   if (result==kTRUE){
1285     result=kTRUE; // just to put break point here
1286   }
1287   return result;
1288 }
1289
1290
1291 void  AliTPCcalibTime::ProcessSame(const AliESDtrack *const track, AliESDfriendTrack *const friendTrack, const AliESDEvent *const event){
1292   //
1293   // Process  TPC tracks crossing CE
1294   //
1295   // 0. Select only track crossing the CE
1296   // 1. Cut on the track length
1297   // 2. Refit the the track on A and C side separatelly
1298   // 3. Fill time histograms
1299   const Int_t kMinNcl=100;
1300   const Int_t kMinNclS=25;  // minimul number of clusters on the sides
1301   const Double_t pimass=TDatabasePDG::Instance()->GetParticle("pi+")->Mass();
1302   const Double_t kMaxDy=1;  // maximal distance in y
1303   const Double_t kMaxDsnp=0.05;  // maximal distance in snp
1304   const Double_t kMaxDtheta=0.05;  // maximal distance in theta
1305   
1306   if (!friendTrack->GetTPCOut()) return;
1307   //
1308   // 0. Select only track crossing the CE
1309   //
1310   if (track->GetInnerParam()->GetZ()*friendTrack->GetTPCOut()->GetZ()>0) return;
1311   //
1312   // 1. cut on track length
1313   //
1314   if (track->GetTPCNcls()<kMinNcl) return;
1315   //
1316   // 2. Refit track sepparatel on A and C side
1317   //
1318   TObject *calibObject;
1319   AliTPCseed *seed = 0;
1320   for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
1321     if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
1322   }
1323   if (!seed) return;
1324   //
1325   AliExternalTrackParam trackIn(*track->GetInnerParam());
1326   AliExternalTrackParam trackOut(*track->GetOuterParam());
1327   Double_t cov[3]={0.01,0.,0.01}; //use the same errors
1328   Double_t xyz[3]={0,0.,0.0};  
1329   Double_t bz   =0;
1330   Int_t nclIn=0,nclOut=0;
1331   trackIn.ResetCovariance(1000.);
1332   trackOut.ResetCovariance(1000.);
1333   //
1334   //2.a Refit inner
1335   // 
1336   Int_t sideIn=0;
1337   for (Int_t irow=0;irow<159;irow++) {
1338     AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
1339     if (!cl) continue;
1340     if (cl->GetX()<80) continue;
1341     if (sideIn==0){
1342       if (cl->GetDetector()%36<18) sideIn=1;
1343       if (cl->GetDetector()%36>=18) sideIn=-1;
1344     }
1345     if (sideIn== -1 && (cl->GetDetector()%36)<18) break;
1346     if (sideIn==  1 &&(cl->GetDetector()%36)>=18) break;
1347     Int_t sector = cl->GetDetector();
1348     Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackIn.GetAlpha();
1349     if (TMath::Abs(dalpha)>0.01){
1350       if (!trackIn.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
1351     }
1352     Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
1353     trackIn.GetXYZ(xyz);
1354     bz = AliTracker::GetBz(xyz);
1355     AliTracker::PropagateTrackToBxByBz(&trackIn,r[0],pimass,1.,kFALSE);
1356     if (!trackIn.PropagateTo(r[0],bz)) break;
1357     nclIn++;
1358     trackIn.Update(&r[1],cov);    
1359   }
1360   //
1361   //2.b Refit outer
1362   //
1363   Int_t sideOut=0;
1364   for (Int_t irow=159;irow>0;irow--) {
1365     AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
1366     if (!cl) continue;
1367     if (cl->GetX()<80) continue;
1368     if (sideOut==0){
1369       if (cl->GetDetector()%36<18) sideOut=1;
1370       if (cl->GetDetector()%36>=18) sideOut=-1;
1371       if (sideIn==sideOut) break;
1372     }
1373     if (sideOut== -1 && (cl->GetDetector()%36)<18) break;
1374     if (sideOut==  1 &&(cl->GetDetector()%36)>=18) break;
1375     //
1376     Int_t sector = cl->GetDetector();
1377     Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackOut.GetAlpha();
1378     if (TMath::Abs(dalpha)>0.01){
1379       if (!trackOut.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
1380     }
1381     Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
1382     trackOut.GetXYZ(xyz);
1383     bz = AliTracker::GetBz(xyz);
1384     AliTracker::PropagateTrackToBxByBz(&trackOut,r[0],pimass,1.,kFALSE);
1385     if (!trackOut.PropagateTo(r[0],bz)) break;
1386     nclOut++;
1387     trackOut.Update(&r[1],cov);    
1388   }
1389   trackOut.Rotate(trackIn.GetAlpha());
1390   Double_t meanX = (trackIn.GetX()+trackOut.GetX())*0.5;
1391   trackIn.PropagateTo(meanX,bz); 
1392   trackOut.PropagateTo(meanX,bz); 
1393   if (TMath::Abs(trackIn.GetY()-trackOut.GetY())>kMaxDy) return;
1394   if (TMath::Abs(trackIn.GetSnp()-trackOut.GetSnp())>kMaxDsnp) return;
1395   if (TMath::Abs(trackIn.GetTgl()-trackOut.GetTgl())>kMaxDtheta) return;
1396   if (TMath::Min(nclIn,nclOut)>kMinNclS){
1397     FillResHistoTPCCE(&trackIn,&trackOut);
1398   }
1399   TTreeSRedirector *cstream = GetDebugStreamer();
1400   if (cstream){
1401     TVectorD gxyz(3);
1402     trackIn.GetXYZ(gxyz.GetMatrixArray());
1403     TTimeStamp tstamp(fTime);
1404     (*cstream)<<"tpctpc"<<
1405       "run="<<fRun<<              //  run number
1406       "event="<<fEvent<<          //  event number
1407       "time="<<fTime<<            //  time stamp of event
1408       "trigger="<<fTrigger<<      //  trigger
1409       "mag="<<fMagF<<             //  magnetic field
1410       //
1411       "sideIn="<<sideIn<<         // side at inner part
1412       "sideOut="<<sideOut<<         // side at puter part
1413       "xyz.="<<&gxyz<<             // global position
1414       "tIn.="<<&trackIn<<         // refitterd track in 
1415       "tOut.="<<&trackOut<<       // refitter track out
1416       "nclIn="<<nclIn<<           // 
1417       "nclOut="<<nclOut<<         //
1418       "\n";  
1419   }
1420   //
1421   // 3. Fill time histograms
1422   // Debug stremaer expression
1423   // chainTPCTPC->Draw("(tIn.fP[1]-tOut.fP[1])*sign(-tIn.fP[3]):tIn.fP[3]","min(nclIn,nclOut)>30","")
1424   if (TMath::Min(nclIn,nclOut)>kMinNclS){
1425     fDz = trackOut.GetZ()-trackIn.GetZ();
1426     if (trackOut.GetTgl()<0) fDz*=-1.;
1427     TTimeStamp tstamp(fTime);
1428     Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
1429     Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
1430     Double_t vecDrift[4]={fTime,(ptrelative0+ptrelative1)/2.0,fDz/500.0,event->GetRunNumber()};
1431     //
1432     // fill histograms per trigger class and itegrated
1433     //
1434     THnSparse* curHist=NULL;
1435     for (Int_t itype=0; itype<2; itype++){
1436       TString name="MEAN_VDRIFT_CROSS_";  
1437       if (itype==0){
1438         name+=event->GetFiredTriggerClasses();
1439         name.ToUpper();
1440       }else{
1441         name+="ALL";
1442       }
1443       curHist=(THnSparseF*)fArrayDz->FindObject(name);
1444       if(!curHist){
1445         curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
1446         fArrayDz->AddLast(curHist);
1447       }
1448       curHist->Fill(vecDrift);
1449     }
1450   }
1451
1452 }
1453
1454 void  AliTPCcalibTime::ProcessAlignITS(AliESDtrack *const track, const AliESDfriendTrack *const friendTrack, const AliESDEvent *const event, AliESDfriend *const esdFriend){
1455   //
1456   // Process track - Update TPC-ITS alignment
1457   // Updates: 
1458   // 0. Apply standartd cuts 
1459   // 1. Recalucluate the current statistic median/RMS
1460   // 2. Apply median+-rms cut
1461   // 3. Update kalman filter
1462   //
1463   const Int_t    kMinTPC  = 80;    // minimal number of TPC cluster
1464   const Int_t    kMinITS  = 3;     // minimal number of ITS cluster
1465   const Double_t kMinZ    = 10;    // maximal dz distance
1466   const Double_t kMaxDy   = 2.;    // maximal dy distance
1467   const Double_t kMaxAngle= 0.07;  // maximal angular distance
1468   const Double_t kSigmaCut= 5;     // maximal sigma distance to median
1469   const Double_t kVdErr   = 0.1;  // initial uncertainty of the vd correction 
1470   const Double_t kT0Err   = 3.;  // initial uncertainty of the T0 time
1471   const Double_t kVdYErr  = 0.05;  // initial uncertainty of the vd correction 
1472   const Double_t kOutCut  = 3.0;   // outlyer cut in AliRelAlgnmentKalman
1473   const Double_t kMinPt   = 0.3;   // minimal pt
1474   const Double_t kMax1Pt=0.5;        //maximal 1/pt distance
1475   const  Int_t     kN=50;         // deepnes of history
1476   static Int_t     kglast=0;
1477   static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1478   //
1479   // 0. Apply standard cuts
1480   // 
1481   Int_t dummycl[1000];
1482   if (track->GetTPCNcls()<kMinTPC) return;  // minimal amount of clusters cut
1483   if (!track->IsOn(AliESDtrack::kTPCrefit)) return;
1484   if (!track->GetInnerParam())   return;
1485   if (!track->GetOuterParam())   return;
1486   if (track->GetInnerParam()->Pt()<kMinPt)  return;
1487   // exclude crossing track
1488   if (track->GetOuterParam()->GetZ()*track->GetInnerParam()->GetZ()<0)   return;
1489   if (TMath::Abs(track->GetInnerParam()->GetZ())<kMinZ/3.)   return;
1490   if (track->GetInnerParam()->GetX()>90)   return;
1491   //
1492   AliExternalTrackParam &pTPC=(AliExternalTrackParam &)(*(track->GetInnerParam()));
1493   //  
1494   AliExternalTrackParam pITS;   // ITS standalone if possible
1495   AliExternalTrackParam pITS2;  //TPC-ITS track
1496   if (friendTrack->GetITSOut()){
1497     pITS2=(*(friendTrack->GetITSOut()));  //TPC-ITS track - snapshot ITS out
1498     pITS2.Rotate(pTPC.GetAlpha());
1499     AliTracker::PropagateTrackToBxByBz(&pITS2,pTPC.GetX(),0.1,0.1,kFALSE);
1500   }
1501
1502   AliESDfriendTrack *itsfriendTrack=0;
1503   //
1504   // try to find standalone ITS track corresponing to the TPC if possible
1505   //
1506   Bool_t hasAlone=kFALSE;
1507   Int_t ntracks=event->GetNumberOfTracks();
1508   for (Int_t i=0; i<ntracks; i++){
1509     AliESDtrack * trackITS = event->GetTrack(i); 
1510     if (!trackITS) continue;
1511     if (trackITS->GetITSclusters(dummycl)<kMinITS) continue;  // minimal amount of clusters
1512     itsfriendTrack = esdFriend->GetTrack(i);
1513     if (!itsfriendTrack) continue;
1514     if (!itsfriendTrack->GetITSOut()) continue;
1515      
1516     if (TMath::Abs(pTPC.GetTgl()-itsfriendTrack->GetITSOut()->GetTgl())> kMaxAngle) continue;
1517     if (TMath::Abs(pTPC.GetSigned1Pt()-itsfriendTrack->GetITSOut()->GetSigned1Pt())> kMax1Pt) continue;
1518     pITS=(*(itsfriendTrack->GetITSOut()));
1519     //
1520     pITS.Rotate(pTPC.GetAlpha());
1521     AliTracker::PropagateTrackToBxByBz(&pITS,pTPC.GetX(),0.1,0.1,kFALSE);
1522     if (TMath::Abs(pTPC.GetY()-pITS.GetY())> kMaxDy) continue;
1523     if (TMath::Abs(pTPC.GetSnp()-pITS.GetSnp())> kMaxAngle) continue;
1524     hasAlone=kTRUE;
1525   }
1526   if (!hasAlone) {
1527     if (track->GetITSclusters(dummycl)<kMinITS) return;
1528     pITS=pITS2;  // use combined track if it has ITS
1529   }
1530   //
1531   if (TMath::Abs(pITS.GetY()-pTPC.GetY())    >kMaxDy)    return;
1532   if (TMath::Abs(pITS.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1533   if (TMath::Abs(pITS.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1534   //
1535   // 1. Update median and RMS info
1536   //
1537   TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
1538   TVectorD vecDeltaN(5);
1539   Double_t sign=(pITS.GetParameter()[1]>0)? 1.:-1.;
1540   vecDelta[4]=0;
1541   for (Int_t i=0;i<4;i++){
1542     vecDelta[i]=(pITS.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1543     kgdP[i][kglast%kN]=vecDelta[i];
1544   }
1545   kglast=(kglast+1);
1546   Int_t entries=(kglast<kN)?kglast:kN;
1547   for (Int_t i=0;i<4;i++){
1548     vecMedian[i] = TMath::Median(entries,kgdP[i]);
1549     vecRMS[i]    = TMath::RMS(entries,kgdP[i]);
1550     vecDeltaN[i] = 0;
1551     if (vecRMS[i]>0.){
1552       vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/vecRMS[i];
1553       vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]);  //sum of abs residuals
1554     }
1555   }
1556   //
1557   // 2. Apply median+-rms cut
1558   //
1559   if (kglast<3)  return;   //median and RMS to be defined
1560   if ( vecDeltaN[4]/4.>kSigmaCut) return;
1561   //
1562   // 3. Update alignment
1563   //
1564   Int_t htime = (fTime-fTimeKalmanBin/2)/fTimeKalmanBin; //time bins number
1565   if (fAlignITSTPC->GetEntriesFast()<htime){
1566     fAlignITSTPC->Expand(htime*2+20);
1567   }
1568   AliRelAlignerKalman* align =  (AliRelAlignerKalman*)fAlignITSTPC->At(htime);
1569   if (!align){
1570     // make Alignment object if doesn't exist
1571     align=new AliRelAlignerKalman(); 
1572     align->SetRunNumber(fRun);
1573     (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1574     (*align->GetStateCov())(7,7)=kT0Err*kT0Err;
1575     (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1576     align->SetOutRejSigma(kOutCut+kOutCut*kN);
1577     //    align->SetRejectOutliers(kFALSE);
1578     align->SetRejectOutliers(kTRUE);
1579     align->SetRejectOutliersSigma2Median(kTRUE);
1580
1581     align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1582     align->SetMagField(fMagF); 
1583     fAlignITSTPC->AddAt(align,htime);
1584   }
1585   align->AddTrackParams(&pITS,&pTPC);
1586   Double_t averageTime =  fTime;
1587   if (align->GetTimeStamp()>0&&align->GetNUpdates()>0){
1588     averageTime=((Double_t(align->GetTimeStamp())*Double_t(align->GetNUpdates())+Double_t(fTime)))/(Double_t(align->GetNUpdates())+1.);
1589   }
1590   align->SetTimeStamp(Int_t(averageTime));
1591
1592   align->SetRunNumber(fRun );
1593   Float_t dca[2],cov[3];
1594   track->GetImpactParameters(dca,cov);
1595   if (TMath::Abs(dca[0])<kMaxDy){
1596     FillResHistoTPCITS(&pTPC,&pITS);
1597     FillResHistoTPC(track);
1598   }
1599   //
1600   Int_t nupdates=align->GetNUpdates();
1601   align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates+1));
1602   //  align->SetRejectOutliers(kFALSE);
1603   align->SetRejectOutliers(kTRUE);
1604   align->SetRejectOutliersSigma2Median(kTRUE);
1605
1606   TTreeSRedirector *cstream = GetDebugStreamer();  
1607   if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1608     TVectorD gpTPC(3), gdTPC(3);
1609     TVectorD gpITS(3), gdITS(3);
1610     pTPC.GetXYZ(gpTPC.GetMatrixArray());
1611     pTPC.GetDirection(gdTPC.GetMatrixArray());
1612     pITS.GetXYZ(gpITS.GetMatrixArray());
1613     pITS.GetDirection(gdITS.GetMatrixArray());
1614     (*cstream)<<"itstpc"<<
1615       "run="<<fRun<<              //  run number
1616       "event="<<fEvent<<          //  event number
1617       "time="<<fTime<<            //  time stamp of event
1618       "trigger="<<fTrigger<<      //  trigger
1619       "mag="<<fMagF<<             //  magnetic field
1620       //
1621       "hasAlone="<<hasAlone<<    // has ITS standalone ?
1622       "track.="<<track<<  // track info
1623       "nmed="<<kglast<<        // number of entries to define median and RMS
1624       "vMed.="<<&vecMedian<<    // median of deltas
1625       "vRMS.="<<&vecRMS<<       // rms of deltas
1626       "vDelta.="<<&vecDelta<<   // delta in respect to median
1627       "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1628       "a.="<<align<<            // current alignment
1629       "pITS.="<<&pITS<<         // track param ITS 
1630       "pITS2.="<<&pITS2<<       // track param ITS+TPC
1631       "pTPC.="<<&pTPC<<         // track param TPC
1632       "gpTPC.="<<&gpTPC<<       // global position  TPC
1633       "gdTPC.="<<&gdTPC<<       // global direction TPC
1634       "gpITS.="<<&gpITS<<       // global position  ITS
1635       "gdITS.="<<&gdITS<<       // global position  ITS
1636       "\n";
1637   }
1638 }
1639
1640
1641
1642
1643 void  AliTPCcalibTime::ProcessAlignTRD(AliESDtrack *const track, const AliESDfriendTrack *const friendTrack){
1644   //
1645   // Process track - Update TPC-TRD alignment
1646   // Updates: 
1647   // 0. Apply standartd cuts 
1648   // 1. Recalucluate the current statistic median/RMS
1649   // 2. Apply median+-rms cut
1650   // 3. Update kalman filter
1651   //
1652   const Int_t    kMinTPC  = 80;    // minimal number of TPC cluster
1653   const Int_t    kMinTRD  = 50;    // minimal number of TRD cluster
1654   const Double_t kMinZ    = 20;    // maximal dz distance
1655   const Double_t kMaxDy   = 5.;    // maximal dy distance
1656   const Double_t kMaxAngle= 0.1;  // maximal angular distance
1657   const Double_t kSigmaCut= 10;     // maximal sigma distance to median
1658   const Double_t kVdErr   = 0.1;  // initial uncertainty of the vd correction 
1659   const Double_t kT0Err   = 3.;  // initial uncertainty of the T0 time
1660   const Double_t kVdYErr  = 0.05;  // initial uncertainty of the vd correction 
1661   const Double_t kOutCut  = 3.0;   // outlyer cut in AliRelAlgnmentKalman
1662   const Double_t kRefX    = 275;   // reference X
1663   const  Int_t     kN=50;         // deepnes of history
1664   static Int_t     kglast=0;
1665   static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1666   //
1667   // 0. Apply standard cuts
1668   //
1669   Int_t dummycl[1000];
1670   if (track->GetTRDclusters(dummycl)<kMinTRD) return;  // minimal amount of clusters
1671   if (track->GetTPCNcls()<kMinTPC) return;  // minimal amount of clusters cut
1672   if (!friendTrack->GetTRDIn()) return;  
1673   if (!track->IsOn(AliESDtrack::kTRDrefit)) return;  
1674   if (!track->IsOn(AliESDtrack::kTRDout)) return;  
1675   if (!track->GetInnerParam())   return;
1676   if (!friendTrack->GetTPCOut())   return;
1677   // exclude crossing track
1678   if (friendTrack->GetTPCOut()->GetZ()*track->GetInnerParam()->GetZ()<0)   return;
1679   if (TMath::Abs(track->GetInnerParam()->GetZ())<kMinZ)   return;
1680   //
1681   AliExternalTrackParam &pTPC=(AliExternalTrackParam &)(*(friendTrack->GetTPCOut()));
1682   AliTracker::PropagateTrackToBxByBz(&pTPC,kRefX,0.1,0.1,kFALSE);
1683   AliExternalTrackParam pTRD(*(friendTrack->GetTRDIn()));
1684   pTRD.Rotate(pTPC.GetAlpha());
1685   //  pTRD.PropagateTo(pTPC.GetX(),fMagF);
1686   AliTracker::PropagateTrackToBxByBz(&pTRD,pTPC.GetX(),0.1,0.1,kFALSE);
1687
1688   ((Double_t*)pTRD.GetCovariance())[2]+=3.*3.;   // increas sys errors
1689   ((Double_t*)pTRD.GetCovariance())[9]+=0.1*0.1; // increse sys errors
1690
1691   if (TMath::Abs(pTRD.GetY()-pTPC.GetY())    >kMaxDy)    return;
1692   if (TMath::Abs(pTRD.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1693   //  if (TMath::Abs(pTRD.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1694   //
1695   // 1. Update median and RMS info
1696   //
1697   TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
1698   TVectorD vecDeltaN(5);
1699   Double_t sign=(pTRD.GetParameter()[1]>0)? 1.:-1.;
1700   vecDelta[4]=0;
1701   for (Int_t i=0;i<4;i++){
1702     vecDelta[i]=(pTRD.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1703     kgdP[i][kglast%kN]=vecDelta[i];
1704   }
1705   kglast=(kglast+1);
1706   Int_t entries=(kglast<kN)?kglast:kN;
1707   for (Int_t i=0;i<4;i++){
1708     vecMedian[i] = TMath::Median(entries,kgdP[i]);
1709
1710     vecRMS[i]    = TMath::RMS(entries,kgdP[i]);
1711     vecDeltaN[i] = 0;
1712     if (vecRMS[i]>0.){
1713       vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/vecRMS[i];
1714       vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]);  //sum of abs residuals
1715     }
1716   }
1717   //
1718   // 2. Apply median+-rms cut
1719   //
1720   if (kglast<3)  return;   //median and RMS to be defined
1721   if ( vecDeltaN[4]/4.>kSigmaCut) return;
1722   //
1723   // 3. Update alignment
1724   //
1725   //Int_t htime = fTime/3600; //time in hours
1726   Int_t htime = (Int_t)(fTime-fTimeKalmanBin/2)/fTimeKalmanBin; //time in half hour
1727   if (fAlignTRDTPC->GetEntriesFast()<htime){
1728     fAlignTRDTPC->Expand(htime*2+20);
1729   }
1730   AliRelAlignerKalman* align =  (AliRelAlignerKalman*)fAlignTRDTPC->At(htime);
1731   if (!align){
1732     // make Alignment object if doesn't exist
1733     align=new AliRelAlignerKalman(); 
1734     align->SetRunNumber(fRun);
1735     (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1736     (*align->GetStateCov())(7,7)=kT0Err*kT0Err;
1737     (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1738     align->SetOutRejSigma(kOutCut+kOutCut*kN);
1739     //    align->SetRejectOutliers(kFALSE);
1740     align->SetRejectOutliers(kTRUE);
1741     align->SetRejectOutliersSigma2Median(kTRUE);
1742
1743     align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1744     align->SetMagField(fMagF); 
1745     fAlignTRDTPC->AddAt(align,htime);
1746   }
1747   align->AddTrackParams(&pTRD,&pTPC);
1748   //align->SetTimeStamp(fTime);
1749   Double_t averageTime =  fTime;
1750   if (align->GetTimeStamp()>0 && align->GetNUpdates()>0) {
1751     averageTime = (((Double_t)fTime) + ((Double_t)align->GetTimeStamp())*align->GetNUpdates()) / (align->GetNUpdates() + 1.);
1752     //printf("align->GetTimeStamp() %d, align->GetNUpdates() %d \n", align->GetTimeStamp(), align->GetNUpdates());
1753   }
1754   align->SetTimeStamp((Int_t)averageTime);
1755
1756   //printf("fTime %d, averageTime %d \n", fTime, (Int_t)averageTime);
1757
1758   align->SetRunNumber(fRun );
1759   Float_t dca[2],cov[3];
1760   track->GetImpactParameters(dca,cov);
1761   if (TMath::Abs(dca[0])<kMaxDy){
1762     FillResHistoTPCTRD(&pTPC,&pTRD);  //only primaries
1763   }
1764   //
1765   Int_t nupdates=align->GetNUpdates();
1766   align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates+1));
1767   //  align->SetRejectOutliers(kFALSE);
1768   align->SetRejectOutliers(kTRUE);
1769   align->SetRejectOutliersSigma2Median(kTRUE);
1770
1771   TTreeSRedirector *cstream = GetDebugStreamer();  
1772   if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1773     TVectorD gpTPC(3), gdTPC(3);
1774     TVectorD gpTRD(3), gdTRD(3);
1775     pTPC.GetXYZ(gpTPC.GetMatrixArray());
1776     pTPC.GetDirection(gdTPC.GetMatrixArray());
1777     pTRD.GetXYZ(gpTRD.GetMatrixArray());
1778     pTRD.GetDirection(gdTRD.GetMatrixArray());
1779     (*cstream)<<"trdtpc"<<
1780       "run="<<fRun<<              //  run number
1781       "event="<<fEvent<<          //  event number
1782       "time="<<fTime<<            //  time stamp of event
1783       "trigger="<<fTrigger<<      //  trigger
1784       "mag="<<fMagF<<             //  magnetic field
1785       //
1786       "nmed="<<kglast<<        // number of entries to define median and RMS
1787       "vMed.="<<&vecMedian<<    // median of deltas
1788       "vRMS.="<<&vecRMS<<       // rms of deltas
1789       "vDelta.="<<&vecDelta<<   // delta in respect to median
1790       "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1791       "t.="<<track<<            // ful track - find proper cuts
1792       "a.="<<align<<            // current alignment
1793       "pTRD.="<<&pTRD<<         // track param TRD
1794       "pTPC.="<<&pTPC<<         // track param TPC
1795       "gpTPC.="<<&gpTPC<<       // global position  TPC
1796       "gdTPC.="<<&gdTPC<<       // global direction TPC
1797       "gpTRD.="<<&gpTRD<<       // global position  TRD
1798       "gdTRD.="<<&gdTRD<<       // global position  TRD
1799       "\n";
1800   }
1801 }
1802
1803
1804 void  AliTPCcalibTime::ProcessAlignTOF(AliESDtrack *const track, const AliESDfriendTrack *const friendTrack){
1805   //
1806   //
1807   // Process track - Update TPC-TOF alignment
1808   // Updates: 
1809   // -1. Make a TOF "track"
1810   // 0. Apply standartd cuts 
1811   // 1. Recalucluate the current statistic median/RMS
1812   // 2. Apply median+-rms cut
1813   // 3. Update kalman filter
1814   //
1815   const Int_t      kMinTPC  = 80;    // minimal number of TPC cluster
1816   //  const Double_t   kMinZ    = 10;    // maximal dz distance
1817   const Double_t   kMaxDy   = 5.;    // maximal dy distance
1818   const Double_t   kMaxAngle= 0.05;  // maximal angular distance
1819   const Double_t   kSigmaCut= 5;     // maximal sigma distance to median
1820   const Double_t   kVdErr   = 0.1;  // initial uncertainty of the vd correction 
1821   const Double_t   kT0Err   = 3.;  // initial uncertainty of the T0 time
1822   const Double_t   kVdYErr  = 0.05;  // initial uncertainty of the vd correction 
1823
1824   const Double_t   kOutCut  = 3.0;   // outlyer cut in AliRelAlgnmentKalman
1825   const  Int_t     kN=50;         // deepnes of history
1826   static Int_t     kglast=0;
1827   static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1828   //
1829   // -1. Make a TOF track-
1830   //     Clusters are not in friends - use alingment points
1831   //
1832   if (track->GetTOFsignal()<=0)  return;
1833   if (!friendTrack->GetTPCOut()) return;
1834   if (!track->GetInnerParam())   return;
1835   if (!friendTrack->GetTPCOut())   return;
1836   const AliTrackPointArray *points=friendTrack->GetTrackPointArray();
1837   if (!points) return;
1838   AliExternalTrackParam pTPC(*(friendTrack->GetTPCOut()));
1839   AliExternalTrackParam pTOF(pTPC);
1840   Double_t mass = TDatabasePDG::Instance()->GetParticle("mu+")->Mass();
1841   Int_t npoints = points->GetNPoints();
1842   AliTrackPoint point;
1843   Int_t naccept=0;
1844   //
1845   for (Int_t ipoint=0;ipoint<npoints;ipoint++){
1846     points->GetPoint(point,ipoint);
1847     Float_t xyz[3];
1848     point.GetXYZ(xyz);
1849     Double_t r=TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]);
1850     if (r<350)  continue;
1851     if (r>400)  continue;
1852     AliTracker::PropagateTrackToBxByBz(&pTPC,r,mass,2.,kTRUE);
1853     AliTracker::PropagateTrackToBxByBz(&pTPC,r,mass,0.1,kTRUE);    
1854     AliTrackPoint lpoint = point.Rotate(pTPC.GetAlpha());
1855     pTPC.PropagateTo(lpoint.GetX(),fMagF);
1856     pTOF=pTPC;
1857     ((Double_t*)pTOF.GetParameter())[0] =lpoint.GetY();
1858     ((Double_t*)pTOF.GetParameter())[1] =lpoint.GetZ();
1859     ((Double_t*)pTOF.GetCovariance())[0]+=3.*3./12.;
1860     ((Double_t*)pTOF.GetCovariance())[2]+=3.*3./12.;
1861     ((Double_t*)pTOF.GetCovariance())[5]+=0.1*0.1;
1862     ((Double_t*)pTOF.GetCovariance())[9]+=0.1*0.1;
1863     naccept++;
1864   }
1865   if (naccept==0) return;  // no tof match clusters
1866   //
1867   // 0. Apply standard cuts
1868   //
1869   if (track->GetTPCNcls()<kMinTPC) return;  // minimal amount of clusters cut
1870   // exclude crossing track
1871   if (friendTrack->GetTPCOut()->GetZ()*track->GetInnerParam()->GetZ()<0)   return;
1872   //
1873   if (TMath::Abs(pTOF.GetY()-pTPC.GetY())    >kMaxDy)    return;
1874   if (TMath::Abs(pTOF.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1875   if (TMath::Abs(pTOF.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1876   //
1877   // 1. Update median and RMS info
1878   //
1879   TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
1880   TVectorD vecDeltaN(5);
1881   Double_t sign=(pTOF.GetParameter()[1]>0)? 1.:-1.;
1882   vecDelta[4]=0;
1883   for (Int_t i=0;i<4;i++){
1884     vecDelta[i]=(pTOF.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1885     kgdP[i][kglast%kN]=vecDelta[i];
1886   }
1887   kglast=(kglast+1);
1888   Int_t entries=(kglast<kN)?kglast:kN;
1889   Bool_t isOK=kTRUE;
1890   for (Int_t i=0;i<4;i++){
1891     vecMedian[i] = TMath::Median(entries,kgdP[i]);
1892     vecRMS[i]    = TMath::RMS(entries,kgdP[i]);
1893     vecDeltaN[i] = 0;
1894     if (vecRMS[i]>0.){
1895       vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/(vecRMS[i]+1.);
1896       vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]);  //sum of abs residuals
1897       if (TMath::Abs(vecDeltaN[i])>kSigmaCut) isOK=kFALSE;
1898     }
1899   }
1900   //
1901   // 2. Apply median+-rms cut
1902   //
1903   if (kglast<10)  return;   //median and RMS to be defined
1904   if (!isOK) return;
1905   //
1906   // 3. Update alignment
1907   //
1908   //Int_t htime = fTime/3600; //time in hours
1909   Int_t htime = (Int_t)(fTime-fTimeKalmanBin)/fTimeKalmanBin; //time bin
1910   if (fAlignTOFTPC->GetEntriesFast()<htime){
1911     fAlignTOFTPC->Expand(htime*2+20);
1912   }
1913   AliRelAlignerKalman* align =  (AliRelAlignerKalman*)fAlignTOFTPC->At(htime);
1914   if (!align){
1915     // make Alignment object if doesn't exist
1916     align=new AliRelAlignerKalman(); 
1917     align->SetRunNumber(fRun);
1918     (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1919     (*align->GetStateCov())(7,7)=kT0Err*kT0Err;
1920     (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1921     align->SetOutRejSigma(kOutCut+kOutCut*kN);
1922     //    align->SetRejectOutliers(kFALSE);
1923     align->SetRejectOutliers(kTRUE);
1924     align->SetRejectOutliersSigma2Median(kTRUE);
1925
1926     align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1927     align->SetMagField(fMagF); 
1928     fAlignTOFTPC->AddAt(align,htime);
1929   }
1930   align->AddTrackParams(&pTOF,&pTPC);
1931   Float_t dca[2],cov[3];
1932   track->GetImpactParameters(dca,cov);
1933   if (TMath::Abs(dca[0])<kMaxDy){
1934     FillResHistoTPCTOF(&pTPC,&pTOF);
1935   }
1936   //align->SetTimeStamp(fTime);
1937   Double_t averageTime =  fTime;
1938   if (align->GetTimeStamp()>0 && align->GetNUpdates()>0) {
1939     averageTime = (((Double_t)fTime) + ((Double_t)align->GetTimeStamp())*align->GetNUpdates()) / (align->GetNUpdates() + 1.);
1940     //printf("align->GetTimeStamp() %d, align->GetNUpdates() %d \n", align->GetTimeStamp(), align->GetNUpdates());
1941   }
1942   align->SetTimeStamp((Int_t)averageTime);
1943
1944   //printf("fTime %d, averageTime %d \n", fTime, (Int_t)averageTime);
1945
1946   align->SetRunNumber(fRun );
1947   //
1948   Int_t nupdates=align->GetNUpdates();
1949   align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates+1));
1950   //  align->SetRejectOutliers(kFALSE);
1951   align->SetRejectOutliers(kTRUE);
1952   align->SetRejectOutliersSigma2Median(kTRUE);
1953
1954   TTreeSRedirector *cstream = GetDebugStreamer();  
1955   if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1956     TVectorD gpTPC(3), gdTPC(3);
1957     TVectorD gpTOF(3), gdTOF(3);
1958     pTPC.GetXYZ(gpTPC.GetMatrixArray());
1959     pTPC.GetDirection(gdTPC.GetMatrixArray());
1960     pTOF.GetXYZ(gpTOF.GetMatrixArray());
1961     pTOF.GetDirection(gdTOF.GetMatrixArray());
1962     (*cstream)<<"toftpc"<<
1963       "run="<<fRun<<              //  run number
1964       "event="<<fEvent<<          //  event number
1965       "time="<<fTime<<            //  time stamp of event
1966       "trigger="<<fTrigger<<      //  trigger
1967       "mag="<<fMagF<<             //  magnetic field
1968       //
1969       "nmed="<<kglast<<        // number of entries to define median and RMS
1970       "vMed.="<<&vecMedian<<    // median of deltas
1971       "vRMS.="<<&vecRMS<<       // rms of deltas
1972       "vDelta.="<<&vecDelta<<   // delta in respect to median
1973       "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1974       "t.="<<track<<            // ful track - find proper cuts
1975       "a.="<<align<<            // current alignment
1976       "pTOF.="<<&pTOF<<         // track param TOF
1977       "pTPC.="<<&pTPC<<         // track param TPC
1978       "gpTPC.="<<&gpTPC<<       // global position  TPC
1979       "gdTPC.="<<&gdTPC<<       // global direction TPC
1980       "gpTOF.="<<&gpTOF<<       // global position  TOF
1981       "gdTOF.="<<&gdTOF<<       // global position  TOF
1982       "\n";
1983   }
1984 }
1985
1986
1987 void  AliTPCcalibTime::BookDistortionMaps(){
1988   //
1989   //   Book ndimensional histograms of distortions/residuals
1990   //   Only primary tracks are selected for analysis
1991   //
1992  
1993   Double_t xminTrack[5], xmaxTrack[5];
1994   Int_t binsTrack[5];
1995   TString axisName[5];
1996   TString axisTitle[5];
1997   //
1998   binsTrack[0]  =50;
1999   axisName[0]   ="#Delta";
2000   axisTitle[0]  ="#Delta";
2001   //
2002   binsTrack[1] =44;
2003   xminTrack[1] =-1.1; xmaxTrack[1]=1.1;
2004   axisName[1]  ="tanTheta";
2005   axisTitle[1]  ="tan(#Theta)";
2006   //
2007   binsTrack[2] =180;
2008   xminTrack[2] =-TMath::Pi(); xmaxTrack[2]=TMath::Pi(); 
2009   axisName[2]  ="phi";
2010   axisTitle[2]  ="#phi";
2011   //
2012   binsTrack[3] =20;
2013   xminTrack[3] =-1.; xmaxTrack[3]=1.;   // 0.33 GeV cut 
2014   axisName[3]  ="snp";
2015   axisTitle[3]  ="snp";
2016   //
2017   binsTrack[4] =10;
2018   xminTrack[4] =120.; xmaxTrack[4]=215.;   // crossing radius for CE only 
2019   axisName[4]  ="r";
2020   axisTitle[4] ="r(cm)";
2021   //
2022   // delta y
2023   xminTrack[0] =-1.5; xmaxTrack[0]=1.5;  // 
2024   fResHistoTPCCE[0] = new THnSparseS("TPCCE#Delta_{Y} (cm)","#Delta_{Y} (cm)",    5, binsTrack,xminTrack, xmaxTrack);
2025   fResHistoTPCITS[0] = new THnSparseS("TPCITS#Delta_{Y} (cm)","#Delta_{Y} (cm)",    4, binsTrack,xminTrack, xmaxTrack);
2026   fResHistoTPCvertex[0]    = new THnSparseS("TPCVertex#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2027   xminTrack[0] =-1.5; xmaxTrack[0]=1.5;  // 
2028   fResHistoTPCTRD[0] = new THnSparseS("TPCTRD#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2029   xminTrack[0] =-5; xmaxTrack[0]=5;  // 
2030   fResHistoTPCTOF[0] = new THnSparseS("TPCTOF#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2031   //
2032   // delta z
2033   xminTrack[0] =-6.; xmaxTrack[0]=6.;  // 
2034   fResHistoTPCCE[1] = new THnSparseS("TPCCE#Delta_{Z} (cm)","#Delta_{Z} (cm)",    5, binsTrack,xminTrack, xmaxTrack);
2035   fResHistoTPCITS[1] = new THnSparseS("TPCITS#Delta_{Z} (cm)","#Delta_{Z} (cm)",    4, binsTrack,xminTrack, xmaxTrack);
2036   fResHistoTPCvertex[1]    = new THnSparseS("TPCVertex#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2037   fResHistoTPCTRD[1] = new THnSparseS("TPCTRD#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2038   xminTrack[0] =-5.; xmaxTrack[0]=5.;  // 
2039   fResHistoTPCTOF[1] = new THnSparseS("TPCTOF#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2040   //
2041   // delta snp-P2
2042   xminTrack[0] =-0.015; xmaxTrack[0]=0.015;  // 
2043   fResHistoTPCCE[2] = new THnSparseS("TPCCE#Delta_{#phi}","#Delta_{#phi}",    5, binsTrack,xminTrack, xmaxTrack);
2044   fResHistoTPCITS[2] = new THnSparseS("TPCITS#Delta_{#phi}","#Delta_{#phi}",    4, binsTrack,xminTrack, xmaxTrack);
2045   fResHistoTPCvertex[2] = new THnSparseS("TPCITSv#Delta_{#phi}","#Delta_{#phi}",    4, binsTrack,xminTrack, xmaxTrack);
2046   fResHistoTPCTRD[2] = new THnSparseS("TPCTRD#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
2047   fResHistoTPCTOF[2] = new THnSparseS("TPCTOF#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
2048   //
2049   // delta theta-P3
2050   xminTrack[0] =-0.05; xmaxTrack[0]=0.05;  // 
2051   fResHistoTPCCE[3] = new THnSparseS("TPCCE#Delta_{#theta}","#Delta_{#theta}",    5, binsTrack,xminTrack, xmaxTrack);
2052   fResHistoTPCITS[3] = new THnSparseS("TPCITS#Delta_{#theta}","#Delta_{#theta}",    4, binsTrack,xminTrack, xmaxTrack);
2053   fResHistoTPCvertex[3] = new THnSparseS("TPCITSv#Delta_{#theta}","#Delta_{#theta}",    4, binsTrack,xminTrack, xmaxTrack);
2054   fResHistoTPCTRD[3] = new THnSparseS("TPCTRD#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
2055   fResHistoTPCTOF[3] = new THnSparseS("TPCTOF#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
2056   //
2057   // delta theta-P4
2058   xminTrack[0] =-0.2; xmaxTrack[0]=0.2;  // 
2059   fResHistoTPCCE[4] = new THnSparseS("TPCCE#Delta_{1/pt}","#Delta_{1/pt}",    5, binsTrack,xminTrack, xmaxTrack);
2060   fResHistoTPCITS[4] = new THnSparseS("TPCITS#Delta_{1/pt}","#Delta_{1/pt}",    4, binsTrack,xminTrack, xmaxTrack);
2061   fResHistoTPCvertex[4] = new THnSparseS("TPCITSv#Delta_{1/pt}","#Delta_{1/pt}",    4, binsTrack,xminTrack, xmaxTrack);
2062   fResHistoTPCTRD[4] = new THnSparseS("TPCTRD#Delta_{1/pt}","#Delta_{1/pt}",    4, binsTrack,xminTrack, xmaxTrack);
2063   fResHistoTPCTOF[4] = new THnSparseS("TPCTOF#Delta_{1/pt}","#Delta_{1/pt}",    4, binsTrack,xminTrack, xmaxTrack);
2064   //
2065   for (Int_t ivar=0;ivar<4;ivar++){
2066     for (Int_t ivar2=0;ivar2<5;ivar2++){      
2067       fResHistoTPCCE[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
2068       fResHistoTPCCE[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
2069       if (ivar2<4){
2070         fResHistoTPCITS[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
2071         fResHistoTPCITS[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
2072         fResHistoTPCTRD[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
2073         fResHistoTPCTRD[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
2074         fResHistoTPCvertex[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
2075         fResHistoTPCvertex[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
2076       }
2077     }
2078   }
2079   //
2080   // Book vertex: time histograms
2081   //
2082   Int_t    binsVertex[2]={500, fTimeBins};
2083   Double_t    aminVertex[2]={-5,fTimeStart};
2084   Double_t    amaxVertex[2]={5, fTimeEnd};
2085   const char* hnames[12]={"TPCXAside", "TPCXCside","TPCXACdiff","TPCXAPCdiff",
2086                           "TPCYAside", "TPCYCside","TPCYACdiff","TPCYAPCdiff",
2087                           "TPCZAPCside", "TPCZAMCside","TPCZACdiff","TPCZAPCdiff"}; 
2088   const char* anames[12]={"x (cm) - A side ", "x (cm) - C side","#Delta_{x} (cm) - TPC-A-C","#Delta_{x} (cm) - TPC-Common",
2089                           "y (cm) - A side ", "y (cm) - C side","#Delta_{x} (cm) - TPC-A-C","#Delta_{y} (cm) - TPC-Common",
2090                           "z (cm)", "#Delta_{Z} (cm) A-C side","#Delta_{x} (cm) - TPC-A-C","#Delta_{Z} (cm) TPC-common"}; 
2091   for (Int_t ihis=0; ihis<12; ihis++) {
2092     if (ihis>=8) aminVertex[0]=-20.;
2093     if (ihis>=8) amaxVertex[0]=20.;
2094     fTPCVertex[ihis]=new THnSparseF(hnames[ihis],hnames[ihis],2,binsVertex,aminVertex,amaxVertex);
2095     fTPCVertex[ihis]->GetAxis(1)->SetTitle("Time");
2096     fTPCVertex[ihis]->GetAxis(0)->SetTitle(anames[ihis]);
2097   }
2098   
2099   Int_t    binsVertexC[2]={40, 300};
2100   Double_t aminVertexC[2]={-20,-30};
2101   Double_t amaxVertexC[2]={20,30};
2102   const char* hnamesC[5]={"TPCA_TPC","TPCC_TPC","TPCA_ITS","TPCC_ITS","TPC_ITS"};
2103   for (Int_t ihis=0; ihis<5; ihis++) {
2104     fTPCVertexCorrelation[ihis]=new THnSparseF(hnamesC[ihis],hnamesC[ihis],2,binsVertexC,aminVertexC,amaxVertexC);
2105     fTPCVertexCorrelation[ihis]->GetAxis(1)->SetTitle("z (cm)");
2106     fTPCVertexCorrelation[ihis]->GetAxis(0)->SetTitle("z (cm)");
2107   }
2108 }
2109
2110
2111 void        AliTPCcalibTime::FillResHistoTPCCE(const AliExternalTrackParam * pTPCIn, const AliExternalTrackParam * pTPCOut ){
2112   //
2113   // fill residual histograms pTPCOut-pTPCin - trac crossing CE
2114   // Histogram 
2115   //
2116   if (fMemoryMode<2) return;
2117   Double_t histoX[5];
2118   Double_t xyz[3];
2119   pTPCIn->GetXYZ(xyz);
2120   Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2121   histoX[1]= pTPCIn->GetTgl();
2122   histoX[2]= phi;
2123   histoX[3]= pTPCIn->GetSnp();
2124   histoX[4]= pTPCIn->GetX();
2125   AliExternalTrackParam lout(*pTPCOut);
2126   lout.Rotate(pTPCIn->GetAlpha());
2127   lout.PropagateTo(pTPCIn->GetX(),fMagF);
2128   //
2129   for (Int_t ihisto=0; ihisto<5; ihisto++){
2130     histoX[0]=lout.GetParameter()[ihisto]-pTPCIn->GetParameter()[ihisto];
2131     fResHistoTPCCE[ihisto]->Fill(histoX);
2132   }
2133 }  
2134 void        AliTPCcalibTime::FillResHistoTPCITS(const AliExternalTrackParam * pTPCIn, const AliExternalTrackParam * pITSOut ){
2135   //
2136   // fill residual histograms pTPCIn-pITSOut
2137   // Histogram is filled only for primary tracks
2138   //
2139   Double_t histoX[4];
2140   Double_t xyz[3];
2141   pTPCIn->GetXYZ(xyz);
2142   Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2143   histoX[1]= pTPCIn->GetTgl();
2144   histoX[2]= phi;
2145   histoX[3]= pTPCIn->GetSnp();
2146   AliExternalTrackParam lits(*pITSOut);
2147   lits.Rotate(pTPCIn->GetAlpha());
2148   lits.PropagateTo(pTPCIn->GetX(),fMagF);
2149   //
2150   for (Int_t ihisto=0; ihisto<5; ihisto++){
2151     histoX[0]=pTPCIn->GetParameter()[ihisto]-lits.GetParameter()[ihisto];
2152     fResHistoTPCITS[ihisto]->Fill(histoX);
2153   }
2154 }  
2155
2156      
2157 void        AliTPCcalibTime::FillResHistoTPC(const AliESDtrack * pTrack){
2158   //
2159   // fill residual histograms pTPC - vertex
2160   // Histogram is filled only for primary tracks
2161   //
2162   if (fMemoryMode<2) return;
2163   Double_t histoX[4];
2164   const AliExternalTrackParam * pTPCIn = pTrack->GetInnerParam();
2165   AliExternalTrackParam pTPCvertex(*(pTrack->GetInnerParam()));
2166   //
2167   AliExternalTrackParam lits(*pTrack);
2168   if (TMath::Abs(pTrack->GetY())>3) return;  // beam pipe
2169   pTPCvertex.Rotate(lits.GetAlpha());
2170   //pTPCvertex.PropagateTo(pTPCvertex->GetX(),fMagF);
2171   AliTracker::PropagateTrackToBxByBz(&pTPCvertex,lits.GetX(),0.1,2,kFALSE);
2172   AliTracker::PropagateTrackToBxByBz(&pTPCvertex,lits.GetX(),0.1,0.1,kFALSE);
2173   Double_t xyz[3];
2174   pTPCIn->GetXYZ(xyz);
2175   Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2176   histoX[1]= pTPCIn->GetTgl();
2177   histoX[2]= phi;
2178   histoX[3]= pTPCIn->GetSnp();
2179   //
2180   Float_t dca[2], cov[3];
2181   pTrack->GetImpactParametersTPC(dca,cov);
2182   for (Int_t ihisto=0; ihisto<5; ihisto++){
2183     histoX[0]=pTPCvertex.GetParameter()[ihisto]-lits.GetParameter()[ihisto];
2184     //    if (ihisto<2) histoX[0]=dca[ihisto];
2185     fResHistoTPCvertex[ihisto]->Fill(histoX);
2186   }
2187 }
2188
2189
2190 void        AliTPCcalibTime::FillResHistoTPCTRD(const AliExternalTrackParam * pTPCOut, const AliExternalTrackParam * pTRDIn ){
2191   //
2192   // fill resuidual histogram TPCout-TRDin
2193   //
2194   if (fMemoryMode<2) return;
2195   Double_t histoX[4];
2196   Double_t xyz[3];
2197   pTPCOut->GetXYZ(xyz);
2198   Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2199   histoX[1]= pTPCOut->GetTgl();
2200   histoX[2]= phi;
2201   histoX[3]= pTPCOut->GetSnp();
2202   //
2203   AliExternalTrackParam ltrd(*pTRDIn);
2204   ltrd.Rotate(pTPCOut->GetAlpha());
2205   //  ltrd.PropagateTo(pTPCOut->GetX(),fMagF);
2206   AliTracker::PropagateTrackToBxByBz(&ltrd,pTPCOut->GetX(),0.1,0.1,kFALSE);
2207
2208   for (Int_t ihisto=0; ihisto<5; ihisto++){
2209     histoX[0]=pTPCOut->GetParameter()[ihisto]-ltrd.GetParameter()[ihisto];
2210     fResHistoTPCTRD[ihisto]->Fill(histoX);
2211   }
2212
2213 }
2214
2215 void        AliTPCcalibTime::FillResHistoTPCTOF(const AliExternalTrackParam * pTPCOut, const AliExternalTrackParam * pTOFIn ){
2216   //
2217   // fill resuidual histogram TPCout-TOFin
2218   // track propagated to the TOF position
2219   if (fMemoryMode<2) return;
2220   Double_t histoX[4];
2221   Double_t xyz[3];
2222
2223   AliExternalTrackParam ltpc(*pTPCOut);
2224   ltpc.Rotate(pTOFIn->GetAlpha());
2225   AliTracker::PropagateTrackToBxByBz(&ltpc,pTOFIn->GetX(),0.1,0.1,kFALSE);
2226   //
2227   ltpc.GetXYZ(xyz);
2228   Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2229   histoX[1]= ltpc.GetTgl();
2230   histoX[2]= phi;
2231   histoX[3]= ltpc.GetSnp();
2232   //
2233   for (Int_t ihisto=0; ihisto<2; ihisto++){
2234     histoX[0]=ltpc.GetParameter()[ihisto]-pTOFIn->GetParameter()[ihisto];
2235     fResHistoTPCTOF[ihisto]->Fill(histoX);
2236   }
2237
2238 }