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::PrintStatus()
399 Int_t status = GetStatus();
401 AliInfo(Form("TOF calibration successful: %s (status=%d)", fgkStatusCodeName[fStatus], status));
403 else if (status > 0) {
404 AliInfo(Form("TOF calibration failed: %s (status=%d)", fgkStatusCodeName[fStatus], status));
406 else if (status < 0) {
407 AliInfo(Form("TOF calibration failed (expected): %s (status=%d)", fgkStatusCodeName[fStatus], status));
412 //_______________________________________________________
415 AliTOFAnalysisTaskCalibPass0::ProcessOutput(const Char_t *filename, AliCDBStorage* db)
421 Int_t ret = DoProcessOutput(filename, db);
426 //_______________________________________________________
429 AliTOFAnalysisTaskCalibPass0::DoProcessOutput(const Char_t *filename, AliCDBStorage *db)
435 /* reset status to OK */
439 TFile *file = TFile::Open(filename);
440 if (!file || !file->IsOpen()) {
441 AliError(Form("cannot open output file %s", filename));
442 fStatus = kInputError;
446 TList *list = (TList *)file->Get("TOFHistos");
447 TH2F *histoVertexTimestamp = NULL;
448 TH2F *histoDeltatTimestamp = NULL;
449 TH2F *histoDeltazEta = NULL;
450 TH2F *histoDeltazCosTheta = NULL;
451 TH2F *histoAcceptedTracksEtaPt = NULL;
452 TH2F *histoMatchedTracksEtaPt = NULL;
454 AliInfo(Form("getting histograms from \"Histos\" list from file %s", filename));
455 histoVertexTimestamp = (TH2F *)list->FindObject("hHistoVertexTimestamp");
456 histoDeltatTimestamp = (TH2F *)list->FindObject("hHistoDeltatTimestamp");
457 histoDeltazEta = (TH2F *)list->FindObject("hHistoDeltazEta");
458 histoDeltazCosTheta = (TH2F *)list->FindObject("hHistoDeltazCosTheta");
459 histoAcceptedTracksEtaPt = (TH2F *)list->FindObject("hHistoAcceptedTracksEtaPt");
460 histoMatchedTracksEtaPt = (TH2F *)list->FindObject("hHistoMatchedTracksEtaPt");
463 AliInfo(Form("getting histograms directly from file %s", filename));
464 histoVertexTimestamp = (TH2F *)file->Get("hHistoVertexTimestamp");
465 histoDeltatTimestamp = (TH2F *)file->Get("hHistoDeltatTimestamp");
466 histoDeltazEta = (TH2F *)file->Get("hHistoDeltazEta");
467 histoDeltazCosTheta = (TH2F *)file->Get("hHistoDeltazCosTheta");
468 histoAcceptedTracksEtaPt = (TH2F *)file->Get("hHistoAcceptedTracksEtaPt");
469 histoMatchedTracksEtaPt = (TH2F *)file->Get("hHistoMatchedTracksEtaPt");
472 if (!histoVertexTimestamp) {
473 AliError(Form("cannot get \"hHistoVertexTimestamp\" object from file %s", filename));
474 fStatus = kInputError;
477 if (!histoDeltatTimestamp) {
478 AliError(Form("cannot get \"hHistoDeltatTimestamp\" object from file %s", filename));
479 fStatus = kInputError;
482 if (!histoDeltazEta) {
483 AliError(Form("cannot get \"hHistoDeltazEta\" object from file %s", filename));
484 fStatus = kInputError;
487 if (!histoDeltazCosTheta) {
488 AliError(Form("cannot get \"hHistoDeltazCosTheta\" object from file %s", filename));
489 fStatus = kInputError;
492 if (!histoAcceptedTracksEtaPt) {
493 AliError(Form("cannot get \"hHistoAccptedTracksEtaPt\" object from file %s", filename));
494 fStatus = kInputError;
497 if (!histoMatchedTracksEtaPt) {
498 AliError(Form("cannot get \"hHistoMatchedTracksEtaPt\" object from file %s", filename));
499 fStatus = kInputError;
503 /* check matching performance */
504 if (!CheckMatchingPerformance(histoDeltazEta, histoAcceptedTracksEtaPt, histoMatchedTracksEtaPt)) {
505 AliError("error while checking matching efficiency");
508 /* calibrate and store */
509 if (!CalibrateAndStore(histoVertexTimestamp, histoDeltatTimestamp, db)) {
510 AliError("error while calibrating and storing");
518 //_______________________________________________________
521 AliTOFAnalysisTaskCalibPass0::CheckMatchingPerformance(const TH2F *histoDeltazCosTheta, const TH2F *histoAcceptedTracksEtaPt, const TH2F *histoMatchedTracksEtaPt) const
524 * check matching performance
528 if (!histoDeltazCosTheta || !histoAcceptedTracksEtaPt || !histoMatchedTracksEtaPt)
530 /* dummy for the time being */
534 //_______________________________________________________
537 AliTOFAnalysisTaskCalibPass0::CalibrateAndStore(TH2F *histoVertexTimestamp, TH2F *histoDeltatTimestamp, AliCDBStorage *db)
540 * calibrate and store
544 if (!histoVertexTimestamp || !histoDeltatTimestamp)
547 /*** GET RUN-NUMBER AND START-TIMESTAMP ***/
550 TObjArray *strarr = NULL;
551 TObjString *ostr = NULL;
552 str = histoVertexTimestamp->GetTitle();
553 strarr = str.Tokenize(",");
555 AliError("problems whith tokenize histogram title");
556 fStatus = kDataError;
561 ostr = (TObjString *)strarr->At(0);
563 AliError("problems while getting run number from histogram title");
564 fStatus = kDataError;
567 str = ostr->GetString();
568 if (!str.BeginsWith("run:")) {
569 AliError("problems while getting run number from histogram title");
570 fStatus = kDataError;
574 Int_t runNb = atoi(str.Data());
576 AliError(Form("bad run number: %d", runNb));
577 fStatus = kDataError;
580 AliInfo(Form("got run number: %d", runNb));
582 /* start timestamp */
583 ostr = (TObjString *)strarr->At(1);
585 AliError("problems while getting start timestamp from histogram title");
586 fStatus = kDataError;
589 str = ostr->GetString();
590 str.Remove(0, 1); /* remove empty space at the beginning */
591 if (!str.BeginsWith("startTimestamp:")) {
592 AliError("problems while getting start timestamp from histogram title");
593 fStatus = kDataError;
597 UInt_t startTimestamp = atoi(str.Data());
598 if (startTimestamp <= 0) {
599 AliError(Form("bad start timestamp: %d", startTimestamp));
600 fStatus = kDataError;
603 TTimeStamp ts = startTimestamp;
604 AliInfo(Form("got start timestamp: %d (%s)", startTimestamp, ts.AsString()));
607 ostr = (TObjString *)strarr->At(2);
609 AliError("problems while getting BPTX from histogram title");
610 fStatus = kDataError;
613 str = ostr->GetString();
614 str.Remove(0, 1); /* remove empty space at the beginning */
615 if (!str.BeginsWith("BPTX:")) {
616 AliError("problems while getting BPTX from histogram title");
617 fStatus = kDataError;
621 Bool_t useBPTX = atoi(str.Data());
622 AliInfo(Form("got BPTX: %d", useBPTX));
624 /*** CALIBRATION STAGE ***/
626 /* get fit function */
627 TF1 *fitFunc = (TF1 *)gROOT->GetFunction("gaus");
630 TH1D *histoVertexTimestamppx = histoVertexTimestamp->ProjectionX("histoVertexTimestamppx");
631 TH1D *histoDeltatTimestamppx = histoDeltatTimestamp->ProjectionX("histoDeltatTimestamppx");
633 /* check statistics */
634 if (histoVertexTimestamppx->Integral() < fgMinVertexIntegral ||
635 histoDeltatTimestamppx->Integral() < fgMinDeltatIntegral) {
636 fStatus = kLowStatistics;
640 /* define mix and max time bin */
641 Int_t minBin = histoVertexTimestamppx->FindFirstBinAbove(0);
642 Int_t maxBin = histoVertexTimestamppx->FindLastBinAbove(0);
643 Float_t minTime = histoVertexTimestamppx->GetBinLowEdge(minBin);
644 Float_t maxTime = histoVertexTimestamppx->GetBinLowEdge(maxBin + 1);
645 AliInfo(Form("min/max time defined: %d < t < %d s [%d, %d]", (Int_t)minTime, (Int_t)maxTime, minBin, maxBin));
647 /* loop over time bins */
649 Float_t time[fgkMaxNumberOfPoints], timeerr[fgkMaxNumberOfPoints];
650 Float_t vertexMean[fgkMaxNumberOfPoints], vertexMeanerr[fgkMaxNumberOfPoints];
651 Float_t vertexSigma[fgkMaxNumberOfPoints], vertexSigmaerr[fgkMaxNumberOfPoints];
652 Float_t timeZeroMean[fgkMaxNumberOfPoints], timeZeroMeanerr[fgkMaxNumberOfPoints];
653 Float_t timeZeroSigma[fgkMaxNumberOfPoints], timeZeroSigmaerr[fgkMaxNumberOfPoints];
654 Float_t averageTimeZero = 0.;
655 Float_t averageTimeSigma = 0.;
656 Float_t averageVertexSpread = 0.;
657 for (Int_t ibin = minBin; ibin <= maxBin; ibin++) {
659 /* define time window */
660 Int_t startBin = ibin;
662 while(histoVertexTimestamppx->Integral(startBin, endBin) < fgMinVertexIntegralSample ||
663 histoDeltatTimestamppx->Integral(startBin, endBin) < fgMinDeltatIntegralSample) {
664 if (endBin < maxBin) endBin++;
665 else if (startBin > minBin) startBin--;
668 if (histoVertexTimestamppx->Integral(startBin, endBin) < fgMinVertexIntegral ||
669 histoDeltatTimestamppx->Integral(startBin, endBin) < fgMinDeltatIntegral) continue;
670 Float_t startTime = histoVertexTimestamppx->GetBinLowEdge(startBin);
671 Float_t endTime = histoVertexTimestamppx->GetBinLowEdge(endBin + 1);
672 Float_t vertexIntegral = histoVertexTimestamppx->Integral(startBin, endBin);
673 Float_t deltatIntegral = histoDeltatTimestamppx->Integral(startBin, endBin);
674 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));
677 TH1D *histoVertexTimestamppy = histoVertexTimestamp->ProjectionY("histoVertexTimestamppy", startBin, endBin);
678 TH1D *histoDeltatTimestamppy = histoDeltatTimestamp->ProjectionY("histoDeltatTimestamppy", startBin, endBin);
681 histoVertexTimestamppx->GetXaxis()->SetRange(startBin, endBin);
682 time[nPoints] = histoVertexTimestamppx->GetMean();
683 timeerr[nPoints] = histoVertexTimestamppx->GetMeanError();
686 if (FitPeak(fitFunc, histoVertexTimestamppy, 10., 3., 3.) != 0) {
687 AliError("troubles fitting vertex, skip");
688 delete histoVertexTimestamppy;
689 delete histoDeltatTimestamppy;
692 vertexMean[nPoints] = fitFunc->GetParameter(1);
693 vertexMeanerr[nPoints] = fitFunc->GetParError(1);
694 vertexSigma[nPoints] = fitFunc->GetParameter(2);
695 vertexSigmaerr[nPoints] = fitFunc->GetParError(2);
696 averageVertexSpread += fitFunc->GetParameter(2);
699 if (FitPeak(fitFunc, histoDeltatTimestamppy, 500., 2., 1.) != 0) {
700 AliError("troubles fitting time-zero TRACKS, skip");
701 delete histoVertexTimestamppy;
702 delete histoDeltatTimestamppy;
705 timeZeroMean[nPoints] = fitFunc->GetParameter(1);
706 timeZeroMeanerr[nPoints] = fitFunc->GetParError(1);
707 timeZeroSigma[nPoints] = fitFunc->GetParameter(2);
708 timeZeroSigmaerr[nPoints] = fitFunc->GetParError(2);
709 averageTimeZero += fitFunc->GetParameter(1);
710 averageTimeSigma += fitFunc->GetParameter(2);
712 /* delete projection-y */
713 delete histoVertexTimestamppy;
714 delete histoDeltatTimestamppy;
716 /* increment n points */
719 /* set current bin */
722 } /* end of loop over time bins */
724 /* delete projection-x */
725 delete histoVertexTimestamppx;
726 delete histoDeltatTimestamppx;
730 AliError("no measurement available, quit");
731 fStatus = kNoMeasurement;
734 AliInfo(Form("average time-zero = %f", averageTimeZero / nPoints));
735 AliInfo(Form("average time-sigma = %f", averageTimeSigma / nPoints));
736 AliInfo(Form("average v-spread = %f", averageVertexSpread / nPoints));
738 /*** CREATE RUN PARAMS OBJECT ***/
741 /* get start time from GRP */
742 TGrid::Connect("alien://", 0, 0, "t");
743 AliCDBManager *cdb = AliCDBManager::Instance();
744 cdb->SetDefaultStorage("raw://");
747 if (!grp.ReadGRPEntry()) {
748 AliError("error while reading GRP entry");
751 UInt_t startTimestamp = grp.GetGRPData()->GetTimeStart();
754 AliInfo(Form("got start time from GRP: %s", ts.AsString()));
758 UInt_t timestamp[fgkMaxNumberOfPoints];
759 Float_t t0[fgkMaxNumberOfPoints];
760 Float_t tofReso[fgkMaxNumberOfPoints];
761 Float_t t0Spread[fgkMaxNumberOfPoints];
762 for (Int_t ipoint = 0; ipoint < nPoints; ipoint++) {
763 timestamp[ipoint] = (UInt_t)time[ipoint] + startTimestamp;
764 t0[ipoint] = timeZeroMean[ipoint];
765 tofReso[ipoint] = timeZeroSigma[ipoint];
766 t0Spread[ipoint] = vertexSigma[ipoint] / 2.99792457999999984e-02;
768 UInt_t run[1] = {runNb};
769 UInt_t runFirstPoint[1] = {0};
770 UInt_t runLastPoint[1] = {nPoints - 1};
772 /* create run params object */
773 AliTOFRunParams obj(nPoints, 1);
774 AliInfo(Form("create run params object for run %d with %d points", runNb, nPoints));
775 obj.SetTimestamp(timestamp);
777 obj.SetTOFResolution(tofReso);
778 obj.SetT0Spread(t0Spread);
780 obj.SetRunFirstPoint(runFirstPoint);
781 obj.SetRunLastPoint(runLastPoint);
782 obj.SetUseLHCClockPhase(useBPTX);
784 /*** CREATE OCDB ENTRY ***/
787 AliError("cannot store object because of NULL storage");
788 fStatus = kStoreError;
792 /* install run params object in OCDB */
793 AliCDBId id("TOF/Calib/RunParams", runNb, runNb);
795 md.SetResponsible("Roberto Preghenella");
796 md.SetComment("offline TOF run parameters");
797 md.SetAliRootVersion(gSystem->Getenv("ARVERSION"));
799 if (!db->Put(&obj, id, &md)) {
800 fStatus = kStoreError;
801 AliError(Form("error while putting object in storage %s", db->GetURI().Data()));
809 //_____________________________________________________________
812 AliTOFAnalysisTaskCalibPass0::FitPeak(TF1 *fitFunc, TH1D *h, Float_t startSigma, Float_t nSigmaMin, Float_t nSigmaMax)
818 Double_t fitCent = h->GetBinCenter(h->GetMaximumBin());
819 Double_t fitMin = fitCent - nSigmaMin * startSigma;
820 Double_t fitMax = fitCent + nSigmaMax * startSigma;
821 if (fitMin < h->GetXaxis()->GetXmin()) fitMin = h->GetXaxis()->GetXmin();
822 if (fitMax > h->GetXaxis()->GetXmax()) fitMax = h->GetXaxis()->GetXmax();
823 fitFunc->SetParameter(1, fitCent);
824 fitFunc->SetParameter(2, startSigma);
825 Int_t fitres = h->Fit(fitFunc, "WWq0", "", fitMin, fitMax);
826 if (fitres != 0) return fitres;
827 /* refit with better range */
828 for (Int_t i = 0; i < 3; i++) {
829 fitCent = fitFunc->GetParameter(1);
830 fitMin = fitCent - nSigmaMin * fitFunc->GetParameter(2);
831 fitMax = fitCent + nSigmaMax * fitFunc->GetParameter(2);
832 if (fitMin < h->GetXaxis()->GetXmin()) fitMin = h->GetXaxis()->GetXmin();
833 if (fitMax > h->GetXaxis()->GetXmax()) fitMax = h->GetXaxis()->GetXmax();
834 fitres = h->Fit(fitFunc, "q0", "", fitMin, fitMax);
835 if (fitres != 0) return fitres;