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