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 ***************************************************************************/
16 ///////////////////////////////////////////////////////////////
18 // This class provides TOF pass0/passX calibration tools //
20 ///////////////////////////////////////////////////////////////
22 #include "AliTOFAnalysisTaskCalibPass0.h"
23 #include "AliESDtrack.h"
24 #include "AliESDEvent.h"
25 #include "AliPhysicsSelection.h"
26 #include "AliESDtrackCuts.h"
27 #include "AliCDBManager.h"
28 #include "AliGRPManager.h"
29 #include "AliGRPObject.h"
30 #include "AliTOFcalib.h"
35 #include "TObjArray.h"
37 #include "TObjString.h"
42 #include "TTimeStamp.h"
43 #include "AliTOFRunParams.h"
44 #include "AliCDBStorage.h"
46 #include "AliCDBMetaData.h"
50 ClassImp(AliTOFAnalysisTaskCalibPass0)
52 //_______________________________________________________
54 const Char_t *AliTOFAnalysisTaskCalibPass0::fgkStatusCodeName[AliTOFAnalysisTaskCalibPass0::kNStatusCodes] = {
64 const Int_t AliTOFAnalysisTaskCalibPass0::fgkMaxNumberOfPoints = 10000; // max number of points
65 Double_t AliTOFAnalysisTaskCalibPass0::fgMinVertexIntegral = 100.;
66 Double_t AliTOFAnalysisTaskCalibPass0::fgMinDeltatIntegral = 2000.;
67 Double_t AliTOFAnalysisTaskCalibPass0::fgMinVertexIntegralSample = 1000.;
68 Double_t AliTOFAnalysisTaskCalibPass0::fgMinDeltatIntegralSample = 20000.;
70 //_______________________________________________________
72 AliTOFAnalysisTaskCalibPass0::AliTOFAnalysisTaskCalibPass0() :
73 AliAnalysisTaskSE("TOFCalib-Pass0"),
76 fEventSelectionFlag(kFALSE),
77 fVertexSelectionFlag(kFALSE),
81 fEventCuts(new AliPhysicsSelection()),
82 fTrackCuts(new AliESDtrackCuts()),
88 fGRPManager(new AliGRPManager()),
90 fTOFcalib(new AliTOFcalib()),
91 fHistoList(new TList()),
92 fHistoVertexTimestamp(NULL),
93 fHistoDeltatTimestamp(NULL),
94 fHistoDeltazEta(NULL),
95 fHistoDeltazCosTheta(NULL),
96 fHistoAcceptedTracksEtaPt(NULL),
97 fHistoMatchedTracksEtaPt(NULL)
100 * default constructor
104 fHistoList->SetOwner(kTRUE);
105 DefineOutput(1, TList::Class());
108 //_______________________________________________________
110 AliTOFAnalysisTaskCalibPass0::~AliTOFAnalysisTaskCalibPass0()
116 if (fHistoList) delete fHistoList;
119 //_______________________________________________________
122 AliTOFAnalysisTaskCalibPass0::InitRun()
129 fESDEvent = dynamic_cast<AliESDEvent *>(InputEvent());
131 AliError("cannot get ESD event");
135 Int_t runNb = fESDEvent->GetRunNumber();
136 /* check run already initialized */
137 if (fInitFlag && fRunNumber == runNb) return kTRUE;
138 if (fInitFlag && fRunNumber != runNb) {
139 AliFatal(Form("run number has changed within same job and this is not allowed: %d -> %d", fRunNumber, runNb));
143 /* init cdb if not yet done */
144 AliCDBManager *cdb = AliCDBManager::Instance();
145 if (!cdb->IsDefaultStorageSet())
146 cdb->SetDefaultStorage("raw://");
147 if (cdb->GetRun() < 0)
150 if (!fTOFcalib->Init()) {
151 AliError("cannot init TOF calib");
155 if (!fGRPManager->ReadGRPEntry()) {
156 AliError("error while reading \"GRPEntry\" from OCDB");
159 fkGRPObject = fGRPManager->GetGRPData();
161 AliError("cannot get \"GRPData\" from GRP manager");
164 fStartTime = fkGRPObject->GetTimeStart();
165 fEndTime = fkGRPObject->GetTimeEnd();
166 AliInfo(Form("got \"GRPData\": startTime=%d, endTime=%d", fStartTime, fEndTime));
168 AliInfo(Form("initialized for run %d", runNb));
172 Int_t useBPTX = fTOFcalib->GetUseLHCClockPhase() ? 1 : 0;
174 /* set histo title with run-number and start-time */
175 fHistoVertexTimestamp->SetTitle(Form("run: %d, startTimestamp: %u, BPTX: %d", fRunNumber, fStartTime, useBPTX));
176 fHistoDeltatTimestamp->SetTitle(Form("run: %d, startTimestamp: %u, BPTX: %d", fRunNumber, fStartTime, useBPTX));
177 fHistoDeltazEta->SetTitle(Form("run: %d, startTimestamp: %u, BPTX: %d", fRunNumber, fStartTime, useBPTX));
178 fHistoDeltazCosTheta->SetTitle(Form("run: %d, startTimestamp: %u, BPTX: %d", fRunNumber, fStartTime, useBPTX));
179 fHistoAcceptedTracksEtaPt->SetTitle(Form("run: %d, startTimestamp: %u, BPTX: %d", fRunNumber, fStartTime, useBPTX));
180 fHistoMatchedTracksEtaPt->SetTitle(Form("run: %d, startTimestamp: %u, BPTX: %d", fRunNumber, fStartTime, useBPTX));
185 //_______________________________________________________
188 AliTOFAnalysisTaskCalibPass0::InitEvent()
195 fESDEvent = dynamic_cast<AliESDEvent *>(InputEvent());
196 if (!fESDEvent) return kFALSE;
197 /* set event time and elapsed time */
198 fEventTime = fESDEvent->GetTimeStamp();
199 fElapsedTime = fESDEvent->GetTimeStamp() - fStartTime;
200 /* event selection */
201 if (fEventSelectionFlag && !fEventCuts->IsCollisionCandidate(fESDEvent)) return kFALSE;
202 /* vertex selection */
203 fkVertex = fESDEvent->GetPrimaryVertexTracks();
204 if (fVertexSelectionFlag && fkVertex->GetNContributors() < 1) {
205 fkVertex = fESDEvent->GetPrimaryVertexSPD();
206 if (fkVertex->GetNContributors() < 1) return kFALSE;
208 if (fVertexSelectionFlag && TMath::Abs(fkVertex->GetZ()) > fVertexCut) return kFALSE;
209 /* calibrate ESD if requested */
210 fTOFcalib->CalibrateESD(fESDEvent);
215 //_______________________________________________________
218 AliTOFAnalysisTaskCalibPass0::HasTOFMeasurement(const AliESDtrack *track) const
221 * has TOF measurement
224 /* check TOF status flags */
225 if (!(track->GetStatus() & AliESDtrack::kTOFout) ||
226 !(track->GetStatus() & AliESDtrack::kTIME)) return kFALSE;
227 /* check integrated length */
228 if (track->GetIntegratedLength() < 350.) return kFALSE;
230 /* TOF measurement ok */
234 //_______________________________________________________
237 AliTOFAnalysisTaskCalibPass0::UserCreateOutputObjects()
240 * user create output objects
244 Int_t timeBins = 288;
245 Float_t timeMin = 0.;
246 Float_t timeMax = 24. * 3600.;
248 Int_t vertexBins = 200;
249 Float_t vertexMin = -50.;
250 Float_t vertexMax = 50.;
252 Int_t deltatBins = 2000;
253 Float_t deltatMin = -24400.;
254 Float_t deltatMax = 24400.;
256 Int_t deltazBins = 200;
257 Float_t deltazMin = -10.;
258 Float_t deltazMax = 10.;
261 Float_t etaMin = -1.;
268 fHistoVertexTimestamp = new TH2F("hHistoVertexTimestamp", "Vertex position;elapsed time (s);z (cm);", timeBins, timeMin, timeMax, vertexBins, vertexMin, vertexMax);
269 fHistoList->Add(fHistoVertexTimestamp);
271 fHistoDeltatTimestamp = new TH2F("hHistoDeltatTimestamp", "Global time shift (T0-fill);elapsed time (s);t - t_{exp}^{(#pi)} (ps);", timeBins, timeMin, timeMax, deltatBins, deltatMin, deltatMax);
272 fHistoList->Add(fHistoDeltatTimestamp);
274 fHistoDeltazEta = new TH2F("hHistoDeltazEta", "Matching residuals (longitudinal);#eta;#Deltaz (cm);", etaBins, etaMin, etaMax, deltazBins, deltazMin, deltazMax);
275 fHistoList->Add(fHistoDeltazEta);
277 fHistoDeltazCosTheta = new TH2F("hHistoDeltazCosTheta", "Matching residuals (longitudinal);cos #theta;#Deltaz (cm);", etaBins, etaMin, etaMax, deltazBins, deltazMin, deltazMax);
278 fHistoList->Add(fHistoDeltazCosTheta);
280 fHistoAcceptedTracksEtaPt = new TH2F("hHistoAcceptedTracksEtaPt", "Accepted tracks;#eta;#p_{T} (GeV/c);", etaBins, etaMin, etaMax, pBins, pMin, pMax);
281 fHistoList->Add(fHistoAcceptedTracksEtaPt);
283 fHistoMatchedTracksEtaPt = new TH2F("hHistoMatchedTracksEtaPt", "Matched tracks;#eta;p_{T} (GeV/c);", etaBins, etaMin, etaMax, pBins, pMin, pMax);
284 fHistoList->Add(fHistoMatchedTracksEtaPt);
287 PostData(1, fHistoList);
290 //_______________________________________________________
293 AliTOFAnalysisTaskCalibPass0::UserExec(Option_t *)
300 if (!InitRun()) return;
302 if (!InitEvent()) return;
304 /*** ACCEPTED EVENT ***/
306 /* fill vertex histo */
307 fHistoVertexTimestamp->Fill(fElapsedTime, fkVertex->GetZ());
309 /* loop over ESD tracks */
310 Int_t nTracks = fESDEvent->GetNumberOfTracks();
312 Double_t eta, costheta, pt, time, timei[AliPID::kSPECIES], deltat, deltaz;
313 for (Int_t itrk = 0; itrk < nTracks; itrk++) {
315 track = fESDEvent->GetTrack(itrk);
316 if (!track) continue;
317 /* check accept track */
318 if (!fTrackCuts->AcceptTrack(track)) continue;
321 costheta = TMath::Cos(track->Theta());
323 /* fill accepted tracks histo */
324 fHistoAcceptedTracksEtaPt->Fill(eta, pt);
325 /* check TOF measurement */
326 if (!HasTOFMeasurement(track)) continue;
328 /*** ACCEPTED TRACK WITH TOF MEASUREMENT ***/
330 /* fill matched tracks histo */
331 fHistoMatchedTracksEtaPt->Fill(eta, pt);
333 time = track->GetTOFsignal();
334 track->GetIntegratedTimes(timei);
335 deltat = time - timei[AliPID::kPion];
336 deltaz = track->GetTOFsignalDz();
339 fHistoDeltatTimestamp->Fill(fElapsedTime, deltat);
340 fHistoDeltazEta->Fill(eta, deltaz);
341 fHistoDeltazCosTheta->Fill(costheta, deltaz);
343 } /* end of loop over ESD tracks */
346 PostData(1, fHistoList);
350 //_______________________________________________________
353 AliTOFAnalysisTaskCalibPass0::GetStatus()
361 /* OK, return zero */
366 /* non-fatal error, return negative status */
373 /* fatal error, return positive status */
380 /* anything else, return negative large number */
386 /* should never arrive here, anyway return negative large number */
390 //_______________________________________________________
393 AliTOFAnalysisTaskCalibPass0::ProcessOutput(const Char_t *filename, const Char_t *dbString)
399 Int_t ret = DoProcessOutput(filename, dbString);
400 Int_t status = GetStatus();
402 AliInfo(Form("TOF calibration successful: %s (status=%d)", fgkStatusCodeName[fStatus], status));
404 else if (status > 0) {
405 AliInfo(Form("TOF calibration failed: %s (status=%d)", fgkStatusCodeName[fStatus], status));
407 else if (status < 0) {
408 AliInfo(Form("TOF calibration failed (expected): %s (status=%d)", fgkStatusCodeName[fStatus], status));
414 //_______________________________________________________
417 AliTOFAnalysisTaskCalibPass0::DoProcessOutput(const Char_t *filename, const Char_t *dbString)
423 /* reset status to OK */
427 TFile *file = TFile::Open(filename);
428 if (!file || !file->IsOpen()) {
429 AliError(Form("cannot open output file %s", filename));
430 fStatus = kInputError;
434 TList *list = (TList *)file->Get("TOFHistos");
435 TH2F *histoVertexTimestamp = NULL;
436 TH2F *histoDeltatTimestamp = NULL;
437 TH2F *histoDeltazEta = NULL;
438 TH2F *histoDeltazCosTheta = NULL;
439 TH2F *histoAcceptedTracksEtaPt = NULL;
440 TH2F *histoMatchedTracksEtaPt = NULL;
442 AliInfo(Form("getting histograms from \"Histos\" list from file %s", filename));
443 histoVertexTimestamp = (TH2F *)list->FindObject("hHistoVertexTimestamp");
444 histoDeltatTimestamp = (TH2F *)list->FindObject("hHistoDeltatTimestamp");
445 histoDeltazEta = (TH2F *)list->FindObject("hHistoDeltazEta");
446 histoDeltazCosTheta = (TH2F *)list->FindObject("hHistoDeltazCosTheta");
447 histoAcceptedTracksEtaPt = (TH2F *)list->FindObject("hHistoAcceptedTracksEtaPt");
448 histoMatchedTracksEtaPt = (TH2F *)list->FindObject("hHistoMatchedTracksEtaPt");
451 AliInfo(Form("getting histograms directly from file %s", filename));
452 histoVertexTimestamp = (TH2F *)file->Get("hHistoVertexTimestamp");
453 histoDeltatTimestamp = (TH2F *)file->Get("hHistoDeltatTimestamp");
454 histoDeltazEta = (TH2F *)file->Get("hHistoDeltazEta");
455 histoDeltazCosTheta = (TH2F *)file->Get("hHistoDeltazCosTheta");
456 histoAcceptedTracksEtaPt = (TH2F *)file->Get("hHistoAcceptedTracksEtaPt");
457 histoMatchedTracksEtaPt = (TH2F *)file->Get("hHistoMatchedTracksEtaPt");
460 if (!histoVertexTimestamp) {
461 AliError(Form("cannot get \"hHistoVertexTimestamp\" object from file %s", filename));
462 fStatus = kInputError;
465 if (!histoDeltatTimestamp) {
466 AliError(Form("cannot get \"hHistoDeltatTimestamp\" object from file %s", filename));
467 fStatus = kInputError;
470 if (!histoDeltazEta) {
471 AliError(Form("cannot get \"hHistoDeltazEta\" object from file %s", filename));
472 fStatus = kInputError;
475 if (!histoDeltazCosTheta) {
476 AliError(Form("cannot get \"hHistoDeltazCosTheta\" object from file %s", filename));
477 fStatus = kInputError;
480 if (!histoAcceptedTracksEtaPt) {
481 AliError(Form("cannot get \"hHistoAccptedTracksEtaPt\" object from file %s", filename));
482 fStatus = kInputError;
485 if (!histoMatchedTracksEtaPt) {
486 AliError(Form("cannot get \"hHistoMatchedTracksEtaPt\" object from file %s", filename));
487 fStatus = kInputError;
491 /* check matching performance */
492 if (!CheckMatchingPerformance(histoDeltazEta, histoAcceptedTracksEtaPt, histoMatchedTracksEtaPt)) {
493 AliError("error while checking matching efficiency");
496 /* calibrate and store */
497 if (!CalibrateAndStore(histoVertexTimestamp, histoDeltatTimestamp, dbString)) {
498 AliError("error while calibrating and storing");
506 //_______________________________________________________
509 AliTOFAnalysisTaskCalibPass0::CheckMatchingPerformance(const TH2F *histoDeltazCosTheta, const TH2F *histoAcceptedTracksEtaPt, const TH2F *histoMatchedTracksEtaPt) const
512 * check matching performance
516 if (!histoDeltazCosTheta || !histoAcceptedTracksEtaPt || !histoMatchedTracksEtaPt)
518 /* dummy for the time being */
522 //_______________________________________________________
525 AliTOFAnalysisTaskCalibPass0::CalibrateAndStore(TH2F *histoVertexTimestamp, TH2F *histoDeltatTimestamp, const Char_t *dbString)
528 * calibrate and store
532 if (!histoVertexTimestamp || !histoDeltatTimestamp)
535 /*** GET RUN-NUMBER AND START-TIMESTAMP ***/
538 TObjArray *strarr = NULL;
539 TObjString *ostr = NULL;
540 str = histoVertexTimestamp->GetTitle();
541 strarr = str.Tokenize(",");
543 AliError("problems whith tokenize histogram title");
544 fStatus = kDataError;
549 ostr = (TObjString *)strarr->At(0);
551 AliError("problems while getting run number from histogram title");
552 fStatus = kDataError;
555 str = ostr->GetString();
556 if (!str.BeginsWith("run:")) {
557 AliError("problems while getting run number from histogram title");
558 fStatus = kDataError;
562 Int_t runNb = atoi(str.Data());
564 AliError(Form("bad run number: %d", runNb));
565 fStatus = kDataError;
568 AliInfo(Form("got run number: %d", runNb));
570 /* start timestamp */
571 ostr = (TObjString *)strarr->At(1);
573 AliError("problems while getting start timestamp from histogram title");
574 fStatus = kDataError;
577 str = ostr->GetString();
578 str.Remove(0, 1); /* remove empty space at the beginning */
579 if (!str.BeginsWith("startTimestamp:")) {
580 AliError("problems while getting start timestamp from histogram title");
581 fStatus = kDataError;
585 UInt_t startTimestamp = atoi(str.Data());
586 if (startTimestamp <= 0) {
587 AliError(Form("bad start timestamp: %d", startTimestamp));
588 fStatus = kDataError;
591 TTimeStamp ts = startTimestamp;
592 AliInfo(Form("got start timestamp: %d (%s)", startTimestamp, ts.AsString()));
595 ostr = (TObjString *)strarr->At(2);
597 AliError("problems while getting BPTX from histogram title");
598 fStatus = kDataError;
601 str = ostr->GetString();
602 str.Remove(0, 1); /* remove empty space at the beginning */
603 if (!str.BeginsWith("BPTX:")) {
604 AliError("problems while getting BPTX from histogram title");
605 fStatus = kDataError;
609 Bool_t useBPTX = atoi(str.Data());
610 AliInfo(Form("got BPTX: %d", useBPTX));
612 /*** CALIBRATION STAGE ***/
614 /* get fit function */
615 TF1 *fitFunc = (TF1 *)gROOT->GetFunction("gaus");
618 TH1D *histoVertexTimestamppx = histoVertexTimestamp->ProjectionX("histoVertexTimestamppx");
619 TH1D *histoDeltatTimestamppx = histoDeltatTimestamp->ProjectionX("histoDeltatTimestamppx");
621 /* check statistics */
622 if (histoVertexTimestamppx->Integral() < fgMinVertexIntegral ||
623 histoDeltatTimestamppx->Integral() < fgMinDeltatIntegral) {
624 fStatus = kLowStatistics;
628 /* define mix and max time bin */
629 Int_t minBin = histoVertexTimestamppx->FindFirstBinAbove(0);
630 Int_t maxBin = histoVertexTimestamppx->FindLastBinAbove(0);
631 Float_t minTime = histoVertexTimestamppx->GetBinLowEdge(minBin);
632 Float_t maxTime = histoVertexTimestamppx->GetBinLowEdge(maxBin + 1);
633 AliInfo(Form("min/max time defined: %d < t < %d s [%d, %d]", (Int_t)minTime, (Int_t)maxTime, minBin, maxBin));
635 /* loop over time bins */
637 Float_t time[fgkMaxNumberOfPoints], timeerr[fgkMaxNumberOfPoints];
638 Float_t vertexMean[fgkMaxNumberOfPoints], vertexMeanerr[fgkMaxNumberOfPoints];
639 Float_t vertexSigma[fgkMaxNumberOfPoints], vertexSigmaerr[fgkMaxNumberOfPoints];
640 Float_t timeZeroMean[fgkMaxNumberOfPoints], timeZeroMeanerr[fgkMaxNumberOfPoints];
641 Float_t timeZeroSigma[fgkMaxNumberOfPoints], timeZeroSigmaerr[fgkMaxNumberOfPoints];
642 Float_t averageTimeZero = 0.;
643 Float_t averageTimeSigma = 0.;
644 Float_t averageVertexSpread = 0.;
645 for (Int_t ibin = minBin; ibin <= maxBin; ibin++) {
647 /* define time window */
648 Int_t startBin = ibin;
650 while(histoVertexTimestamppx->Integral(startBin, endBin) < fgMinVertexIntegralSample ||
651 histoDeltatTimestamppx->Integral(startBin, endBin) < fgMinDeltatIntegralSample) {
652 if (endBin < maxBin) endBin++;
653 else if (startBin > minBin) startBin--;
656 if (histoVertexTimestamppx->Integral(startBin, endBin) < fgMinVertexIntegral ||
657 histoDeltatTimestamppx->Integral(startBin, endBin) < fgMinDeltatIntegral) continue;
658 Float_t startTime = histoVertexTimestamppx->GetBinLowEdge(startBin);
659 Float_t endTime = histoVertexTimestamppx->GetBinLowEdge(endBin + 1);
660 Float_t vertexIntegral = histoVertexTimestamppx->Integral(startBin, endBin);
661 Float_t deltatIntegral = histoDeltatTimestamppx->Integral(startBin, endBin);
662 AliInfo(Form("time window defined: %d < t < %d s [%d, %d]: %d vertices, %d tracks", (Int_t)startTime, (Int_t)endTime, startBin, endBin, (Int_t)vertexIntegral, (Int_t)deltatIntegral));
665 TH1D *histoVertexTimestamppy = histoVertexTimestamp->ProjectionY("histoVertexTimestamppy", startBin, endBin);
666 TH1D *histoDeltatTimestamppy = histoDeltatTimestamp->ProjectionY("histoDeltatTimestamppy", startBin, endBin);
669 histoVertexTimestamppx->GetXaxis()->SetRange(startBin, endBin);
670 time[nPoints] = histoVertexTimestamppx->GetMean();
671 timeerr[nPoints] = histoVertexTimestamppx->GetMeanError();
674 if (FitPeak(fitFunc, histoVertexTimestamppy, 10., 3., 3.) != 0) {
675 AliError("troubles fitting vertex, skip");
676 delete histoVertexTimestamppy;
677 delete histoDeltatTimestamppy;
680 vertexMean[nPoints] = fitFunc->GetParameter(1);
681 vertexMeanerr[nPoints] = fitFunc->GetParError(1);
682 vertexSigma[nPoints] = fitFunc->GetParameter(2);
683 vertexSigmaerr[nPoints] = fitFunc->GetParError(2);
684 averageVertexSpread += fitFunc->GetParameter(2);
687 if (FitPeak(fitFunc, histoDeltatTimestamppy, 500., 2., 1.) != 0) {
688 AliError("troubles fitting time-zero TRACKS, skip");
689 delete histoVertexTimestamppy;
690 delete histoDeltatTimestamppy;
693 timeZeroMean[nPoints] = fitFunc->GetParameter(1);
694 timeZeroMeanerr[nPoints] = fitFunc->GetParError(1);
695 timeZeroSigma[nPoints] = fitFunc->GetParameter(2);
696 timeZeroSigmaerr[nPoints] = fitFunc->GetParError(2);
697 averageTimeZero += fitFunc->GetParameter(1);
698 averageTimeSigma += fitFunc->GetParameter(2);
700 /* delete projection-y */
701 delete histoVertexTimestamppy;
702 delete histoDeltatTimestamppy;
704 /* increment n points */
707 /* set current bin */
710 } /* end of loop over time bins */
712 /* delete projection-x */
713 delete histoVertexTimestamppx;
714 delete histoDeltatTimestamppx;
718 AliError("no measurement available, quit");
719 fStatus = kNoMeasurement;
722 AliInfo(Form("average time-zero = %f", averageTimeZero / nPoints));
723 AliInfo(Form("average time-sigma = %f", averageTimeSigma / nPoints));
724 AliInfo(Form("average v-spread = %f", averageVertexSpread / nPoints));
726 /*** CREATE RUN PARAMS OBJECT ***/
729 /* get start time from GRP */
730 TGrid::Connect("alien://", 0, 0, "t");
731 AliCDBManager *cdb = AliCDBManager::Instance();
732 cdb->SetDefaultStorage("raw://");
735 if (!grp.ReadGRPEntry()) {
736 AliError("error while reading GRP entry");
739 UInt_t startTimestamp = grp.GetGRPData()->GetTimeStart();
742 AliInfo(Form("got start time from GRP: %s", ts.AsString()));
746 UInt_t timestamp[fgkMaxNumberOfPoints];
747 Float_t t0[fgkMaxNumberOfPoints];
748 Float_t tofReso[fgkMaxNumberOfPoints];
749 Float_t t0Spread[fgkMaxNumberOfPoints];
750 for (Int_t ipoint = 0; ipoint < nPoints; ipoint++) {
751 timestamp[ipoint] = (UInt_t)time[ipoint] + startTimestamp;
752 t0[ipoint] = timeZeroMean[ipoint];
753 tofReso[ipoint] = timeZeroSigma[ipoint];
754 t0Spread[ipoint] = vertexSigma[ipoint] / 2.99792457999999984e-02;
756 UInt_t run[1] = {runNb};
757 UInt_t runFirstPoint[1] = {0};
758 UInt_t runLastPoint[1] = {nPoints - 1};
760 /* create run params object */
761 AliTOFRunParams obj(nPoints, 1);
762 AliInfo(Form("create run params object for run %d with %d points", runNb, nPoints));
763 obj.SetTimestamp(timestamp);
765 obj.SetTOFResolution(tofReso);
766 obj.SetT0Spread(t0Spread);
768 obj.SetRunFirstPoint(runFirstPoint);
769 obj.SetRunLastPoint(runLastPoint);
770 obj.SetUseLHCClockPhase(useBPTX);
772 /*** CREATE OCDB ENTRY ***/
775 AliError("cannot store object because of NULL string");
776 fStatus = kStoreError;
780 /* install run params object in OCDB */
781 AliCDBManager *cdb = AliCDBManager::Instance();
782 AliCDBStorage *sto = cdb->GetStorage(dbString);
784 AliError(Form("cannot get storage %s", dbString));
785 fStatus = kStoreError;
788 AliCDBId id("TOF/Calib/RunParams", runNb, runNb);
790 md.SetResponsible("Roberto Preghenella");
791 md.SetComment("offline TOF run parameters");
792 md.SetAliRootVersion(gSystem->Getenv("ARVERSION"));
794 if (!sto->Put(&obj, id, &md)) {
795 fStatus = kStoreError;
796 AliError(Form("error while putting object in storage %s", dbString));
804 //_____________________________________________________________
807 AliTOFAnalysisTaskCalibPass0::FitPeak(TF1 *fitFunc, TH1D *h, Float_t startSigma, Float_t nSigmaMin, Float_t nSigmaMax)
813 Double_t fitCent = h->GetBinCenter(h->GetMaximumBin());
814 Double_t fitMin = fitCent - nSigmaMin * startSigma;
815 Double_t fitMax = fitCent + nSigmaMax * startSigma;
816 if (fitMin < h->GetXaxis()->GetXmin()) fitMin = h->GetXaxis()->GetXmin();
817 if (fitMax > h->GetXaxis()->GetXmax()) fitMax = h->GetXaxis()->GetXmax();
818 fitFunc->SetParameter(1, fitCent);
819 fitFunc->SetParameter(2, startSigma);
820 Int_t fitres = h->Fit(fitFunc, "WWq0", "", fitMin, fitMax);
821 if (fitres != 0) return fitres;
822 /* refit with better range */
823 for (Int_t i = 0; i < 3; i++) {
824 fitCent = fitFunc->GetParameter(1);
825 fitMin = fitCent - nSigmaMin * fitFunc->GetParameter(2);
826 fitMax = fitCent + nSigmaMax * fitFunc->GetParameter(2);
827 if (fitMin < h->GetXaxis()->GetXmin()) fitMin = h->GetXaxis()->GetXmin();
828 if (fitMax > h->GetXaxis()->GetXmax()) fitMax = h->GetXaxis()->GetXmax();
829 fitres = h->Fit(fitFunc, "q0", "", fitMin, fitMax);
830 if (fitres != 0) return fitres;