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