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