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