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