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