9#include <cdc/calibration/CDCdEdx/HadronSaturation.h>
17 m_gamma(0.0019308200),
47 B2INFO(
"\n\t HadronSaturation: filling sample events...\n");
52 std::vector< TString > types;
53 types.push_back(
"pion");
54 types.push_back(
"kaon");
55 types.push_back(
"proton");
56 types.push_back(
"muon");
59 TFile* satFile =
new TFile(infilename);
61 for (
int i = 0; i < int(types.size()); ++i) {
63 TTree* satTree = (TTree*)satFile->Get(types[i]);
64 double satbg, satcosth, satdedx, satdedxerr;
65 satTree->SetBranchAddress(
"bg", &satbg);
66 satTree->SetBranchAddress(
"costh", &satcosth);
67 satTree->SetBranchAddress(
"dedx", &satdedx);
68 satTree->SetBranchAddress(
"dedxerr", &satdedxerr);
71 for (
unsigned int j = 0; j < satTree->GetEntries(); ++j) {
73 if (satdedxerr == 0)
continue;
87 if ((nevents + firstevent) >
int(
m_dedx.size())) nevents =
m_dedx.size();
89 if (firstevent < 0 || (nevents + firstevent) >
int(
m_dedx.size()))
90 B2FATAL(
"HadronSaturation: trying to print events out of range (" <<
m_dedx.size() <<
")");
92 for (
int i = firstevent; i < nevents; ++i)
93 B2INFO(
"\t Event " << i <<
":\t bg = " <<
m_betagamma[i] <<
"\t cos(theta) = " <<
m_costheta[i] <<
"\t dE/dx mean = " <<
101 unsigned int nevt =
m_dedx.size();
102 double chisq = 0, dedxsum = 0;
103 std::vector< double > vdedxavg;
111 for (
unsigned int i = 0; i < nevt; i++) {
113 dedxcor = hadsat.
D2I(
m_costheta[i], hadsat.
I2D(
m_costheta[i], 1.0, alpha, gamma, delta, power, ratio) *
m_dedx[i], alpha, gamma,
114 delta, power, ratio);
127 for (
unsigned int i = 0; i < nevt; i++) {
128 dedxcor = hadsat.
D2I(
m_costheta[i], hadsat.
I2D(
m_costheta[i], 1.0, alpha, gamma, delta, power, ratio) *
m_dedx[i], alpha, gamma,
129 delta, power, ratio);
132 chisq += pow((dedxcor - vdedxavg[j]) /
m_dedxerror[i], 2);
133 B2INFO(
"\t " << i <<
") " << dedxcor <<
"/" << vdedxavg[j] <<
", error was "
135 ": Final " << chisq);
144 result = HC_obj->
myFunction(par[0], par[1], par[2], par[3], par[4]);
151 B2INFO(
"\t Performing the hadron saturation fit...");
154 TFitter* minimizer =
new TFitter(5);
157 minimizer->SetParameter(0,
"alpha",
m_alpha, gRandom->Rndm() * 0.001, -5.0, 5.0);
158 minimizer->SetParameter(1,
"gamma",
m_gamma, gRandom->Rndm() * 0.001, -5.0, 5.0);
159 minimizer->SetParameter(2,
"delta",
m_delta, gRandom->Rndm() * 0.001, 0, 0.50);
160 minimizer->SetParameter(3,
"power",
m_power, gRandom->Rndm() * 0.01, 0.05, 2.5);
161 minimizer->SetParameter(4,
"ratio",
m_ratio, gRandom->Rndm() * 0.01, 0.5, 2);
162 minimizer->FixParameter(2);
163 minimizer->FixParameter(4);
169 minimizer->ExecuteCommand(
"SET START", &strategy, 1);
172 minimizer->ExecuteCommand(
"SET ERR", &up, 1);
176 minimizer->ExecuteCommand(
"SET PRINT", arg, 1);
178 double eps_machine(std::numeric_limits<double>::epsilon());
179 minimizer->ExecuteCommand(
"SET EPS", &eps_machine, 1);
181 double fitpar[5], fiterr[5];
183 for (
int i = 0; i < 30; ++i) {
185 minimizer->SetParameter(0,
"alpha",
m_alpha, gRandom->Rndm() * 0.001, -5.0, 5.0);
186 minimizer->SetParameter(1,
"gamma",
m_gamma, gRandom->Rndm() * 0.001, -5.0, 5.0);
187 minimizer->SetParameter(2,
"delta",
m_delta, gRandom->Rndm() * 0.001, 0, 0.50);
188 minimizer->SetParameter(3,
"power",
m_power, gRandom->Rndm() * 0.01, 0.05, 2.5);
189 minimizer->SetParameter(4,
"ratio",
m_ratio, gRandom->Rndm() * 0.01, 0.5, 2);
190 minimizer->FixParameter(2);
191 minimizer->FixParameter(4);
193 double maxcalls(5000.), tolerance(0.1);
194 double arglist[] = {maxcalls, tolerance};
195 unsigned int nargs(2);
196 minimizer->ExecuteCommand(
"MIGRAD", arglist, nargs);
197 minimizer->ExecuteCommand(
"MIGRAD", arglist, nargs);
198 int status = minimizer->ExecuteCommand(
"HESSE", arglist, nargs);
200 minimizer->PrintResults(1, 0);
201 B2INFO(
"\t iter = " << i <<
": Fit status: " << status);
204 while (status != 0 && counter < 5) {
205 minimizer->SetParameter(0,
"alpha",
m_alpha, gRandom->Rndm() * 0.001, -5.0, 5.0);
206 minimizer->SetParameter(1,
"gamma",
m_gamma, gRandom->Rndm() * 0.001, -5.0, 5.0);
207 minimizer->SetParameter(2,
"delta",
m_delta, gRandom->Rndm() * 0.001, 0, 0.50);
208 minimizer->SetParameter(3,
"power",
m_power, gRandom->Rndm() * 0.01, 0.05, 2.5);
209 minimizer->SetParameter(4,
"ratio",
m_ratio, gRandom->Rndm() * 0.01, 0.5, 2);
210 minimizer->FixParameter(2);
212 minimizer->FixParameter(4);
214 minimizer->ExecuteCommand(
"MIGRAD", arglist, nargs);
215 status = minimizer->ExecuteCommand(
"HESSE", arglist, nargs);
216 B2INFO(
"\t re-Fit iter: " << counter <<
", status code: " << status);
220 B2INFO(
"\t HadronSaturation::ERROR - BAD FIT!");
225 for (
int par = 0; par < 5; ++par) {
226 fitpar[par] = minimizer->GetParameter(par);
227 fiterr[par] = minimizer->GetParError(par);
230 B2INFO(
"\t Final Result: HadronSaturation: fit results");
231 B2INFO(
"\t alpha (" <<
m_alpha <<
"): " << fitpar[0] <<
" +- " << fiterr[0]);
232 B2INFO(
"\t gamma (" <<
m_gamma <<
"): " << fitpar[1] <<
" +- " << fiterr[1]);
233 B2INFO(
"\t delta (" <<
m_delta <<
"): " << fitpar[2] <<
" +- " << fiterr[2]);
234 B2INFO(
"\t power (" <<
m_power <<
"): " << fitpar[3] <<
" +- " << fiterr[3]);
235 B2INFO(
"\t ratio (" <<
m_ratio <<
"): " << fitpar[4] <<
" +- " << fiterr[4]);
239 double chi2, edm, errdef;
241 minimizer->GetStats(chi2, edm, errdef, nvpar, nparx);
242 B2INFO(
"\t\t Fit chi^2: " << chi2);
244 std::ofstream parfile;
245 parfile.open(
"sat-pars.fit.txt");
247 parfile << fitpar[0] << std::endl;
248 parfile << fitpar[1] << std::endl;
249 parfile << fitpar[2] << std::endl;
250 parfile << fitpar[3] << std::endl;
251 parfile << fitpar[4] << std::endl;
Class to hold the hadron saturation functions.
double I2D(double cosTheta, double I) const
hadron saturation parameterization part 2
double D2I(double cosTheta, double D) const
hadron saturation parameterization part 1
Class to perform the hadron saturation calibration.
std::vector< double > m_costheta
a vector to hold cos(theta) values
double m_delta
the delta parameter for the hadron saturation correction
std::vector< double > m_betagamma
a vector to hold beta-gamma values
double m_ratio
the ratio parameter for the hadron saturation correction
std::vector< double > m_dedx
a vector to hold dE/dx measurements
double m_gamma
the gamma parameter for the hadron saturation correction
double myFunction(double alpha, double gamma, double delta, double power, double ratio)
some helper functions for the hadron saturation correction
std::vector< double > m_dedxerror
a vector to hold dE/dx errors
void printEvents(int firstevent, int nevents)
print a sample of events
void fillSample(TString infilename)
fill the vectors below
HadronSaturation()
Constructor: Sets the description, the properties and the parameters of the algorithm.
double m_alpha
the alpha parameter for the hadron saturation correction
void fitSaturation()
perform the hadron saturation fit
void clear()
clear the vectors
double m_power
the power parameter for the hadron saturation correction
static void minuitFunction(int &, double *, double &result, double *para, int)
functions for the hadron saturation correction
int m_cosbins
the number of cosine bins
Abstract base class for different kinds of events.