9#include <top/calibration/TOPLocalCalFitter.h>
19#include <framework/logging/Logger.h>
20#include <top/dbobjects/TOPCalChannelT0.h>
27double crystalball_function(
double x,
double alpha,
double n,
double sigma,
double mean)
30 if (sigma < 0.)
return 0.;
31 double z = (x - mean) / sigma;
32 if (alpha < 0) z = -z;
33 double abs_alpha = std::abs(alpha);
35 return std::exp(- 0.5 * z * z);
37 double nDivAlpha = n / abs_alpha;
38 double AA = std::exp(-0.5 * abs_alpha * abs_alpha);
39 double B = nDivAlpha - abs_alpha;
40 double arg = nDivAlpha / (B - z);
41 return AA * std::pow(arg, n);
45double crystalball_pdf(
double x,
double alpha,
double n,
double sigma,
double mean)
48 if (sigma < 0.)
return 0.;
49 if (n <= 1)
return std::numeric_limits<double>::quiet_NaN();
50 double abs_alpha = std::abs(alpha);
51 double C = n / abs_alpha * 1. / (n - 1.) * std::exp(-alpha * alpha / 2.);
52 double D = std::sqrt(M_PI / 2.) * (1. + erf(abs_alpha / std::sqrt(2.)));
53 double N = 1. / (sigma * (C + D));
54 return N * crystalball_function(x, alpha, n, sigma, mean);
60double TTSPDF(
double x,
double time,
double deltaT,
double sigma1,
double sigma2,
double f1)
62 return f1 * TMath::Gaus(x, time, sigma1, kTRUE) + (1 - f1) * TMath::Gaus(x, time + deltaT, sigma2, kTRUE) ;
71double laserPDF(
double* x,
double* p)
75 double fraction = p[2];
76 double deltaTLaser = p[3];
77 double sigmaRatio = p[4];
78 double deltaTTS = p[5];
82 double deltaTExtra = p[8];
83 double sigmaExtra = p[9];
84 double normExtra = p[10];
86 double deltaTBkg = p[11];
87 double sigmaBkg = p[12];
93 double mainPeak = fraction * TTSPDF(x[0], time, deltaTTS, sigma, TMath::Sqrt(sigma * sigma + sigmaRatio * sigmaRatio), f1);
94 double secondaryPeak = (1. - fraction) * TTSPDF(x[0], time + deltaTLaser, deltaTTS, sigma,
95 TMath::Sqrt(sigma * sigma + sigmaRatio * sigmaRatio), f1);
96 double extraPeak = crystalball_pdf(x[0], alpha, n, sigmaExtra, time + deltaTExtra);
97 double background = TMath::Gaus(x[0], time + deltaTBkg + deltaTTS, sigmaBkg, kTRUE);
99 return norm * (mainPeak + secondaryPeak) + normExtra * extraPeak + bkg * background;
106 "Perform the fit of the laser and pulser runs"
118 B2INFO(
"Getting the TTS parameters from " <<
m_TTSData);
130 std::cout <<
"Running in MC mode, not constraints will be set" << std::endl;
155 m_fitTree =
new TTree(
"fitTree",
"fitTree");
250 m_treeTTS->GetEntry(iChannel + 512 * iSlot);
252 double maxpos = h_profile->GetBinCenter(h_profile->GetMaximumBin());
253 h_profile->GetXaxis()->SetRangeUser(maxpos - 1, maxpos + 2.);
256 double integral = h_profile->Integral();
259 TF1 laser = TF1(
"laser", laserPDF, maxpos - 1, maxpos + 2., 16);
262 laser.SetParameter(0, maxpos);
263 laser.SetParLimits(0, maxpos - 0.06, maxpos + 0.06);
266 laser.SetParameter(1, 0.1);
267 laser.SetParLimits(1, 0.05, 0.25);
269 laser.SetParameter(1, 0.02);
270 laser.SetParLimits(1, 0., 0.04);
274 laser.SetParLimits(2, 0.5, 1.);
282 laser.SetParameter(3, -0.3);
283 laser.SetParLimits(3, -0.4, -0.2);
289 laser.FixParameter(5,
m_mean2);
291 laser.FixParameter(6,
m_f1);
293 laser.FixParameter(6, 0);
296 laser.SetParameter(7, integral * 0.005);
297 laser.SetParLimits(7, 0.2 * integral * 0.005, 2.*integral * 0.005);
300 laser.SetParameter(8, 1.);
301 laser.SetParLimits(8, 0.3, 2.);
302 laser.SetParameter(9, 0.2);
303 laser.SetParLimits(9, 0.08, 1.);
304 laser.SetParameter(10, 0.1 * integral * 0.005);
305 laser.SetParLimits(10, 0., 0.2 * integral * 0.005);
307 laser.SetParameter(14, -2.);
308 laser.SetParameter(15, 2.);
309 laser.SetParLimits(15, 1.01, 20.);
312 laser.SetParameter(11, 1.);
313 laser.SetParLimits(11, 0.1, 5.);
314 laser.SetParameter(12, 0.8);
315 laser.SetParLimits(12, 0., 5.);
316 laser.SetParameter(13, 0.01 * integral * 0.005);
317 laser.SetParLimits(13, 0., 0.2 * integral * 0.005);
334 laser.SetParameter(2, 0.8);
335 laser.SetParLimits(2, 0., 1.);
336 laser.SetParameter(3, -0.1);
337 laser.SetParLimits(3, -0.4, -0.);
339 laser.FixParameter(8, 0);
340 laser.FixParameter(9, 0.1);
341 laser.FixParameter(14, -2.);
342 laser.FixParameter(15, 2);
343 laser.FixParameter(11, 1.);
344 laser.FixParameter(12, 0.1);
345 laser.FixParameter(13, 0.);
346 laser.FixParameter(10, 0.);
356 h_profile->Fit(
"laser",
"R L Q");
360 TF1* peak1 =
new TF1(
"peak1", laserPDF, maxpos - 1, maxpos + 2., 16);
361 TF1* peak2 =
new TF1(
"peak2", laserPDF, maxpos - 1, maxpos + 2., 16);
362 TF1* extra =
new TF1(
"extra", laserPDF, maxpos - 1, maxpos + 2., 16);
363 TF1* background =
new TF1(
"background", laserPDF, maxpos - 1, maxpos + 2., 16);
364 for (
int iPar = 0; iPar < 16; iPar++) {
365 peak1->FixParameter(iPar, laser.GetParameter(iPar));
366 peak2->FixParameter(iPar, laser.GetParameter(iPar));
367 extra->FixParameter(iPar, laser.GetParameter(iPar));
368 background->FixParameter(iPar, laser.GetParameter(iPar));
370 peak1->FixParameter(2, 0.);
371 peak1->FixParameter(7, (1 - laser.GetParameter(2))*laser.GetParameter(7));
372 peak1->FixParameter(10, 0.);
373 peak1->FixParameter(13, 0.);
375 peak2->FixParameter(2, 1.);
376 peak2->FixParameter(7, laser.GetParameter(2)*laser.GetParameter(7));
377 peak2->FixParameter(10, 0.);
378 peak2->FixParameter(13, 0.);
380 extra->FixParameter(7, 0.);
381 extra->FixParameter(13, 0.);
383 background->FixParameter(7, 0.);
384 background->FixParameter(10, 0.);
386 h_profile->GetListOfFunctions()->Add(peak1);
387 h_profile->GetListOfFunctions()->Add(peak2);
388 h_profile->GetListOfFunctions()->Add(extra);
389 h_profile->GetListOfFunctions()->Add(background);
401 m_sigma = laser.GetParameter(1);
415 m_chi2 = laser.GetChisquare() / laser.GetNDF();
428 float maxpos = h_profileFirstPulser->GetBinCenter(h_profileFirstPulser->GetMaximumBin());
429 h_profileFirstPulser->GetXaxis()->SetRangeUser(maxpos - 1, maxpos + 1.);
430 if (h_profileFirstPulser->Integral() > 1000) {
431 TF1 pulser1 = TF1(
"pulser1",
"[0]*TMath::Gaus(x, [1], [2], kTRUE)", maxpos - 1, maxpos + 1.);
432 pulser1.SetParameter(0, 1.);
433 pulser1.SetParameter(1, maxpos);
434 pulser1.SetParameter(2, 0.05);
435 h_profileFirstPulser->Fit(
"pulser1",
"R Q");
438 h_profileFirstPulser->Write();
444 maxpos = h_profileSecondPulser->GetBinCenter(h_profileSecondPulser->GetMaximumBin());
445 h_profileSecondPulser->GetXaxis()->SetRangeUser(maxpos - 1, maxpos + 1.);
446 if (h_profileSecondPulser->Integral() > 1000) {
447 TF1 pulser2 = TF1(
"pulser2",
"[0]*TMath::Gaus(x, [1], [2], kTRUE)", maxpos - 1, maxpos + 1.);
448 pulser2.SetParameter(0, 1.);
449 pulser2.SetParameter(1, maxpos);
450 pulser2.SetParameter(2, 0.05);
451 h_profileSecondPulser->Fit(
"pulser2",
"R Q");
454 h_profileSecondPulser->Write();
466 if (
m_chi2 < 4 && m_sigma < 0.2 && m_yieldLaser > 1000) {
478 if (nEntries != 8192) {
479 B2ERROR(
"fitTree does not contain an entry wit a fit result for each channel. Found " << nEntries <<
480 " instead of 8192. Perhaps you tried to run the commonT0 calculation before finishing the fitting?");
487 short nCal[16] = {0};
488 for (
int i = 0; i < nEntries; i++) {
499 channelT0->suppressAverage();
506 for (
int iSlot = 1; iSlot < 17; iSlot++) {
507 B2INFO(
"--> Number of calibrated channels on Slot " << iSlot <<
" : " << nCal[iSlot - 1] <<
"/512");
508 B2INFO(
"--> Cal on ch 1, 256 and 511: " << channelT0->getT0(iSlot, 0) <<
", " << channelT0->getT0(iSlot,
509 257) <<
", " << channelT0->getT0(iSlot, 511));
510 nCalTot += nCal[iSlot - 1];
513 B2RESULT(
"Channel T0 calibration constants imported to database, calibrated channels: " << nCalTot <<
"/ 8192");
520 for (
int i = 0; i < nEntries; i++) {
524 channelT0Branch->Fill();
525 channelT0ErrBranch->Fill();
543 auto hitTree = getObjectPtr<TTree>(
"hitTree");
544 TH2F* h_hitTime =
new TH2F(
"h_hitTime",
" ", 512 * 16, 0., 512 * 16, 22000, -70, 40.);
546 float amplitude, hitTime;
549 hitTree->SetBranchAddress(
"amplitude", &litude);
550 hitTree->SetBranchAddress(
"hitTime", &hitTime);
551 hitTree->SetBranchAddress(
"channel", &channel);
552 hitTree->SetBranchAddress(
"slot", &slot);
553 hitTree->SetBranchAddress(
"refTimeValid", &refTimeValid);
556 std::vector<TH2F*> h_hitTimeLaserHistos = {};
557 for (
int iLowerEdge = 0; iLowerEdge < (int)
m_binEdges.size() - 1; iLowerEdge++) {
558 TH2F* h_hitTimeLaser =
new TH2F((
"h_hitTimeLaser_" + std::to_string(iLowerEdge + 1)).c_str(),
" ", 512 * 16, 0., 512 * 16, 14000,
560 h_hitTimeLaserHistos.push_back(h_hitTimeLaser);
563 for (
unsigned int i = 0; i < hitTree->GetEntries(); i++) {
564 auto onepc = (
unsigned int)(hitTree->GetEntries() / 100);
566 std::cout <<
"processing hit " << i <<
" of " << hitTree->GetEntries() <<
" (" << i / (onepc * 10) <<
" %)" << std::endl;
567 hitTree->GetEntry(i);
569 int iLowerEdge = std::distance(
m_binEdges.cbegin(), it) - 1;
571 if (iLowerEdge >= 0 && iLowerEdge <
static_cast<int>(
m_binEdges.size()) - 1 && refTimeValid)
572 h_hitTimeLaserHistos[iLowerEdge]->Fill(channel + (slot - 1) * 512, hitTime);
573 if (amplitude > 80. && refTimeValid)
574 h_hitTime->Fill(channel + (slot - 1) * 512, hitTime);
580 for (
short iSlot = 0; iSlot < 16; iSlot++) {
581 std::cout <<
"fitting slot " << iSlot + 1 << std::endl;
582 for (
short iChannel = 0; iChannel < 512; iChannel++) {
583 TH1D* h_profile = h_hitTime->ProjectionY((
"profile_" + std::to_string(iSlot + 1) +
"_" + std::to_string(iChannel)).c_str(),
584 iSlot * 512 + iChannel + 1, iSlot * 512 + iChannel + 1);
586 h_profile->GetXaxis()->SetRangeUser(-10, -10);
588 h_profile->GetXaxis()->SetRangeUser(-65,
595 TH1D* h_profileFirstPulser = h_hitTime->ProjectionY((
"profileFirstPulser_" + std::to_string(iSlot + 1) +
"_" + std::to_string(
596 iChannel)).c_str(), iSlot * 512 + iChannel + 1, iSlot * 512 + iChannel + 1);
597 TH1D* h_profileSecondPulser = h_hitTime->ProjectionY((
"profileSecondPulser_" + std::to_string(iSlot + 1) +
"_" + std::to_string(
598 iChannel)).c_str(), iSlot * 512 + iChannel + 1, iSlot * 512 + iChannel + 1);
599 h_profileFirstPulser->GetXaxis()->SetRangeUser(-10, 10);
600 h_profileSecondPulser->GetXaxis()->SetRangeUser(10, 40);
602 fitPulser(h_profileFirstPulser, h_profileSecondPulser);
607 h_profileFirstPulser->Write();
608 h_profileSecondPulser->Write();
620 std::cout <<
"Fitting in bins" << std::endl;
621 for (
int iLowerEdge = 0; iLowerEdge < (int)
m_binEdges.size() - 1; iLowerEdge++) {
626 for (
short iSlot = 0; iSlot < 16; iSlot++) {
627 std::cout <<
" Fitting slot " << iSlot + 1 << std::endl;
628 for (
short iChannel = 0; iChannel < 512; iChannel++) {
629 TH1D* h_profile = h_hitTimeLaserHistos[iLowerEdge]->ProjectionY((
"profile_" + std::to_string(iSlot + 1) +
"_" + std::to_string(
630 iChannel) +
"_" + std::to_string(iLowerEdge)).c_str(),
631 iSlot * 512 + iChannel + 1, iSlot * 512 + iChannel + 1);
633 h_profile->GetXaxis()->SetRangeUser(-10, -10);
635 h_profile->GetXaxis()->SetRangeUser(-65,
Base class for calibration algorithms.
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
void setDescription(const std::string &description)
Set algorithm description (in constructor)
EResult
The result of calibration.
@ c_OK
Finished successfully =0 in Python.
Channel T0 calibration constants for all 512 channels of 16 modules.
short m_fitStatus
Fit quality flag, propagated to the constants.
float m_yieldLaserErr
Statistical error on yield.
float m_binLowerEdge
Lower edge of the amplitude bin in which this fit is performed.
void fitChannel(short, short, TH1 *)
Fits the laser light on one channel.
float m_nExtraConstraints
parameter n of the tail of the extra peak
float m_chi2
Reduced chi2 of the fit.
float m_timeExtraConstraints
Position of the guassian used to describe the extra peak on the timing distribution tail.
float m_f1
Fraction of the first gaussian on the TTS parametrization.
float m_fraction
Fraction of events in the secondary peak.
float m_sigmaExtra
Gaussian sigma of the extra peak in the timing tail
float m_secondPulserSigma
Time resolution from the fit of the first electronic pulse, from a Gaussian fit.
float m_yieldLaserBackground
Integral of the background gaussian.
TOPLocalCalFitter()
Constructor.
float m_sigma
Gaussian time resolution, fitted.
float m_peakTimeMC
Time of the main peak in teh MC simulation, i.e.
std::string m_output
Name of the output file.
float m_timeBackground
Position of the gaussian used to describe the background, w/ respect to peakTime.
void loadMCInfoTrees()
loads the TTS parameters and the MC truth info
float m_deltaTMC
Time difference between the main peak and the secondary peak in the MC simulation.
float m_peakTimeErr
Statistical error on peakTime.
short m_channel
Channel number (0-511)
float m_channelT0Err
Statistical error on channelT0.
std::string m_TTSData
File with the Fit constraints and MC info.
std::vector< float > m_binEdges
Amplitude bins.
float m_sigmaBackground
Sigma of the gaussian used to describe the background.
float m_peakTime
Fitted time of the main (i.e.
float m_fractionErr
Statistical error on fraction.
bool m_isFitInAmplitudeBins
Enables the fit in amplitude bins.
void calculateChennelT0()
Calculates the commonT0 calibration after the fits have been done.
float m_alphaExtra
alpha parameter of the tail of the extra peak.
float m_rms
RMS of the histogram used for the fit.
float m_mean2
Position of the second gaussian of the TTS parametrization with respect to the first one.
std::string m_fitConstraints
File with the TTS parametrization.
TTree * m_timewalkTree
Output of the fitter.
short m_pixelCol
Pixel column.
float m_yieldLaser
Total number of laser hits from the fitting function integral.
TFile * m_inputTTS
File containing m_treeTTS.
float m_nExtra
parameter n of the tail of the extra peak
float m_alphaExtraConstraints
alpha parameter of the tail of the extra peak.
float m_channelT0
Raw, channelT0 calibration, defined as peakTime-peakTimeMC.
float m_sigmaErr
Statistical error on sigma.
std::string m_fitterMode
Fit mode.
float m_binUpperEdge
Upper edge of the amplitude bin in which this fit is performed.
TFile * m_histFile
Output of the fitter.
TTree * m_treeConstraints
Input to the fitter.
float m_sigmaExtraConstraints
Width of the guassian used to describe the extra peak on the timing distribution tail.
float m_firstPulserSigma
Time resolution from the fit of the first electronic pulse, from a Gaussian fit.
float m_deltaT
Time difference between the main peak and the secondary peak.
TTree * m_fitTree
Output of the fitter.
short m_slot
Slot ID (1-16)
float m_f2
Fraction of the second gaussian on the TTS parametrization.
void determineFitStatus()
determines if the constant obtained by the fit are good or not
short m_pixelRow
Pixel row.
float m_secondPulserTime
Average time of the second electronic pulse respect to the reference pulse, from a gaussian fit.
float m_peakTimeConstraints
Time of the main laser peak in the MC simulation (aka MC correction)
float m_sigma1
Width of the first gaussian on the TTS parametrization.
float m_fractionMC
Fraction of events in the secondary peak form the MC simulation.
float m_sigmaBackgroundConstraints
Sigma of the gaussian used to describe the background.
virtual EResult calibrate() override
Runs the algorithm on events.
TTree * m_treeTTS
Input to the fitter.
void setupOutputTreeAndFile()
prepares the output tree
float m_timeBackgroundConstraints
Position of the gaussian used to describe the background, w/ respect to peakTime.
TFile * m_inputConstraints
File containing m_treeConstraints.
float m_deltaTErr
Statistical error on deltaT.
float m_firstPulserTime
Average time of the first electronic pulse respect to the reference pulse, from a Gaussian fit.
float m_histoIntegral
Integral of the fitted histogram.
float m_yieldLaserExtra
Integral of the extra peak.
float m_sigma2
Width of the second gaussian on the TTS parametrization.
void fitPulser(TH1 *, TH1 *)
Fits the two pulsers.
float m_timeExtra
Position of the extra peak seen in the timing tail, w/ respect to peakTime.
float m_fractionConstraints
Fraction of the main peak.
float m_deltaTConstraints
Distance between the main and the secondary laser peak.
Abstract base class for different kinds of events.