]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibTime.cxx
Add test for sector-by-sector fits from CE calibration. Fail CE entry production
[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(100),
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(100),
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(5);
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   //
809   AliESDfriend *esdFriend=(AliESDfriend*)(((AliESDEvent*)event)->FindListObject("AliESDfriend"));
810   //
811   // Divide tracks to A and C side tracks - using the cluster indexes
812   TObjArray tracksA(ntracks);  
813   TObjArray tracksC(ntracks);  
814   //
815   AliESDVertex *vertexSPD =  (AliESDVertex *)event->GetPrimaryVertexSPD();
816   AliESDVertex *vertex    =  (AliESDVertex *)event->GetPrimaryVertex();
817   AliESDVertex *vertexTracks =  (AliESDVertex *)event->GetPrimaryVertexTracks();
818   Double_t vertexZA[10000], vertexZC[10000];
819   //
820   Int_t ntracksA= 0;
821   Int_t ntracksC= 0;
822   //
823   for (Int_t itrack=0;itrack<ntracks;itrack++) {
824     AliESDtrack *track = event->GetTrack(itrack);
825     AliESDfriendTrack *friendTrack = esdFriend->GetTrack(itrack);
826     if (!friendTrack) continue;
827     if (TMath::Abs(track->GetTgl())>kMaxTgl) continue;
828     if (TMath::Abs(track->Pt())<kMinPt) continue;
829     const AliExternalTrackParam * trackIn  = track->GetInnerParam();
830     TObject *calibObject=0;
831     AliTPCseed *seed = 0;
832     Int_t nA=0, nC=0;
833     for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
834     if (seed) {
835       for (Int_t irow=159;irow>0;irow--) {
836         AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
837         if (!cl) continue;
838         if ((cl->GetDetector()%36)<18) nA++;
839         if ((cl->GetDetector()%36)>=18) nC++;
840       }
841       if ((nA>kMinClusters || nC>kMinClusters) && (nA*nC==0) ){
842         track->GetImpactParameters(dca0[0],dca0[1]);
843         if (TMath::Abs(dca0[0])>kMaxD0) continue;
844         if (TMath::Abs(dca0[1])>kMaxZ0) continue;
845         AliExternalTrackParam pTPCvertex(*trackIn);
846         if (!AliTracker::PropagateTrackToBxByBz(&pTPCvertex,4.+4.*TMath::Abs(dca0[0]),0.1,2,kTRUE)) continue;
847         pTPCvertex.PropagateToDCA(vertex,AliTracker::GetBz(), kMaxD, dcaVertex,0);
848         if (TMath::Abs(dcaVertex[0])>kMaxD) continue;
849         if (nA>kMinClusters &&nC==0) { tracksA.AddLast(pTPCvertex.Clone()); vertexZA[ntracksA++] = pTPCvertex.GetZ();}
850         if (nC>kMinClusters &&nA==0) {tracksC.AddLast(pTPCvertex.Clone());  vertexZC[ntracksC++] = pTPCvertex.GetZ();}
851       }
852     }
853   }
854   Double_t medianZA=TMath::Median(ntracksA, vertexZA);  // tracks median
855   Double_t medianZC=TMath::Median(ntracksC, vertexZC);  // tracks median
856   //
857   ntracksA= tracksA.GetEntriesFast();
858   ntracksC= tracksC.GetEntriesFast();
859   if (ntracksA>kMinTracks && ntracksC>kMinTracks){
860     deltaZ[counterZ%kBuffSize]=medianZA-medianZC;
861     counterZ+=1;
862     Double_t medianDelta=(counterZ>=kBuffSize)? TMath::Median(kBuffSize, deltaZ): TMath::Median(counterZ, deltaZ);
863     if (TMath::Abs(medianDelta-(medianZA-medianZC))>kMaxZ) flags+=16;
864     // increse the error of cumulative vertex at the beginning of event
865     cumulVertexA.Covariance(0,0)+=kCumulCovarXY*kCumulCovarXY;
866     cumulVertexA.Covariance(1,1)+=kCumulCovarXY*kCumulCovarXY;
867     cumulVertexA.Covariance(2,2)+=kCumulCovarZ*kCumulCovarZ;
868     cumulVertexC.Covariance(0,0)+=kCumulCovarXY*kCumulCovarXY;
869     cumulVertexC.Covariance(1,1)+=kCumulCovarXY*kCumulCovarXY;
870     cumulVertexC.Covariance(2,2)+=kCumulCovarZ*kCumulCovarZ;
871     cumulVertexAC.Covariance(0,0)+=kCumulCovarXY*kCumulCovarXY;
872     cumulVertexAC.Covariance(1,1)+=kCumulCovarXY*kCumulCovarXY;
873     cumulVertexAC.Covariance(2,2)+=kCumulCovarZ*kCumulCovarZ;
874     //
875     for (Int_t iA=0; iA<ntracksA; iA++){
876       if (flags!=0) continue;
877       AliExternalTrackParam *aliTrack =  (AliExternalTrackParam *)tracksA.At(iA);
878       if (TMath::Abs(aliTrack->GetZ()-medianZA)>kMaxZ) continue;
879       AliKFParticle part(*aliTrack,211);
880       vertexA+=part;
881     }   
882     for (Int_t iC=0; iC<ntracksC; iC++){
883       if (flags!=0) continue;
884       AliExternalTrackParam *aliTrack =  (AliExternalTrackParam *)tracksC.At(iC);
885       if (TMath::Abs(aliTrack->GetZ()-medianZC)>kMaxZ) continue;
886       AliKFParticle part(*aliTrack,211);
887       vertexC+=part;
888     }    
889     //
890     if (vertexA.GetNDF()<kMinTracks) flags+=32;
891     if (vertexC.GetNDF()<kMinTracks) flags+=32;
892     if (TMath::Abs(vertexA.Z()-medianZA)>kMaxZ) flags+=1;   //apply cuts
893     if (TMath::Abs(vertexC.Z()-medianZC)>kMaxZ) flags+=2;
894     if (TMath::Abs(vertexA.GetChi2()/vertexA.GetNDF()+vertexC.GetChi2()/vertexC.GetNDF())> kMaxChi2) flags+=4;
895     //
896     if (flags==0){
897       for (Int_t iA=0; iA<ntracksA; iA++){
898         if (flags!=0) continue;
899         AliExternalTrackParam *aliTrack =  (AliExternalTrackParam *)tracksA.At(iA);
900         if (TMath::Abs(aliTrack->GetZ()-medianZA)>kMaxZ) continue;
901         AliKFParticle part(*aliTrack,211);
902         cumulVertexA+=part;
903         cumulVertexAC+=part;
904       } 
905       for (Int_t iC=0; iC<ntracksC; iC++){
906         if (flags!=0) continue;
907         AliExternalTrackParam *aliTrack =  (AliExternalTrackParam *)tracksC.At(iC);
908         if (TMath::Abs(aliTrack->GetZ()-medianZC)>kMaxZ) continue;
909         AliKFParticle part(*aliTrack,211);
910         cumulVertexC+=part;
911         cumulVertexAC+=part;
912       }      
913       //
914       if (TMath::Abs(cumulVertexA.X()-vertexA.X())>kMaxDvertex) flags+=64;
915       if (TMath::Abs(cumulVertexA.Y()-vertexA.Y())>kMaxDvertex) flags+=64;
916       if (TMath::Abs(cumulVertexA.Z()-vertexA.Z())>kMaxDvertex) flags+=64;
917       //
918       if (TMath::Abs(cumulVertexC.X()-vertexC.X())>kMaxDvertex) flags+=64;
919       if (TMath::Abs(cumulVertexC.Y()-vertexC.Y())>kMaxDvertex) flags+=64;
920       if (TMath::Abs(cumulVertexC.Z()-vertexC.Z())>kMaxDvertex) flags+=64;
921       
922       
923       if ( flags==0 && cumulVertexC.GetNDF()>kMinTracksVertex&&cumulVertexA.GetNDF()>kMinTracksVertex){
924         Double_t cont[2]={0,fTime};
925         //
926         cont[0]= cumulVertexA.X();
927         fTPCVertex[0]->Fill(cont);
928         cont[0]= cumulVertexC.X();
929         fTPCVertex[1]->Fill(cont);
930         cont[0]= 0.5*(cumulVertexA.X()-cumulVertexC.X());
931         fTPCVertex[2]->Fill(cont);
932         cont[0]= 0.5*(cumulVertexA.X()+cumulVertexC.X())-vertexSPD->GetX();
933         fTPCVertex[3]->Fill(cont);
934         //
935         cont[0]= cumulVertexA.Y();
936         fTPCVertex[4]->Fill(cont);
937         cont[0]= cumulVertexC.Y();
938         fTPCVertex[5]->Fill(cont);
939         cont[0]= 0.5*(cumulVertexA.Y()-cumulVertexC.Y());
940         fTPCVertex[6]->Fill(cont);
941         cont[0]= 0.5*(cumulVertexA.Y()+cumulVertexC.Y())-vertexSPD->GetY();
942         fTPCVertex[7]->Fill(cont);
943         //
944         //
945         cont[0]= 0.5*(cumulVertexA.Z()+cumulVertexC.Z());
946         fTPCVertex[8]->Fill(cont);
947         cont[0]= 0.5*(cumulVertexA.Z()-cumulVertexC.Z());
948         fTPCVertex[9]->Fill(cont);
949         cont[0]= 0.5*(cumulVertexA.Z()-cumulVertexC.Z());
950         fTPCVertex[10]->Fill(cont);      
951         cont[0]= 0.5*(cumulVertexA.Z()+cumulVertexC.Z())-vertexSPD->GetZ();
952         fTPCVertex[11]->Fill(cont);
953         //
954         Double_t correl[2]={0,0};
955         //
956         correl[0]=cumulVertexC.Z();
957         correl[1]=cumulVertexA.Z();
958         fTPCVertexCorrelation[0]->Fill(correl);   // fill A side :TPC
959         correl[0]=cumulVertexA.Z();
960         correl[1]=cumulVertexC.Z(); 
961         fTPCVertexCorrelation[1]->Fill(correl);   // fill C side :TPC
962         //
963         correl[0]=vertexSPD->GetZ();
964         correl[1]=cumulVertexA.Z()-correl[0];
965         fTPCVertexCorrelation[2]->Fill(correl);   // fill A side :ITS
966         correl[1]=cumulVertexC.Z()-correl[0]; 
967         fTPCVertexCorrelation[3]->Fill(correl);   // fill C side :ITS
968         correl[1]=0.5*(cumulVertexA.Z()+cumulVertexC.Z())-correl[0]; 
969         fTPCVertexCorrelation[4]->Fill(correl);   // fill C side :ITS
970       }
971     }        
972     TTreeSRedirector *cstream = GetDebugStreamer();
973     if (cstream){
974       /*
975         TCut cutChi2= "sqrt(vA.fChi2/vA.fNDF+vC.fChi2/vC.fNDF)<10";  // chi2 Cut e.g 10         
976         TCut cutXY= "sqrt((vA.fP[0]-vC.fP[0])^2+(vA.fP[0]-vC.fP[1])^2)<5";   // vertex Cut      
977         TCut cutZ= "abs(vA.fP[2]-mZA)<3&&abs(vC.fP[2]-mZC)<5";           // vertex Cut  
978         tree->Draw("sqrt(vA.fChi2/vA.fNDF)","sqrt(vA.fChi2/vA.fNDF)<100","")
979         
980       */
981       //vertexA.Print();
982       //vertexC.Print();      
983       (*cstream)<<"vertexTPC"<<
984         "flags="<<flags<<        // rejection flags
985         "vSPD.="<<vertexSPD<<    // SPD vertex
986         "vT.="<<vertexTracks<<   // track vertex
987         "v.="<<vertex<<          // esd vertex
988         "mZA="<<medianZA<<       // median Z position at vertex A side
989         "mZC="<<medianZC<<       // median Z position at vertex C side
990         "mDelta="<<medianDelta<< // median delta A side -C side
991         "counter="<<counterZ<<    // counter Z
992         //
993         "vA.="<<&vertexA<<       // vertex A side
994         "vC.="<<&vertexC<<       // vertex C side
995         "cvA.="<<&cumulVertexA<<       // cumulative vertex A side
996         "cvC.="<<&cumulVertexC<<       // cumulative vertex C side
997         "cvAC.="<<&cumulVertexAC<<       // cumulative vertex A+C side
998         "nA="<<ntracksA<<        // contributors
999         "nC="<<ntracksC<<        // contributors
1000         "\n";
1001     }      
1002   }
1003   tracksA.Delete();
1004   tracksC.Delete();
1005 }
1006
1007 void AliTPCcalibTime::Analyze(){
1008   //
1009   // Special macro to analyze result of calibration and extract calibration entries
1010   // Not yet ported to the Analyze function yet
1011   //
1012 }
1013
1014 THnSparse* AliTPCcalibTime::GetHistoDrift(const char* name) const
1015 {
1016   //
1017   // Get histogram for given trigger mask
1018   //
1019   TIterator* iterator = fArrayDz->MakeIterator();
1020   iterator->Reset();
1021   TString newName=name;
1022   newName.ToUpper();
1023   THnSparse* newHist=new THnSparseF(newName,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
1024   THnSparse* addHist=NULL;
1025   while((addHist=(THnSparseF*)iterator->Next())){
1026   if(!addHist) continue;
1027     TString histName=addHist->GetName();
1028     if(!histName.Contains(newName)) continue;
1029     addHist->Print();
1030     newHist->Add(addHist);
1031   }
1032   return newHist;
1033 }
1034
1035 TObjArray* AliTPCcalibTime::GetHistoDrift() const
1036 {
1037   //
1038   // return array of histograms
1039   //
1040   return fArrayDz;
1041 }
1042
1043 TGraphErrors* AliTPCcalibTime::GetGraphDrift(const char* name){
1044   //
1045   // Make a drift velocity (delta Z) graph
1046   //
1047   THnSparse* histoDrift=GetHistoDrift(name);
1048   TGraphErrors* graphDrift=NULL;
1049   if(histoDrift){
1050     graphDrift=FitSlices(histoDrift,2,0,400,100,0.05,0.95, kTRUE);
1051     TString end=histoDrift->GetName();
1052     Int_t pos=end.Index("_");
1053     end=end(pos,end.Capacity()-pos);
1054     TString graphName=graphDrift->ClassName();
1055     graphName+=end;
1056     graphName.ToUpper();
1057     graphDrift->SetName(graphName);
1058   }
1059   return graphDrift;
1060 }
1061
1062 TObjArray* AliTPCcalibTime::GetGraphDrift(){
1063   //
1064   // make a array of drift graphs
1065   //
1066   TObjArray* arrayGraphDrift=new TObjArray();
1067   TIterator* iterator=fArrayDz->MakeIterator();
1068   iterator->Reset();
1069   THnSparse* addHist=NULL;
1070   while((addHist=(THnSparseF*)iterator->Next())) arrayGraphDrift->AddLast(GetGraphDrift(addHist->GetName()));
1071   return arrayGraphDrift;
1072 }
1073
1074 AliSplineFit* AliTPCcalibTime::GetFitDrift(const char* name){
1075   //
1076   // Make a fit AliSplinefit  of drift velocity
1077   //
1078   TGraph* graphDrift=GetGraphDrift(name);
1079   AliSplineFit* fitDrift=NULL;
1080   if(graphDrift && graphDrift->GetN()){
1081     fitDrift=new AliSplineFit();
1082     fitDrift->SetGraph(graphDrift);
1083     fitDrift->SetMinPoints(graphDrift->GetN()+1);
1084     fitDrift->InitKnots(graphDrift,2,0,0.001);
1085     fitDrift->SplineFit(0);
1086     TString end=graphDrift->GetName();
1087     Int_t pos=end.Index("_");
1088     end=end(pos,end.Capacity()-pos);
1089     TString fitName=fitDrift->ClassName();
1090     fitName+=end;
1091     fitName.ToUpper();
1092     //fitDrift->SetName(fitName);
1093     delete graphDrift;
1094     graphDrift=NULL;
1095   }
1096   return fitDrift;
1097 }
1098
1099
1100 Long64_t AliTPCcalibTime::Merge(TCollection *const li) {
1101   //
1102   // Object specific merging procedure
1103   //
1104   TIterator* iter = li->MakeIterator();
1105   AliTPCcalibTime* cal = 0;
1106   //
1107   while ((cal = (AliTPCcalibTime*)iter->Next())) {
1108     if (!cal->InheritsFrom(AliTPCcalibTime::Class())) {
1109       Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName());
1110       return -1;
1111     }
1112     for (Int_t imeas=0; imeas<3; imeas++){
1113       if (cal->GetHistVdriftLaserA(imeas) && cal->GetHistVdriftLaserA(imeas)){
1114         fHistVdriftLaserA[imeas]->Add(cal->GetHistVdriftLaserA(imeas));
1115         fHistVdriftLaserC[imeas]->Add(cal->GetHistVdriftLaserC(imeas));
1116       }
1117     }
1118     //
1119     if (fTPCVertexCorrelation && cal->fTPCVertexCorrelation){
1120       for (Int_t imeas=0; imeas<5; imeas++){
1121         if (fTPCVertexCorrelation[imeas] && cal->fTPCVertexCorrelation[imeas]) fTPCVertexCorrelation[imeas]->Add(cal->fTPCVertexCorrelation[imeas]);
1122       }
1123     }
1124       
1125     if (fTPCVertex && cal->fTPCVertex) 
1126       for (Int_t imeas=0; imeas<12; imeas++){
1127         if (fTPCVertex[imeas] && cal->fTPCVertex[imeas]) fTPCVertex[imeas]->Add(cal->fTPCVertex[imeas]);
1128       }
1129     
1130     if (fMemoryMode>0) for (Int_t imeas=0; imeas<5; imeas++){
1131       if (fMemoryMode>1){
1132         if ( cal->GetResHistoTPCCE(imeas) && cal->GetResHistoTPCCE(imeas)){
1133           fResHistoTPCCE[imeas]->Add(cal->fResHistoTPCCE[imeas]);
1134         }else{
1135           fResHistoTPCCE[imeas]=(THnSparse*)cal->fResHistoTPCCE[imeas]->Clone();
1136         }
1137       }
1138       //
1139       if ((fMemoryMode>0) &&cal->GetResHistoTPCITS(imeas) && cal->GetResHistoTPCITS(imeas)){
1140         if (fMemoryMode>1 || (imeas%2)==1) fResHistoTPCITS[imeas]->Add(cal->fResHistoTPCITS[imeas]);
1141         if (fMemoryMode>1) fResHistoTPCvertex[imeas]->Add(cal->fResHistoTPCvertex[imeas]);
1142       }
1143       //
1144       if ((fMemoryMode>1) && cal->fResHistoTPCTRD[imeas]){
1145         if (fResHistoTPCTRD[imeas])
1146           fResHistoTPCTRD[imeas]->Add(cal->fResHistoTPCTRD[imeas]);
1147         else
1148           fResHistoTPCTRD[imeas]=(THnSparse*)cal->fResHistoTPCTRD[imeas]->Clone();
1149       }
1150       //
1151       if  ((fMemoryMode>1) && cal->fResHistoTPCTOF[imeas]){
1152         if (fResHistoTPCTOF[imeas])
1153           fResHistoTPCTOF[imeas]->Add(cal->fResHistoTPCTOF[imeas]);
1154         else
1155           fResHistoTPCTOF[imeas]=(THnSparse*)cal->fResHistoTPCTOF[imeas]->Clone();      
1156       }
1157       //
1158       if (cal->fArrayLaserA){
1159         fArrayLaserA->Expand(fArrayLaserA->GetEntriesFast()+cal->fArrayLaserA->GetEntriesFast());
1160         fArrayLaserC->Expand(fArrayLaserC->GetEntriesFast()+cal->fArrayLaserC->GetEntriesFast());
1161         for (Int_t ical=0; ical<cal->fArrayLaserA->GetEntriesFast(); ical++){
1162           if (cal->fArrayLaserA->UncheckedAt(ical)) fArrayLaserA->AddLast(cal->fArrayLaserA->UncheckedAt(ical)->Clone());
1163           if (cal->fArrayLaserC->UncheckedAt(ical)) fArrayLaserC->AddLast(cal->fArrayLaserC->UncheckedAt(ical)->Clone());
1164         }
1165       }
1166
1167     }
1168     TObjArray* addArray=cal->GetHistoDrift();
1169     if(!addArray) return 0;
1170     TIterator* iterator = addArray->MakeIterator();
1171     iterator->Reset();
1172     THnSparse* addHist=NULL;
1173     if ((fMemoryMode>1)) while((addHist=(THnSparseF*)iterator->Next())){
1174       if(!addHist) continue;
1175       addHist->Print();
1176       THnSparse* localHist=(THnSparseF*)fArrayDz->FindObject(addHist->GetName());
1177       if(!localHist){
1178         localHist=new THnSparseF(addHist->GetName(),"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
1179         fArrayDz->AddLast(localHist);
1180       }
1181       localHist->Add(addHist);
1182     }
1183
1184     for(Int_t i=0;i<10;i++) if (cal->GetCosmiMatchingHisto(i)) fCosmiMatchingHisto[i]->Add(cal->GetCosmiMatchingHisto(i));
1185     //
1186     // Merge alignment
1187     //
1188     for (Int_t itype=0; itype<3; itype++){
1189       //
1190       //
1191       TObjArray *arr0= 0;
1192       TObjArray *arr1= 0;
1193       if (itype==0) {arr0=fAlignITSTPC; arr1=cal->fAlignITSTPC;}
1194       if (itype==1) {arr0=fAlignTRDTPC; arr1=cal->fAlignTRDTPC;}
1195       if (itype==2) {arr0=fAlignTOFTPC; arr1=cal->fAlignTOFTPC;}
1196       if (!arr1) continue;
1197       if (!arr0) arr0=new TObjArray(arr1->GetEntriesFast());
1198       if (arr1->GetEntriesFast()>arr0->GetEntriesFast()){
1199         arr0->Expand(arr1->GetEntriesFast());
1200       }
1201       for (Int_t i=0;i<arr1->GetEntriesFast(); i++){
1202         AliRelAlignerKalman *kalman1 = (AliRelAlignerKalman *)arr1->UncheckedAt(i);
1203         AliRelAlignerKalman *kalman0 = (AliRelAlignerKalman *)arr0->UncheckedAt(i);
1204         if (!kalman1)  continue;
1205         if (!kalman0) {arr0->AddAt(new AliRelAlignerKalman(*kalman1),i); continue;}
1206         kalman0->SetRejectOutliers(kFALSE);
1207         kalman0->Merge(kalman1);
1208       }
1209     }
1210
1211   }
1212   return 0;
1213 }
1214
1215 Bool_t  AliTPCcalibTime::IsPair(AliExternalTrackParam *tr0, AliExternalTrackParam *tr1){
1216   /*
1217   // 0. Same direction - OPOSITE  - cutDir +cutT    
1218   TCut cutDir("cutDir","dir<-0.99")
1219   // 1. 
1220   TCut cutT("cutT","abs(Tr1.fP[3]+Tr0.fP[3])<0.03")
1221   //
1222   // 2. The same rphi 
1223   TCut cutD("cutD","abs(Tr0.fP[0]+Tr1.fP[0])<5")
1224   //
1225   TCut cutPt("cutPt","abs(Tr1.fP[4]+Tr0.fP[4])<1&&abs(Tr0.fP[4])+abs(Tr1.fP[4])<10");  
1226   // 1/Pt diff cut
1227   */
1228   const Double_t *p0 = tr0->GetParameter();
1229   const Double_t *p1 = tr1->GetParameter();
1230   fCosmiMatchingHisto[0]->Fill(p0[0]+p1[0]);
1231   fCosmiMatchingHisto[1]->Fill(p0[1]-p1[1]);
1232   fCosmiMatchingHisto[2]->Fill(tr1->GetAlpha()-tr0->GetAlpha()+TMath::Pi());
1233   fCosmiMatchingHisto[3]->Fill(p0[3]+p1[3]);
1234   fCosmiMatchingHisto[4]->Fill(p0[4]+p1[4]);
1235   
1236   if (TMath::Abs(p0[3]+p1[3])>fCutTheta) return kFALSE;
1237   if (TMath::Abs(p0[0]+p1[0])>fCutMaxD)  return kFALSE;
1238   if (TMath::Abs(p0[1]-p1[1])>fCutMaxDz)  return kFALSE;
1239   Double_t d0[3], d1[3];
1240   tr0->GetDirection(d0);    
1241   tr1->GetDirection(d1);       
1242   if (d0[0]*d1[0] + d0[1]*d1[1] + d0[2]*d1[2] >fCutMinDir) return kFALSE;
1243
1244   fCosmiMatchingHisto[5]->Fill(p0[0]+p1[0]);
1245   fCosmiMatchingHisto[6]->Fill(p0[1]-p1[1]);
1246   fCosmiMatchingHisto[7]->Fill(tr1->GetAlpha()-tr0->GetAlpha()+TMath::Pi());
1247   fCosmiMatchingHisto[8]->Fill(p0[3]+p1[3]);
1248   fCosmiMatchingHisto[9]->Fill(p0[4]+p1[4]);
1249
1250   return kTRUE;  
1251 }
1252 Bool_t AliTPCcalibTime::IsCross(AliESDtrack *const tr0, AliESDtrack *const tr1){
1253   //
1254   // check if the cosmic pair of tracks crossed A/C side
1255   // 
1256   Bool_t result= tr0->GetOuterParam()->GetZ()*tr1->GetOuterParam()->GetZ()<0;
1257   if (result==kFALSE) return result;
1258   result=kTRUE;
1259   return result;
1260 }
1261
1262 Bool_t AliTPCcalibTime::IsSame(AliESDtrack *const tr0, AliESDtrack *const tr1){
1263   // 
1264   // track crossing the CE
1265   // 0. minimal number of clusters 
1266   // 1. Same sector +-1
1267   // 2. Inner and outer track param on opposite side
1268   // 3. Outer and inner track parameter close each to other
1269   // 3. 
1270   Bool_t result=kTRUE;
1271   //
1272   // inner and outer on opposite sides in z
1273   //
1274   const Int_t knclCut0  = 30;
1275   const Double_t kalphaCut = 0.4;
1276   //
1277   // 0. minimal number of clusters
1278   //
1279   if (tr0->GetTPCNcls()<knclCut0) return kFALSE;
1280   if (tr1->GetTPCNcls()<knclCut0) return kFALSE;
1281   //
1282   // 1. alpha cut - sector+-1
1283   //
1284   if (TMath::Abs(tr0->GetOuterParam()->GetAlpha()-tr1->GetOuterParam()->GetAlpha())>kalphaCut) return kFALSE;
1285   //
1286   // 2. Z crossing
1287   //
1288   if (tr0->GetOuterParam()->GetZ()*tr0->GetInnerParam()->GetZ()>0) result&=kFALSE;
1289   if (tr1->GetOuterParam()->GetZ()*tr1->GetInnerParam()->GetZ()>0) result&=kFALSE;
1290   if (result==kFALSE){
1291     return result;
1292   }
1293   //
1294   //
1295   const Double_t *p0I = tr0->GetInnerParam()->GetParameter();
1296   const Double_t *p1I = tr1->GetInnerParam()->GetParameter();
1297   const Double_t *p0O = tr0->GetOuterParam()->GetParameter();
1298   const Double_t *p1O = tr1->GetOuterParam()->GetParameter();
1299   //
1300   if (TMath::Abs(p0I[0]-p1I[0])>fCutMaxD)  result&=kFALSE;
1301   if (TMath::Abs(p0I[1]-p1I[1])>fCutMaxDz) result&=kFALSE;
1302   if (TMath::Abs(p0I[2]-p1I[2])>fCutTheta) result&=kFALSE;
1303   if (TMath::Abs(p0I[3]-p1I[3])>fCutTheta) result&=kFALSE;
1304   if (TMath::Abs(p0O[0]-p1O[0])>fCutMaxD)  result&=kFALSE;
1305   if (TMath::Abs(p0O[1]-p1O[1])>fCutMaxDz) result&=kFALSE;
1306   if (TMath::Abs(p0O[2]-p1O[2])>fCutTheta) result&=kFALSE;
1307   if (TMath::Abs(p0O[3]-p1O[3])>fCutTheta) result&=kFALSE;
1308   if (result==kTRUE){
1309     result=kTRUE; // just to put break point here
1310   }
1311   return result;
1312 }
1313
1314
1315 void  AliTPCcalibTime::ProcessSame(AliESDtrack *const track, AliESDfriendTrack *const friendTrack, const AliESDEvent *const event){
1316   //
1317   // Process  TPC tracks crossing CE
1318   //
1319   // 0. Select only track crossing the CE
1320   // 1. Cut on the track length
1321   // 2. Refit the the track on A and C side separatelly
1322   // 3. Fill time histograms
1323   const Int_t kMinNcl=100;
1324   const Int_t kMinNclS=25;  // minimul number of clusters on the sides
1325   const Double_t pimass=TDatabasePDG::Instance()->GetParticle("pi+")->Mass();
1326   const Double_t kMaxDy=1;  // maximal distance in y
1327   const Double_t kMaxDsnp=0.05;  // maximal distance in snp
1328   const Double_t kMaxDtheta=0.05;  // maximal distance in theta
1329   
1330   if (!friendTrack->GetTPCOut()) return;
1331   //
1332   // 0. Select only track crossing the CE
1333   //
1334   if (track->GetInnerParam()->GetZ()*friendTrack->GetTPCOut()->GetZ()>0) return;
1335   //
1336   // 1. cut on track length
1337   //
1338   if (track->GetTPCNcls()<kMinNcl) return;
1339   //
1340   // 2. Refit track sepparatel on A and C side
1341   //
1342   TObject *calibObject;
1343   AliTPCseed *seed = 0;
1344   for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
1345     if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
1346   }
1347   if (!seed) return;
1348   //
1349   AliExternalTrackParam trackIn(*track->GetInnerParam());
1350   AliExternalTrackParam trackOut(*track->GetOuterParam());
1351   Double_t cov[3]={0.01,0.,0.01}; //use the same errors
1352   Double_t xyz[3]={0,0.,0.0};  
1353   Double_t bz   =0;
1354   Int_t nclIn=0,nclOut=0;
1355   trackIn.ResetCovariance(1000.);
1356   trackOut.ResetCovariance(1000.);
1357   //
1358   //2.a Refit inner
1359   // 
1360   Int_t sideIn=0;
1361   for (Int_t irow=0;irow<159;irow++) {
1362     AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
1363     if (!cl) continue;
1364     if (cl->GetX()<80) continue;
1365     if (sideIn==0){
1366       if (cl->GetDetector()%36<18) sideIn=1;
1367       if (cl->GetDetector()%36>=18) sideIn=-1;
1368     }
1369     if (sideIn== -1 && (cl->GetDetector()%36)<18) break;
1370     if (sideIn==  1 &&(cl->GetDetector()%36)>=18) break;
1371     Int_t sector = cl->GetDetector();
1372     Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackIn.GetAlpha();
1373     if (TMath::Abs(dalpha)>0.01){
1374       if (!trackIn.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
1375     }
1376     Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
1377     trackIn.GetXYZ(xyz);
1378     bz = AliTracker::GetBz(xyz);
1379     AliTracker::PropagateTrackToBxByBz(&trackIn,r[0],1.,pimass,kFALSE);
1380     if (!trackIn.PropagateTo(r[0],bz)) break;
1381     nclIn++;
1382     trackIn.Update(&r[1],cov);    
1383   }
1384   //
1385   //2.b Refit outer
1386   //
1387   Int_t sideOut=0;
1388   for (Int_t irow=159;irow>0;irow--) {
1389     AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
1390     if (!cl) continue;
1391     if (cl->GetX()<80) continue;
1392     if (sideOut==0){
1393       if (cl->GetDetector()%36<18) sideOut=1;
1394       if (cl->GetDetector()%36>=18) sideOut=-1;
1395       if (sideIn==sideOut) break;
1396     }
1397     if (sideOut== -1 && (cl->GetDetector()%36)<18) break;
1398     if (sideOut==  1 &&(cl->GetDetector()%36)>=18) break;
1399     //
1400     Int_t sector = cl->GetDetector();
1401     Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackOut.GetAlpha();
1402     if (TMath::Abs(dalpha)>0.01){
1403       if (!trackOut.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
1404     }
1405     Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
1406     trackOut.GetXYZ(xyz);
1407     bz = AliTracker::GetBz(xyz);
1408     AliTracker::PropagateTrackToBxByBz(&trackOut,r[0],1.,pimass,kFALSE);
1409     if (!trackOut.PropagateTo(r[0],bz)) break;
1410     nclOut++;
1411     trackOut.Update(&r[1],cov);    
1412   }
1413   trackOut.Rotate(trackIn.GetAlpha());
1414   Double_t meanX = (trackIn.GetX()+trackOut.GetX())*0.5;
1415   trackIn.PropagateTo(meanX,bz); 
1416   trackOut.PropagateTo(meanX,bz); 
1417   if (TMath::Abs(trackIn.GetY()-trackOut.GetY())>kMaxDy) return;
1418   if (TMath::Abs(trackIn.GetSnp()-trackOut.GetSnp())>kMaxDsnp) return;
1419   if (TMath::Abs(trackIn.GetTgl()-trackOut.GetTgl())>kMaxDtheta) return;
1420   if (TMath::Min(nclIn,nclOut)>kMinNclS){
1421     FillResHistoTPCCE(&trackIn,&trackOut);
1422   }
1423   TTreeSRedirector *cstream = GetDebugStreamer();
1424   if (cstream){
1425     TVectorD gxyz(3);
1426     trackIn.GetXYZ(gxyz.GetMatrixArray());
1427     TTimeStamp tstamp(fTime);
1428     (*cstream)<<"tpctpc"<<
1429       "run="<<fRun<<              //  run number
1430       "event="<<fEvent<<          //  event number
1431       "time="<<fTime<<            //  time stamp of event
1432       "trigger="<<fTrigger<<      //  trigger
1433       "mag="<<fMagF<<             //  magnetic field
1434       //
1435       "sideIn="<<sideIn<<         // side at inner part
1436       "sideOut="<<sideOut<<         // side at puter part
1437       "xyz.="<<&gxyz<<             // global position
1438       "tIn.="<<&trackIn<<         // refitterd track in 
1439       "tOut.="<<&trackOut<<       // refitter track out
1440       "nclIn="<<nclIn<<           // 
1441       "nclOut="<<nclOut<<         //
1442       "\n";  
1443   }
1444   //
1445   // 3. Fill time histograms
1446   // Debug stremaer expression
1447   // chainTPCTPC->Draw("(tIn.fP[1]-tOut.fP[1])*sign(-tIn.fP[3]):tIn.fP[3]","min(nclIn,nclOut)>30","")
1448   if (TMath::Min(nclIn,nclOut)>kMinNclS){
1449     fDz = trackOut.GetZ()-trackIn.GetZ();
1450     if (trackOut.GetTgl()<0) fDz*=-1.;
1451     TTimeStamp tstamp(fTime);
1452     Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
1453     Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
1454     Double_t vecDrift[4]={fTime,(ptrelative0+ptrelative1)/2.0,fDz/500.0,event->GetRunNumber()};
1455     //
1456     // fill histograms per trigger class and itegrated
1457     //
1458     THnSparse* curHist=NULL;
1459     for (Int_t itype=0; itype<2; itype++){
1460       TString name="MEAN_VDRIFT_CROSS_";  
1461       if (itype==0){
1462         name+=event->GetFiredTriggerClasses();
1463         name.ToUpper();
1464       }else{
1465         name+="ALL";
1466       }
1467       curHist=(THnSparseF*)fArrayDz->FindObject(name);
1468       if(!curHist){
1469         curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
1470         fArrayDz->AddLast(curHist);
1471       }
1472       curHist->Fill(vecDrift);
1473     }
1474   }
1475
1476 }
1477
1478 void  AliTPCcalibTime::ProcessAlignITS(AliESDtrack *const track, AliESDfriendTrack *const friendTrack, const AliESDEvent *const event, AliESDfriend *const esdFriend){
1479   //
1480   // Process track - Update TPC-ITS alignment
1481   // Updates: 
1482   // 0. Apply standartd cuts 
1483   // 1. Recalucluate the current statistic median/RMS
1484   // 2. Apply median+-rms cut
1485   // 3. Update kalman filter
1486   //
1487   const Int_t    kMinTPC  = 80;    // minimal number of TPC cluster
1488   const Int_t    kMinITS  = 3;     // minimal number of ITS cluster
1489   const Double_t kMinZ    = 10;    // maximal dz distance
1490   const Double_t kMaxDy   = 2.;    // maximal dy distance
1491   const Double_t kMaxAngle= 0.07;  // maximal angular distance
1492   const Double_t kSigmaCut= 5;     // maximal sigma distance to median
1493   const Double_t kVdErr   = 0.1;  // initial uncertainty of the vd correction 
1494   const Double_t kT0Err   = 3.;  // initial uncertainty of the T0 time
1495   const Double_t kVdYErr  = 0.05;  // initial uncertainty of the vd correction 
1496   const Double_t kOutCut  = 1.0;   // outlyer cut in AliRelAlgnmentKalman
1497   const Double_t kMinPt   = 0.3;   // minimal pt
1498   const Double_t kMax1Pt=0.5;        //maximal 1/pt distance
1499   const  Int_t     kN=50;         // deepnes of history
1500   static Int_t     kglast=0;
1501   static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1502   //
1503   // 0. Apply standard cuts
1504   // 
1505   Int_t dummycl[1000];
1506   if (track->GetTPCNcls()<kMinTPC) return;  // minimal amount of clusters cut
1507   if (!track->IsOn(AliESDtrack::kTPCrefit)) return;
1508   if (!track->GetInnerParam())   return;
1509   if (!track->GetOuterParam())   return;
1510   if (track->GetInnerParam()->Pt()<kMinPt)  return;
1511   // exclude crossing track
1512   if (track->GetOuterParam()->GetZ()*track->GetInnerParam()->GetZ()<0)   return;
1513   if (TMath::Abs(track->GetInnerParam()->GetZ())<kMinZ/3.)   return;
1514   if (track->GetInnerParam()->GetX()>90)   return;
1515   //
1516   AliExternalTrackParam &pTPC=(AliExternalTrackParam &)(*(track->GetInnerParam()));
1517   //  
1518   AliExternalTrackParam pITS;   // ITS standalone if possible
1519   AliExternalTrackParam pITS2;  //TPC-ITS track
1520   if (friendTrack->GetITSOut()){
1521     pITS2=(*(friendTrack->GetITSOut()));  //TPC-ITS track - snapshot ITS out
1522     pITS2.Rotate(pTPC.GetAlpha());
1523     AliTracker::PropagateTrackToBxByBz(&pITS2,pTPC.GetX(),0.1,0.1,kFALSE);
1524   }
1525
1526   AliESDfriendTrack *itsfriendTrack=0;
1527   //
1528   // try to find standalone ITS track corresponing to the TPC if possible
1529   //
1530   Bool_t hasAlone=kFALSE;
1531   Int_t ntracks=event->GetNumberOfTracks();
1532   for (Int_t i=0; i<ntracks; i++){
1533     AliESDtrack * trackITS = event->GetTrack(i); 
1534     if (!trackITS) continue;
1535     if (trackITS->GetITSclusters(dummycl)<kMinITS) continue;  // minimal amount of clusters
1536     itsfriendTrack = esdFriend->GetTrack(i);
1537     if (!itsfriendTrack) continue;
1538     if (!itsfriendTrack->GetITSOut()) continue;
1539      
1540     if (TMath::Abs(pTPC.GetTgl()-itsfriendTrack->GetITSOut()->GetTgl())> kMaxAngle) continue;
1541     if (TMath::Abs(pTPC.GetSigned1Pt()-itsfriendTrack->GetITSOut()->GetSigned1Pt())> kMax1Pt) continue;
1542     pITS=(*(itsfriendTrack->GetITSOut()));
1543     //
1544     pITS.Rotate(pTPC.GetAlpha());
1545     AliTracker::PropagateTrackToBxByBz(&pITS,pTPC.GetX(),0.1,0.1,kFALSE);
1546     if (TMath::Abs(pTPC.GetY()-pITS.GetY())> kMaxDy) continue;
1547     if (TMath::Abs(pTPC.GetSnp()-pITS.GetSnp())> kMaxAngle) continue;
1548     hasAlone=kTRUE;
1549   }
1550   if (!hasAlone) {
1551     if (track->GetITSclusters(dummycl)<kMinITS) return;
1552     pITS=pITS2;  // use combined track if it has ITS
1553   }
1554   //
1555   if (TMath::Abs(pITS.GetY()-pTPC.GetY())    >kMaxDy)    return;
1556   if (TMath::Abs(pITS.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1557   if (TMath::Abs(pITS.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1558   //
1559   // 1. Update median and RMS info
1560   //
1561   TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
1562   TVectorD vecDeltaN(5);
1563   Double_t sign=(pITS.GetParameter()[1]>0)? 1.:-1.;
1564   vecDelta[4]=0;
1565   for (Int_t i=0;i<4;i++){
1566     vecDelta[i]=(pITS.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1567     kgdP[i][kglast%kN]=vecDelta[i];
1568   }
1569   kglast=(kglast+1);
1570   Int_t entries=(kglast<kN)?kglast:kN;
1571   for (Int_t i=0;i<4;i++){
1572     vecMedian[i] = TMath::Median(entries,kgdP[i]);
1573     vecRMS[i]    = TMath::RMS(entries,kgdP[i]);
1574     vecDeltaN[i] = 0;
1575     if (vecRMS[i]>0.){
1576       vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/vecRMS[i];
1577       vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]);  //sum of abs residuals
1578     }
1579   }
1580   //
1581   // 2. Apply median+-rms cut
1582   //
1583   if (kglast<3)  return;   //median and RMS to be defined
1584   if ( vecDeltaN[4]/4.>kSigmaCut) return;
1585   //
1586   // 3. Update alignment
1587   //
1588   Int_t htime = (fTime-fTimeKalmanBin/2)/fTimeKalmanBin; //time bins number
1589   if (fAlignITSTPC->GetEntriesFast()<htime){
1590     fAlignITSTPC->Expand(htime*2+20);
1591   }
1592   AliRelAlignerKalman* align =  (AliRelAlignerKalman*)fAlignITSTPC->At(htime);
1593   if (!align){
1594     // make Alignment object if doesn't exist
1595     align=new AliRelAlignerKalman(); 
1596     align->SetRunNumber(fRun);
1597     (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1598     (*align->GetStateCov())(7,7)=kT0Err*kT0Err;
1599     (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1600     align->SetOutRejSigma(kOutCut+kOutCut*kN);
1601     align->SetRejectOutliers(kFALSE);
1602
1603     align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1604     align->SetMagField(fMagF); 
1605     fAlignITSTPC->AddAt(align,htime);
1606   }
1607   align->AddTrackParams(&pITS,&pTPC);
1608   Double_t averageTime =  fTime;
1609   if (align->GetTimeStamp()>0&&align->GetNUpdates()>0){
1610     averageTime=((Double_t(align->GetTimeStamp())*Double_t(align->GetNUpdates())+Double_t(fTime)))/(Double_t(align->GetNUpdates())+1.);
1611   }
1612   align->SetTimeStamp(Int_t(averageTime));
1613
1614   align->SetRunNumber(fRun );
1615   Float_t dca[2],cov[3];
1616   track->GetImpactParameters(dca,cov);
1617   if (TMath::Abs(dca[0])<kMaxDy){
1618     FillResHistoTPCITS(&pTPC,&pITS);
1619     FillResHistoTPC(track);
1620   }
1621   //
1622   Int_t nupdates=align->GetNUpdates();
1623   align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1624   align->SetRejectOutliers(kFALSE);
1625   TTreeSRedirector *cstream = GetDebugStreamer();  
1626   if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1627     TVectorD gpTPC(3), gdTPC(3);
1628     TVectorD gpITS(3), gdITS(3);
1629     pTPC.GetXYZ(gpTPC.GetMatrixArray());
1630     pTPC.GetDirection(gdTPC.GetMatrixArray());
1631     pITS.GetXYZ(gpITS.GetMatrixArray());
1632     pITS.GetDirection(gdITS.GetMatrixArray());
1633     (*cstream)<<"itstpc"<<
1634       "run="<<fRun<<              //  run number
1635       "event="<<fEvent<<          //  event number
1636       "time="<<fTime<<            //  time stamp of event
1637       "trigger="<<fTrigger<<      //  trigger
1638       "mag="<<fMagF<<             //  magnetic field
1639       //
1640       "hasAlone="<<hasAlone<<    // has ITS standalone ?
1641       "track.="<<track<<  // track info
1642       "nmed="<<kglast<<        // number of entries to define median and RMS
1643       "vMed.="<<&vecMedian<<    // median of deltas
1644       "vRMS.="<<&vecRMS<<       // rms of deltas
1645       "vDelta.="<<&vecDelta<<   // delta in respect to median
1646       "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1647       "t.="<<track<<            // ful track - find proper cuts
1648       "a.="<<align<<            // current alignment
1649       "pITS.="<<&pITS<<         // track param ITS 
1650       "pITS2.="<<&pITS2<<       // track param ITS+TPC
1651       "pTPC.="<<&pTPC<<         // track param TPC
1652       "gpTPC.="<<&gpTPC<<       // global position  TPC
1653       "gdTPC.="<<&gdTPC<<       // global direction TPC
1654       "gpITS.="<<&gpITS<<       // global position  ITS
1655       "gdITS.="<<&gdITS<<       // global position  ITS
1656       "\n";
1657   }
1658 }
1659
1660
1661
1662
1663 void  AliTPCcalibTime::ProcessAlignTRD(AliESDtrack *const track, AliESDfriendTrack *const friendTrack){
1664   //
1665   // Process track - Update TPC-TRD alignment
1666   // Updates: 
1667   // 0. Apply standartd cuts 
1668   // 1. Recalucluate the current statistic median/RMS
1669   // 2. Apply median+-rms cut
1670   // 3. Update kalman filter
1671   //
1672   const Int_t    kMinTPC  = 80;    // minimal number of TPC cluster
1673   const Int_t    kMinTRD  = 50;    // minimal number of TRD cluster
1674   const Double_t kMinZ    = 20;    // maximal dz distance
1675   const Double_t kMaxDy   = 5.;    // maximal dy distance
1676   const Double_t kMaxAngle= 0.1;  // maximal angular distance
1677   const Double_t kSigmaCut= 10;     // maximal sigma distance to median
1678   const Double_t kVdErr   = 0.1;  // initial uncertainty of the vd correction 
1679   const Double_t kT0Err   = 3.;  // initial uncertainty of the T0 time
1680   const Double_t kVdYErr  = 0.05;  // initial uncertainty of the vd correction 
1681   const Double_t kOutCut  = 1.0;   // outlyer cut in AliRelAlgnmentKalman
1682   const Double_t kRefX    = 275;   // reference X
1683   const  Int_t     kN=50;         // deepnes of history
1684   static Int_t     kglast=0;
1685   static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1686   //
1687   // 0. Apply standard cuts
1688   //
1689   Int_t dummycl[1000];
1690   if (track->GetTRDclusters(dummycl)<kMinTRD) return;  // minimal amount of clusters
1691   if (track->GetTPCNcls()<kMinTPC) return;  // minimal amount of clusters cut
1692   if (!friendTrack->GetTRDIn()) return;  
1693   if (!track->IsOn(AliESDtrack::kTRDrefit)) return;  
1694   if (!track->IsOn(AliESDtrack::kTRDout)) return;  
1695   if (!track->GetInnerParam())   return;
1696   if (!friendTrack->GetTPCOut())   return;
1697   // exclude crossing track
1698   if (friendTrack->GetTPCOut()->GetZ()*track->GetInnerParam()->GetZ()<0)   return;
1699   if (TMath::Abs(track->GetInnerParam()->GetZ())<kMinZ)   return;
1700   //
1701   AliExternalTrackParam &pTPC=(AliExternalTrackParam &)(*(friendTrack->GetTPCOut()));
1702   AliTracker::PropagateTrackToBxByBz(&pTPC,kRefX,0.1,0.1,kFALSE);
1703   AliExternalTrackParam pTRD(*(friendTrack->GetTRDIn()));
1704   pTRD.Rotate(pTPC.GetAlpha());
1705   //  pTRD.PropagateTo(pTPC.GetX(),fMagF);
1706   AliTracker::PropagateTrackToBxByBz(&pTRD,pTPC.GetX(),0.1,0.1,kFALSE);
1707
1708   ((Double_t*)pTRD.GetCovariance())[2]+=3.*3.;   // increas sys errors
1709   ((Double_t*)pTRD.GetCovariance())[9]+=0.1*0.1; // increse sys errors
1710
1711   if (TMath::Abs(pTRD.GetY()-pTPC.GetY())    >kMaxDy)    return;
1712   if (TMath::Abs(pTRD.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1713   //  if (TMath::Abs(pTRD.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1714   //
1715   // 1. Update median and RMS info
1716   //
1717   TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
1718   TVectorD vecDeltaN(5);
1719   Double_t sign=(pTRD.GetParameter()[1]>0)? 1.:-1.;
1720   vecDelta[4]=0;
1721   for (Int_t i=0;i<4;i++){
1722     vecDelta[i]=(pTRD.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1723     kgdP[i][kglast%kN]=vecDelta[i];
1724   }
1725   kglast=(kglast+1);
1726   Int_t entries=(kglast<kN)?kglast:kN;
1727   for (Int_t i=0;i<4;i++){
1728     vecMedian[i] = TMath::Median(entries,kgdP[i]);
1729
1730     vecRMS[i]    = TMath::RMS(entries,kgdP[i]);
1731     vecDeltaN[i] = 0;
1732     if (vecRMS[i]>0.){
1733       vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/vecRMS[i];
1734       vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]);  //sum of abs residuals
1735     }
1736   }
1737   //
1738   // 2. Apply median+-rms cut
1739   //
1740   if (kglast<3)  return;   //median and RMS to be defined
1741   if ( vecDeltaN[4]/4.>kSigmaCut) return;
1742   //
1743   // 3. Update alignment
1744   //
1745   //Int_t htime = fTime/3600; //time in hours
1746   Int_t htime = (Int_t)(fTime-fTimeKalmanBin/2)/fTimeKalmanBin; //time in half hour
1747   if (fAlignTRDTPC->GetEntriesFast()<htime){
1748     fAlignTRDTPC->Expand(htime*2+20);
1749   }
1750   AliRelAlignerKalman* align =  (AliRelAlignerKalman*)fAlignTRDTPC->At(htime);
1751   if (!align){
1752     // make Alignment object if doesn't exist
1753     align=new AliRelAlignerKalman(); 
1754     align->SetRunNumber(fRun);
1755     (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1756     (*align->GetStateCov())(7,7)=kT0Err*kT0Err;
1757     (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1758     align->SetOutRejSigma(kOutCut+kOutCut*kN);
1759     align->SetRejectOutliers(kFALSE);
1760     align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1761     align->SetMagField(fMagF); 
1762     fAlignTRDTPC->AddAt(align,htime);
1763   }
1764   align->AddTrackParams(&pTRD,&pTPC);
1765   //align->SetTimeStamp(fTime);
1766   Double_t averageTime =  fTime;
1767   if (align->GetTimeStamp()>0 && align->GetNUpdates()>0) {
1768     averageTime = (((Double_t)fTime) + ((Double_t)align->GetTimeStamp())*align->GetNUpdates()) / (align->GetNUpdates() + 1.);
1769     //printf("align->GetTimeStamp() %d, align->GetNUpdates() %d \n", align->GetTimeStamp(), align->GetNUpdates());
1770   }
1771   align->SetTimeStamp((Int_t)averageTime);
1772
1773   //printf("fTime %d, averageTime %d \n", fTime, (Int_t)averageTime);
1774
1775   align->SetRunNumber(fRun );
1776   Float_t dca[2],cov[3];
1777   track->GetImpactParameters(dca,cov);
1778   if (TMath::Abs(dca[0])<kMaxDy){
1779     FillResHistoTPCTRD(&pTPC,&pTRD);  //only primaries
1780   }
1781   //
1782   Int_t nupdates=align->GetNUpdates();
1783   align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1784   align->SetRejectOutliers(kFALSE);
1785   TTreeSRedirector *cstream = GetDebugStreamer();  
1786   if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1787     TVectorD gpTPC(3), gdTPC(3);
1788     TVectorD gpTRD(3), gdTRD(3);
1789     pTPC.GetXYZ(gpTPC.GetMatrixArray());
1790     pTPC.GetDirection(gdTPC.GetMatrixArray());
1791     pTRD.GetXYZ(gpTRD.GetMatrixArray());
1792     pTRD.GetDirection(gdTRD.GetMatrixArray());
1793     (*cstream)<<"trdtpc"<<
1794       "run="<<fRun<<              //  run number
1795       "event="<<fEvent<<          //  event number
1796       "time="<<fTime<<            //  time stamp of event
1797       "trigger="<<fTrigger<<      //  trigger
1798       "mag="<<fMagF<<             //  magnetic field
1799       //
1800       "nmed="<<kglast<<        // number of entries to define median and RMS
1801       "vMed.="<<&vecMedian<<    // median of deltas
1802       "vRMS.="<<&vecRMS<<       // rms of deltas
1803       "vDelta.="<<&vecDelta<<   // delta in respect to median
1804       "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1805       "t.="<<track<<            // ful track - find proper cuts
1806       "a.="<<align<<            // current alignment
1807       "pTRD.="<<&pTRD<<         // track param TRD
1808       "pTPC.="<<&pTPC<<         // track param TPC
1809       "gpTPC.="<<&gpTPC<<       // global position  TPC
1810       "gdTPC.="<<&gdTPC<<       // global direction TPC
1811       "gpTRD.="<<&gpTRD<<       // global position  TRD
1812       "gdTRD.="<<&gdTRD<<       // global position  TRD
1813       "\n";
1814   }
1815 }
1816
1817
1818 void  AliTPCcalibTime::ProcessAlignTOF(AliESDtrack *const track, AliESDfriendTrack *const friendTrack){
1819   //
1820   //
1821   // Process track - Update TPC-TOF alignment
1822   // Updates: 
1823   // -1. Make a TOF "track"
1824   // 0. Apply standartd cuts 
1825   // 1. Recalucluate the current statistic median/RMS
1826   // 2. Apply median+-rms cut
1827   // 3. Update kalman filter
1828   //
1829   const Int_t      kMinTPC  = 80;    // minimal number of TPC cluster
1830   //  const Double_t   kMinZ    = 10;    // maximal dz distance
1831   const Double_t   kMaxDy   = 5.;    // maximal dy distance
1832   const Double_t   kMaxAngle= 0.05;  // maximal angular distance
1833   const Double_t   kSigmaCut= 5;     // maximal sigma distance to median
1834   const Double_t   kVdErr   = 0.1;  // initial uncertainty of the vd correction 
1835   const Double_t   kT0Err   = 3.;  // initial uncertainty of the T0 time
1836   const Double_t   kVdYErr  = 0.05;  // initial uncertainty of the vd correction 
1837
1838   const Double_t   kOutCut  = 1.0;   // outlyer cut in AliRelAlgnmentKalman
1839   const  Int_t     kN=50;         // deepnes of history
1840   static Int_t     kglast=0;
1841   static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1842   //
1843   // -1. Make a TOF track-
1844   //     Clusters are not in friends - use alingment points
1845   //
1846   if (track->GetTOFsignal()<=0)  return;
1847   if (!friendTrack->GetTPCOut()) return;
1848   if (!track->GetInnerParam())   return;
1849   if (!friendTrack->GetTPCOut())   return;
1850   const AliTrackPointArray *points=friendTrack->GetTrackPointArray();
1851   if (!points) return;
1852   AliExternalTrackParam pTPC(*(friendTrack->GetTPCOut()));
1853   AliExternalTrackParam pTOF(pTPC);
1854   Double_t mass = TDatabasePDG::Instance()->GetParticle("mu+")->Mass();
1855   Int_t npoints = points->GetNPoints();
1856   AliTrackPoint point;
1857   Int_t naccept=0;
1858   //
1859   for (Int_t ipoint=0;ipoint<npoints;ipoint++){
1860     points->GetPoint(point,ipoint);
1861     Float_t xyz[3];
1862     point.GetXYZ(xyz);
1863     Double_t r=TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]);
1864     if (r<350)  continue;
1865     if (r>400)  continue;
1866     AliTracker::PropagateTrackToBxByBz(&pTPC,r,mass,2.,kTRUE);
1867     AliTracker::PropagateTrackToBxByBz(&pTPC,r,mass,0.1,kTRUE);    
1868     AliTrackPoint lpoint = point.Rotate(pTPC.GetAlpha());
1869     pTPC.PropagateTo(lpoint.GetX(),fMagF);
1870     pTOF=pTPC;
1871     ((Double_t*)pTOF.GetParameter())[0] =lpoint.GetY();
1872     ((Double_t*)pTOF.GetParameter())[1] =lpoint.GetZ();
1873     ((Double_t*)pTOF.GetCovariance())[0]+=3.*3./12.;
1874     ((Double_t*)pTOF.GetCovariance())[2]+=3.*3./12.;
1875     ((Double_t*)pTOF.GetCovariance())[5]+=0.1*0.1;
1876     ((Double_t*)pTOF.GetCovariance())[9]+=0.1*0.1;
1877     naccept++;
1878   }
1879   if (naccept==0) return;  // no tof match clusters
1880   //
1881   // 0. Apply standard cuts
1882   //
1883   if (track->GetTPCNcls()<kMinTPC) return;  // minimal amount of clusters cut
1884   // exclude crossing track
1885   if (friendTrack->GetTPCOut()->GetZ()*track->GetInnerParam()->GetZ()<0)   return;
1886   //
1887   if (TMath::Abs(pTOF.GetY()-pTPC.GetY())    >kMaxDy)    return;
1888   if (TMath::Abs(pTOF.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1889   if (TMath::Abs(pTOF.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1890   //
1891   // 1. Update median and RMS info
1892   //
1893   TVectorD vecDelta(5),vecMedian(5), vecRMS(5);
1894   TVectorD vecDeltaN(5);
1895   Double_t sign=(pTOF.GetParameter()[1]>0)? 1.:-1.;
1896   vecDelta[4]=0;
1897   for (Int_t i=0;i<4;i++){
1898     vecDelta[i]=(pTOF.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1899     kgdP[i][kglast%kN]=vecDelta[i];
1900   }
1901   kglast=(kglast+1);
1902   Int_t entries=(kglast<kN)?kglast:kN;
1903   Bool_t isOK=kTRUE;
1904   for (Int_t i=0;i<4;i++){
1905     vecMedian[i] = TMath::Median(entries,kgdP[i]);
1906     vecRMS[i]    = TMath::RMS(entries,kgdP[i]);
1907     vecDeltaN[i] = 0;
1908     if (vecRMS[i]>0.){
1909       vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/(vecRMS[i]+1.);
1910       vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]);  //sum of abs residuals
1911       if (TMath::Abs(vecDeltaN[i])>kSigmaCut) isOK=kFALSE;
1912     }
1913   }
1914   //
1915   // 2. Apply median+-rms cut
1916   //
1917   if (kglast<10)  return;   //median and RMS to be defined
1918   if (!isOK) return;
1919   //
1920   // 3. Update alignment
1921   //
1922   //Int_t htime = fTime/3600; //time in hours
1923   Int_t htime = (Int_t)(fTime-fTimeKalmanBin)/fTimeKalmanBin; //time bin
1924   if (fAlignTOFTPC->GetEntriesFast()<htime){
1925     fAlignTOFTPC->Expand(htime*2+20);
1926   }
1927   AliRelAlignerKalman* align =  (AliRelAlignerKalman*)fAlignTOFTPC->At(htime);
1928   if (!align){
1929     // make Alignment object if doesn't exist
1930     align=new AliRelAlignerKalman(); 
1931     align->SetRunNumber(fRun);
1932     (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1933     (*align->GetStateCov())(7,7)=kT0Err*kT0Err;
1934     (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1935     align->SetOutRejSigma(kOutCut+kOutCut*kN);
1936     align->SetRejectOutliers(kFALSE);
1937     align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1938     align->SetMagField(fMagF); 
1939     fAlignTOFTPC->AddAt(align,htime);
1940   }
1941   align->AddTrackParams(&pTOF,&pTPC);
1942   Float_t dca[2],cov[3];
1943   track->GetImpactParameters(dca,cov);
1944   if (TMath::Abs(dca[0])<kMaxDy){
1945     FillResHistoTPCTOF(&pTPC,&pTOF);
1946   }
1947   //align->SetTimeStamp(fTime);
1948   Double_t averageTime =  fTime;
1949   if (align->GetTimeStamp()>0 && align->GetNUpdates()>0) {
1950     averageTime = (((Double_t)fTime) + ((Double_t)align->GetTimeStamp())*align->GetNUpdates()) / (align->GetNUpdates() + 1.);
1951     //printf("align->GetTimeStamp() %d, align->GetNUpdates() %d \n", align->GetTimeStamp(), align->GetNUpdates());
1952   }
1953   align->SetTimeStamp((Int_t)averageTime);
1954
1955   //printf("fTime %d, averageTime %d \n", fTime, (Int_t)averageTime);
1956
1957   align->SetRunNumber(fRun );
1958   //
1959   Int_t nupdates=align->GetNUpdates();
1960   align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1961   align->SetRejectOutliers(kFALSE);
1962   TTreeSRedirector *cstream = GetDebugStreamer();  
1963   if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1964     TVectorD gpTPC(3), gdTPC(3);
1965     TVectorD gpTOF(3), gdTOF(3);
1966     pTPC.GetXYZ(gpTPC.GetMatrixArray());
1967     pTPC.GetDirection(gdTPC.GetMatrixArray());
1968     pTOF.GetXYZ(gpTOF.GetMatrixArray());
1969     pTOF.GetDirection(gdTOF.GetMatrixArray());
1970     (*cstream)<<"toftpc"<<
1971       "run="<<fRun<<              //  run number
1972       "event="<<fEvent<<          //  event number
1973       "time="<<fTime<<            //  time stamp of event
1974       "trigger="<<fTrigger<<      //  trigger
1975       "mag="<<fMagF<<             //  magnetic field
1976       //
1977       "nmed="<<kglast<<        // number of entries to define median and RMS
1978       "vMed.="<<&vecMedian<<    // median of deltas
1979       "vRMS.="<<&vecRMS<<       // rms of deltas
1980       "vDelta.="<<&vecDelta<<   // delta in respect to median
1981       "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1982       "t.="<<track<<            // ful track - find proper cuts
1983       "a.="<<align<<            // current alignment
1984       "pTOF.="<<&pTOF<<         // track param TOF
1985       "pTPC.="<<&pTPC<<         // track param TPC
1986       "gpTPC.="<<&gpTPC<<       // global position  TPC
1987       "gdTPC.="<<&gdTPC<<       // global direction TPC
1988       "gpTOF.="<<&gpTOF<<       // global position  TOF
1989       "gdTOF.="<<&gdTOF<<       // global position  TOF
1990       "\n";
1991   }
1992 }
1993
1994
1995 void  AliTPCcalibTime::BookDistortionMaps(){
1996   //
1997   //   Book ndimensional histograms of distortions/residuals
1998   //   Only primary tracks are selected for analysis
1999   //
2000  
2001   Double_t xminTrack[5], xmaxTrack[5];
2002   Int_t binsTrack[5];
2003   TString axisName[5];
2004   TString axisTitle[5];
2005   //
2006   binsTrack[0]  =50;
2007   axisName[0]   ="#Delta";
2008   axisTitle[0]  ="#Delta";
2009   //
2010   binsTrack[1] =44;
2011   xminTrack[1] =-1.1; xmaxTrack[1]=1.1;
2012   axisName[1]  ="tanTheta";
2013   axisTitle[1]  ="tan(#Theta)";
2014   //
2015   binsTrack[2] =180;
2016   xminTrack[2] =-TMath::Pi(); xmaxTrack[2]=TMath::Pi(); 
2017   axisName[2]  ="phi";
2018   axisTitle[2]  ="#phi";
2019   //
2020   binsTrack[3] =20;
2021   xminTrack[3] =-1.; xmaxTrack[3]=1.;   // 0.33 GeV cut 
2022   axisName[3]  ="snp";
2023   axisTitle[3]  ="snp";
2024   //
2025   binsTrack[4] =10;
2026   xminTrack[4] =120.; xmaxTrack[4]=215.;   // crossing radius for CE only 
2027   axisName[4]  ="r";
2028   axisTitle[4] ="r(cm)";
2029   //
2030   // delta y
2031   xminTrack[0] =-1.5; xmaxTrack[0]=1.5;  // 
2032   fResHistoTPCCE[0] = new THnSparseS("TPCCE#Delta_{Y} (cm)","#Delta_{Y} (cm)",    5, binsTrack,xminTrack, xmaxTrack);
2033   fResHistoTPCITS[0] = new THnSparseS("TPCITS#Delta_{Y} (cm)","#Delta_{Y} (cm)",    4, binsTrack,xminTrack, xmaxTrack);
2034   fResHistoTPCvertex[0]    = new THnSparseS("TPCVertex#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2035   xminTrack[0] =-1.5; xmaxTrack[0]=1.5;  // 
2036   fResHistoTPCTRD[0] = new THnSparseS("TPCTRD#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2037   xminTrack[0] =-5; xmaxTrack[0]=5;  // 
2038   fResHistoTPCTOF[0] = new THnSparseS("TPCTOF#Delta_{Y} (cm)","#Delta_{Y} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2039   //
2040   // delta z
2041   xminTrack[0] =-3.; xmaxTrack[0]=3.;  // 
2042   fResHistoTPCCE[1] = new THnSparseS("TPCCE#Delta_{Z} (cm)","#Delta_{Z} (cm)",    5, binsTrack,xminTrack, xmaxTrack);
2043   fResHistoTPCITS[1] = new THnSparseS("TPCITS#Delta_{Z} (cm)","#Delta_{Z} (cm)",    4, binsTrack,xminTrack, xmaxTrack);
2044   fResHistoTPCvertex[1]    = new THnSparseS("TPCVertex#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2045   fResHistoTPCTRD[1] = new THnSparseS("TPCTRD#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2046   xminTrack[0] =-5.; xmaxTrack[0]=5.;  // 
2047   fResHistoTPCTOF[1] = new THnSparseS("TPCTOF#Delta_{Z} (cm)","#Delta_{Z} (cm)", 4, binsTrack,xminTrack, xmaxTrack);
2048   //
2049   // delta snp-P2
2050   xminTrack[0] =-0.015; xmaxTrack[0]=0.015;  // 
2051   fResHistoTPCCE[2] = new THnSparseS("TPCCE#Delta_{#phi}","#Delta_{#phi}",    5, binsTrack,xminTrack, xmaxTrack);
2052   fResHistoTPCITS[2] = new THnSparseS("TPCITS#Delta_{#phi}","#Delta_{#phi}",    4, binsTrack,xminTrack, xmaxTrack);
2053   fResHistoTPCvertex[2] = new THnSparseS("TPCITSv#Delta_{#phi}","#Delta_{#phi}",    4, binsTrack,xminTrack, xmaxTrack);
2054   fResHistoTPCTRD[2] = new THnSparseS("TPCTRD#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
2055   fResHistoTPCTOF[2] = new THnSparseS("TPCTOF#Delta_{#phi}","#Delta_{#phi}", 4, binsTrack,xminTrack, xmaxTrack);
2056   //
2057   // delta theta-P3
2058   xminTrack[0] =-0.025; xmaxTrack[0]=0.025;  // 
2059   fResHistoTPCCE[3] = new THnSparseS("TPCCE#Delta_{#theta}","#Delta_{#theta}",    5, binsTrack,xminTrack, xmaxTrack);
2060   fResHistoTPCITS[3] = new THnSparseS("TPCITS#Delta_{#theta}","#Delta_{#theta}",    4, binsTrack,xminTrack, xmaxTrack);
2061   fResHistoTPCvertex[3] = new THnSparseS("TPCITSv#Delta_{#theta}","#Delta_{#theta}",    4, binsTrack,xminTrack, xmaxTrack);
2062   fResHistoTPCTRD[3] = new THnSparseS("TPCTRD#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
2063   fResHistoTPCTOF[3] = new THnSparseS("TPCTOF#Delta_{#theta}","#Delta_{#theta}", 4, binsTrack,xminTrack, xmaxTrack);
2064   //
2065   // delta theta-P4
2066   xminTrack[0] =-0.2; xmaxTrack[0]=0.2;  // 
2067   fResHistoTPCCE[4] = new THnSparseS("TPCCE#Delta_{1/pt}","#Delta_{1/pt}",    5, binsTrack,xminTrack, xmaxTrack);
2068   fResHistoTPCITS[4] = new THnSparseS("TPCITS#Delta_{1/pt}","#Delta_{1/pt}",    4, binsTrack,xminTrack, xmaxTrack);
2069   fResHistoTPCvertex[4] = new THnSparseS("TPCITSv#Delta_{1/pt}","#Delta_{1/pt}",    4, binsTrack,xminTrack, xmaxTrack);
2070   fResHistoTPCTRD[4] = new THnSparseS("TPCTRD#Delta_{1/pt}","#Delta_{1/pt}",    4, binsTrack,xminTrack, xmaxTrack);
2071   fResHistoTPCTOF[4] = new THnSparseS("TPCTOF#Delta_{1/pt}","#Delta_{1/pt}",    4, binsTrack,xminTrack, xmaxTrack);
2072   //
2073   for (Int_t ivar=0;ivar<4;ivar++){
2074     for (Int_t ivar2=0;ivar2<5;ivar2++){      
2075       fResHistoTPCCE[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
2076       fResHistoTPCCE[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
2077       if (ivar2<4){
2078         fResHistoTPCITS[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
2079         fResHistoTPCITS[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
2080         fResHistoTPCTRD[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
2081         fResHistoTPCTRD[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
2082         fResHistoTPCvertex[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
2083         fResHistoTPCvertex[ivar]->GetAxis(ivar2)->SetTitle(axisTitle[ivar2].Data());
2084       }
2085     }
2086   }
2087   //
2088   // Book vertex: time histograms
2089   //
2090   Int_t    binsVertex[2]={500, fTimeBins};
2091   Double_t    aminVertex[2]={-5,fTimeStart};
2092   Double_t    amaxVertex[2]={5, fTimeEnd};
2093   const char* hnames[12]={"TPCXAside", "TPCXCside","TPCXACdiff","TPCXAPCdiff",
2094                           "TPCYAside", "TPCYCside","TPCYACdiff","TPCYAPCdiff",
2095                           "TPCZAPCside", "TPCZAMCside","TPCZACdiff","TPCZAPCdiff"}; 
2096   const char* anames[12]={"x (cm) - A side ", "x (cm) - C side","#Delta_{x} (cm) - TPC-A-C","#Delta_{x} (cm) - TPC-Common",
2097                           "y (cm) - A side ", "y (cm) - C side","#Delta_{x} (cm) - TPC-A-C","#Delta_{y} (cm) - TPC-Common",
2098                           "z (cm)", "#Delta_{Z} (cm) A-C side","#Delta_{x} (cm) - TPC-A-C","#Delta_{Z} (cm) TPC-common"}; 
2099   for (Int_t ihis=0; ihis<12; ihis++) {
2100     if (ihis>=8) aminVertex[0]=-20.;
2101     if (ihis>=8) amaxVertex[0]=20.;
2102     fTPCVertex[ihis]=new THnSparseF(hnames[ihis],hnames[ihis],2,binsVertex,aminVertex,amaxVertex);
2103     fTPCVertex[ihis]->GetAxis(1)->SetTitle("Time");
2104     fTPCVertex[ihis]->GetAxis(0)->SetTitle(anames[ihis]);
2105   }
2106   
2107   Int_t    binsVertexC[2]={40, 300};
2108   Double_t aminVertexC[2]={-20,-30};
2109   Double_t amaxVertexC[2]={20,30};
2110   const char* hnamesC[5]={"TPCA_TPC","TPCC_TPC","TPCA_ITS","TPCC_ITS","TPC_ITS"};
2111   for (Int_t ihis=0; ihis<5; ihis++) {
2112     fTPCVertexCorrelation[ihis]=new THnSparseF(hnamesC[ihis],hnamesC[ihis],2,binsVertexC,aminVertexC,amaxVertexC);
2113     fTPCVertexCorrelation[ihis]->GetAxis(1)->SetTitle("z (cm)");
2114     fTPCVertexCorrelation[ihis]->GetAxis(0)->SetTitle("z (cm)");
2115   }
2116 }
2117
2118
2119 void        AliTPCcalibTime::FillResHistoTPCCE(const AliExternalTrackParam * pTPCIn, const AliExternalTrackParam * pTPCOut ){
2120   //
2121   // fill residual histograms pTPCOut-pTPCin - trac crossing CE
2122   // Histogram 
2123   //
2124   if (fMemoryMode<2) return;
2125   Double_t histoX[5];
2126   Double_t xyz[3];
2127   pTPCIn->GetXYZ(xyz);
2128   Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2129   histoX[1]= pTPCIn->GetTgl();
2130   histoX[2]= phi;
2131   histoX[3]= pTPCIn->GetSnp();
2132   histoX[4]= pTPCIn->GetX();
2133   AliExternalTrackParam lout(*pTPCOut);
2134   lout.Rotate(pTPCIn->GetAlpha());
2135   lout.PropagateTo(pTPCIn->GetX(),fMagF);
2136   //
2137   for (Int_t ihisto=0; ihisto<5; ihisto++){
2138     histoX[0]=lout.GetParameter()[ihisto]-pTPCIn->GetParameter()[ihisto];
2139     fResHistoTPCCE[ihisto]->Fill(histoX);
2140   }
2141 }  
2142 void        AliTPCcalibTime::FillResHistoTPCITS(const AliExternalTrackParam * pTPCIn, const AliExternalTrackParam * pITSOut ){
2143   //
2144   // fill residual histograms pTPCIn-pITSOut
2145   // Histogram is filled only for primary tracks
2146   //
2147   Double_t histoX[4];
2148   Double_t xyz[3];
2149   pTPCIn->GetXYZ(xyz);
2150   Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2151   histoX[1]= pTPCIn->GetTgl();
2152   histoX[2]= phi;
2153   histoX[3]= pTPCIn->GetSnp();
2154   AliExternalTrackParam lits(*pITSOut);
2155   lits.Rotate(pTPCIn->GetAlpha());
2156   lits.PropagateTo(pTPCIn->GetX(),fMagF);
2157   //
2158   for (Int_t ihisto=0; ihisto<5; ihisto++){
2159     histoX[0]=pTPCIn->GetParameter()[ihisto]-lits.GetParameter()[ihisto];
2160     fResHistoTPCITS[ihisto]->Fill(histoX);
2161   }
2162 }  
2163
2164      
2165 void        AliTPCcalibTime::FillResHistoTPC(const AliESDtrack * pTrack){
2166   //
2167   // fill residual histograms pTPC - vertex
2168   // Histogram is filled only for primary tracks
2169   //
2170   if (fMemoryMode<2) return;
2171   Double_t histoX[4];
2172   const AliExternalTrackParam * pTPCIn = pTrack->GetInnerParam();
2173   AliExternalTrackParam pTPCvertex(*(pTrack->GetInnerParam()));
2174   //
2175   AliExternalTrackParam lits(*pTrack);
2176   if (TMath::Abs(pTrack->GetY())>3) return;  // beam pipe
2177   pTPCvertex.Rotate(lits.GetAlpha());
2178   //pTPCvertex.PropagateTo(pTPCvertex->GetX(),fMagF);
2179   AliTracker::PropagateTrackToBxByBz(&pTPCvertex,lits.GetX(),0.1,2,kFALSE);
2180   AliTracker::PropagateTrackToBxByBz(&pTPCvertex,lits.GetX(),0.1,0.1,kFALSE);
2181   Double_t xyz[3];
2182   pTPCIn->GetXYZ(xyz);
2183   Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2184   histoX[1]= pTPCIn->GetTgl();
2185   histoX[2]= phi;
2186   histoX[3]= pTPCIn->GetSnp();
2187   //
2188   Float_t dca[2], cov[3];
2189   pTrack->GetImpactParametersTPC(dca,cov);
2190   for (Int_t ihisto=0; ihisto<5; ihisto++){
2191     histoX[0]=pTPCvertex.GetParameter()[ihisto]-lits.GetParameter()[ihisto];
2192     //    if (ihisto<2) histoX[0]=dca[ihisto];
2193     fResHistoTPCvertex[ihisto]->Fill(histoX);
2194   }
2195 }
2196
2197
2198 void        AliTPCcalibTime::FillResHistoTPCTRD(const AliExternalTrackParam * pTPCOut, const AliExternalTrackParam * pTRDIn ){
2199   //
2200   // fill resuidual histogram TPCout-TRDin
2201   //
2202   if (fMemoryMode<2) return;
2203   Double_t histoX[4];
2204   Double_t xyz[3];
2205   pTPCOut->GetXYZ(xyz);
2206   Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2207   histoX[1]= pTPCOut->GetTgl();
2208   histoX[2]= phi;
2209   histoX[3]= pTPCOut->GetSnp();
2210   //
2211   AliExternalTrackParam ltrd(*pTRDIn);
2212   ltrd.Rotate(pTPCOut->GetAlpha());
2213   //  ltrd.PropagateTo(pTPCOut->GetX(),fMagF);
2214   AliTracker::PropagateTrackToBxByBz(&ltrd,pTPCOut->GetX(),0.1,0.1,kFALSE);
2215
2216   for (Int_t ihisto=0; ihisto<5; ihisto++){
2217     histoX[0]=pTPCOut->GetParameter()[ihisto]-ltrd.GetParameter()[ihisto];
2218     fResHistoTPCTRD[ihisto]->Fill(histoX);
2219   }
2220
2221 }
2222
2223 void        AliTPCcalibTime::FillResHistoTPCTOF(const AliExternalTrackParam * pTPCOut, const AliExternalTrackParam * pTOFIn ){
2224   //
2225   // fill resuidual histogram TPCout-TOFin
2226   // track propagated to the TOF position
2227   if (fMemoryMode<2) return;
2228   Double_t histoX[4];
2229   Double_t xyz[3];
2230
2231   AliExternalTrackParam ltpc(*pTPCOut);
2232   ltpc.Rotate(pTOFIn->GetAlpha());
2233   AliTracker::PropagateTrackToBxByBz(&ltpc,pTOFIn->GetX(),0.1,0.1,kFALSE);
2234   //
2235   ltpc.GetXYZ(xyz);
2236   Double_t phi= TMath::ATan2(xyz[1],xyz[0]);
2237   histoX[1]= ltpc.GetTgl();
2238   histoX[2]= phi;
2239   histoX[3]= ltpc.GetSnp();
2240   //
2241   for (Int_t ihisto=0; ihisto<2; ihisto++){
2242     histoX[0]=ltpc.GetParameter()[ihisto]-pTOFIn->GetParameter()[ihisto];
2243     fResHistoTPCTOF[ihisto]->Fill(histoX);
2244   }
2245
2246 }