1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 Comments to be written here:
19 1. What do we calibrate.
21 Time dependence of gain and drift velocity in order to account for changes in: temperature, pressure, gas composition.
23 AliTPCcalibTime *calibTime = new AliTPCcalibTime("cosmicTime","cosmicTime",0, 1213.9e+06, 1213.96e+06, 0.04e+04, 0.04e+04);
25 2. How to interpret results
29 a) determine the required time range:
31 AliXRDPROOFtoolkit tool;
32 TChain * chain = tool.MakeChain("pass2.txt","esdTree",0,6000);
33 chain->Draw("GetTimeStamp()")
35 b) analyse calibration object on Proof in calibration train
37 AliTPCcalibTime *calibTime = new AliTPCcalibTime("cosmicTime","cosmicTime", StartTimeStamp, EndTimeStamp, IntegrationTimeVdrift);
41 gSystem->Load("libANALYSIS");
42 gSystem->Load("libTPCcalib");
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")
49 4. Analysis using debug streamers.
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);
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);
68 #include "Riostream.h"
74 #include "THnSparse.h"
82 #include "TGraphErrors.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"
92 #include "AliTracker.h"
94 #include "AliTPCCalROC.h"
95 #include "AliTPCParam.h"
99 #include "AliTPCcalibTime.h"
100 #include "AliRelAlignerKalman.h"
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"
110 #include "TDatabasePDG.h"
111 #include "AliTrackPointArray.h"
113 ClassImp(AliTPCcalibTime)
116 AliTPCcalibTime::AliTPCcalibTime()
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
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
141 // fBinsVdrift(fTimeBins,fPtBins,fVdriftBins),
142 // fXminVdrift(fTimeStart,fPtStart,fVdriftStart),
143 // fXmaxVdrift(fTimeEnd,fPtEnd,fVdriftEnd)
145 AliInfo("Default Constructor");
146 for (Int_t i=0;i<3;i++) {
147 fHistVdriftLaserA[i]=0;
148 fHistVdriftLaserC[i]=0;
150 for (Int_t i=0;i<10;i++) {
151 fCosmiMatchingHisto[i]=0;
155 AliTPCcalibTime::AliTPCcalibTime(const Text_t *name, const Text_t *title, UInt_t StartTime, UInt_t EndTime, Int_t deltaIntegrationTimeVdrift)
157 fLaser(0), // pointer to laser calibration
158 fDz(0), // current delta z
159 fCutMaxD(5*0.5356), // maximal distance in rfi ditection
160 fCutMaxDz(40), // maximal distance in rfi ditection
161 fCutTheta(5*0.004644),// maximal distan theta
162 fCutMinDir(-0.99), // direction vector products
164 fArrayDz(0), //Tmap of V drifts for different triggers
165 fAlignITSTPC(0), //alignemnt array ITS TPC match
166 fAlignTRDTPC(0), //alignemnt array TRD TPC match
167 fAlignTOFTPC(0), //alignemnt array TOF TPC match
183 for (Int_t i=0;i<3;i++) {
184 fHistVdriftLaserA[i]=0;
185 fHistVdriftLaserC[i]=0;
188 AliInfo("Non Default Constructor");
189 fTimeBins =(EndTime-StartTime)/deltaIntegrationTimeVdrift;
190 fTimeStart =StartTime; //(((TObjString*)(mapGRP->GetValue("fAliceStartTime")))->GetString()).Atoi();
191 fTimeEnd =EndTime; //(((TObjString*)(mapGRP->GetValue("fAliceStopTime")))->GetString()).Atoi();
202 Int_t binsVdriftLaser[4] = {fTimeBins , fPtBins , fVdriftBins*20, fRunBins };
203 Double_t xminVdriftLaser[4] = {fTimeStart, fPtStart, fVdriftStart , fRunStart};
204 Double_t xmaxVdriftLaser[4] = {fTimeEnd , fPtEnd , fVdriftEnd , fRunEnd };
205 TString axisTitle[4]={
211 TString histoName[3]={
218 for (Int_t i=0;i<3;i++) {
219 fHistVdriftLaserA[i] = new THnSparseF("HistVdriftLaser","HistVdriftLaser;time;p/T ratio;Vdrift;run",4,binsVdriftLaser,xminVdriftLaser,xmaxVdriftLaser);
220 fHistVdriftLaserC[i] = new THnSparseF("HistVdriftLaser","HistVdriftLaser;time;p/T ratio;Vdrift;run",4,binsVdriftLaser,xminVdriftLaser,xmaxVdriftLaser);
221 fHistVdriftLaserA[i]->SetName(histoName[i]);
222 fHistVdriftLaserC[i]->SetName(histoName[i]);
223 for (Int_t iaxis=0; iaxis<4;iaxis++){
224 fHistVdriftLaserA[i]->GetAxis(iaxis)->SetName(axisTitle[iaxis]);
225 fHistVdriftLaserC[i]->GetAxis(iaxis)->SetName(axisTitle[iaxis]);
228 fBinsVdrift[0] = fTimeBins;
229 fBinsVdrift[1] = fPtBins;
230 fBinsVdrift[2] = fVdriftBins;
231 fBinsVdrift[3] = fRunBins;
232 fXminVdrift[0] = fTimeStart;
233 fXminVdrift[1] = fPtStart;
234 fXminVdrift[2] = fVdriftStart;
235 fXminVdrift[3] = fRunStart;
236 fXmaxVdrift[0] = fTimeEnd;
237 fXmaxVdrift[1] = fPtEnd;
238 fXmaxVdrift[2] = fVdriftEnd;
239 fXmaxVdrift[3] = fRunEnd;
241 fArrayDz=new TObjArray();
242 fAlignITSTPC = new TObjArray; //alignemnt array ITS TPC match
243 fAlignTRDTPC = new TObjArray; //alignemnt array ITS TPC match
244 fAlignTOFTPC = new TObjArray; //alignemnt array ITS TPC match
245 fAlignITSTPC->SetOwner(kTRUE);
246 fAlignTRDTPC->SetOwner(kTRUE);
247 fAlignTOFTPC->SetOwner(kTRUE);
249 // fArrayDz->AddLast(fHistVdriftLaserA[0]);
250 // fArrayDz->AddLast(fHistVdriftLaserA[1]);
251 // fArrayDz->AddLast(fHistVdriftLaserA[2]);
252 // fArrayDz->AddLast(fHistVdriftLaserC[0]);
253 // fArrayDz->AddLast(fHistVdriftLaserC[1]);
254 // fArrayDz->AddLast(fHistVdriftLaserC[2]);
256 fCosmiMatchingHisto[0]=new TH1F("Cosmics matching","p0-all" ,100,-10*0.5356 ,10*0.5356 );
257 fCosmiMatchingHisto[1]=new TH1F("Cosmics matching","p1-all" ,100,-10*4.541 ,10*4.541 );
258 fCosmiMatchingHisto[2]=new TH1F("Cosmics matching","p2-all" ,100,-10*0.01134 ,10*0.01134 );
259 fCosmiMatchingHisto[3]=new TH1F("Cosmics matching","p3-all" ,100,-10*0.004644,10*0.004644);
260 fCosmiMatchingHisto[4]=new TH1F("Cosmics matching","p4-all" ,100,-10*0.03773 ,10*0.03773 );
261 fCosmiMatchingHisto[5]=new TH1F("Cosmics matching","p0-isPair",100,-10*0.5356 ,10*0.5356 );
262 fCosmiMatchingHisto[6]=new TH1F("Cosmics matching","p1-isPair",100,-10*4.541 ,10*4.541 );
263 fCosmiMatchingHisto[7]=new TH1F("Cosmics matching","p2-isPair",100,-10*0.01134 ,10*0.01134 );
264 fCosmiMatchingHisto[8]=new TH1F("Cosmics matching","p3-isPair",100,-10*0.004644,10*0.004644);
265 fCosmiMatchingHisto[9]=new TH1F("Cosmics matching","p4-isPair",100,-10*0.03773 ,10*0.03773 );
266 // Char_t nameHisto[3]={'p','0','\n'};
267 // for (Int_t i=0;i<10;i++){
268 // fCosmiMatchingHisto[i]=new TH1F("Cosmics matching",nameHisto,8192,0,0);
270 // if(i==4) nameHisto[1]='0';
274 AliTPCcalibTime::~AliTPCcalibTime(){
278 for(Int_t i=0;i<3;i++){
279 if(fHistVdriftLaserA[i]){
280 delete fHistVdriftLaserA[i];
281 fHistVdriftLaserA[i]=NULL;
283 if(fHistVdriftLaserC[i]){
284 delete fHistVdriftLaserC[i];
285 fHistVdriftLaserC[i]=NULL;
289 fArrayDz->SetOwner();
294 for(Int_t i=0;i<5;i++){
295 if(fCosmiMatchingHisto[i]){
296 delete fCosmiMatchingHisto[i];
297 fCosmiMatchingHisto[i]=NULL;
300 fAlignITSTPC->SetOwner(kTRUE);
301 fAlignTRDTPC->SetOwner(kTRUE);
302 fAlignTOFTPC->SetOwner(kTRUE);
304 fAlignITSTPC->Delete();
305 fAlignTRDTPC->Delete();
306 fAlignTOFTPC->Delete();
312 Bool_t AliTPCcalibTime::IsLaser(AliESDEvent */*event*/){
313 return kTRUE; //More accurate creteria to be added
315 Bool_t AliTPCcalibTime::IsCosmics(AliESDEvent */*event*/){
316 return kTRUE; //More accurate creteria to be added
318 Bool_t AliTPCcalibTime::IsBeam(AliESDEvent */*event*/){
319 return kTRUE; //More accurate creteria to be added
321 void AliTPCcalibTime::ResetCurrent(){
322 fDz=0; //Reset current dz
324 void AliTPCcalibTime::Process(AliESDEvent *event){
326 if (event->GetNumberOfTracks()<2) return;
328 if(IsLaser (event)) ProcessLaser (event);
329 if(IsCosmics(event)) ProcessCosmic(event);
330 if(IsBeam (event)) ProcessBeam (event);
333 void AliTPCcalibTime::ProcessLaser(AliESDEvent *event){
335 // Fit drift velocity using laser
338 const Int_t kMinTracks = 40; // minimal number of laser tracks
339 const Int_t kMinTracksSide = 20; // minimal number of tracks per side
340 const Float_t kMaxDeltaZ = 30.; // maximal trigger delay
341 const Float_t kMaxDeltaV = 0.05; // maximal deltaV
342 const Float_t kMaxRMS = 0.1; // maximal RMS of tracks
345 TCut cutRMS("sqrt(laserA.fElements[4])<0.1&&sqrt(laserC.fElements[4])<0.1");
346 TCut cutZ("abs(laserA.fElements[0]-laserC.fElements[0])<3");
347 TCut cutV("abs(laserA.fElements[1]-laserC.fElements[1])<0.01");
348 TCut cutY("abs(laserA.fElements[2]-laserC.fElements[2])<2");
349 TCut cutAll = cutRMS+cutZ+cutV+cutY;
351 if (event->GetNumberOfTracks()<kMinTracks) return;
353 if(!fLaser) fLaser = new AliTPCcalibLaser("laserTPC","laserTPC",kFALSE);
354 fLaser->Process(event);
355 if (fLaser->GetNtracks()<kMinTracks) return; // small amount of tracks cut
356 if (fLaser->fFitAside->GetNrows()==0 && fLaser->fFitCside->GetNrows()==0) return; // no fit neither a or C side
358 // debug streamer - activate stream level
359 // Use it for tuning of the cuts
361 // cuts to be applied
363 Int_t isReject[2]={0,0};
366 if (TMath::Abs((*fLaser->fFitAside)[3]) < kMinTracksSide) isReject[0]|=1;
367 if (TMath::Abs((*fLaser->fFitCside)[3]) < kMinTracksSide) isReject[1]|=1;
368 // unreasonable z offset
369 if (TMath::Abs((*fLaser->fFitAside)[0])>kMaxDeltaZ) isReject[0]|=2;
370 if (TMath::Abs((*fLaser->fFitCside)[0])>kMaxDeltaZ) isReject[1]|=2;
371 // unreasonable drift velocity
372 if (TMath::Abs((*fLaser->fFitAside)[1]-1)>kMaxDeltaV) isReject[0]|=4;
373 if (TMath::Abs((*fLaser->fFitCside)[1]-1)>kMaxDeltaV) isReject[1]|=4;
375 if (TMath::Sqrt(TMath::Abs((*fLaser->fFitAside)[4]))>kMaxRMS ) isReject[0]|=8;
376 if (TMath::Sqrt(TMath::Abs((*fLaser->fFitCside)[4]))>kMaxRMS ) isReject[1]|=8;
382 printf("Trigger: %s\n",event->GetFiredTriggerClasses().Data());
384 TTreeSRedirector *cstream = GetDebugStreamer();
386 TTimeStamp tstamp(fTime);
387 Float_t valuePressure0 = AliTPCcalibDB::GetPressure(tstamp,fRun,0);
388 Float_t valuePressure1 = AliTPCcalibDB::GetPressure(tstamp,fRun,1);
389 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
390 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
391 Double_t temp0 = AliTPCcalibDB::GetTemperature(tstamp,fRun,0);
392 Double_t temp1 = AliTPCcalibDB::GetTemperature(tstamp,fRun,1);
393 Double_t vdcorr = AliTPCcalibDB::Instance()->GetVDriftCorrectionTime(tstamp,fRun,0,1);
394 TVectorD vecGoofie(20);
395 AliDCSSensorArray* goofieArray = AliTPCcalibDB::Instance()->GetGoofieSensors(fRun);
397 for (Int_t isensor=0; isensor<goofieArray->NumSensors();isensor++){
398 AliDCSSensor *gsensor = goofieArray->GetSensor(isensor);
399 if (gsensor) vecGoofie[isensor]=gsensor->GetValue(tstamp);
402 (*cstream)<<"laserInfo"<<
403 "run="<<fRun<< // run number
404 "event="<<fEvent<< // event number
405 "time="<<fTime<< // time stamp of event
406 "trigger="<<fTrigger<< // trigger
407 "mag="<<fMagF<< // magnetic field
408 // Environment values
409 "press0="<<valuePressure0<<
410 "press1="<<valuePressure1<<
411 "pt0="<<ptrelative0<<
412 "pt1="<<ptrelative1<<
415 "vecGoofie.="<<&vecGoofie<<
418 "rejectA="<<isReject[0]<<
419 "rejectC="<<isReject[1]<<
420 "laserA.="<<fLaser->fFitAside<<
421 "laserC.="<<fLaser->fFitCside<<
422 "laserAC.="<<fLaser->fFitACside<<
423 "trigger="<<event->GetFiredTriggerClasses()<<
430 TVectorD vdriftA(5), vdriftC(5),vdriftAC(5);
431 vdriftA=*(fLaser->fFitAside);
432 vdriftC=*(fLaser->fFitCside);
433 vdriftAC=*(fLaser->fFitACside);
434 Int_t npointsA=0, npointsC=0;
435 Float_t chi2A=0, chi2C=0;
436 npointsA= TMath::Nint(vdriftA[3]);
438 npointsC= TMath::Nint(vdriftC[3]);
441 TTimeStamp tstamp(fTime);
442 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
443 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
444 Double_t driftA=0, driftC=0;
445 if (vdriftA[1]>1.-kMaxDeltaV) driftA = 1./vdriftA[1]-1.;
446 if (vdriftC[1]>1.-kMaxDeltaV) driftC = 1./vdriftC[1]-1.;
448 Double_t vecDriftLaserA[4]={fTime,(ptrelative0+ptrelative1)/2.0,driftA,event->GetRunNumber()};
449 Double_t vecDriftLaserC[4]={fTime,(ptrelative0+ptrelative1)/2.0,driftC,event->GetRunNumber()};
450 // Double_t vecDrift[4] ={fTime,(ptrelative0+ptrelative1)/2.0,1./((*(fLaser->fFitACside))[1])-1,event->GetRunNumber()};
452 for (Int_t icalib=0;icalib<3;icalib++){
453 if (icalib==0){ //z0 shift
454 vecDriftLaserA[2]=vdriftA[0]/250.;
455 vecDriftLaserC[2]=vdriftC[0]/250.;
457 if (icalib==1){ //vdrel shift
458 vecDriftLaserA[2]=driftA;
459 vecDriftLaserC[2]=driftC;
461 if (icalib==2){ //gy shift - full gy - full drift
462 vecDriftLaserA[2]=vdriftA[2]/250.;
463 vecDriftLaserC[2]=vdriftC[2]/250.;
465 if (isReject[0]==0) fHistVdriftLaserA[icalib]->Fill(vecDriftLaserA);
466 if (isReject[1]==0) fHistVdriftLaserC[icalib]->Fill(vecDriftLaserC);
469 // THnSparse* curHist=new THnSparseF("","HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
470 // TString shortName=curHist->ClassName();
471 // shortName+="_MEAN_DRIFT_LASER_";
477 // name+=event->GetFiredTriggerClasses();
479 // curHist=(THnSparseF*)fArrayDz->FindObject(name);
481 // curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
482 // fArrayDz->AddLast(curHist);
484 // curHist->Fill(vecDrift);
489 // curHist=(THnSparseF*)fArrayDz->FindObject(name);
491 // curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
492 // fArrayDz->AddLast(curHist);
494 // curHist->Fill(vecDrift);
497 void AliTPCcalibTime::ProcessCosmic(AliESDEvent *event){
499 Printf("ERROR: ESD not available");
502 if (event->GetTimeStamp() == 0 ) {
503 Printf("no time stamp!");
510 // Track0 is choosen in upper TPC part
511 // Track1 is choosen in lower TPC part
513 const Int_t kMinClustersCross =30;
514 const Int_t kMinClusters =80;
515 Int_t ntracks=event->GetNumberOfTracks();
516 if (ntracks==0) return;
517 if (ntracks > fCutTracks) return;
519 if (GetDebugLevel()>20) printf("Hallo world: Im here\n");
520 AliESDfriend *ESDfriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
522 TObjArray tpcSeeds(ntracks);
523 Double_t vtxx[3]={0,0,0};
524 Double_t svtxx[3]={0.000001,0.000001,100.};
525 AliESDVertex vtx(vtxx,svtxx);
529 TArrayI clusterSideA(ntracks);
530 TArrayI clusterSideC(ntracks);
531 for (Int_t i=0;i<ntracks;++i) {
534 AliESDtrack *track = event->GetTrack(i);
536 const AliExternalTrackParam * trackIn = track->GetInnerParam();
537 const AliExternalTrackParam * trackOut = track->GetOuterParam();
538 if (!trackIn) continue;
539 if (!trackOut) continue;
541 AliESDfriendTrack *friendTrack = ESDfriend->GetTrack(i);
542 if (friendTrack) ProcessSame(track,friendTrack,event);
543 if (friendTrack) ProcessAlignITS(track,friendTrack,event,ESDfriend);
544 if (friendTrack) ProcessAlignTRD(track,friendTrack);
545 if (friendTrack) ProcessAlignTOF(track,friendTrack);
546 TObject *calibObject;
547 AliTPCseed *seed = 0;
548 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
550 tpcSeeds.AddAt(seed,i);
552 for (Int_t irow=159;irow>0;irow--) {
553 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
555 if ((cl->GetDetector()%36)<18) nA++;
556 if ((cl->GetDetector()%36)>=18) nC++;
562 if (ntracks<2) return;
567 for (Int_t i=0;i<ntracks;++i) {
568 AliESDtrack *track0 = event->GetTrack(i);
569 // track0 - choosen upper part
570 if (!track0) continue;
571 if (!track0->GetOuterParam()) continue;
572 if (track0->GetOuterParam()->GetAlpha()<0) continue;
574 track0->GetDirection(d1);
575 for (Int_t j=0;j<ntracks;++j) {
577 AliESDtrack *track1 = event->GetTrack(j);
579 if (!track1) continue;
580 if (!track1->GetOuterParam()) continue;
581 if (track0->GetTPCNcls()+ track1->GetTPCNcls()< kMinClusters) continue;
582 Int_t nAC = TMath::Max( TMath::Min(clusterSideA[i], clusterSideC[j]),
583 TMath::Min(clusterSideC[i], clusterSideA[j]));
584 if (nAC<kMinClustersCross) continue;
585 Int_t nA0=clusterSideA[i];
586 Int_t nC0=clusterSideC[i];
587 Int_t nA1=clusterSideA[j];
588 Int_t nC1=clusterSideC[j];
589 // if (track1->GetOuterParam()->GetAlpha()>0) continue;
592 track1->GetDirection(d2);
594 AliTPCseed * seed0 = (AliTPCseed*) tpcSeeds.At(i);
595 AliTPCseed * seed1 = (AliTPCseed*) tpcSeeds.At(j);
596 if (! seed0) continue;
597 if (! seed1) continue;
598 Float_t dir = (d1[0]*d2[0] + d1[1]*d2[1] + d1[2]*d2[2]);
599 Float_t dist0 = track0->GetLinearD(0,0);
600 Float_t dist1 = track1->GetLinearD(0,0);
602 // conservative cuts - convergence to be guarantied
603 // applying before track propagation
604 if (TMath::Abs(TMath::Abs(dist0)-TMath::Abs(dist1))>fCutMaxD) continue; // distance to the 0,0
605 if (TMath::Abs(dir)<TMath::Abs(fCutMinDir)) continue; // direction vector product
606 Float_t bz = AliTracker::GetBz();
607 Float_t dvertex0[2]; //distance to 0,0
608 Float_t dvertex1[2]; //distance to 0,0
609 track0->GetDZ(0,0,0,bz,dvertex0);
610 track1->GetDZ(0,0,0,bz,dvertex1);
611 if (TMath::Abs(dvertex0[1])>250) continue;
612 if (TMath::Abs(dvertex1[1])>250) continue;
616 Float_t dmax = TMath::Max(TMath::Abs(dist0),TMath::Abs(dist1));
617 AliExternalTrackParam param0(*track0);
618 AliExternalTrackParam param1(*track1);
620 // Propagate using Magnetic field and correct fo material budget
622 AliTracker::PropagateTrackTo(¶m0,dmax+1,0.0005,3,kTRUE);
623 AliTracker::PropagateTrackTo(¶m1,dmax+1,0.0005,3,kTRUE);
625 // Propagate rest to the 0,0 DCA - z should be ignored
628 param0.PropagateToDCA(&vtx,bz,1000);
630 param1.PropagateToDCA(&vtx,bz,1000);
631 param0.GetDZ(0,0,0,bz,dvertex0);
632 param1.GetDZ(0,0,0,bz,dvertex1);
637 Bool_t isPair = IsPair(¶m0,¶m1);
638 Bool_t isCross = IsCross(track0, track1);
639 Bool_t isSame = IsSame(track0, track1);
641 THnSparse* hist=new THnSparseF("","HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
642 TString shortName=hist->ClassName();
643 shortName+="_MEAN_VDRIFT_COSMICS_";
647 if((isSame) || (isCross && isPair)){
648 if (track0->GetTPCNcls()+ track1->GetTPCNcls()> 80) {
649 fDz = param0.GetZ() - param1.GetZ();
650 Double_t sign=(nA0>nA1)? 1:-1;
652 TTimeStamp tstamp(fTime);
653 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
654 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
655 Double_t vecDrift[4]={fTime,(ptrelative0+ptrelative1)/2.0,fDz/500.0,event->GetRunNumber()};
656 THnSparse* curHist=NULL;
660 name+=event->GetFiredTriggerClasses();
662 curHist=(THnSparseF*)fArrayDz->FindObject(name);
664 curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
665 fArrayDz->AddLast(curHist);
667 // curHist=(THnSparseF*)(fMapDz->GetValue(event->GetFiredTriggerClasses()));
669 // curHist=new THnSparseF(event->GetFiredTriggerClasses(),"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
670 // fMapDz->Add(new TObjString(event->GetFiredTriggerClasses()),curHist);
672 curHist->Fill(vecDrift);
677 curHist=(THnSparseF*)fArrayDz->FindObject(name);
679 curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
680 fArrayDz->AddLast(curHist);
682 // curHist=(THnSparseF*)(fMapDz->GetValue("all"));
684 // curHist=new THnSparseF("all","HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
685 // fMapDz->Add(new TObjString("all"),curHist);
687 curHist->Fill(vecDrift);
690 TTreeSRedirector *cstream = GetDebugStreamer();
693 (*cstream)<<"trackInfo"<<
704 "isCross="<<isCross<<
712 } // end 2nd order loop
713 } // end 1st order loop
716 TTreeSRedirector *cstream = GetDebugStreamer();
718 TTimeStamp tstamp(fTime);
719 Float_t valuePressure0 = AliTPCcalibDB::GetPressure(tstamp,fRun,0);
720 Float_t valuePressure1 = AliTPCcalibDB::GetPressure(tstamp,fRun,1);
721 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
722 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
723 Double_t temp0 = AliTPCcalibDB::GetTemperature(tstamp,fRun,0);
724 Double_t temp1 = AliTPCcalibDB::GetTemperature(tstamp,fRun,1);
725 Double_t vdcorr = AliTPCcalibDB::Instance()->GetVDriftCorrectionTime(tstamp,fRun,0,1);
726 TVectorD vecGoofie(20);
727 AliDCSSensorArray* goofieArray = AliTPCcalibDB::Instance()->GetGoofieSensors(fRun);
729 for (Int_t isensor=0; isensor<goofieArray->NumSensors();isensor++){
730 AliDCSSensor *gsensor = goofieArray->GetSensor(isensor);
731 if (gsensor) vecGoofie[isensor]=gsensor->GetValue(tstamp);
734 (*cstream)<<"timeInfo"<<
735 "run="<<fRun<< // run number
736 "event="<<fEvent<< // event number
737 "time="<<fTime<< // time stamp of event
738 "trigger="<<fTrigger<< // trigger
739 "mag="<<fMagF<< // magnetic field
740 // Environment values
741 "press0="<<valuePressure0<<
742 "press1="<<valuePressure1<<
743 "pt0="<<ptrelative0<<
744 "pt1="<<ptrelative1<<
747 "vecGoofie.=<<"<<&vecGoofie<<
750 // accumulated values
752 "fDz="<<fDz<< //! current delta z
753 "trigger="<<event->GetFiredTriggerClasses()<<
757 if (GetDebugLevel()>20) printf("Trigger: %s\n",event->GetFiredTriggerClasses().Data());
760 void AliTPCcalibTime::ProcessBeam(AliESDEvent */*event*/){
763 void AliTPCcalibTime::Analyze(){}
765 THnSparse* AliTPCcalibTime::GetHistoDrift(const char* name){
766 TIterator* iterator = fArrayDz->MakeIterator();
768 TString newName=name;
770 THnSparse* newHist=new THnSparseF(newName,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
771 THnSparse* addHist=NULL;
772 while((addHist=(THnSparseF*)iterator->Next())){
773 if(!addHist) continue;
774 TString histName=addHist->GetName();
775 if(!histName.Contains(newName)) continue;
777 newHist->Add(addHist);
782 TObjArray* AliTPCcalibTime::GetHistoDrift(){
786 TGraphErrors* AliTPCcalibTime::GetGraphDrift(const char* name){
787 THnSparse* histoDrift=GetHistoDrift(name);
788 TGraphErrors* graphDrift=NULL;
790 graphDrift=FitSlices(histoDrift,2,0,400,100,0.05,0.95, kTRUE);
791 TString end=histoDrift->GetName();
792 Int_t pos=end.Index("_");
793 end=end(pos,end.Capacity()-pos);
794 TString graphName=graphDrift->ClassName();
797 graphDrift->SetName(graphName);
802 TObjArray* AliTPCcalibTime::GetGraphDrift(){
803 TObjArray* arrayGraphDrift=new TObjArray();
804 TIterator* iterator=fArrayDz->MakeIterator();
806 THnSparse* addHist=NULL;
807 while((addHist=(THnSparseF*)iterator->Next())) arrayGraphDrift->AddLast(GetGraphDrift(addHist->GetName()));
808 return arrayGraphDrift;
811 AliSplineFit* AliTPCcalibTime::GetFitDrift(const char* name){
812 TGraph* graphDrift=GetGraphDrift(name);
813 AliSplineFit* fitDrift=NULL;
814 if(graphDrift && graphDrift->GetN()){
815 fitDrift=new AliSplineFit();
816 fitDrift->SetGraph(graphDrift);
817 fitDrift->SetMinPoints(graphDrift->GetN()+1);
818 fitDrift->InitKnots(graphDrift,2,0,0.001);
819 fitDrift->SplineFit(0);
820 TString end=graphDrift->GetName();
821 Int_t pos=end.Index("_");
822 end=end(pos,end.Capacity()-pos);
823 TString fitName=fitDrift->ClassName();
826 //fitDrift->SetName(fitName);
833 //TObjArray* AliTPCcalibTime::GetFitDrift(){
834 // TObjArray* arrayFitDrift=new TObjArray();
835 // TIterator* iterator = fArrayDz->MakeIterator();
836 // iterator->Reset();
837 // THnSparse* addHist=NULL;
838 // while((addHist=(THnSparseF*)iterator->Next())) arrayFitDrift->AddLast(GetFitDrift(addHist->GetName()));
839 // return arrayFitDrift;
842 Long64_t AliTPCcalibTime::Merge(TCollection *li) {
843 TIterator* iter = li->MakeIterator();
844 AliTPCcalibTime* cal = 0;
846 while ((cal = (AliTPCcalibTime*)iter->Next())) {
847 if (!cal->InheritsFrom(AliTPCcalibTime::Class())) {
848 Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName());
851 for (Int_t imeas=0; imeas<3; imeas++){
852 if (cal->GetHistVdriftLaserA(imeas) && cal->GetHistVdriftLaserA(imeas)){
853 fHistVdriftLaserA[imeas]->Add(cal->GetHistVdriftLaserA(imeas));
854 fHistVdriftLaserC[imeas]->Add(cal->GetHistVdriftLaserC(imeas));
857 TObjArray* addArray=cal->GetHistoDrift();
858 if(!addArray) return 0;
859 TIterator* iterator = addArray->MakeIterator();
861 THnSparse* addHist=NULL;
862 while((addHist=(THnSparseF*)iterator->Next())){
863 if(!addHist) continue;
865 THnSparse* localHist=(THnSparseF*)fArrayDz->FindObject(addHist->GetName());
867 localHist=new THnSparseF(addHist->GetName(),"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
868 fArrayDz->AddLast(localHist);
870 localHist->Add(addHist);
872 // TMap * addMap=cal->GetHistoDrift();
873 // if(!addMap) return 0;
874 // TIterator* iterator = addMap->MakeIterator();
875 // iterator->Reset();
877 // while((addPair=(TPair *)(addMap->FindObject(iterator->Next())))){
878 // THnSparse* addHist=dynamic_cast<THnSparseF*>(addPair->Value());
879 // if (!addHist) continue;
881 // THnSparse* localHist=dynamic_cast<THnSparseF*>(fMapDz->GetValue(addHist->GetName()));
883 // localHist=new THnSparseF(addHist->GetName(),"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
884 // fMapDz->Add(new TObjString(addHist->GetName()),localHist);
886 // localHist->Add(addHist);
888 for(Int_t i=0;i<10;i++) if (cal->GetCosmiMatchingHisto(i)) fCosmiMatchingHisto[i]->Add(cal->GetCosmiMatchingHisto(i));
892 for (Int_t itype=0; itype<3; itype++){
897 if (itype==0) {arr0=fAlignITSTPC; arr1=cal->fAlignITSTPC;}
898 if (itype==1) {arr0=fAlignTRDTPC; arr1=cal->fAlignTRDTPC;}
899 if (itype==2) {arr0=fAlignTOFTPC; arr1=cal->fAlignTOFTPC;}
901 if (!arr0) arr0=new TObjArray(arr1->GetEntriesFast());
902 if (arr1->GetEntriesFast()>arr0->GetEntriesFast()){
903 arr0->Expand(arr1->GetEntriesFast());
905 for (Int_t i=0;i<arr1->GetEntriesFast(); i++){
906 AliRelAlignerKalman *kalman1 = (AliRelAlignerKalman *)arr1->UncheckedAt(i);
907 AliRelAlignerKalman *kalman0 = (AliRelAlignerKalman *)arr0->UncheckedAt(i);
908 if (!kalman1) continue;
909 if (!kalman0) {arr0->AddAt(new AliRelAlignerKalman(*kalman1),i); continue;}
910 kalman0->SetRejectOutliers(kFALSE);
911 kalman0->Merge(kalman1);
919 Bool_t AliTPCcalibTime::IsPair(AliExternalTrackParam *tr0, AliExternalTrackParam *tr1){
921 // 0. Same direction - OPOSITE - cutDir +cutT
922 TCut cutDir("cutDir","dir<-0.99")
924 TCut cutT("cutT","abs(Tr1.fP[3]+Tr0.fP[3])<0.03")
927 TCut cutD("cutD","abs(Tr0.fP[0]+Tr1.fP[0])<5")
929 TCut cutPt("cutPt","abs(Tr1.fP[4]+Tr0.fP[4])<1&&abs(Tr0.fP[4])+abs(Tr1.fP[4])<10");
932 const Double_t *p0 = tr0->GetParameter();
933 const Double_t *p1 = tr1->GetParameter();
934 fCosmiMatchingHisto[0]->Fill(p0[0]+p1[0]);
935 fCosmiMatchingHisto[1]->Fill(p0[1]-p1[1]);
936 fCosmiMatchingHisto[2]->Fill(tr1->GetAlpha()-tr0->GetAlpha()+TMath::Pi());
937 fCosmiMatchingHisto[3]->Fill(p0[3]+p1[3]);
938 fCosmiMatchingHisto[4]->Fill(p0[4]+p1[4]);
940 if (TMath::Abs(p0[3]+p1[3])>fCutTheta) return kFALSE;
941 if (TMath::Abs(p0[0]+p1[0])>fCutMaxD) return kFALSE;
942 if (TMath::Abs(p0[1]-p1[1])>fCutMaxDz) return kFALSE;
943 Double_t d0[3], d1[3];
944 tr0->GetDirection(d0);
945 tr1->GetDirection(d1);
946 if (d0[0]*d1[0] + d0[1]*d1[1] + d0[2]*d1[2] >fCutMinDir) return kFALSE;
948 fCosmiMatchingHisto[5]->Fill(p0[0]+p1[0]);
949 fCosmiMatchingHisto[6]->Fill(p0[1]-p1[1]);
950 fCosmiMatchingHisto[7]->Fill(tr1->GetAlpha()-tr0->GetAlpha()+TMath::Pi());
951 fCosmiMatchingHisto[8]->Fill(p0[3]+p1[3]);
952 fCosmiMatchingHisto[9]->Fill(p0[4]+p1[4]);
956 Bool_t AliTPCcalibTime::IsCross(AliESDtrack *tr0, AliESDtrack *tr1){
958 // check if the cosmic pair of tracks crossed A/C side
960 Bool_t result= tr0->GetOuterParam()->GetZ()*tr1->GetOuterParam()->GetZ()<0;
961 if (result==kFALSE) return result;
966 Bool_t AliTPCcalibTime::IsSame(AliESDtrack *tr0, AliESDtrack *tr1){
968 // track crossing the CE
969 // 0. minimal number of clusters
970 // 1. Same sector +-1
971 // 2. Inner and outer track param on opposite side
972 // 3. Outer and inner track parameter close each to other
976 // inner and outer on opposite sides in z
978 const Int_t knclCut0 = 30;
979 const Double_t kalphaCut = 0.4;
981 // 0. minimal number of clusters
983 if (tr0->GetTPCNcls()<knclCut0) return kFALSE;
984 if (tr1->GetTPCNcls()<knclCut0) return kFALSE;
986 // 1. alpha cut - sector+-1
988 if (TMath::Abs(tr0->GetOuterParam()->GetAlpha()-tr1->GetOuterParam()->GetAlpha())>kalphaCut) return kFALSE;
992 if (tr0->GetOuterParam()->GetZ()*tr0->GetInnerParam()->GetZ()>0) result&=kFALSE;
993 if (tr1->GetOuterParam()->GetZ()*tr1->GetInnerParam()->GetZ()>0) result&=kFALSE;
999 const Double_t *p0I = tr0->GetInnerParam()->GetParameter();
1000 const Double_t *p1I = tr1->GetInnerParam()->GetParameter();
1001 const Double_t *p0O = tr0->GetOuterParam()->GetParameter();
1002 const Double_t *p1O = tr1->GetOuterParam()->GetParameter();
1004 if (TMath::Abs(p0I[0]-p1I[0])>fCutMaxD) result&=kFALSE;
1005 if (TMath::Abs(p0I[1]-p1I[1])>fCutMaxDz) result&=kFALSE;
1006 if (TMath::Abs(p0I[2]-p1I[2])>fCutTheta) result&=kFALSE;
1007 if (TMath::Abs(p0I[3]-p1I[3])>fCutTheta) result&=kFALSE;
1008 if (TMath::Abs(p0O[0]-p1O[0])>fCutMaxD) result&=kFALSE;
1009 if (TMath::Abs(p0O[1]-p1O[1])>fCutMaxDz) result&=kFALSE;
1010 if (TMath::Abs(p0O[2]-p1O[2])>fCutTheta) result&=kFALSE;
1011 if (TMath::Abs(p0O[3]-p1O[3])>fCutTheta) result&=kFALSE;
1013 result=kTRUE; // just to put break point here
1019 void AliTPCcalibTime::ProcessSame(AliESDtrack* track, AliESDfriendTrack *friendTrack,AliESDEvent *event){
1021 // Process TPC tracks crossing CE
1023 // 0. Select only track crossing the CE
1024 // 1. Cut on the track length
1025 // 2. Refit the terack on A and C side separatelly
1026 // 3. Fill time histograms
1027 const Int_t kMinNcl=100;
1028 const Int_t kMinNclS=25; // minimul number of clusters on the sides
1029 if (!friendTrack->GetTPCOut()) return;
1031 // 0. Select only track crossing the CE
1033 if (track->GetInnerParam()->GetZ()*friendTrack->GetTPCOut()->GetZ()>0) return;
1035 // 1. cut on track length
1037 if (track->GetTPCNcls()<kMinNcl) return;
1039 // 2. Refit track sepparatel on A and C side
1041 TObject *calibObject;
1042 AliTPCseed *seed = 0;
1043 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
1044 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
1048 AliExternalTrackParam trackIn(*track->GetInnerParam());
1049 AliExternalTrackParam trackOut(*track->GetOuterParam());
1050 Double_t cov[3]={0.01,0.,0.01}; //use the same errors
1051 Double_t xyz[3]={0,0.,0.0};
1053 Int_t nclIn=0,nclOut=0;
1054 trackIn.ResetCovariance(30.);
1055 trackOut.ResetCovariance(30.);
1059 for (Int_t irow=0;irow<159;irow++) {
1060 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
1062 if (cl->GetX()<80) continue;
1063 if (track->GetInnerParam()->GetZ()<0 &&(cl->GetDetector()%36)<18) break;
1064 if (track->GetInnerParam()->GetZ()>0 &&(cl->GetDetector()%36)>=18) break;
1065 Int_t sector = cl->GetDetector();
1066 Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackIn.GetAlpha();
1067 if (TMath::Abs(dalpha)>0.01){
1068 if (!trackIn.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
1070 Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
1071 trackIn.GetXYZ(xyz);
1072 bz = AliTracker::GetBz(xyz);
1073 if (!trackIn.PropagateTo(r[0],bz)) break;
1075 trackIn.Update(&r[1],cov);
1080 for (Int_t irow=159;irow>0;irow--) {
1081 AliTPCclusterMI *cl=seed->GetClusterPointer(irow);
1083 if (cl->GetX()<80) continue;
1084 if (cl->GetZ()*track->GetOuterParam()->GetZ()<0) break;
1085 if (friendTrack->GetTPCOut()->GetZ()<0 &&(cl->GetDetector()%36)<18) break;
1086 if (friendTrack->GetTPCOut()->GetZ()>0 &&(cl->GetDetector()%36)>=18) break;
1087 Int_t sector = cl->GetDetector();
1088 Float_t dalpha = TMath::DegToRad()*(sector%18*20.+10.)-trackOut.GetAlpha();
1089 if (TMath::Abs(dalpha)>0.01){
1090 if (!trackOut.Rotate(TMath::DegToRad()*(sector%18*20.+10.))) break;
1092 Double_t r[3]={cl->GetX(),cl->GetY(),cl->GetZ()};
1093 trackOut.GetXYZ(xyz);
1094 bz = AliTracker::GetBz(xyz);
1095 if (!trackOut.PropagateTo(r[0],bz)) break;
1097 trackOut.Update(&r[1],cov);
1099 trackOut.Rotate(trackIn.GetAlpha());
1100 Double_t meanX = (trackIn.GetX()+trackOut.GetX())*0.5;
1101 trackIn.PropagateTo(meanX,bz);
1102 trackOut.PropagateTo(meanX,bz);
1103 TTreeSRedirector *cstream = GetDebugStreamer();
1106 trackIn.GetXYZ(gxyz.GetMatrixArray());
1107 TTimeStamp tstamp(fTime);
1108 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
1109 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
1110 Double_t vdcorr = AliTPCcalibDB::Instance()->GetVDriftCorrectionTime(tstamp,fRun,0,1);
1111 (*cstream)<<"tpctpc"<<
1112 "run="<<fRun<< // run number
1113 "event="<<fEvent<< // event number
1114 "time="<<fTime<< // time stamp of event
1115 "trigger="<<fTrigger<< // trigger
1116 "mag="<<fMagF<< // magnetic field
1117 "ptrel0.="<<ptrelative0<<
1118 "ptrel1.="<<ptrelative1<<
1119 "vdcorr="<<vdcorr<< // drift correction applied
1121 "xyz.="<<&gxyz<< // global position
1122 "tIn.="<<&trackIn<< // refitterd track in
1123 "tOut.="<<&trackOut<< // refitter track out
1124 "nclIn="<<nclIn<< //
1125 "nclOut="<<nclOut<< //
1129 // 3. Fill time histograms
1130 // Debug stremaer expression
1131 // chainTPCTPC->Draw("(tIn.fP[1]-tOut.fP[1])*sign(-tIn.fP[3]):tIn.fP[3]","min(nclIn,nclOut)>30","")
1132 if (TMath::Min(nclIn,nclOut)>kMinNclS){
1133 fDz = trackOut.GetZ()-trackIn.GetZ();
1134 if (trackOut.GetTgl()<0) fDz*=-1.;
1135 TTimeStamp tstamp(fTime);
1136 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
1137 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
1138 Double_t vecDrift[4]={fTime,(ptrelative0+ptrelative1)/2.0,fDz/500.0,event->GetRunNumber()};
1140 // fill histograms per trigger class and itegrated
1142 THnSparse* curHist=NULL;
1143 for (Int_t itype=0; itype<2; itype++){
1144 TString name="MEAN_VDRIFT_CROSS_";
1146 name+=event->GetFiredTriggerClasses();
1151 curHist=(THnSparseF*)fArrayDz->FindObject(name);
1153 curHist=new THnSparseF(name,"HistVdrift;time;p/T ratio;Vdrift;run",4,fBinsVdrift,fXminVdrift,fXmaxVdrift);
1154 fArrayDz->AddLast(curHist);
1156 curHist->Fill(vecDrift);
1162 void AliTPCcalibTime::ProcessAlignITS(AliESDtrack* track, AliESDfriendTrack *friendTrack, AliESDEvent *event,AliESDfriend *ESDfriend){
1164 // Process track - Update TPC-ITS alignment
1166 // 0. Apply standartd cuts
1167 // 1. Recalucluate the current statistic median/RMS
1168 // 2. Apply median+-rms cut
1169 // 3. Update kalman filter
1171 const Int_t kMinTPC = 80; // minimal number of TPC cluster
1172 const Int_t kMinITS = 3; // minimal number of ITS cluster
1173 const Double_t kMinZ = 10; // maximal dz distance
1174 const Double_t kMaxDy = 2.; // maximal dy distance
1175 const Double_t kMaxAngle= 0.015; // maximal angular distance
1176 const Double_t kSigmaCut= 5; // maximal sigma distance to median
1177 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
1178 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1179 const Double_t kOutCut = 1.0; // outlyer cut in AliRelAlgnmentKalman
1180 const Double_t kMinPt = 0.3; // minimal pt
1181 const Int_t kN=500; // deepnes of history
1182 static Int_t kglast=0;
1183 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1186 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";
1189 // 0. Apply standard cuts
1191 Int_t dummycl[1000];
1192 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
1193 if (track->GetITSclusters(dummycl)<kMinITS) return; // minimal amount of clusters
1194 if (!track->IsOn(AliESDtrack::kTPCrefit)) return;
1195 if (!friendTrack->GetITSOut()) return;
1196 if (!track->GetInnerParam()) return;
1197 if (!track->GetOuterParam()) return;
1198 if (track->GetInnerParam()->Pt()<kMinPt) return;
1199 // exclude crossing track
1200 if (track->GetOuterParam()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
1201 if (TMath::Abs(track->GetInnerParam()->GetZ())<kMinZ) return;
1202 if (track->GetInnerParam()->GetX()>90) return;
1204 AliExternalTrackParam &pTPC=(AliExternalTrackParam &)(*(track->GetInnerParam()));
1205 AliExternalTrackParam pITS(*(friendTrack->GetITSOut())); // ITS standalone if possible
1206 AliExternalTrackParam pITS2(*(friendTrack->GetITSOut())); //TPC-ITS track
1207 pITS2.Rotate(pTPC.GetAlpha());
1208 pITS2.PropagateTo(pTPC.GetX(),fMagF);
1209 AliESDfriendTrack *itsfriendTrack=0;
1211 // try to find standalone ITS track corresponing to the TPC if possible
1213 Bool_t hasAlone=kFALSE;
1214 Int_t ntracks=event->GetNumberOfTracks();
1215 for (Int_t i=0; i<ntracks; i++){
1216 AliESDtrack *trackS = event->GetTrack(i);
1217 if (trackS->GetTPCNcls()>0) continue; //continue if has TPC info
1218 itsfriendTrack = ESDfriend->GetTrack(i);
1219 if (!itsfriendTrack) continue;
1220 if (!itsfriendTrack->GetITSOut()) continue;
1221 if (TMath::Abs(pITS2.GetTgl()-itsfriendTrack->GetITSOut()->GetTgl())> kMaxAngle) continue;
1222 pITS=(*(itsfriendTrack->GetITSOut()));
1224 pITS.Rotate(pTPC.GetAlpha());
1225 pITS.PropagateTo(pTPC.GetX(),fMagF);
1226 if (TMath::Abs(pITS2.GetY()-pITS.GetY())> kMaxDy) continue;
1229 if (!hasAlone) pITS=pITS2;
1231 if (TMath::Abs(pITS.GetY()-pTPC.GetY()) >kMaxDy) return;
1232 if (TMath::Abs(pITS.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1233 if (TMath::Abs(pITS.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1235 // 1. Update median and RMS info
1237 TVectorD vecDelta(4),vecMedian(4), vecRMS(4);
1238 TVectorD vecDeltaN(5);
1239 Double_t sign=(pITS.GetParameter()[1]>0)? 1.:-1.;
1241 for (Int_t i=0;i<4;i++){
1242 vecDelta[i]=(pITS.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1243 kgdP[i][kglast%kN]=vecDelta[i];
1246 Int_t entries=(kglast<kN)?kglast:kN;
1247 for (Int_t i=0;i<4;i++){
1248 vecMedian[i] = TMath::Median(entries,kgdP[i]);
1249 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1252 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/vecRMS[i];
1253 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1257 // 2. Apply median+-rms cut
1259 if (kglast<3) return; //median and RMS to be defined
1260 if ( vecDeltaN[4]/4.>kSigmaCut) return;
1262 // 3. Update alignment
1264 Int_t htime = fTime/3600; //time in hours
1265 if (fAlignITSTPC->GetEntries()<htime){
1266 fAlignITSTPC->Expand(htime*2+20);
1268 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignITSTPC->At(htime);
1270 // make Alignment object if doesn't exist
1271 align=new AliRelAlignerKalman();
1272 align->SetRunNumber(fRun);
1273 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1274 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1275 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1276 align->SetRejectOutliers(kFALSE);
1278 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1279 align->SetMagField(fMagF);
1280 fAlignITSTPC->AddAt(align,htime);
1282 align->AddTrackParams(&pITS,&pTPC);
1283 align->SetTimeStamp(fTime);
1284 align->SetRunNumber(fRun );
1286 Int_t nupdates=align->GetNUpdates();
1287 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1288 align->SetRejectOutliers(kFALSE);
1289 TTreeSRedirector *cstream = GetDebugStreamer();
1290 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1291 TTimeStamp tstamp(fTime);
1292 Float_t valuePressure0 = AliTPCcalibDB::GetPressure(tstamp,fRun,0);
1293 Float_t valuePressure1 = AliTPCcalibDB::GetPressure(tstamp,fRun,1);
1294 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
1295 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
1296 Double_t temp0 = AliTPCcalibDB::GetTemperature(tstamp,fRun,0);
1297 Double_t temp1 = AliTPCcalibDB::GetTemperature(tstamp,fRun,1);
1298 TVectorD vecGoofie(20);
1299 AliDCSSensorArray* goofieArray = AliTPCcalibDB::Instance()->GetGoofieSensors(fRun);
1301 for (Int_t isensor=0; isensor<goofieArray->NumSensors();isensor++){
1302 AliDCSSensor *gsensor = goofieArray->GetSensor(isensor);
1303 if (gsensor) vecGoofie[isensor]=gsensor->GetValue(tstamp);
1306 TVectorD gpTPC(3), gdTPC(3);
1307 TVectorD gpITS(3), gdITS(3);
1308 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1309 pTPC.GetDirection(gdTPC.GetMatrixArray());
1310 pITS.GetXYZ(gpITS.GetMatrixArray());
1311 pITS.GetDirection(gdITS.GetMatrixArray());
1312 Double_t vdcorr = AliTPCcalibDB::Instance()->GetVDriftCorrectionTime(tstamp,fRun,0,1);
1313 (*cstream)<<"itstpc"<<
1314 "run="<<fRun<< // run number
1315 "event="<<fEvent<< // event number
1316 "time="<<fTime<< // time stamp of event
1317 "trigger="<<fTrigger<< // trigger
1318 "mag="<<fMagF<< // magnetic field
1319 // Environment values
1320 "press0="<<valuePressure0<<
1321 "press1="<<valuePressure1<<
1322 "pt0="<<ptrelative0<<
1323 "pt1="<<ptrelative1<<
1326 "vecGoofie.="<<&vecGoofie<<
1327 "vdcorr="<<vdcorr<< // drift correction applied
1329 "hasAlone="<<hasAlone<< // has ITS standalone ?
1330 "track.="<<track<< // track info
1331 "nmed="<<kglast<< // number of entries to define median and RMS
1332 "vMed.="<<&vecMedian<< // median of deltas
1333 "vRMS.="<<&vecRMS<< // rms of deltas
1334 "vDelta.="<<&vecDelta<< // delta in respect to median
1335 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1336 "t.="<<track<< // ful track - find proper cuts
1337 "a.="<<align<< // current alignment
1338 "pITS.="<<&pITS<< // track param ITS
1339 "pITS2.="<<&pITS2<< // track param ITS+TPC
1340 "pTPC.="<<&pTPC<< // track param TPC
1341 "gpTPC.="<<&gpTPC<< // global position TPC
1342 "gdTPC.="<<&gdTPC<< // global direction TPC
1343 "gpITS.="<<&gpITS<< // global position ITS
1344 "gdITS.="<<&gdITS<< // global position ITS
1352 void AliTPCcalibTime::ProcessAlignTRD(AliESDtrack* track, AliESDfriendTrack *friendTrack){
1354 // Process track - Update TPC-TRD alignment
1356 // 0. Apply standartd cuts
1357 // 1. Recalucluate the current statistic median/RMS
1358 // 2. Apply median+-rms cut
1359 // 3. Update kalman filter
1361 const Int_t kMinTPC = 80; // minimal number of TPC cluster
1362 const Int_t kMinTRD = 50; // minimal number of TRD cluster
1363 const Double_t kMinZ = 20; // maximal dz distance
1364 const Double_t kMaxDy = 2.; // maximal dy distance
1365 const Double_t kMaxAngle= 0.015; // maximal angular distance
1366 const Double_t kSigmaCut= 5; // maximal sigma distance to median
1367 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
1368 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1369 const Double_t kOutCut = 1.0; // outlyer cut in AliRelAlgnmentKalman
1370 const Int_t kN=500; // deepnes of history
1371 static Int_t kglast=0;
1372 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1374 // 0. Apply standard cuts
1376 Int_t dummycl[1000];
1377 if (track->GetTRDclusters(dummycl)<kMinTRD) return; // minimal amount of clusters
1378 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
1379 if (!friendTrack->GetTRDIn()) return;
1380 if (!track->GetInnerParam()) return;
1381 if (!track->GetOuterParam()) return;
1382 // exclude crossing track
1383 if (track->GetOuterParam()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
1384 if (TMath::Abs(track->GetInnerParam()->GetZ())<kMinZ) return;
1386 AliExternalTrackParam &pTPC=(AliExternalTrackParam &)(*(track->GetOuterParam()));
1387 AliExternalTrackParam pTRD(*(friendTrack->GetTRDIn()));
1388 pTRD.Rotate(pTPC.GetAlpha());
1389 pTRD.PropagateTo(pTPC.GetX(),fMagF);
1390 ((Double_t*)pTRD.GetCovariance())[2]+=3.*3.; // increas sys errors
1391 ((Double_t*)pTRD.GetCovariance())[9]+=0.1*0.1; // increse sys errors
1393 if (TMath::Abs(pTRD.GetY()-pTPC.GetY()) >kMaxDy) return;
1394 if (TMath::Abs(pTRD.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1395 if (TMath::Abs(pTRD.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1397 // 1. Update median and RMS info
1399 TVectorD vecDelta(4),vecMedian(4), vecRMS(4);
1400 TVectorD vecDeltaN(5);
1401 Double_t sign=(pTRD.GetParameter()[1]>0)? 1.:-1.;
1403 for (Int_t i=0;i<4;i++){
1404 vecDelta[i]=(pTRD.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1405 kgdP[i][kglast%kN]=vecDelta[i];
1408 Int_t entries=(kglast<kN)?kglast:kN;
1409 for (Int_t i=0;i<4;i++){
1410 vecMedian[i] = TMath::Median(entries,kgdP[i]);
1411 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1414 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/vecRMS[i];
1415 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1419 // 2. Apply median+-rms cut
1421 if (kglast<3) return; //median and RMS to be defined
1422 if ( vecDeltaN[4]/4.>kSigmaCut) return;
1424 // 3. Update alignment
1426 Int_t htime = fTime/3600; //time in hours
1427 if (fAlignTRDTPC->GetEntries()<htime){
1428 fAlignTRDTPC->Expand(htime*2+20);
1430 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignTRDTPC->At(htime);
1432 // make Alignment object if doesn't exist
1433 align=new AliRelAlignerKalman();
1434 align->SetRunNumber(fRun);
1435 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1436 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1437 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1438 align->SetRejectOutliers(kFALSE);
1439 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1440 align->SetMagField(fMagF);
1441 fAlignTRDTPC->AddAt(align,htime);
1443 align->AddTrackParams(&pTRD,&pTPC);
1444 align->SetTimeStamp(fTime);
1445 align->SetRunNumber(fRun );
1447 Int_t nupdates=align->GetNUpdates();
1448 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1449 align->SetRejectOutliers(kFALSE);
1450 TTreeSRedirector *cstream = GetDebugStreamer();
1451 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1452 TTimeStamp tstamp(fTime);
1453 Float_t valuePressure0 = AliTPCcalibDB::GetPressure(tstamp,fRun,0);
1454 Float_t valuePressure1 = AliTPCcalibDB::GetPressure(tstamp,fRun,1);
1455 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
1456 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
1457 Double_t temp0 = AliTPCcalibDB::GetTemperature(tstamp,fRun,0);
1458 Double_t temp1 = AliTPCcalibDB::GetTemperature(tstamp,fRun,1);
1459 TVectorD vecGoofie(20);
1460 AliDCSSensorArray* goofieArray = AliTPCcalibDB::Instance()->GetGoofieSensors(fRun);
1462 for (Int_t isensor=0; isensor<goofieArray->NumSensors();isensor++){
1463 AliDCSSensor *gsensor = goofieArray->GetSensor(isensor);
1464 if (gsensor) vecGoofie[isensor]=gsensor->GetValue(tstamp);
1467 TVectorD gpTPC(3), gdTPC(3);
1468 TVectorD gpTRD(3), gdTRD(3);
1469 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1470 pTPC.GetDirection(gdTPC.GetMatrixArray());
1471 pTRD.GetXYZ(gpTRD.GetMatrixArray());
1472 pTRD.GetDirection(gdTRD.GetMatrixArray());
1473 Double_t vdcorr = AliTPCcalibDB::Instance()->GetVDriftCorrectionTime(tstamp,fRun,0,1);
1474 (*cstream)<<"trdtpc"<<
1475 "run="<<fRun<< // run number
1476 "event="<<fEvent<< // event number
1477 "time="<<fTime<< // time stamp of event
1478 "trigger="<<fTrigger<< // trigger
1479 "mag="<<fMagF<< // magnetic field
1480 // Environment values
1481 "press0="<<valuePressure0<<
1482 "press1="<<valuePressure1<<
1483 "pt0="<<ptrelative0<<
1484 "pt1="<<ptrelative1<<
1487 "vecGoofie.="<<&vecGoofie<<
1488 "vdcorr="<<vdcorr<< // drift correction applied
1490 "nmed="<<kglast<< // number of entries to define median and RMS
1491 "vMed.="<<&vecMedian<< // median of deltas
1492 "vRMS.="<<&vecRMS<< // rms of deltas
1493 "vDelta.="<<&vecDelta<< // delta in respect to median
1494 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1495 "t.="<<track<< // ful track - find proper cuts
1496 "a.="<<align<< // current alignment
1497 "pTRD.="<<&pTRD<< // track param TRD
1498 "pTPC.="<<&pTPC<< // track param TPC
1499 "gpTPC.="<<&gpTPC<< // global position TPC
1500 "gdTPC.="<<&gdTPC<< // global direction TPC
1501 "gpTRD.="<<&gpTRD<< // global position TRD
1502 "gdTRD.="<<&gdTRD<< // global position TRD
1508 void AliTPCcalibTime::ProcessAlignTOF(AliESDtrack* track, AliESDfriendTrack *friendTrack){
1511 // Process track - Update TPC-TOF alignment
1513 // -1. Make a TOF "track"
1514 // 0. Apply standartd cuts
1515 // 1. Recalucluate the current statistic median/RMS
1516 // 2. Apply median+-rms cut
1517 // 3. Update kalman filter
1519 const Int_t kMinTPC = 80; // minimal number of TPC cluster
1520 // const Double_t kMinZ = 10; // maximal dz distance
1521 const Double_t kMaxDy = 5.; // maximal dy distance
1522 const Double_t kMaxAngle= 0.015; // maximal angular distance
1523 const Double_t kSigmaCut= 5; // maximal sigma distance to median
1524 const Double_t kVdErr = 0.1; // initial uncertainty of the vd correction
1525 const Double_t kVdYErr = 0.05; // initial uncertainty of the vd correction
1527 const Double_t kOutCut = 1.0; // outlyer cut in AliRelAlgnmentKalman
1528 const Int_t kN=1000; // deepnes of history
1529 static Int_t kglast=0;
1530 static Double_t* kgdP[4]={new Double_t[kN], new Double_t[kN], new Double_t[kN], new Double_t[kN]};
1532 // -1. Make a TOF track-
1533 // Clusters are not in friends - use alingment points
1535 if (track->GetTOFsignal()<=0) return;
1536 if (!friendTrack->GetTPCOut()) return;
1537 if (!track->GetInnerParam()) return;
1538 if (!track->GetOuterParam()) return;
1539 const AliTrackPointArray *points=friendTrack->GetTrackPointArray();
1540 if (!points) return;
1541 AliExternalTrackParam pTPC(*(track->GetOuterParam()));
1542 AliExternalTrackParam pTOF(pTPC);
1543 Double_t mass = TDatabasePDG::Instance()->GetParticle("mu+")->Mass();
1544 Int_t npoints = points->GetNPoints();
1545 AliTrackPoint point;
1548 for (Int_t ipoint=0;ipoint<npoints;ipoint++){
1549 points->GetPoint(point,ipoint);
1552 Double_t r=TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]);
1553 if (r<350) continue;
1554 if (r>400) continue;
1555 AliTracker::PropagateTrackToBxByBz(&pTPC,r,mass,2.,kTRUE);
1556 AliTracker::PropagateTrackToBxByBz(&pTPC,r,mass,0.1,kTRUE);
1557 AliTrackPoint lpoint = point.Rotate(pTPC.GetAlpha());
1558 pTPC.PropagateTo(lpoint.GetX(),fMagF);
1560 ((Double_t*)pTOF.GetParameter())[0] =lpoint.GetY();
1561 ((Double_t*)pTOF.GetParameter())[1] =lpoint.GetZ();
1562 ((Double_t*)pTOF.GetCovariance())[0]+=3.*3./12.;
1563 ((Double_t*)pTOF.GetCovariance())[2]+=3.*3./12.;
1564 ((Double_t*)pTOF.GetCovariance())[5]+=0.1*0.1;
1565 ((Double_t*)pTOF.GetCovariance())[9]+=0.1*0.1;
1568 if (naccept==0) return; // no tof match clusters
1570 // 0. Apply standard cuts
1572 if (track->GetTPCNcls()<kMinTPC) return; // minimal amount of clusters cut
1573 // exclude crossing track
1574 if (track->GetOuterParam()->GetZ()*track->GetInnerParam()->GetZ()<0) return;
1576 if (TMath::Abs(pTOF.GetY()-pTPC.GetY()) >kMaxDy) return;
1577 if (TMath::Abs(pTOF.GetSnp()-pTPC.GetSnp())>kMaxAngle) return;
1578 if (TMath::Abs(pTOF.GetTgl()-pTPC.GetTgl())>kMaxAngle) return;
1580 // 1. Update median and RMS info
1582 TVectorD vecDelta(4),vecMedian(4), vecRMS(4);
1583 TVectorD vecDeltaN(5);
1584 Double_t sign=(pTOF.GetParameter()[1]>0)? 1.:-1.;
1586 for (Int_t i=0;i<4;i++){
1587 vecDelta[i]=(pTOF.GetParameter()[i]-pTPC.GetParameter()[i])*sign;
1588 kgdP[i][kglast%kN]=vecDelta[i];
1591 Int_t entries=(kglast<kN)?kglast:kN;
1593 for (Int_t i=0;i<4;i++){
1594 vecMedian[i] = TMath::Median(entries,kgdP[i]);
1595 vecRMS[i] = TMath::RMS(entries,kgdP[i]);
1598 vecDeltaN[i] = (vecDelta[i]-vecMedian[i])/(vecRMS[i]+1.);
1599 vecDeltaN[4]+= TMath::Abs(vecDeltaN[i]); //sum of abs residuals
1600 if (TMath::Abs(vecDeltaN[i])>kSigmaCut) isOK=kFALSE;
1604 // 2. Apply median+-rms cut
1606 if (kglast<10) return; //median and RMS to be defined
1609 // 3. Update alignment
1611 Int_t htime = fTime/3600; //time in hours
1612 if (fAlignTOFTPC->GetEntries()<htime){
1613 fAlignTOFTPC->Expand(htime*2+20);
1615 AliRelAlignerKalman* align = (AliRelAlignerKalman*)fAlignTOFTPC->At(htime);
1617 // make Alignment object if doesn't exist
1618 align=new AliRelAlignerKalman();
1619 align->SetRunNumber(fRun);
1620 (*align->GetStateCov())(6,6)=kVdErr*kVdErr;
1621 (*align->GetStateCov())(8,8)=kVdYErr*kVdYErr;
1622 align->SetOutRejSigma(kOutCut+kOutCut*kN);
1623 align->SetRejectOutliers(kFALSE);
1624 align->SetTPCvd(AliTPCcalibDB::Instance()->GetParameters()->GetDriftV()/1000000.);
1625 align->SetMagField(fMagF);
1626 fAlignTOFTPC->AddAt(align,htime);
1628 align->AddTrackParams(&pTOF,&pTPC);
1629 align->SetTimeStamp(fTime);
1630 align->SetRunNumber(fRun );
1632 Int_t nupdates=align->GetNUpdates();
1633 align->SetOutRejSigma(kOutCut+kOutCut*kN/Double_t(nupdates));
1634 align->SetRejectOutliers(kFALSE);
1635 TTreeSRedirector *cstream = GetDebugStreamer();
1636 if (cstream && align->GetState() && align->GetState()->GetNrows()>2 ){
1637 TTimeStamp tstamp(fTime);
1638 Float_t valuePressure0 = AliTPCcalibDB::GetPressure(tstamp,fRun,0);
1639 Float_t valuePressure1 = AliTPCcalibDB::GetPressure(tstamp,fRun,1);
1640 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,0);
1641 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,fRun,1);
1642 Double_t temp0 = AliTPCcalibDB::GetTemperature(tstamp,fRun,0);
1643 Double_t temp1 = AliTPCcalibDB::GetTemperature(tstamp,fRun,1);
1644 TVectorD vecGoofie(20);
1645 AliDCSSensorArray* goofieArray = AliTPCcalibDB::Instance()->GetGoofieSensors(fRun);
1647 for (Int_t isensor=0; isensor<goofieArray->NumSensors();isensor++){
1648 AliDCSSensor *gsensor = goofieArray->GetSensor(isensor);
1649 if (gsensor) vecGoofie[isensor]=gsensor->GetValue(tstamp);
1652 TVectorD gpTPC(3), gdTPC(3);
1653 TVectorD gpTOF(3), gdTOF(3);
1654 pTPC.GetXYZ(gpTPC.GetMatrixArray());
1655 pTPC.GetDirection(gdTPC.GetMatrixArray());
1656 pTOF.GetXYZ(gpTOF.GetMatrixArray());
1657 pTOF.GetDirection(gdTOF.GetMatrixArray());
1658 Double_t vdcorr = AliTPCcalibDB::Instance()->GetVDriftCorrectionTime(tstamp,fRun,0,1);
1659 (*cstream)<<"toftpc"<<
1660 "run="<<fRun<< // run number
1661 "event="<<fEvent<< // event number
1662 "time="<<fTime<< // time stamp of event
1663 "trigger="<<fTrigger<< // trigger
1664 "mag="<<fMagF<< // magnetic field
1665 // Environment values
1666 "press0="<<valuePressure0<<
1667 "press1="<<valuePressure1<<
1668 "pt0="<<ptrelative0<<
1669 "pt1="<<ptrelative1<<
1672 "vecGoofie.="<<&vecGoofie<<
1673 "vdcorr="<<vdcorr<< // drift correction applied
1675 "nmed="<<kglast<< // number of entries to define median and RMS
1676 "vMed.="<<&vecMedian<< // median of deltas
1677 "vRMS.="<<&vecRMS<< // rms of deltas
1678 "vDelta.="<<&vecDelta<< // delta in respect to median
1679 "vDeltaN.="<<&vecDeltaN<< // normalized delta in respect to median
1680 "t.="<<track<< // ful track - find proper cuts
1681 "a.="<<align<< // current alignment
1682 "pTOF.="<<&pTOF<< // track param TOF
1683 "pTPC.="<<&pTPC<< // track param TPC
1684 "gpTPC.="<<&gpTPC<< // global position TPC
1685 "gdTPC.="<<&gdTPC<< // global direction TPC
1686 "gpTOF.="<<&gpTOF<< // global position TOF
1687 "gdTOF.="<<&gdTOF<< // global position TOF