ATLAS Offline Software
Loading...
Searching...
No Matches
ADCMTHistos.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7#include <cmath>
8
9#include "TDirectory.h"
10#include "TF1.h"
11#include "TH1.h"
12
13namespace MuonCalib {
14
17 inline Double_t adcfitf_skewnormal(Double_t *x, Double_t *par) {
18 // par[0] = skew gauss norm
19 // par[1] = skew gauss mean (i.e. mu)
20 // par[2] = skew gauss sigma (i.e sigma)
21 // par[3] = skew factor (i.e. alpha)
22 // Numeric constants
23 Double_t invsq2pi = 1.0 / (std::sqrt(2 * M_PI));
24 Double_t twoPi = 2 * M_PI;
25
26 Double_t delta_value = par[3] / (std::sqrt(1. + par[3] * par[3]));
27 Double_t omega_square = (par[2] * par[2]) / (1. - 4. * delta_value * delta_value / (twoPi));
28 Double_t omega_value = std::sqrt(omega_square);
29 Double_t xi_value = par[1] - delta_value * omega_value * 2. * invsq2pi;
30 Double_t Gauss_part =
31 (invsq2pi / omega_value) * std::exp(-((x[0] - xi_value) * (x[0] - xi_value)) / (2.0 * omega_square)); // phi(x)
32
33 Double_t Erf_part = 0.5 * (1 + std::erf(par[3] * (x[0] - xi_value) / (omega_value)));
34
35 Double_t SkewNormal_Value = par[0] * 2. * Gauss_part * Erf_part;
36
37 Double_t MyGaussFuncValue = SkewNormal_Value;
38
39 return (MyGaussFuncValue);
40 }
41
42 void ADCMTHistos::FillA(double a) { m_adc->Fill(static_cast<Axis_t>(a)); }
43
45 // Initialize(int id, const T0MTSettings & settings) //
47
48 void ADCMTHistos::Initialize(int id, const T0MTSettings *settings, const char *hname) {
49 m_settings = settings;
50 char buf[100];
51 if (hname == nullptr)
52 snprintf(buf, 100, "adc_spec_%d", id);
53 else
54 snprintf(buf, 100, "adc_spec_%s", hname);
55 m_adc = new TH1F(buf, "", settings->NBinsADC(), settings->ADCMin(), settings->ADCMax());
56 m_id = id;
57 m_adc_ok = false;
58 }
59
61 m_adc_fit = std::make_unique<TF1>("adc_fun", adcfitf_skewnormal, 50, 320, 4);
62 Double_t average = m_adc->GetMean();
63 Double_t max = m_adc->GetMaximum();
64 m_adc_fit->SetParameters(max, average, 40, 1); // initialize value
65 m_adc_fit->SetLineColor(kRed);
66 m_adc_fit->SetParNames("Max", "Mean", "Sigma", "Skew_factor");
67 m_adc_fit->SetParLimits(0, 0, 1000000);
68 m_adc_fit->SetParLimits(1, 100, 200);
69 m_adc_fit->SetParLimits(2, 0, 200);
70 m_adc->Fit("adc_fun", "Q+", "", 50, 320);
71 m_adc_ok = true;
72 return true;
73 }
74
75} // namespace MuonCalib
#define M_PI
static Double_t a
#define x
#define max(a, b)
Definition cfImp.cxx:41
bool m_adc_ok
is true if the adc-fit is ok;
Definition ADCMTHistos.h:96
void Initialize(int id, const T0MTSettings *settings, const char *hname=NULL)
Initialize class.
std::unique_ptr< TF1 > m_adc_fit
function which is fitted to the adc-spectrum
Definition ADCMTHistos.h:94
void FillA(double a)
fill adc value
TH1F * m_adc
pulse height spectrum
Definition ADCMTHistos.h:90
const T0MTSettings * m_settings
Pointer to settings class.
Definition ADCMTHistos.h:98
Settings for the T0 calibration (histogram booking and fitting parameters) Parameters for pattern rec...
double ADCMax() const
const int & NBinsADC() const
Number of bins for ADC histogram and range.
double ADCMin() const
CscCalcPed - algorithm that finds the Cathode Strip Chamber pedestals from an RDO.
Double_t adcfitf_skewnormal(Double_t *x, Double_t *par)
skew normal ADC fitting