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 fHistoDeltatEta(NULL),
96 fHistoDeltazCosTheta(NULL),
97 fHistoAcceptedTracksEtaPt(NULL),
98 fHistoMatchedTracksEtaPt(NULL)
101 * default constructor
105 fHistoList->SetOwner(kTRUE);
106 DefineOutput(1, TList::Class());
109 //_______________________________________________________
111 AliTOFAnalysisTaskCalibPass0::~AliTOFAnalysisTaskCalibPass0()
117 if (fHistoList) delete fHistoList;
120 //_______________________________________________________
123 AliTOFAnalysisTaskCalibPass0::InitRun()
130 fESDEvent = dynamic_cast<AliESDEvent *>(InputEvent());
132 AliError("cannot get ESD event");
136 Int_t runNb = fESDEvent->GetRunNumber();
137 /* check run already initialized */
138 if (fInitFlag && fRunNumber == runNb) return kTRUE;
139 if (fInitFlag && fRunNumber != runNb) {
140 AliFatal(Form("run number has changed within same job and this is not allowed: %d -> %d", fRunNumber, runNb));
144 /* init cdb if not yet done */
145 AliCDBManager *cdb = AliCDBManager::Instance();
146 if (!cdb->IsDefaultStorageSet())
147 cdb->SetDefaultStorage("raw://");
148 if (cdb->GetRun() < 0)
151 if (!fTOFcalib->Init()) {
152 AliError("cannot init TOF calib");
156 if (!fGRPManager->ReadGRPEntry()) {
157 AliError("error while reading \"GRPEntry\" from OCDB");
160 fkGRPObject = fGRPManager->GetGRPData();
162 AliError("cannot get \"GRPData\" from GRP manager");
165 fStartTime = fkGRPObject->GetTimeStart();
166 fEndTime = fkGRPObject->GetTimeEnd();
167 AliInfo(Form("got \"GRPData\": startTime=%d, endTime=%d", fStartTime, fEndTime));
169 AliInfo(Form("initialized for run %d", runNb));
173 Int_t useBPTX = fTOFcalib->GetUseLHCClockPhase() ? 1 : 0;
175 /* set histo title with run-number and start-time */
176 fHistoVertexTimestamp->SetTitle(Form("run: %d, startTimestamp: %u, BPTX: %d", fRunNumber, fStartTime, useBPTX));
177 fHistoDeltatTimestamp->SetTitle(Form("run: %d, startTimestamp: %u, BPTX: %d", fRunNumber, fStartTime, useBPTX));
178 fHistoDeltazEta->SetTitle(Form("run: %d, startTimestamp: %u, BPTX: %d", fRunNumber, fStartTime, useBPTX));
179 fHistoDeltatEta->SetTitle(Form("run: %d, startTimestamp: %u, BPTX: %d", fRunNumber, fStartTime, useBPTX));
180 fHistoDeltazCosTheta->SetTitle(Form("run: %d, startTimestamp: %u, BPTX: %d", fRunNumber, fStartTime, useBPTX));
181 fHistoAcceptedTracksEtaPt->SetTitle(Form("run: %d, startTimestamp: %u, BPTX: %d", fRunNumber, fStartTime, useBPTX));
182 fHistoMatchedTracksEtaPt->SetTitle(Form("run: %d, startTimestamp: %u, BPTX: %d", fRunNumber, fStartTime, useBPTX));
187 //_______________________________________________________
190 AliTOFAnalysisTaskCalibPass0::InitEvent()
197 fESDEvent = dynamic_cast<AliESDEvent *>(InputEvent());
198 if (!fESDEvent) return kFALSE;
199 /* set event time and elapsed time */
200 fEventTime = fESDEvent->GetTimeStamp();
201 fElapsedTime = fESDEvent->GetTimeStamp() - fStartTime;
202 /* event selection */
203 if (fEventSelectionFlag && !fEventCuts->IsCollisionCandidate(fESDEvent)) return kFALSE;
204 /* vertex selection */
205 fkVertex = fESDEvent->GetPrimaryVertexTracks();
206 if (fVertexSelectionFlag && fkVertex->GetNContributors() < 1) {
207 fkVertex = fESDEvent->GetPrimaryVertexSPD();
208 if (fkVertex->GetNContributors() < 1) return kFALSE;
210 if (fVertexSelectionFlag && TMath::Abs(fkVertex->GetZ()) > fVertexCut) return kFALSE;
211 /* calibrate ESD if requested */
212 fTOFcalib->CalibrateESD(fESDEvent);
217 //_______________________________________________________
220 AliTOFAnalysisTaskCalibPass0::HasTOFMeasurement(const AliESDtrack *track) const
223 * has TOF measurement
226 /* check TOF status flags */
227 if (!(track->GetStatus() & AliESDtrack::kTOFout) ||
228 !(track->GetStatus() & AliESDtrack::kTIME)) return kFALSE;
229 /* check integrated length */
230 if (track->GetIntegratedLength() < 350.) return kFALSE;
232 /* TOF measurement ok */
236 //_______________________________________________________
239 AliTOFAnalysisTaskCalibPass0::UserCreateOutputObjects()
242 * user create output objects
246 Int_t timeBins = 288;
247 Float_t timeMin = 0.;
248 Float_t timeMax = 24. * 3600.;
250 Int_t vertexBins = 200;
251 Float_t vertexMin = -50.;
252 Float_t vertexMax = 50.;
254 Int_t deltatBins = 2000;
255 Float_t deltatMin = -24400.;
256 Float_t deltatMax = 24400.;
258 Int_t deltazBins = 200;
259 Float_t deltazMin = -10.;
260 Float_t deltazMax = 10.;
263 Float_t etaMin = -1.;
270 fHistoVertexTimestamp = new TH2F("hHistoVertexTimestamp", "Vertex position;elapsed time (s);z (cm);", timeBins, timeMin, timeMax, vertexBins, vertexMin, vertexMax);
271 fHistoList->Add(fHistoVertexTimestamp);
273 fHistoDeltatTimestamp = new TH2F("hHistoDeltatTimestamp", "Global time shift (T0-fill);elapsed time (s);t - t_{exp}^{(#pi)} (ps);", timeBins, timeMin, timeMax, deltatBins, deltatMin, deltatMax);
274 fHistoList->Add(fHistoDeltatTimestamp);
276 fHistoDeltazEta = new TH2F("hHistoDeltazEta", "Matching residuals (longitudinal);#eta;#Deltaz (cm);", etaBins, etaMin, etaMax, deltazBins, deltazMin, deltazMax);
277 fHistoList->Add(fHistoDeltazEta);
279 fHistoDeltatEta = new TH2F("hHistoDeltatEta", "Global time shift (T0-fill) vs #eta; #eta; t - t_{exp}^{(#pi)} (ps);", etaBins, etaMin, etaMax, deltatBins, deltatMin, deltatMax);
280 fHistoList->Add(fHistoDeltatEta);
282 fHistoDeltazCosTheta = new TH2F("hHistoDeltazCosTheta", "Matching residuals (longitudinal);cos #theta;#Deltaz (cm);", etaBins, etaMin, etaMax, deltazBins, deltazMin, deltazMax);
283 fHistoList->Add(fHistoDeltazCosTheta);
285 fHistoAcceptedTracksEtaPt = new TH2F("hHistoAcceptedTracksEtaPt", "Accepted tracks;#eta;#p_{T} (GeV/c);", etaBins, etaMin, etaMax, pBins, pMin, pMax);
286 fHistoList->Add(fHistoAcceptedTracksEtaPt);
288 fHistoMatchedTracksEtaPt = new TH2F("hHistoMatchedTracksEtaPt", "Matched tracks;#eta;p_{T} (GeV/c);", etaBins, etaMin, etaMax, pBins, pMin, pMax);
289 fHistoList->Add(fHistoMatchedTracksEtaPt);
292 PostData(1, fHistoList);
295 //_______________________________________________________
298 AliTOFAnalysisTaskCalibPass0::UserExec(Option_t *)
305 if (!InitRun()) return;
307 if (!InitEvent()) return;
309 /*** ACCEPTED EVENT ***/
311 /* fill vertex histo */
312 fHistoVertexTimestamp->Fill(fElapsedTime, fkVertex->GetZ());
314 /* loop over ESD tracks */
315 Int_t nTracks = fESDEvent->GetNumberOfTracks();
317 Double_t eta, costheta, pt, time, timei[AliPID::kSPECIES], deltat, deltaz;
318 for (Int_t itrk = 0; itrk < nTracks; itrk++) {
320 track = fESDEvent->GetTrack(itrk);
321 if (!track) continue;
322 /* check accept track */
323 if (!fTrackCuts->AcceptTrack(track)) continue;
326 costheta = TMath::Cos(track->Theta());
328 /* fill accepted tracks histo */
329 fHistoAcceptedTracksEtaPt->Fill(eta, pt);
330 /* check TOF measurement */
331 if (!HasTOFMeasurement(track)) continue;
333 /*** ACCEPTED TRACK WITH TOF MEASUREMENT ***/
335 /* fill matched tracks histo */
336 fHistoMatchedTracksEtaPt->Fill(eta, pt);
338 time = track->GetTOFsignal();
339 track->GetIntegratedTimes(timei);
340 deltat = time - timei[AliPID::kPion];
341 deltaz = track->GetTOFsignalDz();
344 fHistoDeltatTimestamp->Fill(fElapsedTime, deltat);
345 fHistoDeltazEta->Fill(eta, deltaz);
346 fHistoDeltatEta->Fill(eta, deltat);
347 fHistoDeltazCosTheta->Fill(costheta, deltaz);
349 } /* end of loop over ESD tracks */
352 PostData(1, fHistoList);
356 //_______________________________________________________
359 AliTOFAnalysisTaskCalibPass0::GetStatus()
367 /* OK, return zero */
372 /* non-fatal error, return negative status */
379 /* fatal error, return positive status */
386 /* anything else, return negative large number */
392 /* should never arrive here, anyway return negative large number */
396 //_______________________________________________________
399 AliTOFAnalysisTaskCalibPass0::PrintStatus()
405 Int_t status = GetStatus();
407 AliInfo(Form("TOF calibration successful: %s (status=%d)", fgkStatusCodeName[fStatus], status));
409 else if (status > 0) {
410 AliInfo(Form("TOF calibration failed: %s (status=%d)", fgkStatusCodeName[fStatus], status));
412 else if (status < 0) {
413 AliInfo(Form("TOF calibration failed (expected): %s (status=%d)", fgkStatusCodeName[fStatus], status));
418 //_______________________________________________________
421 AliTOFAnalysisTaskCalibPass0::ProcessOutput(const Char_t *filename, AliCDBStorage* db)
427 Int_t ret = DoProcessOutput(filename, db);
432 //_______________________________________________________
435 AliTOFAnalysisTaskCalibPass0::DoProcessOutput(const Char_t *filename, AliCDBStorage *db)
441 /* reset status to OK */
445 TFile *file = TFile::Open(filename);
446 if (!file || !file->IsOpen()) {
447 AliError(Form("cannot open output file %s", filename));
448 fStatus = kInputError;
452 TList *list = (TList *)file->Get("TOFHistos");
453 TH2F *histoVertexTimestamp = NULL;
454 TH2F *histoDeltatTimestamp = NULL;
455 TH2F *histoDeltazEta = NULL;
456 TH2F *histoDeltatEta = NULL;
457 TH2F *histoDeltazCosTheta = NULL;
458 TH2F *histoAcceptedTracksEtaPt = NULL;
459 TH2F *histoMatchedTracksEtaPt = NULL;
461 AliInfo(Form("getting histograms from \"Histos\" list from file %s", filename));
462 histoVertexTimestamp = (TH2F *)list->FindObject("hHistoVertexTimestamp");
463 histoDeltatTimestamp = (TH2F *)list->FindObject("hHistoDeltatTimestamp");
464 histoDeltazEta = (TH2F *)list->FindObject("hHistoDeltazEta");
465 histoDeltatEta = (TH2F *)list->FindObject("hHistoDeltatEta");
466 histoDeltazCosTheta = (TH2F *)list->FindObject("hHistoDeltazCosTheta");
467 histoAcceptedTracksEtaPt = (TH2F *)list->FindObject("hHistoAcceptedTracksEtaPt");
468 histoMatchedTracksEtaPt = (TH2F *)list->FindObject("hHistoMatchedTracksEtaPt");
471 AliInfo(Form("getting histograms directly from file %s", filename));
472 histoVertexTimestamp = (TH2F *)file->Get("hHistoVertexTimestamp");
473 histoDeltatTimestamp = (TH2F *)file->Get("hHistoDeltatTimestamp");
474 histoDeltazEta = (TH2F *)file->Get("hHistoDeltazEta");
475 histoDeltatEta = (TH2F *)file->Get("hHistoDeltatEta");
476 histoDeltazCosTheta = (TH2F *)file->Get("hHistoDeltazCosTheta");
477 histoAcceptedTracksEtaPt = (TH2F *)file->Get("hHistoAcceptedTracksEtaPt");
478 histoMatchedTracksEtaPt = (TH2F *)file->Get("hHistoMatchedTracksEtaPt");
481 if (!histoVertexTimestamp) {
482 AliError(Form("cannot get \"hHistoVertexTimestamp\" object from file %s", filename));
483 fStatus = kInputError;
486 if (!histoDeltatTimestamp) {
487 AliError(Form("cannot get \"hHistoDeltatTimestamp\" object from file %s", filename));
488 fStatus = kInputError;
491 if (!histoDeltazEta) {
492 AliError(Form("cannot get \"hHistoDeltazEta\" object from file %s", filename));
493 fStatus = kInputError;
496 if (!histoDeltatEta) {
497 AliError(Form("cannot get \"hHistoDeltatEta\" object from file %s", filename));
498 fStatus = kInputError;
501 if (!histoDeltazCosTheta) {
502 AliError(Form("cannot get \"hHistoDeltazCosTheta\" object from file %s", filename));
503 fStatus = kInputError;
506 if (!histoAcceptedTracksEtaPt) {
507 AliError(Form("cannot get \"hHistoAccptedTracksEtaPt\" object from file %s", filename));
508 fStatus = kInputError;
511 if (!histoMatchedTracksEtaPt) {
512 AliError(Form("cannot get \"hHistoMatchedTracksEtaPt\" object from file %s", filename));
513 fStatus = kInputError;
517 /* check matching performance */
518 if (!CheckMatchingPerformance(histoDeltazEta, histoAcceptedTracksEtaPt, histoMatchedTracksEtaPt)) {
519 AliError("error while checking matching efficiency");
522 /* calibrate and store */
523 if (!CalibrateAndStore(histoVertexTimestamp, histoDeltatTimestamp, db)) {
524 AliError("error while calibrating and storing");
532 //_______________________________________________________
535 AliTOFAnalysisTaskCalibPass0::CheckMatchingPerformance(const TH2F *histoDeltazCosTheta, const TH2F *histoAcceptedTracksEtaPt, const TH2F *histoMatchedTracksEtaPt) const
538 * check matching performance
542 if (!histoDeltazCosTheta || !histoAcceptedTracksEtaPt || !histoMatchedTracksEtaPt)
544 /* dummy for the time being */
548 //_______________________________________________________
551 AliTOFAnalysisTaskCalibPass0::CalibrateAndStore(TH2F *histoVertexTimestamp, TH2F *histoDeltatTimestamp, AliCDBStorage *db)
554 * calibrate and store
558 if (!histoVertexTimestamp || !histoDeltatTimestamp)
561 /*** GET RUN-NUMBER AND START-TIMESTAMP ***/
564 TObjArray *strarr = NULL;
565 TObjString *ostr = NULL;
566 str = histoVertexTimestamp->GetTitle();
567 strarr = str.Tokenize(",");
569 AliError("problems whith tokenize histogram title");
570 fStatus = kDataError;
575 ostr = (TObjString *)strarr->At(0);
577 AliError("problems while getting run number from histogram title");
578 fStatus = kDataError;
582 str = ostr->GetString();
583 if (!str.BeginsWith("run:")) {
584 AliError("problems while getting run number from histogram title");
585 fStatus = kDataError;
590 Int_t runNb = atoi(str.Data());
592 AliError(Form("bad run number: %d", runNb));
593 fStatus = kDataError;
597 AliInfo(Form("got run number: %d", runNb));
599 /* start timestamp */
600 ostr = (TObjString *)strarr->At(1);
602 AliError("problems while getting start timestamp from histogram title");
603 fStatus = kDataError;
607 str = ostr->GetString();
608 str.Remove(0, 1); /* remove empty space at the beginning */
609 if (!str.BeginsWith("startTimestamp:")) {
610 AliError("problems while getting start timestamp from histogram title");
611 fStatus = kDataError;
616 UInt_t startTimestamp = atoi(str.Data());
617 if (startTimestamp <= 0) {
618 AliError(Form("bad start timestamp: %d", startTimestamp));
619 fStatus = kDataError;
623 TTimeStamp ts = startTimestamp;
624 AliInfo(Form("got start timestamp: %d (%s)", startTimestamp, ts.AsString()));
627 ostr = (TObjString *)strarr->At(2);
629 AliError("problems while getting BPTX from histogram title");
630 fStatus = kDataError;
634 str = ostr->GetString();
635 str.Remove(0, 1); /* remove empty space at the beginning */
636 if (!str.BeginsWith("BPTX:")) {
637 AliError("problems while getting BPTX from histogram title");
638 fStatus = kDataError;
643 Bool_t useBPTX = atoi(str.Data());
644 AliInfo(Form("got BPTX: %d", useBPTX));
648 /*** CALIBRATION STAGE ***/
650 /* get fit function */
651 TF1 *fitFunc = (TF1 *)gROOT->GetFunction("gaus");
654 TH1D *histoVertexTimestamppx = histoVertexTimestamp->ProjectionX("histoVertexTimestamppx");
655 TH1D *histoDeltatTimestamppx = histoDeltatTimestamp->ProjectionX("histoDeltatTimestamppx");
657 /* check statistics */
658 if (histoVertexTimestamppx->Integral() < fgMinVertexIntegral ||
659 histoDeltatTimestamppx->Integral() < fgMinDeltatIntegral) {
660 fStatus = kLowStatistics;
664 /* define mix and max time bin */
665 Int_t minBin = histoVertexTimestamppx->FindFirstBinAbove(0);
666 Int_t maxBin = histoVertexTimestamppx->FindLastBinAbove(0);
667 Float_t minTime = histoVertexTimestamppx->GetBinLowEdge(minBin);
668 Float_t maxTime = histoVertexTimestamppx->GetBinLowEdge(maxBin + 1);
669 AliInfo(Form("min/max time defined: %d < t < %d s [%d, %d]", (Int_t)minTime, (Int_t)maxTime, minBin, maxBin));
671 /* loop over time bins */
673 Float_t time[fgkMaxNumberOfPoints], timeerr[fgkMaxNumberOfPoints];
674 Float_t vertexMean[fgkMaxNumberOfPoints], vertexMeanerr[fgkMaxNumberOfPoints];
675 Float_t vertexSigma[fgkMaxNumberOfPoints], vertexSigmaerr[fgkMaxNumberOfPoints];
676 Float_t timeZeroMean[fgkMaxNumberOfPoints], timeZeroMeanerr[fgkMaxNumberOfPoints];
677 Float_t timeZeroSigma[fgkMaxNumberOfPoints], timeZeroSigmaerr[fgkMaxNumberOfPoints];
678 Float_t averageTimeZero = 0.;
679 Float_t averageTimeSigma = 0.;
680 Float_t averageVertexSpread = 0.;
681 for (Int_t ibin = minBin; ibin <= maxBin; ibin++) {
683 /* define time window */
684 Int_t startBin = ibin;
686 while(histoVertexTimestamppx->Integral(startBin, endBin) < fgMinVertexIntegralSample ||
687 histoDeltatTimestamppx->Integral(startBin, endBin) < fgMinDeltatIntegralSample) {
688 if (endBin < maxBin) endBin++;
689 else if (startBin > minBin) startBin--;
692 if (histoVertexTimestamppx->Integral(startBin, endBin) < fgMinVertexIntegral ||
693 histoDeltatTimestamppx->Integral(startBin, endBin) < fgMinDeltatIntegral) continue;
694 Float_t startTime = histoVertexTimestamppx->GetBinLowEdge(startBin);
695 Float_t endTime = histoVertexTimestamppx->GetBinLowEdge(endBin + 1);
696 Float_t vertexIntegral = histoVertexTimestamppx->Integral(startBin, endBin);
697 Float_t deltatIntegral = histoDeltatTimestamppx->Integral(startBin, endBin);
698 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));
701 TH1D *histoVertexTimestamppy = histoVertexTimestamp->ProjectionY("histoVertexTimestamppy", startBin, endBin);
702 TH1D *histoDeltatTimestamppy = histoDeltatTimestamp->ProjectionY("histoDeltatTimestamppy", startBin, endBin);
705 histoVertexTimestamppx->GetXaxis()->SetRange(startBin, endBin);
706 time[nPoints] = histoVertexTimestamppx->GetMean();
707 timeerr[nPoints] = histoVertexTimestamppx->GetMeanError();
710 if (FitPeak(fitFunc, histoVertexTimestamppy, 10., 3., 3.) != 0) {
711 AliError("troubles fitting vertex, skip");
712 delete histoVertexTimestamppy;
713 delete histoDeltatTimestamppy;
716 vertexMean[nPoints] = fitFunc->GetParameter(1);
717 vertexMeanerr[nPoints] = fitFunc->GetParError(1);
718 vertexSigma[nPoints] = fitFunc->GetParameter(2);
719 vertexSigmaerr[nPoints] = fitFunc->GetParError(2);
720 averageVertexSpread += fitFunc->GetParameter(2);
723 if (FitPeak(fitFunc, histoDeltatTimestamppy, 500., 2., 1.) != 0) {
724 AliError("troubles fitting time-zero TRACKS, skip");
725 delete histoVertexTimestamppy;
726 delete histoDeltatTimestamppy;
729 timeZeroMean[nPoints] = fitFunc->GetParameter(1);
730 timeZeroMeanerr[nPoints] = fitFunc->GetParError(1);
731 timeZeroSigma[nPoints] = fitFunc->GetParameter(2);
732 timeZeroSigmaerr[nPoints] = fitFunc->GetParError(2);
733 averageTimeZero += fitFunc->GetParameter(1);
734 averageTimeSigma += fitFunc->GetParameter(2);
736 /* delete projection-y */
737 delete histoVertexTimestamppy;
738 delete histoDeltatTimestamppy;
740 /* increment n points */
743 /* set current bin */
746 } /* end of loop over time bins */
748 /* delete projection-x */
749 delete histoVertexTimestamppx;
750 delete histoDeltatTimestamppx;
754 AliError("no measurement available, quit");
755 fStatus = kNoMeasurement;
758 AliInfo(Form("average time-zero = %f", averageTimeZero / nPoints));
759 AliInfo(Form("average time-sigma = %f", averageTimeSigma / nPoints));
760 AliInfo(Form("average v-spread = %f", averageVertexSpread / nPoints));
762 /*** CREATE RUN PARAMS OBJECT ***/
765 /* get start time from GRP */
766 TGrid::Connect("alien://", 0, 0, "t");
767 AliCDBManager *cdb = AliCDBManager::Instance();
768 cdb->SetDefaultStorage("raw://");
771 if (!grp.ReadGRPEntry()) {
772 AliError("error while reading GRP entry");
775 UInt_t startTimestamp = grp.GetGRPData()->GetTimeStart();
778 AliInfo(Form("got start time from GRP: %s", ts.AsString()));
782 UInt_t timestamp[fgkMaxNumberOfPoints];
783 Float_t t0[fgkMaxNumberOfPoints];
784 Float_t tofReso[fgkMaxNumberOfPoints];
785 Float_t t0Spread[fgkMaxNumberOfPoints];
786 for (Int_t ipoint = 0; ipoint < nPoints; ipoint++) {
787 timestamp[ipoint] = (UInt_t)time[ipoint] + startTimestamp;
788 t0[ipoint] = timeZeroMean[ipoint];
789 tofReso[ipoint] = timeZeroSigma[ipoint];
790 t0Spread[ipoint] = vertexSigma[ipoint] / 2.99792457999999984e-02;
792 UInt_t run[1] = {runNb};
793 UInt_t runFirstPoint[1] = {0};
794 UInt_t runLastPoint[1] = {nPoints - 1};
796 /* create run params object */
797 AliTOFRunParams obj(nPoints, 1);
798 AliInfo(Form("create run params object for run %d with %d points", runNb, nPoints));
799 obj.SetTimestamp(timestamp);
801 obj.SetTOFResolution(tofReso);
802 obj.SetT0Spread(t0Spread);
804 obj.SetRunFirstPoint(runFirstPoint);
805 obj.SetRunLastPoint(runLastPoint);
806 obj.SetUseLHCClockPhase(useBPTX);
808 /*** CREATE OCDB ENTRY ***/
811 AliError("cannot store object because of NULL storage");
812 fStatus = kStoreError;
816 /* install run params object in OCDB */
817 AliCDBId id("TOF/Calib/RunParams", runNb, runNb);
819 md.SetResponsible("Roberto Preghenella");
820 md.SetComment("offline TOF run parameters");
821 md.SetAliRootVersion(gSystem->Getenv("ARVERSION"));
823 if (!db->Put(&obj, id, &md)) {
824 fStatus = kStoreError;
825 AliError(Form("error while putting object in storage %s", db->GetURI().Data()));
833 //_____________________________________________________________
836 AliTOFAnalysisTaskCalibPass0::FitPeak(TF1 *fitFunc, TH1D *h, Float_t startSigma, Float_t nSigmaMin, Float_t nSigmaMax)
842 Double_t fitCent = h->GetBinCenter(h->GetMaximumBin());
843 Double_t fitMin = fitCent - nSigmaMin * startSigma;
844 Double_t fitMax = fitCent + nSigmaMax * startSigma;
845 if (fitMin < h->GetXaxis()->GetXmin()) fitMin = h->GetXaxis()->GetXmin();
846 if (fitMax > h->GetXaxis()->GetXmax()) fitMax = h->GetXaxis()->GetXmax();
847 fitFunc->SetParameter(1, fitCent);
848 fitFunc->SetParameter(2, startSigma);
849 Int_t fitres = h->Fit(fitFunc, "WWq0", "", fitMin, fitMax);
850 if (fitres != 0) return fitres;
851 /* refit with better range */
852 for (Int_t i = 0; i < 3; i++) {
853 fitCent = fitFunc->GetParameter(1);
854 fitMin = fitCent - nSigmaMin * fitFunc->GetParameter(2);
855 fitMax = fitCent + nSigmaMax * fitFunc->GetParameter(2);
856 if (fitMin < h->GetXaxis()->GetXmin()) fitMin = h->GetXaxis()->GetXmin();
857 if (fitMax > h->GetXaxis()->GetXmax()) fitMax = h->GetXaxis()->GetXmax();
858 fitres = h->Fit(fitFunc, "q0", "", fitMin, fitMax);
859 if (fitres != 0) return fitres;