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 Int_t AliTOFAnalysisTaskCalibPass0::fgkMaxNumberOfPoints = 10000; // max number of points
55 const Double_t AliTOFAnalysisTaskCalibPass0::fgkMinVertexIntegral = 1000.;
56 const Double_t AliTOFAnalysisTaskCalibPass0::fgkMinDeltatIntegral = 20000.;
58 //_______________________________________________________
60 AliTOFAnalysisTaskCalibPass0::AliTOFAnalysisTaskCalibPass0() :
61 AliAnalysisTaskSE("TOFCalib-Pass0"),
63 fEventSelectionFlag(kFALSE),
64 fVertexSelectionFlag(kFALSE),
68 fEventCuts(new AliPhysicsSelection()),
69 fTrackCuts(new AliESDtrackCuts()),
75 fGRPManager(new AliGRPManager()),
77 fTOFcalib(new AliTOFcalib()),
78 fHistoList(new TList()),
79 fHistoVertexTimestamp(NULL),
80 fHistoDeltatTimestamp(NULL),
81 fHistoDeltazEta(NULL),
82 fHistoDeltazCosTheta(NULL),
83 fHistoAcceptedTracksEtaPt(NULL),
84 fHistoMatchedTracksEtaPt(NULL)
91 fHistoList->SetOwner(kTRUE);
92 DefineOutput(1, TList::Class());
95 //_______________________________________________________
97 AliTOFAnalysisTaskCalibPass0::~AliTOFAnalysisTaskCalibPass0()
103 if (fHistoList) delete fHistoList;
106 //_______________________________________________________
109 AliTOFAnalysisTaskCalibPass0::InitRun()
116 fESDEvent = dynamic_cast<AliESDEvent *>(InputEvent());
118 AliError("cannot get ESD event");
122 Int_t runNb = fESDEvent->GetRunNumber();
123 /* check run already initialized */
124 if (fInitFlag && fRunNumber == runNb) return kTRUE;
125 if (fInitFlag && fRunNumber != runNb) {
126 AliFatal(Form("run number has changed within same job and this is not allowed: %d -> %d", fRunNumber, runNb));
130 /* init cdb if not yet done */
131 AliCDBManager *cdb = AliCDBManager::Instance();
132 if (!cdb->IsDefaultStorageSet())
133 cdb->SetDefaultStorage("raw://");
134 if (cdb->GetRun() < 0)
137 if (!fTOFcalib->Init()) {
138 AliError("cannot init TOF calib");
142 if (!fGRPManager->ReadGRPEntry()) {
143 AliError("error while reading \"GRPEntry\" from OCDB");
146 fkGRPObject = fGRPManager->GetGRPData();
148 AliError("cannot get \"GRPData\" from GRP manager");
151 fStartTime = fkGRPObject->GetTimeStart();
152 fEndTime = fkGRPObject->GetTimeEnd();
153 AliInfo(Form("got \"GRPData\": startTime=%d, endTime=%d", fStartTime, fEndTime));
155 AliInfo(Form("initialized for run %d", runNb));
159 /* set histo title with run-number and start-time */
160 fHistoVertexTimestamp->SetTitle(Form("run: %d, startTimestamp: %u", fRunNumber, fStartTime));
161 fHistoDeltatTimestamp->SetTitle(Form("run: %d, startTimestamp: %u", fRunNumber, fStartTime));
162 fHistoDeltazEta->SetTitle(Form("run: %d, startTimestamp: %u", fRunNumber, fStartTime));
163 fHistoDeltazCosTheta->SetTitle(Form("run: %d, startTimestamp: %u", fRunNumber, fStartTime));
164 fHistoAcceptedTracksEtaPt->SetTitle(Form("run: %d, startTimestamp: %u", fRunNumber, fStartTime));
165 fHistoMatchedTracksEtaPt->SetTitle(Form("run: %d, startTimestamp: %u", fRunNumber, fStartTime));
170 //_______________________________________________________
173 AliTOFAnalysisTaskCalibPass0::InitEvent()
180 fESDEvent = dynamic_cast<AliESDEvent *>(InputEvent());
181 if (!fESDEvent) return kFALSE;
182 /* set event time and elapsed time */
183 fEventTime = fESDEvent->GetTimeStamp();
184 fElapsedTime = fESDEvent->GetTimeStamp() - fStartTime;
185 /* event selection */
186 if (fEventSelectionFlag && !fEventCuts->IsCollisionCandidate(fESDEvent)) return kFALSE;
187 /* vertex selection */
188 fkVertex = fESDEvent->GetPrimaryVertexTracks();
189 if (fVertexSelectionFlag && fkVertex->GetNContributors() < 1) {
190 fkVertex = fESDEvent->GetPrimaryVertexSPD();
191 if (fkVertex->GetNContributors() < 1) return kFALSE;
193 if (fVertexSelectionFlag && TMath::Abs(fkVertex->GetZ()) > fVertexCut) return kFALSE;
194 /* calibrate ESD if requested */
195 fTOFcalib->CalibrateESD(fESDEvent);
200 //_______________________________________________________
203 AliTOFAnalysisTaskCalibPass0::HasTOFMeasurement(const AliESDtrack *track) const
206 * has TOF measurement
209 /* check TOF status flags */
210 if (!(track->GetStatus() & AliESDtrack::kTOFout) ||
211 !(track->GetStatus() & AliESDtrack::kTIME)) return kFALSE;
212 /* check integrated length */
213 if (track->GetIntegratedLength() < 350.) return kFALSE;
215 /* TOF measurement ok */
219 //_______________________________________________________
222 AliTOFAnalysisTaskCalibPass0::UserCreateOutputObjects()
225 * user create output objects
229 Int_t timeBins = 288;
230 Float_t timeMin = 0.;
231 Float_t timeMax = 24. * 3600.;
233 Int_t vertexBins = 200;
234 Float_t vertexMin = -50.;
235 Float_t vertexMax = 50.;
237 Int_t deltatBins = 2000;
238 Float_t deltatMin = -24400.;
239 Float_t deltatMax = 24400.;
241 Int_t deltazBins = 200;
242 Float_t deltazMin = -10.;
243 Float_t deltazMax = 10.;
246 Float_t etaMin = -1.;
253 fHistoVertexTimestamp = new TH2F("hHistoVertexTimestamp", "Vertex position;elapsed time (s);z (cm);", timeBins, timeMin, timeMax, vertexBins, vertexMin, vertexMax);
254 fHistoList->Add(fHistoVertexTimestamp);
256 fHistoDeltatTimestamp = new TH2F("hHistoDeltatTimestamp", "Global time shift (T0-fill);elapsed time (s);t - t_{exp}^{(#pi)} (ps);", timeBins, timeMin, timeMax, deltatBins, deltatMin, deltatMax);
257 fHistoList->Add(fHistoDeltatTimestamp);
259 fHistoDeltazEta = new TH2F("hHistoDeltazEta", "Matching residuals (longitudinal);#eta;#Deltaz (cm);", etaBins, etaMin, etaMax, deltazBins, deltazMin, deltazMax);
260 fHistoList->Add(fHistoDeltazEta);
262 fHistoDeltazCosTheta = new TH2F("hHistoDeltazCosTheta", "Matching residuals (longitudinal);cos #theta;#Deltaz (cm);", etaBins, etaMin, etaMax, deltazBins, deltazMin, deltazMax);
263 fHistoList->Add(fHistoDeltazCosTheta);
265 fHistoAcceptedTracksEtaPt = new TH2F("hHistoAcceptedTracksEtaPt", "Accepted tracks;#eta;#p_{T} (GeV/c);", etaBins, etaMin, etaMax, pBins, pMin, pMax);
266 fHistoList->Add(fHistoAcceptedTracksEtaPt);
268 fHistoMatchedTracksEtaPt = new TH2F("hHistoMatchedTracksEtaPt", "Matched tracks;#eta;p_{T} (GeV/c);", etaBins, etaMin, etaMax, pBins, pMin, pMax);
269 fHistoList->Add(fHistoMatchedTracksEtaPt);
272 PostData(1, fHistoList);
275 //_______________________________________________________
278 AliTOFAnalysisTaskCalibPass0::UserExec(Option_t *)
285 if (!InitRun()) return;
287 if (!InitEvent()) return;
289 /*** ACCEPTED EVENT ***/
291 /* fill vertex histo */
292 fHistoVertexTimestamp->Fill(fElapsedTime, fkVertex->GetZ());
294 /* loop over ESD tracks */
295 Int_t nTracks = fESDEvent->GetNumberOfTracks();
297 Double_t eta, costheta, pt, time, timei[AliPID::kSPECIES], deltat, deltaz;
298 for (Int_t itrk = 0; itrk < nTracks; itrk++) {
300 track = fESDEvent->GetTrack(itrk);
301 if (!track) continue;
302 /* check accept track */
303 if (!fTrackCuts->AcceptTrack(track)) continue;
306 costheta = TMath::Cos(track->Theta());
308 /* fill accepted tracks histo */
309 fHistoAcceptedTracksEtaPt->Fill(eta, pt);
310 /* check TOF measurement */
311 if (!HasTOFMeasurement(track)) continue;
313 /*** ACCEPTED TRACK WITH TOF MEASUREMENT ***/
315 /* fill matched tracks histo */
316 fHistoMatchedTracksEtaPt->Fill(eta, pt);
318 time = track->GetTOFsignal();
319 track->GetIntegratedTimes(timei);
320 deltat = time - timei[AliPID::kPion];
321 deltaz = track->GetTOFsignalDz();
324 fHistoDeltatTimestamp->Fill(fElapsedTime, deltat);
325 fHistoDeltazEta->Fill(eta, deltaz);
326 fHistoDeltazCosTheta->Fill(costheta, deltaz);
328 } /* end of loop over ESD tracks */
331 PostData(1, fHistoList);
335 //_______________________________________________________
338 AliTOFAnalysisTaskCalibPass0::ProcessOutput(const Char_t *filename, const Char_t *dbString)
345 TFile *file = TFile::Open(filename);
346 if (!file || !file->IsOpen()) {
347 AliError(Form("cannot open output file %s", filename));
351 TList *list = (TList *)file->Get("Histos");
352 TH2F *histoVertexTimestamp = NULL;
353 TH2F *histoDeltatTimestamp = NULL;
354 TH2F *histoDeltazEta = NULL;
355 TH2F *histoDeltazCosTheta = NULL;
356 TH2F *histoAcceptedTracksEtaPt = NULL;
357 TH2F *histoMatchedTracksEtaPt = NULL;
359 AliInfo(Form("getting histograms from \"Histos\" list from file %s", filename));
360 histoVertexTimestamp = (TH2F *)list->FindObject("hHistoVertexTimestamp");
361 histoDeltatTimestamp = (TH2F *)list->FindObject("hHistoDeltatTimestamp");
362 histoDeltazEta = (TH2F *)list->FindObject("hHistoDeltazEta");
363 histoDeltazCosTheta = (TH2F *)list->FindObject("hHistoDeltazCosTheta");
364 histoAcceptedTracksEtaPt = (TH2F *)list->FindObject("hHistoAcceptedTracksEtaPt");
365 histoMatchedTracksEtaPt = (TH2F *)list->FindObject("hHistoMatchedTracksEtaPt");
368 AliInfo(Form("getting histograms directly from file %s", filename));
369 histoVertexTimestamp = (TH2F *)file->Get("hHistoVertexTimestamp");
370 histoDeltatTimestamp = (TH2F *)file->Get("hHistoDeltatTimestamp");
371 histoDeltazEta = (TH2F *)file->Get("hHistoDeltazEta");
372 histoDeltazCosTheta = (TH2F *)file->Get("hHistoDeltazCosTheta");
373 histoAcceptedTracksEtaPt = (TH2F *)file->Get("hHistoAcceptedTracksEtaPt");
374 histoMatchedTracksEtaPt = (TH2F *)file->Get("hHistoMatchedTracksEtaPt");
377 if (!histoVertexTimestamp) {
378 AliError(Form("cannot get \"hHistoVertexTimestamp\" object from file %s", filename));
381 if (!histoDeltatTimestamp) {
382 AliError(Form("cannot get \"hHistoDeltatTimestamp\" object from file %s", filename));
385 if (!histoDeltazEta) {
386 AliError(Form("cannot get \"hHistoDeltazEta\" object from file %s", filename));
389 if (!histoDeltazCosTheta) {
390 AliError(Form("cannot get \"hHistoDeltazCosTheta\" object from file %s", filename));
393 if (!histoAcceptedTracksEtaPt) {
394 AliError(Form("cannot get \"hHistoAccptedTracksEtaPt\" object from file %s", filename));
397 if (!histoMatchedTracksEtaPt) {
398 AliError(Form("cannot get \"hHistoMatchedTracksEtaPt\" object from file %s", filename));
401 /* check matching performance */
402 if (!CheckMatchingPerformance(histoDeltazEta, histoAcceptedTracksEtaPt, histoMatchedTracksEtaPt)) {
403 AliError("error while checking matching efficiency");
406 /* calibrate and store */
407 if (!CalibrateAndStore(histoVertexTimestamp, histoDeltatTimestamp, dbString)) {
408 AliError("error while calibrating and storing");
416 //_______________________________________________________
419 AliTOFAnalysisTaskCalibPass0::CheckMatchingPerformance(const TH2F *histoDeltazCosTheta, const TH2F *histoAcceptedTracksEtaPt, const TH2F *histoMatchedTracksEtaPt) const
422 * check matching performance
426 if (!histoDeltazCosTheta || !histoAcceptedTracksEtaPt || !histoMatchedTracksEtaPt)
428 /* dummy for the time being */
432 //_______________________________________________________
435 AliTOFAnalysisTaskCalibPass0::CalibrateAndStore(TH2F *histoVertexTimestamp, TH2F *histoDeltatTimestamp, const Char_t *dbString)
438 * calibrate and store
442 if (!histoVertexTimestamp || !histoDeltatTimestamp)
445 /*** GET RUN-NUMBER AND START-TIMESTAMP ***/
448 TObjArray *strarr = NULL;
449 TObjString *ostr = NULL;
450 str = histoVertexTimestamp->GetTitle();
451 strarr = str.Tokenize(",");
453 AliError("problems whith tokenize histogram title");
458 ostr = (TObjString *)strarr->At(0);
460 AliError("problems while getting run number from histogram title");
463 str = ostr->GetString();
464 if (!str.BeginsWith("run:")) {
465 AliError("problems while getting run number from histogram title");
469 Int_t runNb = atoi(str.Data());
471 AliError(Form("bad run number: %d", runNb));
474 AliInfo(Form("got run number: %d", runNb));
476 /* start timestamp */
477 ostr = (TObjString *)strarr->At(1);
479 AliError("problems while getting start timestamp from histogram title");
482 str = ostr->GetString();
483 str.Remove(0, 1); /* remove empty space at the beginning */
484 if (!str.BeginsWith("startTimestamp:")) {
485 AliError("problems while getting start timestamp from histogram title");
489 UInt_t startTimestamp = atoi(str.Data());
490 if (startTimestamp <= 0) {
491 AliError(Form("bad start timestamp: %d", startTimestamp));
494 TTimeStamp ts = startTimestamp;
495 AliInfo(Form("got start timestamp: %d (%s)", startTimestamp, ts.AsString()));
497 /*** CALIBRATION STAGE ***/
499 /* get fit function */
500 TF1 *fitFunc = (TF1 *)gROOT->GetFunction("gaus");
503 TH1D *histoVertexTimestamppx = histoVertexTimestamp->ProjectionX("histoVertexTimestamppx");
504 TH1D *histoDeltatTimestamppx = histoDeltatTimestamp->ProjectionX("histoDeltatTimestamppx");
506 /* define mix and max time bin */
507 Int_t minBin = histoVertexTimestamppx->FindFirstBinAbove(0);
508 Int_t maxBin = histoVertexTimestamppx->FindLastBinAbove(0);
509 Float_t minTime = histoVertexTimestamppx->GetBinLowEdge(minBin);
510 Float_t maxTime = histoVertexTimestamppx->GetBinLowEdge(maxBin + 1);
511 AliInfo(Form("min/max time defined: %d < t < %d s [%d, %d]", (Int_t)minTime, (Int_t)maxTime, minBin, maxBin));
513 /* loop over time bins */
515 Float_t time[fgkMaxNumberOfPoints], timeerr[fgkMaxNumberOfPoints];
516 Float_t vertexMean[fgkMaxNumberOfPoints], vertexMeanerr[fgkMaxNumberOfPoints];
517 Float_t vertexSigma[fgkMaxNumberOfPoints], vertexSigmaerr[fgkMaxNumberOfPoints];
518 Float_t timeZeroMean[fgkMaxNumberOfPoints], timeZeroMeanerr[fgkMaxNumberOfPoints];
519 Float_t timeZeroSigma[fgkMaxNumberOfPoints], timeZeroSigmaerr[fgkMaxNumberOfPoints];
520 Float_t averageTimeZero = 0.;
521 Float_t averageTimeSigma = 0.;
522 Float_t averageVertexSpread = 0.;
523 for (Int_t ibin = minBin; ibin <= maxBin; ibin++) {
525 /* define time window */
526 Int_t startBin = ibin;
528 while(histoVertexTimestamppx->Integral(startBin, endBin) < fgkMinVertexIntegral ||
529 histoDeltatTimestamppx->Integral(startBin, endBin) < fgkMinDeltatIntegral) {
530 if (endBin < maxBin) endBin++;
531 else if (startBin > minBin) startBin--;
534 if (histoVertexTimestamppx->Integral(startBin, endBin) <= 0 ||
535 histoDeltatTimestamppx->Integral(startBin, endBin) <= 0) continue;
536 Float_t startTime = histoVertexTimestamppx->GetBinLowEdge(startBin);
537 Float_t endTime = histoVertexTimestamppx->GetBinLowEdge(endBin + 1);
538 Float_t vertexIntegral = histoVertexTimestamppx->Integral(startBin, endBin);
539 Float_t deltatIntegral = histoDeltatTimestamppx->Integral(startBin, endBin);
540 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));
543 TH1D *histoVertexTimestamppy = histoVertexTimestamp->ProjectionY("histoVertexTimestamppy", startBin, endBin);
544 TH1D *histoDeltatTimestamppy = histoDeltatTimestamp->ProjectionY("histoDeltatTimestamppy", startBin, endBin);
547 histoVertexTimestamppx->GetXaxis()->SetRange(startBin, endBin);
548 time[nPoints] = histoVertexTimestamppx->GetMean();
549 timeerr[nPoints] = histoVertexTimestamppx->GetMeanError();
552 if (FitPeak(fitFunc, histoVertexTimestamppy, 10., 3., 3.) != 0) {
553 AliError("troubles fitting vertex, skip");
554 delete histoVertexTimestamppy;
555 delete histoDeltatTimestamppy;
558 vertexMean[nPoints] = fitFunc->GetParameter(1);
559 vertexMeanerr[nPoints] = fitFunc->GetParError(1);
560 vertexSigma[nPoints] = fitFunc->GetParameter(2);
561 vertexSigmaerr[nPoints] = fitFunc->GetParError(2);
562 averageVertexSpread += fitFunc->GetParameter(2);
565 if (FitPeak(fitFunc, histoDeltatTimestamppy, 500., 2., 1.) != 0) {
566 AliError("troubles fitting time-zero TRACKS, skip");
567 delete histoVertexTimestamppy;
568 delete histoDeltatTimestamppy;
571 timeZeroMean[nPoints] = fitFunc->GetParameter(1);
572 timeZeroMeanerr[nPoints] = fitFunc->GetParError(1);
573 timeZeroSigma[nPoints] = fitFunc->GetParameter(2);
574 timeZeroSigmaerr[nPoints] = fitFunc->GetParError(2);
575 averageTimeZero += fitFunc->GetParameter(1);
576 averageTimeSigma += fitFunc->GetParameter(2);
578 /* delete projection-y */
579 delete histoVertexTimestamppy;
580 delete histoDeltatTimestamppy;
582 /* increment n points */
585 /* set current bin */
588 } /* end of loop over time bins */
590 /* delete projection-x */
591 delete histoVertexTimestamppx;
592 delete histoDeltatTimestamppx;
596 AliError("no measurement available, quit");
599 AliInfo(Form("average time-zero = %f", averageTimeZero / nPoints));
600 AliInfo(Form("average time-sigma = %f", averageTimeSigma / nPoints));
601 AliInfo(Form("average v-spread = %f", averageVertexSpread / nPoints));
603 /*** CREATE RUN PARAMS OBJECT ***/
606 /* get start time from GRP */
607 TGrid::Connect("alien://", 0, 0, "t");
608 AliCDBManager *cdb = AliCDBManager::Instance();
609 cdb->SetDefaultStorage("raw://");
612 if (!grp.ReadGRPEntry()) {
613 AliError("error while reading GRP entry");
616 UInt_t startTimestamp = grp.GetGRPData()->GetTimeStart();
619 AliInfo(Form("got start time from GRP: %s", ts.AsString()));
623 UInt_t timestamp[fgkMaxNumberOfPoints];
624 Float_t t0[fgkMaxNumberOfPoints];
625 Float_t tofReso[fgkMaxNumberOfPoints];
626 Float_t t0Spread[fgkMaxNumberOfPoints];
627 for (Int_t ipoint = 0; ipoint < nPoints; ipoint++) {
628 timestamp[ipoint] = (UInt_t)time[ipoint] + startTimestamp;
629 t0[ipoint] = timeZeroMean[ipoint];
630 tofReso[ipoint] = timeZeroSigma[ipoint];
631 t0Spread[ipoint] = vertexSigma[ipoint] / 2.99792457999999984e-02;
633 UInt_t run[1] = {runNb};
634 UInt_t runFirstPoint[1] = {0};
635 UInt_t runLastPoint[1] = {nPoints - 1};
637 /* create run params object */
638 AliTOFRunParams obj(nPoints, 1);
639 AliInfo(Form("create run params object for run %d with %d points", runNb, nPoints));
640 obj.SetTimestamp(timestamp);
642 obj.SetTOFResolution(tofReso);
643 obj.SetT0Spread(t0Spread);
645 obj.SetRunFirstPoint(runFirstPoint);
646 obj.SetRunLastPoint(runLastPoint);
648 /*** CREATE OCDB ENTRY ***/
651 AliError("cannot store object because of NULL string");
655 /* install run params object in OCDB */
656 AliCDBManager *cdb = AliCDBManager::Instance();
657 AliCDBStorage *sto = cdb->GetStorage(dbString);
659 AliError(Form("cannot get storage %s", dbString));
662 AliCDBId id("TOF/Calib/RunParams", runNb, runNb);
664 md.SetResponsible("Roberto Preghenella");
665 md.SetComment("offline TOF run parameters");
666 md.SetAliRootVersion(gSystem->Getenv("ARVERSION"));
668 if (!sto->Put(&obj, id, &md)) {
669 AliError(Form("error while putting object in storage %s", dbString));
676 //_____________________________________________________________
679 AliTOFAnalysisTaskCalibPass0::FitPeak(TF1 *fitFunc, TH1D *h, Float_t startSigma, Float_t nSigmaMin, Float_t nSigmaMax)
685 Double_t fitCent = h->GetBinCenter(h->GetMaximumBin());
686 Double_t fitMin = fitCent - nSigmaMin * startSigma;
687 Double_t fitMax = fitCent + nSigmaMax * startSigma;
688 if (fitMin < h->GetXaxis()->GetXmin()) fitMin = h->GetXaxis()->GetXmin();
689 if (fitMax > h->GetXaxis()->GetXmax()) fitMax = h->GetXaxis()->GetXmax();
690 fitFunc->SetParameter(1, fitCent);
691 fitFunc->SetParameter(2, startSigma);
692 Int_t fitres = h->Fit(fitFunc, "WWq0", "", fitMin, fitMax);
693 if (fitres != 0) return fitres;
694 /* refit with better range */
695 for (Int_t i = 0; i < 3; i++) {
696 fitCent = fitFunc->GetParameter(1);
697 fitMin = fitCent - nSigmaMin * fitFunc->GetParameter(2);
698 fitMax = fitCent + nSigmaMax * fitFunc->GetParameter(2);
699 if (fitMin < h->GetXaxis()->GetXmin()) fitMin = h->GetXaxis()->GetXmin();
700 if (fitMax > h->GetXaxis()->GetXmax()) fitMax = h->GetXaxis()->GetXmax();
701 fitres = h->Fit(fitFunc, "q0", "", fitMin, fitMax);
702 if (fitres != 0) return fitres;