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