]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/Calib/AliTPCcalibCosmic.cxx
419e833b52238a9eb0ad59a6b63295a02d096715
[u/mrichter/AliRoot.git] / TPC / Calib / AliTPCcalibCosmic.cxx
1
2
3 /**************************************************************************
4  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  *                                                                        *
6  * Author: The ALICE Off-line Project.                                    *
7  * Contributors are mentioned in the code where appropriate.              *
8  *                                                                        *
9  * Permission to use, copy, modify and distribute this software and its   *
10  * documentation strictly for non-commercial purposes is hereby granted   *
11  * without fee, provided that the above copyright notice appears in all   *
12  * copies and that both the copyright notice and this permission notice   *
13  * appear in the supporting documentation. The authors make no claims     *
14  * about the suitability of this software for any purpose. It is          *
15  * provided "as is" without express or implied warranty.                  *
16  **************************************************************************/
17
18 /*
19     Comments to be written here: 
20     1. What do we calibrate.
21     2. How to interpret results
22     3. Simple example
23     4. Analysis using debug streamers.
24
25
26
27     3.Simple example
28     // To make cosmic scan the user interaction neccessary
29     //
30      
31   */
32
33
34
35 #include "Riostream.h"
36 #include "TChain.h"
37 #include "TTree.h"
38 #include "TH1F.h"
39 #include "TH2F.h"
40 #include "TList.h"
41 #include "TMath.h"
42 #include "TCanvas.h"
43 #include "TFile.h"
44 #include "TF1.h"
45 #include "THnSparse.h"
46 #include "TDatabasePDG.h"
47
48 #include "AliTPCclusterMI.h"
49 #include "AliTPCseed.h"
50 #include "AliESDVertex.h"
51 #include "AliESDEvent.h"
52 #include "AliESDtrack.h"
53 #include "AliESDfriend.h"
54 #include "AliESDInputHandler.h"
55 #include "AliAnalysisManager.h"
56
57 #include "AliVEvent.h"
58 #include "AliVTrack.h"
59 #include "AliVfriendEvent.h"
60 #include "AliVfriendTrack.h"
61
62 #include "AliTracker.h"
63 #include "AliMagF.h"
64 #include "AliTPCCalROC.h"
65 #include "AliTPCParam.h"
66 #include "AliLog.h"
67
68 #include "AliTPCcalibCosmic.h"
69 #include "TTreeStream.h"
70 #include "AliTPCTracklet.h"
71 //#include "AliESDcosmic.h"
72 #include "AliRecoParam.h"
73 #include "AliMultiplicity.h"
74 #include "AliTPCTransform.h"
75 #include "AliTPCcalibDB.h"
76 #include "AliTPCseed.h"
77 #include "AliGRPObject.h"
78 #include "AliTPCCorrection.h"
79 ClassImp(AliTPCcalibCosmic)
80
81
82 AliTPCcalibCosmic::AliTPCcalibCosmic() 
83   :AliTPCcalibBase(),
84    fHistNTracks(0),
85    fClusters(0),
86    fModules(0),
87    fHistPt(0),
88    fDeDx(0),
89    fDeDxMIP(0),
90    fMIPvalue(1), 
91    fCutMaxD(5),        // maximal distance in rfi ditection
92    fCutMaxDz(40),      // maximal distance in z ditection
93    fCutTheta(0.03),    // maximal distan theta
94    fCutMinDir(-0.99),   // direction vector products
95    fCosmicTree(0)      // tree with cosmic data
96 {  
97   //
98   // CONSTRUCTOR - SEE COMMENTS ABOVE
99   //
100   AliInfo("Default Constructor");    
101   for (Int_t ihis=0; ihis<6;ihis++){
102     fHistoDelta[ihis]=0;
103     fHistoPull[ihis]=0;
104   }
105   for (Int_t ihis=0; ihis<4;ihis++){
106     fHistodEdxMax[ihis]    =0;
107     fHistodEdxTot[ihis]    =0;
108   }
109 }
110
111
112 AliTPCcalibCosmic::AliTPCcalibCosmic(const Text_t *name, const Text_t *title) 
113   :AliTPCcalibBase(),
114    fHistNTracks(0),
115    fClusters(0),
116    fModules(0),
117    fHistPt(0),
118    fDeDx(0),
119    fDeDxMIP(0),
120    fMIPvalue(1),
121    fCutMaxD(5),        // maximal distance in rfi ditection 
122    fCutMaxDz(40),      // maximal distance in z ditection
123    fCutTheta(0.03),    // maximal distan theta
124    fCutMinDir(-0.99),  // direction vector products
125    fCosmicTree(0)      // tree with cosmic data
126 {  
127   //
128   // cONSTRUCTOR - SEE COMENTS ABOVE
129   //
130   SetName(name);
131   SetTitle(title);
132
133   fHistNTracks = new TH1F("ntracks","Number of Tracks per Event; number of tracks per event; number of tracks",501,-0.5,500.5);
134   fClusters = new TH1F("signal","Number of Clusters per track; number of clusters per track n_{cl}; counts",160,0,160);
135   fModules = new TH2F("sector","Acorde hits; z (cm); x(cm)",1200,-650,650,600,-700,700);
136   fHistPt = new TH1F("Pt","Pt distribution; p_{T} (GeV); counts",2000,0,50);
137   fDeDx = new TH2F("DeDx","dEdx; momentum p (GeV); TPC signal (a.u.)",500,0.01,100.,500,2.,1000);
138   BinLogX(fDeDx);
139   fDeDxMIP =  new TH1F("DeDxMIP","MIP region; TPC signal (a.u.);counts ",500,2.,1000);
140   Init();
141   AliInfo("Non Default Constructor");  
142   //
143 }
144
145 AliTPCcalibCosmic::~AliTPCcalibCosmic(){
146   //
147   // destructor
148   //
149   for (Int_t ihis=0; ihis<6;ihis++){
150     delete fHistoDelta[ihis];
151     delete fHistoPull[ihis];
152   }
153   for (Int_t ihis=0; ihis<4;ihis++){
154     delete fHistodEdxTot[ihis];
155     delete fHistodEdxMax[ihis];
156   }
157
158   delete fHistNTracks;            //  histogram showing number of ESD tracks per event
159   delete fClusters;               //  histogram showing the number of clusters per track
160   delete fModules;                //  2d histogram of tracks which are propagated to the ACORDE scintillator array
161   delete fHistPt;                 //  Pt histogram of reconstructed tracks
162   delete fDeDx;                   //  dEdx spectrum showing the different particle types
163   delete fDeDxMIP;                //  TPC signal close to the MIP region of muons 0.4 < p < 0.45 GeV
164 }
165
166
167 void AliTPCcalibCosmic::Init(){
168   //
169   // init component
170   // Make performance histograms
171   //
172
173   // tracking performance bins
174   // 0 - delta of interest
175   // 1 - min (track0, track1) number of clusters
176   // 2 - R  - vertex radius
177   // 3 - P1 - mean z
178   // 4 - P2 - snp(phi)    at inner wall of TPC
179   // 5 - P3 - tan(theta)  at inner wall of TPC
180   // 6 - P4 - 1/pt mean
181   // 7 - pt - pt mean
182   // 8 - alpha
183   // 9 - is corss indicator
184   Int_t ndim=10;
185   Double_t xminTrack[10], xmaxTrack[10];
186   Int_t binsTrack[10];
187   TString axisName[10];
188   //
189   binsTrack[0] =100;
190   axisName[0]  ="#Delta";
191   //
192   binsTrack[1] =8;
193   xminTrack[1] =80; xmaxTrack[1]=160;
194   axisName[1]  ="N_{cl}";
195   //
196   binsTrack[2] =10;
197   xminTrack[2] =0; xmaxTrack[2]=90;  // 
198   axisName[2]  ="dca_{r} (cm)";
199   //
200   binsTrack[3] =25;
201   xminTrack[3] =-250; xmaxTrack[3]=250;  // 
202   axisName[3]  ="z (cm)";
203   //
204   binsTrack[4] =10;
205   xminTrack[4] =-0.8; xmaxTrack[4]=0.8;  // 
206   axisName[4]  ="sin(#phi)";
207   //
208   binsTrack[5] =10;
209   xminTrack[5] =-1; xmaxTrack[5]=1;  // 
210   axisName[5]  ="tan(#theta)";
211   //
212   binsTrack[6] =40;
213   xminTrack[6] =-2; xmaxTrack[6]=2;  // 
214   axisName[6]  ="1/pt (1/GeV)";
215   //
216   binsTrack[7] =50;
217   xminTrack[7] =1; xmaxTrack[7]=1000;  // 
218   axisName[7]  ="pt (GeV)";
219   //
220   binsTrack[8] =18;
221   xminTrack[8] =0; xmaxTrack[8]=TMath::Pi();  // 
222   axisName[8]  ="alpha";
223   //
224   binsTrack[9] =3;
225   xminTrack[9] =-0.1; xmaxTrack[9]=2.1;  // 
226   axisName[9]  ="cross";
227   //
228   // delta y
229   xminTrack[0] =-1; xmaxTrack[0]=1;  // 
230   fHistoDelta[0] = new THnSparseS("#Delta_{Y} (cm)","#Delta_{Y} (cm)", ndim, binsTrack,xminTrack, xmaxTrack);
231   xminTrack[0] =-5; xmaxTrack[0]=5;  // 
232   fHistoPull[0] = new THnSparseS("#Delta_{Y} (unit)","#Delta_{Y} (unit)", ndim, binsTrack,xminTrack, xmaxTrack);
233   //
234   // delta z
235   xminTrack[0] =-1; xmaxTrack[0]=1;  // 
236   fHistoDelta[1] = new THnSparseS("#Delta_{Z} (cm)","#Delta_{Z} (cm)", ndim, binsTrack,xminTrack, xmaxTrack);
237   xminTrack[0] =-5; xmaxTrack[0]=5;  // 
238   fHistoPull[1] = new THnSparseS("#Delta_{Z} (unit)","#Delta_{Z} (unit)", ndim, binsTrack,xminTrack, xmaxTrack);
239   //
240   // delta P2
241   xminTrack[0] =-10; xmaxTrack[0]=10;  // 
242   fHistoDelta[2] = new THnSparseS("#Delta_{#phi} (mrad)","#Delta_{#phi} (mrad)", ndim, binsTrack,xminTrack, xmaxTrack);
243   xminTrack[0] =-5; xmaxTrack[0]=5;  // 
244   fHistoPull[2] = new THnSparseS("#Delta_{#phi} (unit)","#Delta_{#phi} (unit)", ndim, binsTrack,xminTrack, xmaxTrack);
245   //
246   // delta P3
247   xminTrack[0] =-10; xmaxTrack[0]=10;  // 
248   fHistoDelta[3] = new THnSparseS("#Delta_{#theta} (mrad)","#Delta_{#theta} (mrad)", ndim, binsTrack,xminTrack, xmaxTrack);
249   xminTrack[0] =-5; xmaxTrack[0]=5;  // 
250   fHistoPull[3] = new THnSparseS("#Delta_{#theta} (unit)","#Delta_{#theta} (unit)", ndim, binsTrack,xminTrack, xmaxTrack);
251   //
252   // delta P4
253   xminTrack[0] =-0.2; xmaxTrack[0]=0.2;  // 
254   fHistoDelta[4] = new THnSparseS("#Delta_{1/pt} (1/GeV)","#Delta_{1/pt} (1/GeV)", ndim, binsTrack,xminTrack, xmaxTrack);
255   xminTrack[0] =-5; xmaxTrack[0]=5;  // 
256   fHistoPull[4] = new THnSparseS("#Delta_{1/pt} (unit)","#Delta_{1/pt} (unit)", ndim, binsTrack,xminTrack, xmaxTrack);
257   
258   //
259   // delta Pt
260   xminTrack[0] =-0.5; xmaxTrack[0]=0.5;  // 
261   fHistoDelta[5] = new THnSparseS("#Delta_{pt}/p_{t}","#Delta_{pt}/p_{t}", ndim, binsTrack,xminTrack, xmaxTrack);
262   xminTrack[0] =-5; xmaxTrack[0]=5;  // 
263   fHistoPull[5] = new THnSparseS("#Delta_{pt}/p_{t} (unit)","#Delta_{pt}/p_{t} (unit)", ndim, binsTrack,xminTrack, xmaxTrack);
264   //
265
266   for (Int_t idedx=0;idedx<4;idedx++){
267     xminTrack[0] =0.5; xmaxTrack[0]=1.5;  // 
268     binsTrack[1] =40;
269     xminTrack[1] =10; xmaxTrack[1]=160;
270
271     fHistodEdxMax[idedx] = new THnSparseS(Form("dEdx_{MaxUp}/dEdx_{MaxDown}_Pad%d",idedx),
272                                           Form("dEdx_{MaxUp}/dEdx_{MaxDown}_Pad%d",idedx), 
273                                           ndim, binsTrack,xminTrack, xmaxTrack);
274     fHistodEdxTot[idedx] = new THnSparseS(Form("dEdx_{TotUp}/dEdx_{TotDown}_Pad%d",idedx),
275                                           Form("dEdx_{TotUp}/dEdx_{TotDown}_Pad%d",idedx), 
276                                           ndim, binsTrack,xminTrack, xmaxTrack);
277   }
278   
279
280
281   for (Int_t ivar=0;ivar<6;ivar++){
282     for (Int_t ivar2=0;ivar2<ndim;ivar2++){      
283       fHistoDelta[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
284       fHistoDelta[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
285       fHistoPull[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
286       fHistoPull[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
287       BinLogX(fHistoDelta[ivar],7);
288       BinLogX(fHistoPull[ivar],7);
289       if (ivar<4){
290         fHistodEdxMax[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
291         fHistodEdxMax[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
292         fHistodEdxTot[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
293         fHistodEdxTot[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
294         BinLogX(fHistodEdxMax[ivar],7);
295         BinLogX(fHistodEdxTot[ivar],7);
296       }
297     }
298   }
299 }
300
301
302 void AliTPCcalibCosmic::Add(const AliTPCcalibCosmic* cosmic){
303   //
304   // merge the content of the cosmic componentnts
305   //
306   for (Int_t ivar=0; ivar<6;ivar++){
307     if (fHistoDelta[ivar] && cosmic->fHistoDelta[ivar]){
308       fHistoDelta[ivar]->Add(cosmic->fHistoDelta[ivar]);
309     }
310     if (fHistoPull[ivar] && cosmic->fHistoPull[ivar]){
311       fHistoPull[ivar]->Add(cosmic->fHistoPull[ivar]);
312     }
313   }
314   for (Int_t ivar=0; ivar<4;ivar++){
315     if (fHistodEdxMax[ivar] && cosmic->fHistodEdxMax[ivar]){
316       fHistodEdxMax[ivar]->Add(cosmic->fHistodEdxMax[ivar]);
317     }
318     if (fHistodEdxTot[ivar] && cosmic->fHistodEdxTot[ivar]){
319       fHistodEdxTot[ivar]->Add(cosmic->fHistodEdxTot[ivar]);
320     }
321   }
322   if (cosmic->fCosmicTree){
323     if (!fCosmicTree) {
324       fCosmicTree = new TTree("pairs","pairs");
325       fCosmicTree->SetDirectory(0);
326     }
327     AliTPCcalibCosmic::AddTree(fCosmicTree,cosmic->fCosmicTree);
328   }
329 }
330
331
332
333
334 void AliTPCcalibCosmic::Process(AliVEvent *event) {
335   //
336   // Process of the ESD event  - fill calibration components
337   //
338   if (!event) {
339     Printf("ERROR: event not available");
340     return;
341   }  
342    
343   //
344   //Int_t isOK=kTRUE;
345   // COSMIC not signed properly
346   //  UInt_t specie = event->GetEventSpecie();  // select only cosmic events
347   //if (specie==AliRecoParam::kCosmic || specie==AliRecoParam::kCalib) {
348   //  isOK = kTRUE;
349   //}
350   //if (!isOK) return;
351   // Work around
352   FindCosmicPairs(event);
353   //const AliMultiplicity *multiplicity = event->GetMultiplicity();
354   //  Int_t ntracklets = multiplicity->GetNumberOfTracklets();
355   //if (ntracklets>6) return; // filter out "normal" event with high multiplicity
356   //const TString &trigger = event->GetFiredTriggerClasses();
357   //if (trigger.Contains("C0OB0")==0) return;
358
359
360   FindPairs(event); // nearly everything takes place in find pairs...
361
362   if (GetDebugLevel()>20) printf("Hallo world: Im here and processing an event\n");
363   Int_t ntracks=event->GetNumberOfTracks(); 
364   fHistNTracks->Fill(ntracks);
365   
366 }
367
368
369 void AliTPCcalibCosmic::FillHistoPerformance(const AliExternalTrackParam *par0, const AliExternalTrackParam *par1, const AliExternalTrackParam *inner0, const AliExternalTrackParam */*inner1*/, AliTPCseed *seed0,  AliTPCseed *seed1, const AliExternalTrackParam *param0Combined , Int_t cross){
370   //
371   // par0,par1       - parameter of tracks at DCA 0
372   // inner0,inner1   - parameter of tracks at the TPC entrance
373   // seed0, seed1    - detailed track information
374   // param0Combined  - Use combined track parameters for binning
375   //
376   Int_t kMinCldEdx =20;
377   Int_t ncl0 = seed0->GetNumberOfClusters();
378   Int_t ncl1 = seed1->GetNumberOfClusters();
379   const Double_t kpullCut    = 10;
380   Double_t x[10];
381   Double_t xyz0[3];
382   Double_t xyz1[3];
383   par0->GetXYZ(xyz0);
384   par1->GetXYZ(xyz1);
385   Double_t radius0 = TMath::Sqrt(xyz0[0]*xyz0[0]+xyz0[1]*xyz0[1]);
386   Double_t radius1 = TMath::Sqrt(xyz1[0]*xyz1[0]+xyz1[1]*xyz1[1]);
387   inner0->GetXYZ(xyz0);
388   Double_t alpha = TMath::ATan2(xyz0[1],xyz0[0]);
389   // bin parameters
390   x[1] = TMath::Min(ncl0,ncl1);
391   x[2] = (radius0+radius1)*0.5;
392   x[3] = param0Combined->GetZ();
393   x[4] = inner0->GetSnp();
394   x[5] = param0Combined->GetTgl();
395   x[6] = param0Combined->GetSigned1Pt();
396   x[7] = param0Combined->Pt();
397   x[8] = alpha;
398   x[9] = cross;
399   // deltas
400   Double_t delta[6];
401   Double_t sigma[6];
402   delta[0] = (par0->GetY()+par1->GetY());
403   delta[1] = (par0->GetZ()-par1->GetZ());
404   delta[2] = (par0->GetAlpha()-par1->GetAlpha()-TMath::Pi());
405   delta[3] = (par0->GetTgl()+par1->GetTgl());
406   delta[4] = (par0->GetParameter()[4]+par1->GetParameter()[4]);
407   delta[5] = (par0->Pt()-par1->Pt())/((par0->Pt()+par1->Pt())*0.5);
408   //
409   sigma[0] = TMath::Sqrt(par0->GetSigmaY2()+par1->GetSigmaY2());
410   sigma[1] = TMath::Sqrt(par0->GetSigmaZ2()+par1->GetSigmaZ2());
411   sigma[2] = TMath::Sqrt(par0->GetSigmaSnp2()+par1->GetSigmaSnp2());
412   sigma[3] = TMath::Sqrt(par0->GetSigmaTgl2()+par1->GetSigmaTgl2());
413   sigma[4] = TMath::Sqrt(par0->GetSigma1Pt2()+par1->GetSigma1Pt2());
414   sigma[5] = sigma[4]*((par0->Pt()+par1->Pt())*0.5);
415   //
416   Bool_t isOK = kTRUE;
417   for (Int_t ivar=0;ivar<6;ivar++){
418     if (sigma[ivar]==0) isOK=kFALSE;
419     x[0]= delta[ivar]/sigma[ivar];
420     if (TMath::Abs(x[0])>kpullCut) isOK = kFALSE;
421   }
422   //
423
424   if (isOK) for (Int_t ivar=0;ivar<6;ivar++){
425     x[0]= delta[ivar];    // Modifiation 10.10 use not normalized deltas
426     if (ivar==2 || ivar ==3) x[0]*=1000;  // angles in radian
427     fHistoDelta[ivar]->Fill(x);
428     if (sigma[ivar]>0){
429       x[0]= delta[ivar]/sigma[ivar];
430       fHistoPull[ivar]->Fill(x);
431     }
432   }
433
434   //                                            
435   // Fill dedx performance
436   //
437   for (Int_t ipad=0; ipad<4;ipad++){
438     //
439     //
440     //
441     Int_t row0=0;
442     Int_t row1=160;
443     if (ipad==0) row1=63;
444     if (ipad==1) {row0=63; row1=63+64;}
445     if (ipad==2) {row0=128;}
446     Int_t   nclUp       = TMath::Nint(seed0->CookdEdxAnalytical(0.01,0.7,0,row0,row1,2));
447     Int_t   nclDown     = TMath::Nint(seed1->CookdEdxAnalytical(0.01,0.7,0,row0,row1,2));
448     Int_t   minCl       = TMath::Min(nclUp,nclDown);
449     if (minCl<kMinCldEdx) continue;
450     x[1] = minCl;
451     //
452     Float_t dEdxTotUp   = seed0->CookdEdxAnalytical(0.01,0.7,0,row0,row1);
453     Float_t dEdxTotDown = seed1->CookdEdxAnalytical(0.01,0.7,0,row0,row1);
454     Float_t dEdxMaxUp   = seed0->CookdEdxAnalytical(0.01,0.7,1,row0,row1);
455     Float_t dEdxMaxDown = seed1->CookdEdxAnalytical(0.01,0.7,1,row0,row1);
456     //
457     if (dEdxTotDown<=0) continue;
458     if (dEdxMaxDown<=0) continue;
459     x[0]=dEdxTotUp/dEdxTotDown;
460     fHistodEdxTot[ipad]->Fill(x);
461     x[0]=dEdxMaxUp/dEdxMaxDown;
462     fHistodEdxMax[ipad]->Fill(x);
463   }
464
465
466   
467 }
468
469 void AliTPCcalibCosmic::FindPairs(const AliVEvent *event){
470   //
471   // Find cosmic pairs
472   // 
473   // Track0 is choosen in upper TPC part
474   // Track1 is choosen in lower TPC part
475   //
476   if (GetDebugLevel()>20) printf("Hallo world: Im here\n");
477   AliVfriendEvent *friendEvent=event->FindFriend();
478   Int_t ntracks=event->GetNumberOfTracks(); 
479   TObjArray  tpcSeeds(ntracks);
480   if (ntracks==0) return;
481   Double_t vtxx[3]={0,0,0};
482   Double_t svtxx[3]={0.000001,0.000001,100.};
483   AliESDVertex vtx(vtxx,svtxx);
484   //
485   //track loop
486   //
487   for (Int_t i=0;i<ntracks;++i) {
488    AliVTrack *track = event->GetVTrack(i);
489    if(!track) continue;
490    fClusters->Fill(track->GetTPCNcls()); 
491   
492    const AliExternalTrackParam * trackIn = track->GetInnerParam();
493    const AliExternalTrackParam * trackOut = track->GetOuterParam();
494    if (!trackIn) continue;
495    if (!trackOut) continue;
496    if (ntracks>4 && TMath::Abs(trackIn->GetTgl())<0.0015) continue;  // filter laser 
497
498
499    const AliVfriendTrack *friendTrack = friendEvent->GetTrack(i);
500    if (!friendTrack) continue;
501    TObject *calibObject;
502    AliTPCseed *seed = 0;   
503    for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
504      if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
505    }
506    if (seed) tpcSeeds.AddAt(seed,i);
507
508    Double_t meanP = 0.5*(trackIn->GetP() + trackOut->GetP());
509    if (seed && track->GetTPCNcls() > 80 + 60/(1+TMath::Exp(-meanP+5))) {
510      fDeDx->Fill(meanP, seed->CookdEdxNorm(0.0,0.45,0,0,159));
511      //
512      if (meanP > 0.4 && meanP < 0.45) fDeDxMIP->Fill(seed->CookdEdxNorm(0.0,0.45,0,0,159));
513      //
514      // if (GetDebugLevel()>0&&meanP>0.2&&seed->CookdEdxNorm(0.0,0.45,0,0,159)>300) {
515 //        //TFile *curfile = AliAnalysisManager::GetAnalysisManager()->GetTree()->GetCurrentFile();
516 //        //if (curfile) printf(">>> p+ in file: %s \t event: %i \t Number of ESD tracks: %i \n", curfile->GetName(), (int)event->GetEventNumberInFile(), (int)ntracks);
517 //        // if (track->GetOuterParam()->GetAlpha()<0) cout << " Polartiy: " << track->GetSign() << endl;
518 //      }
519
520    }
521
522   }
523
524   if (ntracks<2) return;
525   //
526   // Find pairs
527   //
528   for (Int_t i=0;i<ntracks;++i) {
529       AliVTrack *track0 = event->GetVTrack(i);
530     // track0 - choosen upper part
531     if (!track0) continue;
532     if (!track0->GetOuterParam()) continue;
533     if (track0->GetOuterParam()->GetAlpha()<0) continue;
534     Double_t dir0[3];
535     track0->GetDirection(dir0);    
536     for (Int_t j=0;j<ntracks;++j) {
537       if (i==j) continue;
538       AliVTrack *track1 = event->GetVTrack(j);
539       //track 1 lower part
540       if (!track1) continue;
541       if (!track1->GetOuterParam()) continue;
542       if (track1->GetOuterParam()->GetAlpha()>0) continue;
543       //
544       Double_t dir1[3];
545       track1->GetDirection(dir1);
546       
547       AliTPCseed * seed0 = (AliTPCseed*) tpcSeeds.At(i);
548       AliTPCseed * seed1 = (AliTPCseed*) tpcSeeds.At(j);
549       if (! seed0) continue;
550       if (! seed1) continue;
551       Float_t dedx0 = seed0->CookdEdxNorm(0.05,0.55,0,0,159);
552       Float_t dedx1 = seed1->CookdEdxNorm(0.05,0.55,0,0,159);
553       //
554       Float_t dedx0I = seed0->CookdEdxNorm(0.05,0.55,0,0,63);
555       Float_t dedx1I = seed1->CookdEdxNorm(0.05,0.55,0,0,63);
556       //
557       Float_t dedx0O = seed0->CookdEdxNorm(0.05,0.55,0,64,159);
558       Float_t dedx1O = seed1->CookdEdxNorm(0.05,0.55,0,64,159);
559       //
560       Float_t dir = (dir0[0]*dir1[0] + dir0[1]*dir1[1] + dir0[2]*dir1[2]);
561       Float_t d0  = track0->GetLinearD(0,0);
562       Float_t d1  = track1->GetLinearD(0,0);
563       //
564       // conservative cuts - convergence to be guarantied
565       // applying before track propagation
566       if (TMath::Abs(d0+d1)>fCutMaxD) continue;   // distance to the 0,0
567       if (dir>fCutMinDir) continue;               // direction vector product
568       Float_t bz = AliTracker::GetBz();
569       Float_t dvertex0[2];   //distance to 0,0
570       Float_t dvertex1[2];   //distance to 0,0 
571       track0->GetDZ(0,0,0,bz,dvertex0);
572       track1->GetDZ(0,0,0,bz,dvertex1);
573       if (TMath::Abs(dvertex0[1])>250) continue;
574       if (TMath::Abs(dvertex1[1])>250) continue;
575       //
576       //
577       //
578       Float_t dmax = TMath::Max(TMath::Abs(d0),TMath::Abs(d1));
579       AliExternalTrackParam param0;
580       param0.CopyFromVTrack(track0);
581
582       AliExternalTrackParam param1;
583       param1.CopyFromVTrack(track1);
584       //
585       // Propagate using Magnetic field and correct fo material budget
586       //
587       Double_t sign0=-1;
588       Double_t sign1=1;
589       Double_t maxsnp=0.90;
590       AliTracker::PropagateTrackToBxByBz(&param0,dmax+1,TDatabasePDG::Instance()->GetParticle("e-")->Mass(),3,kTRUE,maxsnp,sign0);
591       AliTracker::PropagateTrackToBxByBz(&param1,dmax+1,TDatabasePDG::Instance()->GetParticle("e-")->Mass(),3,kTRUE,maxsnp,sign1);
592       //
593       // Propagate rest to the 0,0 DCA - z should be ignored
594       //
595       Bool_t b0 = param0.PropagateToDCA(&vtx,bz,1000);
596       Bool_t b1 = param1.PropagateToDCA(&vtx,bz,1000);
597       //      
598       param0.GetDZ(0,0,0,bz,dvertex0);
599       param1.GetDZ(0,0,0,bz,dvertex1);
600       if (TMath::Abs(param0.GetZ()-param1.GetZ())>fCutMaxDz) continue;
601       //
602       Double_t xyz0[3];//,pxyz0[3];
603       Double_t xyz1[3];//,pxyz1[3];
604       param0.GetXYZ(xyz0);
605       param1.GetXYZ(xyz1);
606       Bool_t isPair = IsPair(&param0,&param1);
607       //
608       if (isPair) FillAcordeHist(track0);
609       if (isPair &&param0.Pt()>1) {
610         const TString &trigger = event->GetFiredTriggerClasses(); 
611         UInt_t specie = event->GetEventSpecie();
612         printf("COSMIC ?\t%s\t%d\t%f\t%f\n", trigger.Data(),specie, param0.GetZ(), param1.GetZ()); 
613       }
614       //
615       // combined track params 
616       //
617       AliExternalTrackParam *par0U=MakeCombinedTrack(&param0,&param1);
618       AliExternalTrackParam *par1U=MakeCombinedTrack(&param1,&param0);
619       
620       //
621       if (fStreamLevel>0){
622         TTreeSRedirector * cstream =  GetDebugStreamer();
623         //printf("My stream=%p\n",(void*)cstream);
624         AliExternalTrackParam *ip0 = (AliExternalTrackParam *)track0->GetInnerParam();
625         AliExternalTrackParam *ip1 = (AliExternalTrackParam *)track1->GetInnerParam();
626         AliExternalTrackParam *op0 = (AliExternalTrackParam *)track0->GetOuterParam();
627         AliExternalTrackParam *op1 = (AliExternalTrackParam *)track1->GetOuterParam();
628         Bool_t isCrossI = ip0->GetZ()*ip1->GetZ()<0;
629         Bool_t isCrossO = op0->GetZ()*op1->GetZ()<0;
630         Double_t alpha0 = TMath::ATan2(dir0[1],dir0[0]);
631         Double_t alpha1 = TMath::ATan2(dir1[1],dir1[0]);
632         //
633         //
634         //
635         Int_t cross =0;  // 0 no cross, 2 cross on both sides
636         if (isCrossI) cross+=1; 
637         if (isCrossO) cross+=1; 
638         FillHistoPerformance(&param0, &param1, ip0, ip1, seed0, seed1,par0U, cross);
639         if (cstream) {
640           (*cstream) << "Track0" <<
641             "run="<<fRun<<              //  run number
642             "event="<<fEvent<<          //  event number
643             "time="<<fTime<<            //  time stamp of event
644             "trigger="<<fTrigger<<      //  trigger
645             "triggerClass="<<&fTriggerClass<<      //  trigger
646             "mag="<<fMagF<<             //  magnetic field
647             "dir="<<dir<<               //  direction
648             "OK="<<isPair<<             //  will be accepted
649             "b0="<<b0<<                 //  propagate status
650             "b1="<<b1<<                 //  propagate status
651             "crossI="<<isCrossI<<       //  cross inner
652             "crossO="<<isCrossO<<       //  cross outer
653             //
654             "Orig0.=" << track0 <<      //  original track  0
655             "Orig1.=" << track1 <<      //  original track  1
656             "Tr0.="<<&param0<<          //  track propagated to the DCA 0,0
657             "Tr1.="<<&param1<<          //  track propagated to the DCA 0,0        
658             "Ip0.="<<ip0<<              //  inner param - upper
659             "Ip1.="<<ip1<<              //  inner param - lower
660             "Op0.="<<op0<<              //  outer param - upper
661             "Op1.="<<op1<<              //  outer param - lower
662             "Up0.="<<par0U<<           //  combined track 0
663             "Up1.="<<par1U<<           //  combined track 1
664             //
665             "v00="<<dvertex0[0]<<       //  distance using kalman
666             "v01="<<dvertex0[1]<<       // 
667             "v10="<<dvertex1[0]<<       //
668             "v11="<<dvertex1[1]<<       // 
669             "d0="<<d0<<                 //  linear distance to 0,0
670             "d1="<<d1<<                 //  linear distance to 0,0
671             //
672             //
673             //
674             "x00="<<xyz0[0]<<           // global position close to vertex
675             "x01="<<xyz0[1]<<
676             "x02="<<xyz0[2]<<
677             //
678             "x10="<<xyz1[0]<<           // global position close to vertex
679             "x11="<<xyz1[1]<<
680             "x12="<<xyz1[2]<<
681             //
682             "alpha0="<<alpha0<<
683             "alpha1="<<alpha1<<
684             "dir00="<<dir0[0]<<           // direction upper
685             "dir01="<<dir0[1]<<
686             "dir02="<<dir0[2]<<
687             //
688             "dir10="<<dir1[0]<<           // direction lower
689             "dir11="<<dir1[1]<<
690             "dir12="<<dir1[2]<<
691             //
692             //
693             "Seed0.=" << seed0 <<       //  original seed 0
694             "Seed1.=" << seed1 <<       //  original seed 1
695             //
696             "dedx0="<<dedx0<<           //  dedx0 - all
697             "dedx1="<<dedx1<<           //  dedx1 - all
698             //
699             "dedx0I="<<dedx0I<<         //  dedx0 - inner ROC
700             "dedx1I="<<dedx1I<<         //  dedx1 - inner ROC
701             //
702             "dedx0O="<<dedx0O<<         //  dedx0 - outer ROC
703             "dedx1O="<<dedx1O<<         //  dedx1 - outer ROC
704             "\n";
705         }
706       }
707       delete par0U;
708       delete par1U;
709     }
710   }  
711 }    
712
713
714
715
716 void  AliTPCcalibCosmic::FillAcordeHist(AliVTrack *upperTrack) {
717
718   // Pt cut to select straight tracks which can be easily propagated to ACORDE which is outside the magnetic field
719   if (upperTrack->Pt() < 10 || upperTrack->GetTPCNcls() < 80) return;
720     
721   const Double_t acordePlane = 850.; // distance of the central Acorde detectors to the beam line at y =0
722   const Double_t roof = 210.5;       // distance from x =0 to end of magnet roof
723
724   Double_t r[3];
725   upperTrack->GetXYZ(r);
726   Double_t d[3];
727   upperTrack->GetDirection(d);
728   Double_t x,z;
729   z = r[2] + (d[2]/d[1])*(acordePlane - r[1]);
730   x = r[0] + (d[0]/d[1])*(acordePlane - r[1]);
731   
732   if (x > roof) {
733     x = r[0] + (d[0]/(d[0]+d[1]))*(acordePlane+roof-r[0]-r[1]);
734     z = r[2] + (d[2]/(d[0]+d[1]))*(acordePlane+roof-r[0]-r[1]);
735   }
736   if (x < -roof) {
737     x = r[0] + (d[0]/(d[1]-d[0]))*(acordePlane+roof+r[0]-r[1]);       
738     z = r[2] + (d[2]/(d[1]-d[0]))*(acordePlane+roof+r[0]-r[1]);
739   } 
740
741   fModules->Fill(z, x);
742  
743 }
744
745
746
747 Long64_t AliTPCcalibCosmic::Merge(TCollection *const li) {
748   //
749   // component merging
750   //
751
752   TIterator* iter = li->MakeIterator();
753   AliTPCcalibCosmic* cal = 0;
754
755   while ((cal = (AliTPCcalibCosmic*)iter->Next())) {
756     if (!cal->InheritsFrom(AliTPCcalibCosmic::Class())) {
757       //Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName());
758       return -1;
759     }
760     
761     fHistNTracks->Add(cal->GetHistNTracks());
762     fClusters->Add(cal-> GetHistClusters());
763     fModules->Add(cal->GetHistAcorde());
764     fHistPt->Add(cal->GetHistPt());
765     fDeDx->Add(cal->GetHistDeDx());
766     fDeDxMIP->Add(cal->GetHistMIP());
767     Add(cal);
768   }
769   return 0;
770   
771 }
772
773
774 Bool_t  AliTPCcalibCosmic::IsPair(AliExternalTrackParam *tr0, AliExternalTrackParam *tr1) const{
775   //
776   //
777   /*
778   // 0. Same direction - OPOSITE  - cutDir +cutT    
779   TCut cutDir("cutDir","dir<-0.99")
780   // 1. 
781   TCut cutT("cutT","abs(Tr1.fP[3]+Tr0.fP[3])<0.03")
782   //
783   // 2. The same rphi 
784   TCut cutD("cutD","abs(Tr0.fP[0]+Tr1.fP[0])<5")
785   //
786   //
787   //
788   TCut cutPt("cutPt","abs(Tr1.fP[4]+Tr0.fP[4])<1&&abs(Tr0.fP[4])+abs(Tr1.fP[4])<10");  
789   // 1/Pt diff cut
790   */
791   const Double_t *p0 = tr0->GetParameter();
792   const Double_t *p1 = tr1->GetParameter();
793   if (TMath::Abs(p0[3]+p1[3])>fCutTheta) return kFALSE;
794   if (TMath::Abs(p0[1]-p1[1])>fCutMaxDz) return kFALSE;
795   if (TMath::Abs(p0[0]+p1[0])>fCutMaxD)  return kFALSE;
796   
797   Double_t d0[3], d1[3];
798   tr0->GetDirection(d0);    
799   tr1->GetDirection(d1);       
800   if (d0[0]*d1[0] + d0[1]*d1[1] + d0[2]*d1[2] >fCutMinDir) return kFALSE;
801   //
802   return kTRUE;  
803 }
804  
805
806
807 Double_t AliTPCcalibCosmic::CalculateMIPvalue(TH1F * hist) {
808   //
809   // Calculate the MIP value - gaussian fit used
810   //
811
812   TF1 * funcDoubleGaus = new TF1("funcDoubleGaus", "gaus(0)+gaus(3)",0,1000);
813   funcDoubleGaus->SetParameters(hist->GetEntries()*0.75,hist->GetMean()/1.3,hist->GetMean()*0.10,
814                                 hist->GetEntries()*0.25,hist->GetMean()*1.3,hist->GetMean()*0.10);
815   hist->Fit(funcDoubleGaus);
816   Double_t aMIPvalue = TMath::Min(funcDoubleGaus->GetParameter(1),funcDoubleGaus->GetParameter(4));
817
818   delete funcDoubleGaus;
819
820   return aMIPvalue;
821
822 }
823
824
825
826
827 void AliTPCcalibCosmic::CalculateBetheParams(TH2F */*hist*/, Double_t * /*initialParam*/) {
828   //
829   // Not implemented yet
830   //
831   return;
832
833 }
834
835
836 void AliTPCcalibCosmic::BinLogX(THnSparse *const h, Int_t axisDim) {
837
838   // Method for the correct logarithmic binning of histograms
839
840   TAxis *axis = h->GetAxis(axisDim);
841   int bins = axis->GetNbins();
842
843   Double_t from = axis->GetXmin();
844   Double_t to = axis->GetXmax();
845   Double_t *newBins = new Double_t[bins + 1];
846
847   newBins[0] = from;
848   Double_t factor = pow(to/from, 1./bins);
849
850   for (int i = 1; i <= bins; i++) {
851    newBins[i] = factor * newBins[i-1];
852   }
853   axis->Set(bins, newBins);
854   delete [] newBins;
855
856 }
857
858
859 void AliTPCcalibCosmic::BinLogX(TH1 *const h) {
860
861   // Method for the correct logarithmic binning of histograms
862
863   TAxis *axis = h->GetXaxis();
864   int bins = axis->GetNbins();
865
866   Double_t from = axis->GetXmin();
867   Double_t to = axis->GetXmax();
868   Double_t *newBins = new Double_t[bins + 1];
869    
870   newBins[0] = from;
871   Double_t factor = pow(to/from, 1./bins);
872   
873   for (int i = 1; i <= bins; i++) {
874    newBins[i] = factor * newBins[i-1];
875   }
876   axis->Set(bins, newBins);
877   delete [] newBins;
878   
879 }
880
881
882 AliExternalTrackParam *AliTPCcalibCosmic::MakeTrack(const AliExternalTrackParam *track0, const AliExternalTrackParam *track1){
883   //
884   // Make a atrack using the kalman update of track0 and track1
885   //
886   AliExternalTrackParam *par1R= new AliExternalTrackParam(*track1);
887   par1R->Rotate(track0->GetAlpha());
888   par1R->PropagateTo(track0->GetX(),AliTracker::GetBz()); 
889   //
890   //
891   Double_t * param = (Double_t*)par1R->GetParameter();
892   Double_t * covar = (Double_t*)par1R->GetCovariance();
893
894   param[0]*=1;  //OK
895   param[1]*=1;  //OK
896   param[2]*=1;  //?
897   param[3]*=-1; //OK
898   param[4]*=-1; //OK
899   //
900   covar[6] *=-1.; covar[7] *=-1.; covar[8] *=-1.;
901   //covar[10]*=-1.; covar[11]*=-1.; covar[12]*=-1.;
902   covar[13]*=-1.;
903   return par1R;
904 }
905
906 AliExternalTrackParam *AliTPCcalibCosmic::MakeCombinedTrack(const AliExternalTrackParam *track0, const AliExternalTrackParam *track1){
907   //
908   // Make combined track
909   //
910   //
911   AliExternalTrackParam * par1T = MakeTrack(track0,track1);
912   AliExternalTrackParam * par0U = new AliExternalTrackParam(*track0);
913   //
914   UpdateTrack(*par0U,*par1T);
915   delete par1T;
916   return par0U;
917 }
918
919
920 void AliTPCcalibCosmic::UpdateTrack(AliExternalTrackParam &track1, const AliExternalTrackParam &track2){
921   //
922   // Update track 1 with track 2
923   //
924   //
925   //
926   TMatrixD vecXk(5,1);    // X vector
927   TMatrixD covXk(5,5);    // X covariance 
928   TMatrixD matHk(5,5);    // vector to mesurement
929   TMatrixD measR(5,5);    // measurement error 
930   TMatrixD vecZk(5,1);    // measurement
931   //
932   TMatrixD vecYk(5,1);    // Innovation or measurement residual
933   TMatrixD matHkT(5,5);
934   TMatrixD matSk(5,5);    // Innovation (or residual) covariance
935   TMatrixD matKk(5,5);    // Optimal Kalman gain
936   TMatrixD mat1(5,5);     // update covariance matrix
937   TMatrixD covXk2(5,5);   // 
938   TMatrixD covOut(5,5);
939   //
940   Double_t *param1=(Double_t*) track1.GetParameter();
941   Double_t *covar1=(Double_t*) track1.GetCovariance();
942   Double_t *param2=(Double_t*) track2.GetParameter();
943   Double_t *covar2=(Double_t*) track2.GetCovariance();
944   //
945   // copy data to the matrix
946   for (Int_t ipar=0; ipar<5; ipar++){
947     for (Int_t jpar=0; jpar<5; jpar++){
948       covXk(ipar,jpar) = covar1[track1.GetIndex(ipar, jpar)];
949       measR(ipar,jpar) = covar2[track2.GetIndex(ipar, jpar)];
950       matHk(ipar,jpar)=0;
951       mat1(ipar,jpar)=0;
952     }
953     vecXk(ipar,0) = param1[ipar];
954     vecZk(ipar,0) = param2[ipar];
955     matHk(ipar,ipar)=1;
956     mat1(ipar,ipar)=0;
957   }
958   //
959   //
960   //
961   //
962   //
963   vecYk = vecZk-matHk*vecXk;                 // Innovation or measurement residual
964   matHkT=matHk.T(); matHk.T();
965   matSk = (matHk*(covXk*matHkT))+measR;      // Innovation (or residual) covariance
966   matSk.Invert();
967   matKk = (covXk*matHkT)*matSk;              //  Optimal Kalman gain
968   vecXk += matKk*vecYk;                      //  updated vector 
969   covXk2 = (mat1-(matKk*matHk));
970   covOut =  covXk2*covXk; 
971   //
972   //
973   //
974   // copy from matrix to parameters
975   if (0) {
976     vecXk.Print();
977     vecZk.Print();
978     //
979     measR.Print();
980     covXk.Print();
981     covOut.Print();
982     //
983     track1.Print();
984     track2.Print();
985   }
986
987   for (Int_t ipar=0; ipar<5; ipar++){
988     param1[ipar]= vecXk(ipar,0) ;
989     for (Int_t jpar=0; jpar<5; jpar++){
990       covar1[track1.GetIndex(ipar, jpar)]=covOut(ipar,jpar);
991     }
992   }
993 }
994
995
996
997 void AliTPCcalibCosmic::FindCosmicPairs(const AliVEvent *event) {
998   //
999   // find cosmic pairs trigger by random trigger
1000   //
1001   //
1002   AliESDVertex vtxSPD;
1003   event->GetPrimaryVertexSPD(vtxSPD);
1004   AliESDVertex *vertexSPD=&vtxSPD;
1005
1006   AliESDVertex vtxTPC;
1007   event->GetPrimaryVertexTPC(vtxTPC);
1008   AliESDVertex *vertexTPC=&vtxTPC;
1009
1010   AliVfriendEvent *friendEvent=event->FindFriend();
1011   const Double_t kMinPt=1;
1012   const Double_t kMinPtMax=0.8;
1013   const Double_t kMinNcl=50;
1014   const Double_t kMaxDelta[5]={2,600,0.02,0.02,0.1};
1015   Int_t ntracks=event->GetNumberOfTracks(); 
1016   //  Float_t dcaTPC[2]={0,0};
1017   // Float_t covTPC[3]={0,0,0};
1018
1019   UInt_t specie = event->GetEventSpecie();  // skip laser events
1020   if (specie==AliRecoParam::kCalib) return;
1021   
1022
1023
1024   for (Int_t itrack0=0;itrack0<ntracks;itrack0++) {
1025     AliVTrack *track0 = event->GetVTrack(itrack0);
1026     if (!track0) continue;
1027     if (!track0->IsOn(AliVTrack::kTPCrefit)) continue;
1028
1029     if (TMath::Abs(AliTracker::GetBz())>1&&track0->Pt()<kMinPt) continue;
1030     if (track0->GetTPCncls()<kMinNcl) continue;
1031     if (TMath::Abs(track0->GetY())<kMaxDelta[0]) continue; 
1032     if (track0->GetKinkIndex(0)>0) continue;
1033     const Double_t * par0=track0->GetParameter(); //track param at rhe DCA
1034     //rm primaries
1035     //
1036     //track0->GetImpactParametersTPC(dcaTPC,covTPC);
1037     //if (TMath::Abs(dcaTPC[0])<kMaxDelta[0]) continue;
1038     //if (TMath::Abs(dcaTPC[1])<kMaxDelta[0]*2) continue;
1039     //    const AliExternalTrackParam * trackIn0 = track0->GetInnerParam();
1040     for (Int_t itrack1=itrack0+1;itrack1<ntracks;itrack1++) {
1041       AliVTrack *track1 = event->GetVTrack(itrack1);
1042       if (!track1) continue;  
1043       if (!track1->IsOn(AliVTrack::kTPCrefit)) continue;
1044       if (track1->GetKinkIndex(0)>0) continue;
1045       if (TMath::Abs(AliTracker::GetBz())>1&&track1->Pt()<kMinPt) continue;
1046       if (track1->GetTPCncls()<kMinNcl) continue;
1047       if (TMath::Abs(AliTracker::GetBz())>1&&TMath::Max(track1->Pt(), track0->Pt())<kMinPtMax) continue;
1048       if (TMath::Abs(track1->GetY())<kMaxDelta[0]) continue;
1049       //track1->GetImpactParametersTPC(dcaTPC,covTPC);
1050       //      if (TMath::Abs(dcaTPC[0])<kMaxDelta[0]) continue;
1051       //if (TMath::Abs(dcaTPC[1])<kMaxDelta[0]*2) continue;
1052       //
1053       const Double_t* par1=track1->GetParameter(); //track param at rhe DCA
1054       //
1055       Bool_t isPair=kTRUE;
1056       for (Int_t ipar=0; ipar<5; ipar++){
1057         if (ipar==4&&TMath::Abs(AliTracker::GetBz())<1) continue; // 1/pt not defined for B field off
1058         if (TMath::Abs(TMath::Abs(par0[ipar])-TMath::Abs(par1[ipar]))>kMaxDelta[ipar]) isPair=kFALSE;
1059       }
1060       if (!isPair) continue;
1061       if (TMath::Abs(TMath::Abs(track0->GetAlpha()-track1->GetAlpha())-TMath::Pi())>kMaxDelta[2]) isPair=kFALSE;
1062       //delta with correct sign
1063       /*
1064         TCut cut0="abs(t1.fP[0]+t0.fP[0])<2"
1065         TCut cut3="abs(t1.fP[3]+t0.fP[3])<0.02"
1066         TCut cut4="abs(t1.fP[4]+t0.fP[4])<0.2"
1067       */
1068       if  (TMath::Abs(par0[0]+par1[0])>kMaxDelta[0]) isPair=kFALSE; //delta y   opposite sign
1069       if  (TMath::Abs(par0[3]+par1[3])>kMaxDelta[3]) isPair=kFALSE; //delta tgl opposite sign
1070       if  (TMath::Abs(AliTracker::GetBz())>1 && TMath::Abs(par0[4]+par1[4])>kMaxDelta[4]) isPair=kFALSE; //delta 1/pt opposite sign
1071       if (!isPair) continue;
1072       TString filename(AliAnalysisManager::GetAnalysisManager()->GetTree()->GetCurrentFile()->GetName());
1073       Int_t eventNumber = event->GetEventNumberInFile(); 
1074       Bool_t hasFriend=(friendEvent) ? (friendEvent->GetTrack(itrack0)!=0):0;
1075       Bool_t hasITS=(track0->GetNcls(0)+track1->GetNcls(0)>4);
1076       printf("DUMPHPTCosmic:%s|%f|%d|%d|%d\n",filename.Data(),(TMath::Min(track0->Pt(),track1->Pt())), eventNumber,hasFriend,hasITS);
1077
1078
1079       //      const AliExternalTrackParam * trackIn1 = track1->GetInnerParam();      
1080       //
1081       //       
1082       TTreeSRedirector * pcstream =  GetDebugStreamer();
1083       Int_t ntracksSPD = vertexSPD->GetNContributors();
1084       Int_t ntracksTPC = vertexTPC->GetNContributors();
1085       //
1086       const AliVfriendTrack *friendTrack0 = friendEvent->GetTrack(itrack0);
1087       if (!friendTrack0) continue;
1088       const AliVfriendTrack *friendTrack1 = friendEvent->GetTrack(itrack1);
1089       if (!friendTrack1) continue;
1090       TObject *calibObject;
1091       AliTPCseed *seed0 = 0;   
1092       AliTPCseed *seed1 = 0;
1093       //
1094       for (Int_t l=0;(calibObject=friendTrack0->GetCalibObject(l));++l) {
1095         if ((seed0=dynamic_cast<AliTPCseed*>(calibObject))) break;
1096       }
1097       for (Int_t l=0;(calibObject=friendTrack1->GetCalibObject(l));++l) {
1098         if ((seed1=dynamic_cast<AliTPCseed*>(calibObject))) break;
1099       }
1100
1101       //
1102       if (pcstream){
1103         (*pcstream)<<"pairs"<<
1104           "run="<<fRun<<              //  run number
1105           "event="<<fEvent<<          //  event number
1106           "time="<<fTime<<            //  time stamp of event
1107           "trigger="<<fTrigger<<      //  trigger
1108           "triggerClass="<<&fTriggerClass<<      //  trigger
1109           "mag="<<fMagF<<             //  magnetic field
1110           //
1111           "nSPD="<<ntracksSPD<<
1112           "nTPC="<<ntracksTPC<<
1113           "vSPD.="<<vertexSPD<<         //primary vertex -SPD
1114           "vTPC.="<<vertexTPC<<         //primary vertex -TPC
1115           "t0.="<<track0<<              //track0
1116           "t1.="<<track1<<              //track1
1117       //"ft0.="<<dummyfriendTrack0<<       //track0
1118       //"ft1.="<<dummyfriendTrack1<<       //track1
1119           "s0.="<<seed0<<               //track0
1120           "s1.="<<seed1<<               //track1
1121           "\n";      
1122       }
1123
1124       //**********************TEMPORARY!!*******************************************
1125       // more investigation is needed with Tree ///!!!
1126       //all dummy stuff here is just for code to compile and work with ESD
1127
1128       AliESDfriendTrack *dummyfriendTrack0 = (AliESDfriendTrack*)friendTrack0;
1129       AliESDfriendTrack *dummyfriendTrack1 = (AliESDfriendTrack*)friendTrack1;
1130
1131       AliESDtrack *dummytrack0 = (AliESDtrack*)track0;
1132       AliESDtrack *dummytrack1 = (AliESDtrack*)track1;
1133
1134       if ((pcstream)&&(dummyfriendTrack0)){
1135           (*pcstream)<<"ft0.="<<dummyfriendTrack0<<"\n";
1136       }
1137       if ((pcstream)&&(dummyfriendTrack1)){
1138           (*pcstream)<<"ft1.="<<dummyfriendTrack1<<"\n";
1139       }
1140
1141       if (!fCosmicTree) {
1142         fCosmicTree = new TTree("pairs","pairs");
1143         fCosmicTree->SetDirectory(0);
1144       }
1145       if (fCosmicTree->GetEntries()==0){
1146         //
1147         fCosmicTree->SetDirectory(0);
1148     fCosmicTree->Branch("t0.",&dummytrack0);
1149     fCosmicTree->Branch("t1.",&dummytrack1);
1150     fCosmicTree->Branch("ft0.",&dummyfriendTrack0);
1151     fCosmicTree->Branch("ft1.",&dummyfriendTrack1);
1152       }else{
1153     fCosmicTree->SetBranchAddress("t0.",&dummytrack0);
1154     fCosmicTree->SetBranchAddress("t1.",&dummytrack1);
1155     fCosmicTree->SetBranchAddress("ft0.",&dummyfriendTrack0);
1156     fCosmicTree->SetBranchAddress("ft1.",&dummyfriendTrack1);
1157       }
1158       fCosmicTree->Fill();
1159     }
1160   }
1161 }
1162
1163
1164 void  AliTPCcalibCosmic::Terminate(){
1165   //
1166   // copy the cosmic tree to memory resident tree
1167   //
1168   static Int_t counter=0;
1169   printf("AliTPCcalibCosmic::Terminate\t%d\n",counter);
1170   counter++;
1171   AliTPCcalibBase::Terminate();
1172 }
1173
1174
1175 void AliTPCcalibCosmic::AddTree(TTree* treeOutput, TTree * treeInput){
1176   //
1177   // Add the content of tree: 
1178   // Notice automatic copy of tree in ROOT does not work for such complicated tree
1179   //  
1180   return;
1181   //if (TMath::Abs(fMagF)<0.1) return; // work around - otherwise crashes 
1182   AliESDtrack *track0=new AliESDtrack;     ///!!!
1183   AliESDtrack *track1=new AliESDtrack;
1184   AliESDfriendTrack *ftrack0=new AliESDfriendTrack;
1185   AliESDfriendTrack *ftrack1=new AliESDfriendTrack;
1186   treeInput->SetBranchAddress("t0.",&track0);   
1187   treeInput->SetBranchAddress("t1.",&track1);
1188   treeInput->SetBranchAddress("ft0.",&ftrack0); 
1189   treeInput->SetBranchAddress("ft1.",&ftrack1);
1190   treeOutput->SetDirectory(0);
1191   //
1192   Int_t entries= treeInput->GetEntries();
1193   Int_t step=1+Int_t(TMath::Log(1+treeOutput->GetEntries()/10.));
1194   for (Int_t i=0; i<entries; i+=step){
1195     treeInput->SetBranchAddress("t0.",&track0); 
1196     treeInput->SetBranchAddress("t1.",&track1);
1197     treeInput->SetBranchAddress("ft0.",&ftrack0);       
1198     treeInput->SetBranchAddress("ft1.",&ftrack1);
1199     treeInput->GetEntry(i);
1200     if (!track0) continue;
1201     if (!track1) continue;
1202     if (!ftrack0) continue;
1203     if (!ftrack1) continue;
1204     if (track0->GetTPCncls()<=0) continue;
1205     if (track1->GetTPCncls()<=0) continue;
1206     if (!track0->GetInnerParam()) continue;
1207     if (!track1->GetInnerParam()) continue;
1208     if (!track0->GetTPCInnerParam()) continue;
1209     if (!track1->GetTPCInnerParam()) continue;
1210     //track0
1211     treeOutput->SetBranchAddress("t0.",&track0);        
1212     treeOutput->SetBranchAddress("t1.",&track1);
1213     treeOutput->SetBranchAddress("ft0.",&ftrack0);      
1214     treeOutput->SetBranchAddress("ft1.",&ftrack1);    
1215     treeOutput->Fill();
1216     delete track0;
1217     delete track1;
1218     delete ftrack0;
1219     delete ftrack1;
1220     track0=0;
1221     track1=0;
1222     ftrack0=0;
1223     ftrack1=0;
1224   }
1225 }
1226
1227
1228
1229 void AliTPCcalibCosmic::MakeFitTree(TTree * treeInput, TTreeSRedirector *pcstream, const TObjArray * corrArray, Int_t step, Int_t run){
1230   //
1231   // Make fit tree
1232   // refit the tracks with original points + corrected points for each correction
1233   // Input:
1234   //   treeInput - tree with cosmic tracks
1235   //   pcstream  - debug output
1236
1237   // Algorithm:
1238   // Loop over pair of cosmic tracks:
1239   //   1. Find trigger offset between cosmic event and "physic" trigger
1240   //   2. Refit tracks with current transformation
1241   //   3. Refit tracks using additional "primitive" distortion on top
1242   // Best correction estimated as linear combination of corrections 
1243   // minimizing the observed distortions
1244   // Observed distortions - matching in the y,z, snp, theta and 1/pt
1245   //
1246   const Double_t kResetCov=20.;
1247   const Double_t kMaxDelta[5]={1,40,0.03,0.01,0.2};
1248   const Double_t kSigma=2.;    
1249   const Double_t kMaxTime=1050;
1250   const Double_t kMaxSnp=0.8;
1251   Int_t ncorr=corrArray->GetEntries();
1252   AliTPCTransform *transform = AliTPCcalibDB::Instance()->GetTransform() ;
1253   AliTPCParam     *param     = AliTPCcalibDB::Instance()->GetParameters();
1254   AliGRPObject*  grp = AliTPCcalibDB::Instance()->GetGRP(run);
1255   Double_t time=0.5*(grp->GetTimeStart() +grp->GetTimeEnd()); 
1256   transform->SetCurrentRun(run);
1257   transform->SetCurrentTimeStamp(TMath::Nint(time));
1258   Double_t covar[15];
1259   for (Int_t i=0;i<15;i++) covar[i]=0;
1260   covar[0]=kSigma*kSigma;
1261   covar[2]=kSigma*kSigma;
1262   covar[5]=kSigma*kSigma/Float_t(150*150);
1263   covar[9]=kSigma*kSigma/Float_t(150*150);
1264   covar[14]=0.2*0.2;
1265   Double_t *distortions = new Double_t[ncorr+1];
1266
1267   AliESDtrack *track0=new AliESDtrack;
1268   AliESDtrack *track1=new AliESDtrack;
1269   AliESDfriendTrack *ftrack0=new AliESDfriendTrack;
1270   AliESDfriendTrack *ftrack1=new AliESDfriendTrack;
1271   treeInput->SetBranchAddress("t0.",&track0);   
1272   treeInput->SetBranchAddress("t1.",&track1);
1273   treeInput->SetBranchAddress("ft0.",&ftrack0); 
1274   treeInput->SetBranchAddress("ft1.",&ftrack1);
1275   Int_t entries= treeInput->GetEntries();
1276   for (Int_t i=0; i<entries; i+=step){    
1277     treeInput->GetEntry(i);
1278     if (i%20==0) printf("%d\n",i);
1279     if (!ftrack0->GetTPCOut()) continue;
1280     if (!ftrack1->GetTPCOut()) continue;
1281     AliTPCseed *seed0=0;
1282     AliTPCseed *seed1=0;
1283     TObject *calibObject;
1284     for (Int_t l=0;(calibObject=ftrack0->GetCalibObject(l));++l) {
1285       if ((seed0=dynamic_cast<AliTPCseed*>(calibObject))) break;
1286     }
1287     for (Int_t l=0;(calibObject=ftrack1->GetCalibObject(l));++l) {
1288       if ((seed1=dynamic_cast<AliTPCseed*>(calibObject))) break;
1289     }
1290     if (!seed0) continue;
1291     if (!seed1) continue;
1292     if (TMath::Abs(seed0->GetSnp())>kMaxSnp) continue;
1293     if (TMath::Abs(seed1->GetSnp())>kMaxSnp) continue;
1294     //
1295     //
1296     Int_t nclA0=0, nclC0=0;     // number of clusters
1297     Int_t nclA1=0, nclC1=0;     // number of clusters
1298     Int_t ncl0=0,ncl1=0;
1299     Double_t rmin0=300, rmax0=-300;  // variables to estimate the time0 - trigger offset
1300     Double_t rmin1=300, rmax1=-300;
1301     Double_t tmin0=2000, tmax0=-2000;
1302     Double_t tmin1=2000, tmax1=-2000;
1303     //
1304     //
1305     // calculate trigger offset usig "missing clusters"
1306     for (Int_t irow=0; irow<159; irow++){
1307       AliTPCclusterMI *cluster0=seed0->GetClusterPointer(irow);
1308       if (cluster0 &&cluster0->GetX()>10){
1309         if (cluster0->GetX()<rmin0) { rmin0=cluster0->GetX(); tmin0=cluster0->GetTimeBin();}
1310         if (cluster0->GetX()>rmax0) { rmax0=cluster0->GetX(); tmax0=cluster0->GetTimeBin();}
1311         ncl0++;
1312         if (cluster0->GetDetector()%36<18) nclA0++;
1313         if (cluster0->GetDetector()%36>=18) nclC0++;
1314       }  
1315       AliTPCclusterMI *cluster1=seed1->GetClusterPointer(irow);
1316       if (cluster1&&cluster1->GetX()>10){
1317         if (cluster1->GetX()<rmin1) { rmin1=cluster1->GetX();  tmin1=cluster1->GetTimeBin();}
1318         if (cluster1->GetX()>rmax1) { rmax1=cluster1->GetX(); tmax1=cluster1->GetTimeBin();}
1319         ncl1++;
1320         if (cluster1->GetDetector()%36<18) nclA1++;
1321         if (cluster1->GetDetector()%36>=18) nclC1++;
1322       }
1323     }
1324     Int_t cosmicType=0;  // types of cosmic topology
1325     if ((nclA0>nclC0) && (nclA1>nclC1)) cosmicType=0; // AA side
1326     if ((nclA0<nclC0) && (nclA1<nclC1)) cosmicType=1; // CC side
1327     if ((nclA0>nclC0) && (nclA1<nclC1)) cosmicType=2; // AC side
1328     if ((nclA0<nclC0) && (nclA1>nclC1)) cosmicType=3; // CA side
1329     //if ((nclA0>nclC0) && (nclA1<nclC1)) cosmicType=6; // AC side out of time
1330     //if ((nclA0>nclC0) && (nclA1<nclC1)) cosmicType=7; // CA side out of time
1331     //
1332     Double_t deltaTime = 0;   // correction for trigger not in time - equalizing the time arival
1333     if ((cosmicType>1)&&TMath::Abs(track1->GetZ()-track0->GetZ())>4){
1334       cosmicType+=4;
1335       deltaTime=0.5*(track1->GetZ()-track0->GetZ())/param->GetZWidth();
1336       if (nclA0>nclC0) deltaTime*=-1; // if A side track
1337     }
1338     //
1339     TVectorD vectorDT(8);
1340     Int_t crossCounter=0;
1341     Double_t deltaTimeCross = AliTPCcalibCosmic::GetDeltaTime(rmin0, rmax0, rmin1, rmax1, tmin0, tmax0, tmin1, tmax1, TMath::Abs(track0->GetY()),vectorDT);
1342     Bool_t isOKTrigger=kTRUE;
1343     for (Int_t ic=0; ic<6;ic++) {
1344       if (TMath::Abs(vectorDT[ic])>0) {
1345         if (vectorDT[ic]+vectorDT[6]<0) isOKTrigger=kFALSE;
1346         if (vectorDT[ic]+vectorDT[7]>kMaxTime) isOKTrigger=kFALSE;
1347         if (isOKTrigger){
1348           crossCounter++; 
1349         }
1350       }
1351     }
1352     Double_t deltaTimeCluster=deltaTime;
1353     if ((cosmicType==0 || cosmicType==1) && crossCounter>0){
1354       deltaTimeCluster=deltaTimeCross;
1355       cosmicType+=8;
1356     }
1357     if (nclA0*nclC0>0 || nclA1*nclC1>0) cosmicType+=16;  // mixed A side C side - bad for visualization
1358     //
1359     // Apply current transformation
1360     //
1361     //
1362     for (Int_t irow=0; irow<159; irow++){
1363       AliTPCclusterMI *cluster0=seed0->GetClusterPointer(irow);
1364       if (cluster0 &&cluster0->GetX()>10){
1365         Double_t x0[3]={ static_cast<Double_t>(cluster0->GetRow()),cluster0->GetPad(),cluster0->GetTimeBin()+deltaTimeCluster};
1366         Int_t index0[1]={cluster0->GetDetector()};
1367         transform->Transform(x0,index0,0,1);  
1368         cluster0->SetX(x0[0]);
1369         cluster0->SetY(x0[1]);
1370         cluster0->SetZ(x0[2]);
1371         //
1372       }
1373       AliTPCclusterMI *cluster1=seed1->GetClusterPointer(irow);
1374       if (cluster1&&cluster1->GetX()>10){
1375         Double_t x1[3]={ static_cast<Double_t>(cluster1->GetRow()),cluster1->GetPad(),cluster1->GetTimeBin()+deltaTimeCluster};
1376         Int_t index1[1]={cluster1->GetDetector()};
1377         transform->Transform(x1,index1,0,1);  
1378         cluster1->SetX(x1[0]);
1379         cluster1->SetY(x1[1]);
1380         cluster1->SetZ(x1[2]);
1381       }
1382     }
1383     //
1384     //
1385     Double_t alpha=track0->GetAlpha();   // rotation frame
1386     Double_t cos = TMath::Cos(alpha);
1387     Double_t sin = TMath::Sin(alpha);
1388     Double_t mass =  TDatabasePDG::Instance()->GetParticle("mu+")->Mass();
1389     AliExternalTrackParam  btrack0=*(ftrack0->GetTPCOut());
1390     AliExternalTrackParam  btrack1=*(ftrack1->GetTPCOut());
1391     btrack0.Rotate(alpha);
1392     btrack1.Rotate(alpha);
1393     // change the sign for track 1
1394     Double_t* par1=(Double_t*)btrack0.GetParameter();
1395     par1[3]*=-1;
1396     par1[4]*=-1;
1397     btrack0.AddCovariance(covar);
1398     btrack1.AddCovariance(covar);
1399     btrack0.ResetCovariance(kResetCov);
1400     btrack1.ResetCovariance(kResetCov);
1401     Bool_t isOK=kTRUE;
1402     Bool_t isOKT=kTRUE;
1403     TObjArray tracks0(ncorr+1);
1404     TObjArray tracks1(ncorr+1);
1405     //    
1406     Double_t dEdx0Tot=seed0->CookdEdxAnalytical(0.02,0.6,kTRUE);
1407     Double_t dEdx1Tot=seed1->CookdEdxAnalytical(0.02,0.6,kTRUE);
1408     Double_t dEdx0Max=seed0->CookdEdxAnalytical(0.02,0.6,kFALSE);
1409     Double_t dEdx1Max=seed1->CookdEdxAnalytical(0.02,0.6,kFALSE);
1410     //if (TMath::Abs((dEdx0Max+1)/(dEdx0Tot+1)-1.)>0.1) isOK=kFALSE;
1411     //if (TMath::Abs((dEdx1Max+1)/(dEdx1Tot+1)-1.)>0.1) isOK=kFALSE;
1412     ncl0=0; ncl1=0;
1413     for (Int_t icorr=-1; icorr<ncorr; icorr++){
1414       AliExternalTrackParam  rtrack0=btrack0;
1415       AliExternalTrackParam  rtrack1=btrack1;
1416       AliTPCCorrection *corr = 0;
1417       if (icorr>=0) corr = (AliTPCCorrection*)corrArray->At(icorr);
1418       //
1419       for (Int_t irow=159; irow>0; irow--){ 
1420         AliTPCclusterMI *cluster=seed0->GetClusterPointer(irow);
1421         if (!cluster) continue;
1422         if (!isOKT) break;
1423         Double_t rD[3]={cluster->GetX(),cluster->GetY(),cluster->GetZ()};
1424         transform->RotatedGlobal2Global(cluster->GetDetector()%36,rD);  // transform to global
1425         Float_t  r[3]={static_cast<Float_t>(rD[0]),static_cast<Float_t>(rD[1]),static_cast<Float_t>(rD[2])};
1426         if (corr){
1427           corr->DistortPoint(r, cluster->GetDetector());
1428         }
1429         //
1430         Double_t cov[3]={0.01,0.,0.01}; 
1431         Double_t lx =cos*r[0]+sin*r[1];      
1432         Double_t ly =-sin*r[0]+cos*r[1];
1433         rD[1]=ly; rD[0]=lx; rD[2]=r[2];  //transform to track local
1434         if (!AliTracker::PropagateTrackToBxByBz(&rtrack0, lx,mass,1.,kFALSE)) isOKT=kFALSE;;
1435         if (!rtrack0.Update(&rD[1],cov)) isOKT =kFALSE;
1436         if (icorr<0) ncl0++;
1437       }
1438       //
1439       for (Int_t irow=159; irow>0; irow--){ 
1440         AliTPCclusterMI *cluster=seed1->GetClusterPointer(irow);
1441         if (!cluster) continue;
1442         if (!isOKT) break;
1443         Double_t rD[3]={cluster->GetX(),cluster->GetY(),cluster->GetZ()};
1444         transform->RotatedGlobal2Global(cluster->GetDetector()%36,rD);
1445         Float_t  r[3]={static_cast<Float_t>(rD[0]),static_cast<Float_t>(rD[1]),static_cast<Float_t>(rD[2])};
1446         if (corr){
1447           corr->DistortPoint(r, cluster->GetDetector());
1448         }
1449         //
1450         Double_t cov[3]={0.01,0.,0.01}; 
1451         Double_t lx =cos*r[0]+sin*r[1];      
1452         Double_t ly =-sin*r[0]+cos*r[1];
1453         rD[1]=ly; rD[0]=lx; rD[2]=r[2];
1454         if (!AliTracker::PropagateTrackToBxByBz(&rtrack1, lx,mass,1.,kFALSE)) isOKT=kFALSE;
1455         if (!rtrack1.Update(&rD[1],cov)) isOKT=kFALSE;
1456         if (icorr<0) ncl1++;
1457       }
1458       if (!AliTracker::PropagateTrackToBxByBz(&rtrack0, 0,mass,10.,kFALSE)) isOKT=kFALSE;
1459       if (!AliTracker::PropagateTrackToBxByBz(&rtrack1, 0,mass,10.,kFALSE)) isOKT=kFALSE;
1460       if (!AliTracker::PropagateTrackToBxByBz(&rtrack0, 0,mass,1.,kFALSE))  isOKT=kFALSE;
1461       if (!AliTracker::PropagateTrackToBxByBz(&rtrack1, 0,mass,1.,kFALSE))  isOKT=kFALSE;
1462       const Double_t *param0=rtrack0.GetParameter();
1463       const Double_t *param1=rtrack1.GetParameter();
1464       for (Int_t ipar=0; ipar<4; ipar++){
1465         if (TMath::Abs(param1[ipar]-param0[ipar])>kMaxDelta[ipar]) isOK=kFALSE;
1466       }
1467       tracks0.AddAt(rtrack0.Clone(), icorr+1);
1468       tracks1.AddAt(rtrack1.Clone(), icorr+1);
1469       AliExternalTrackParam out0=*(ftrack0->GetTPCOut());
1470       AliExternalTrackParam out1=*(ftrack1->GetTPCOut());
1471       Int_t nentries=TMath::Min(ncl0,ncl1);
1472
1473       if (icorr<0) {
1474         (*pcstream)<<"cosmic"<<
1475           "isOK="<<isOK<<              // correct all propagation update and also residuals
1476           "isOKT="<<isOKT<<            // correct all propagation update
1477           "isOKTrigger="<<isOKTrigger<< // correct? estimate of trigger offset
1478           "id="<<cosmicType<<
1479           //
1480           //
1481           "cross="<<crossCounter<<
1482           "vDT.="<<&vectorDT<<
1483           //
1484           "dTime="<<deltaTime<<        // delta time using the A-c side cross
1485           "dTimeCross="<<deltaTimeCross<< // delta time using missing clusters
1486           //
1487           "dEdx0Max="<<dEdx0Max<<
1488           "dEdx0Tot="<<dEdx0Tot<<
1489           "dEdx1Max="<<dEdx1Max<<
1490           "dEdx1Tot="<<dEdx1Tot<<
1491           //
1492           "track0.="<<track0<<         // original track 0
1493           "track1.="<<track1<<         // original track 1
1494           "out0.="<<&out0<<             // outer track  0
1495           "out1.="<<&out1<<             // outer track  1
1496           "rtrack0.="<<&rtrack0<<      // refitted track with current transform
1497           "rtrack1.="<<&rtrack1<<     //          
1498           "ncl0="<<ncl0<<
1499           "ncl1="<<ncl1<<
1500           "entries="<<nentries<<       // number of clusters
1501           "\n";
1502       }
1503     }
1504     //
1505
1506     if (isOK){        
1507       Int_t nentries=TMath::Min(ncl0,ncl1);    
1508       for (Int_t ipar=0; ipar<5; ipar++){
1509         for (Int_t icorr=-1; icorr<ncorr; icorr++){
1510           AliTPCCorrection *corr = 0;
1511           if (icorr>=0) corr = (AliTPCCorrection*)corrArray->At(icorr);
1512           //
1513           AliExternalTrackParam *param0=(AliExternalTrackParam *) tracks0.At(icorr+1);
1514           AliExternalTrackParam *param1=(AliExternalTrackParam *) tracks1.At(icorr+1);
1515           distortions[icorr+1]=param1->GetParameter()[ipar]-param0->GetParameter()[ipar];
1516           if (icorr>=0){
1517             distortions[icorr+1]-=distortions[0];
1518           }
1519           //
1520           if (icorr<0){
1521             Double_t bz=AliTrackerBase::GetBz();
1522             Double_t gxyz[3];
1523             param0->GetXYZ(gxyz);
1524             Int_t dtype=20;
1525             Double_t theta=param0->GetParameter()[3];
1526             Double_t phi = alpha;
1527             Double_t snp = track0->GetInnerParam()->GetSnp();
1528             Double_t mean= distortions[0];
1529             Int_t index = param0->GetIndex(ipar,ipar);
1530             Double_t rms=TMath::Sqrt(param1->GetCovariance()[index]+param1->GetCovariance()[index]);
1531             if (crossCounter<1) rms*=1;
1532             Double_t sector=9*phi/TMath::Pi();
1533             Double_t dsec   = sector-TMath::Nint(sector+0.5);
1534             Double_t gx=gxyz[0],gy=gxyz[1],gz=gxyz[2];
1535             Double_t refX=TMath::Sqrt(gx*gx+gy*gy);
1536             Double_t dRrec=0;
1537             //      Double_t pt=(param0->GetSigned1Pt()+param1->GetSigned1Pt())*0.5;
1538             Double_t pt=(param0->GetSigned1Pt()+param1->GetSigned1Pt())*0.5;
1539
1540             (*pcstream)<<"fit"<<  // dump valus for fit
1541               "run="<<run<<       //run number
1542               "bz="<<bz<<         // magnetic filed used
1543               "dtype="<<dtype<<   // detector match type 20
1544               "ptype="<<ipar<<   // parameter type
1545               "theta="<<theta<<   // theta
1546               "phi="<<phi<<       // phi 
1547               "snp="<<snp<<       // snp
1548               "mean="<<mean<<     // mean dist value
1549               "rms="<<rms<<       // rms
1550               "sector="<<sector<<
1551               "dsec="<<dsec<<
1552               //
1553               "refX="<<refX<<      // reference radius
1554               "gx="<<gx<<         // global position
1555               "gy="<<gy<<         // global position
1556               "gz="<<gz<<         // global position
1557               "dRrec="<<dRrec<<      // delta Radius in reconstruction
1558               "pt="<<pt<<            //1/pt
1559               "id="<<cosmicType<<     //type of the cosmic used      
1560               "entries="<<nentries;// number of entries in bin
1561             (*pcstream)<<"cosmicDebug"<<
1562               "p0.="<<param0<<    // dump distorted track 0
1563               "p1.="<<param1;    // dump distorted track 1
1564           }
1565           if (icorr>=0){
1566             (*pcstream)<<"fit"<<
1567               Form("%s=",corr->GetName())<<distortions[icorr+1];    // dump correction value
1568             (*pcstream)<<"cosmicDebug"<<
1569               Form("%s=",corr->GetName())<<distortions[icorr+1]<<    // dump correction value
1570               Form("%sp0.=",corr->GetName())<<param0<<    // dump distorted track 0
1571               Form("%sp1.=",corr->GetName())<<param1;    // dump distorted track 1
1572           }
1573         } //loop corrections      
1574         (*pcstream)<<"fit"<<"isOK="<<isOK<<"\n";
1575         (*pcstream)<<"cosmicDebug"<<"isOK="<<isOK<<"\n";
1576       } //loop over parameters
1577     } // dump results
1578   }//loop tracks
1579   delete [] distortions;
1580 }
1581
1582
1583
1584 Double_t AliTPCcalibCosmic::GetDeltaTime(Double_t rmin0, Double_t rmax0, Double_t rmin1, Double_t rmax1, Double_t tmin0, Double_t tmax0, Double_t tmin1, Double_t tmax1, Double_t dcaR, TVectorD &vectorDT)
1585 {
1586   //
1587   // Estimate trigger offset between random cosmic event and "physics" trigger
1588   // Efficiency about 50 % of cases:
1589   // Cases:
1590   // 0. Tracks crossing A side and C side - no match in z - 30 % of cases
1591   // 1. Track only on one side and  dissapear at small or at high radiuses - 50 % of cases
1592   //    1.a) rmax<Rc    && tmax<Tcmax && tmax>tmin    => deltaT=Tcmax-tmax 
1593   //    1.b) rmin>Rcmin && tmin<Tcmax && tmin>tmax    => deltaT=Tcmax-tmin  
1594   //                      // case the z matching gives proper time
1595   //    1.c) rmax<Rc    && tmax>Tcmin && tmax<tmin    => deltaT=-tmax
1596   //
1597   // check algorithm:
1598   //    TCut cutStop = "(min(rmax0,rmax1)<235||abs(rmin0-rmin1)>10)"; // tracks not registered for full time
1599
1600   // Combinations:
1601   // 0-1 - forbidden
1602   // 0-2 - forbidden
1603   // 0-3 - occur - wrong correlation
1604   // 1-2 - occur - wrong correlation
1605   // 1-3 - forbidden
1606   // 2-3 - occur - small number of outlyers -20%
1607   // Frequency:
1608   //   0 - 106
1609   //   1 - 265
1610   //   2 - 206
1611   //   3 - 367
1612   //
1613   const Double_t kMaxRCut=235;  // max radius
1614   const Double_t kMinRCut=TMath::Max(dcaR,90.);  // min radius
1615   const Double_t kMaxDCut=30;   // max distance for minimal radius
1616   const Double_t kMinTime=110;
1617   const Double_t kMaxTime=950;  
1618   Double_t deltaT=0;
1619   Int_t counter=0;
1620   vectorDT[6]=TMath::Min(TMath::Min(tmin0,tmin1),TMath::Min(tmax0,tmax1));
1621   vectorDT[7]=TMath::Max(TMath::Max(tmin0,tmin1),TMath::Max(tmax0,tmax1));
1622   if (TMath::Min(rmax0,rmax1)<kMaxRCut){
1623     // max cross - deltaT>0
1624     if (rmax0<kMaxRCut && tmax0 <kMaxTime && tmax0>tmin0) vectorDT[0]=kMaxTime-tmax0; // disapear at CE
1625     if (rmax1<kMaxRCut && tmax1 <kMaxTime && tmax1>tmin1) vectorDT[1]=kMaxTime-tmax1; // disapear at CE
1626     // min cross - deltaT<0 - OK they are correlated
1627     if (rmax0<kMaxRCut && tmax0 >kMinTime && tmax0<tmin0) vectorDT[2]=-tmax0;         // disapear at ROC
1628     if (rmax1<kMaxRCut && tmax1 >kMinTime && tmax1<tmin1) vectorDT[3]=-tmax1;         // disapear at ROC
1629   }  
1630   if (rmin0> kMinRCut+kMaxDCut && tmin0 <kMaxTime && tmin0>tmax0) vectorDT[4]=kMaxTime-tmin0;
1631   if (rmin1> kMinRCut+kMaxDCut && tmin1 <kMaxTime && tmin1>tmax1) vectorDT[5]=kMaxTime-tmin1;
1632   Bool_t isOK=kTRUE;
1633   for (Int_t i=0; i<6;i++) {
1634     if (TMath::Abs(vectorDT[i])>0) {
1635       counter++; 
1636       if (vectorDT[i]+vectorDT[6]<0) isOK=kFALSE;
1637       if (vectorDT[i]+vectorDT[7]>kMaxTime) isOK=kFALSE;
1638       if (isOK) deltaT=vectorDT[i];
1639     }
1640   }
1641   return deltaT;  
1642 }
1643