10#include <ecl/calibration/eclWaveformTemplateCalibrationC2Algorithm.h>
13#include <ecl/dataobjects/ECLElementNumbers.h>
14#include <ecl/digitization/OfflineFitFunction.h>
15#include <ecl/dbobjects/ECLDigitWaveformParameters.h>
16#include <ecl/dbobjects/ECLCrystalCalib.h>
19#include <framework/database/DBObjPtr.h>
33using namespace Calibration;
40 "Perform the photon template shape calibration using waveforms from high energy crystals from e+e- --> gamma gamma events"
49 std::vector<TF1*> FitFunctions;
50 const double numberofADCPoints = 31.0;
52 double fitf(
double* x,
double* par)
55 double xtoeval = std::fmod(x[0], numberofADCPoints);
56 int whichFitFunctions = x[0] / numberofADCPoints;
58 for (
int i = 0; i < (int)FitFunctions.size(); i++) {
59 FitFunctions[i]->SetParameter(0, par[i]);
60 FitFunctions[i]->SetParameter(1, par[FitFunctions.size() + i]);
61 FitFunctions[i]->SetParameter(2, par[(2 * FitFunctions.size()) + i]);
62 FitFunctions[i]->FixParameter(3, 0);
63 for (
int k = 0; k < 10; k++) {
64 FitFunctions[i]->SetParameter(4 + k, par[(3 * FitFunctions.size()) + k]);
68 return FitFunctions[whichFitFunctions]->Eval(xtoeval * 0.5);
76 B2INFO(
"Reading ECLCrystalCalib payload: eclWaveformTemplateCalibrationC1MaxResLimit");
77 DBObjPtr<ECLCrystalCalib> existingeclWaveformTemplateCalibrationC1MaxResLimit(
"eclWaveformTemplateCalibrationC1MaxResLimit");
79 ExpRun chosenRun = runs.front();
86 std::vector<double> cellIDArray;
87 std::vector<double> maxResidualArray;
88 std::vector<double> limitResidualArray;
89 std::vector<double> parLimitFactorArray;
92 TFile* histfile =
new TFile(
m_outputName.c_str(),
"recreate");
95 TFile* f_PhotonTemplateOutput =
new TFile(Form(
"PhotonShapes_Low%d_High%d.root",
m_firstCellID,
m_lastCellID),
"RECREATE");
96 TTree* mtree =
new TTree(
"mtree",
"");
97 std::vector<double> PhotonWaveformArray(100000);
98 mtree->Branch(
"PhotonArray", PhotonWaveformArray.data(),
"PhotonWaveformArray[100000]/D");
101 auto tree = getObjectPtr<TTree>(
"tree");
103 tree->SetBranchAddress(
"CellID", &CellID);
107 tree->SetBranchAddress(Form(
"ADC%d", i), &Waveform[i]);
111 std::time_t t = std::time(0);
114 int AttemptCounter = 0;
121 double ParMin11t[11];
127 if (CellID_i > 7776 || CellID_i < 1153) {
128 ParMin11t[0] = 20.3216;
129 ParMin11t[1] = -0.0206266;
130 ParMin11t[2] = 0.313928;
131 ParMin11t[3] = 0.589646;
132 ParMin11t[4] = 0.455526;
133 ParMin11t[5] = 1.03656;
134 ParMin11t[6] = 0.000822467;
135 ParMin11t[7] = 45.1574;
136 ParMin11t[8] = 0.716034;
137 ParMin11t[9] = 0.616753;
138 ParMin11t[10] = 0.0851222;
140 ParMin11t[0] = 24.6176;
141 ParMin11t[1] = 0.00725002;
142 ParMin11t[2] = 0.601578;
143 ParMin11t[3] = 0.491976;
144 ParMin11t[4] = 0.601034;
145 ParMin11t[5] = 0.601684;
146 ParMin11t[6] = -0.0103788;
147 ParMin11t[7] = 2.22615;
148 ParMin11t[8] = 0.671294;
149 ParMin11t[9] = 0.529878;
150 ParMin11t[10] = 0.0757927;
153 double resLimit = 2 * existingeclWaveformTemplateCalibrationC1MaxResLimit->getCalibVector()[CellID_i -
155 double resLimitOriginal = resLimit;
158 std::vector<int> EntriesToSkip;
161 double maxResidual = 1000.0;
165 while (PASS ==
false) {
168 std::vector<double> xValuesToFit;
169 std::vector<double> yValuesToFit;
172 std::vector<double> guessBaseline;
173 std::vector<double> guessAmp;
174 std::vector<double> guessTime;
177 std::vector<int> NtupleEntries;
180 int counterWaveforms = 0;
182 for (
int i = 0; i < tree->GetEntries(); i++) {
185 bool skipEvent =
false;
186 for (
int k = 0; k < (int)EntriesToSkip.size(); k++) {
187 if (EntriesToSkip[k] == i) skipEvent =
true;
189 if (skipEvent)
continue;
193 if (CellID != CellID_i)
continue;
198 xValuesToFit.push_back(counter);
199 yValuesToFit.push_back(Waveform[j]);
200 if (Waveform[j] > maxval) {
201 maxval = Waveform[j];
208 guessBaseline.push_back(Waveform[0]);
209 guessAmp.push_back(maxval);
210 guessTime.push_back((maxIndex - 4.5) * 0.5);
212 NtupleEntries.push_back(i);
213 B2INFO(
"Entry: " << i);
222 B2INFO(
"CellID " << CellID_i <<
" counterWaveforms = " << counterWaveforms);
226 B2INFO(
"eclWaveformTemplateCalibrationC2Algorithm: warning total entries for cell ID " << CellID_i <<
" is only: " <<
228 EntriesToSkip.clear();
230 B2INFO(
"eclWaveformTemplateCalibrationC2Algorithm: warning " << CellID_i <<
" resLimit is doubled to" << resLimit <<
" start was "
231 << resLimitOriginal);
236 auto gWaveformToFit =
new TGraph(xValuesToFit.size(), xValuesToFit.data(), yValuesToFit.data());
237 gWaveformToFit->SetName(Form(
"gWaveformToFit_%d",
int(CellID_i)));
242 FitFunctions.clear();
243 for (
int i = 0; i < counterWaveforms; i++) {
244 FitFunctions.push_back(
new TF1(Form(
"Shp_%d", i), Belle2::ECL::WaveFuncTwoComponent, 0, 30.5, 26));
245 FitFunctions[i]->SetNpx(10000);
246 FitFunctions[i]->FixParameter(3, 0);
247 for (
int k = 0; k < 10; k++) {
248 FitFunctions[i]->SetParameter(4 + k, ParMin11t[k + 1]);
249 FitFunctions[i]->FixParameter(10 + 4 + k, ParMin11t[k + 1]);
251 FitFunctions[i]->FixParameter(24, ParMin11t[0]);
252 FitFunctions[i]->FixParameter(25, 1);
256 TF1* TotalFitFunction =
new TF1(
"TotalFitFunction", fitf, 0, counterWaveforms *
m_NumberofADCPoints,
257 (3 * FitFunctions.size()) + 10);
260 int FFsize = FitFunctions.size();
261 for (
int i = 0; i < FFsize; i++) {
262 TotalFitFunction->SetParameter(i, guessTime[i]);
263 TotalFitFunction->SetParameter(FFsize + i, guessBaseline[i]);
264 TotalFitFunction->SetParameter((2 * FFsize) + i, guessAmp[i]);
265 for (
int k = 0; k < 10; k++) {
266 TotalFitFunction->SetParameter((3 * FFsize) + k, ParMin11t[k + 1]);
268 TotalFitFunction->SetParLimits((3 * FFsize) + k, ParMin11t[k + 1] -
m_ParamLimitFactor * fabs(ParMin11t[k + 1]),
271 TotalFitFunction->ReleaseParameter((3 * FFsize) + k);
277 gWaveformToFit->Fit(
"TotalFitFunction",
"Q M W N 0 R",
"", 0, counterWaveforms *
m_NumberofADCPoints);
280 std::vector<int> FitResultY;
281 std::vector<int> FitResultX;
282 int maxResidualWaveformID = 0;
286 double npts = xValuesToFit.size();
287 double maxResidualOld = 0.0;
288 for (
int k = 0; k < npts; k++) {
289 double xVal = xValuesToFit[k];
290 double yVal = TotalFitFunction->Eval(xVal);
291 FitResultX.push_back(xVal);
292 FitResultY.push_back(yVal);
293 double diff = fabs(yValuesToFit[k] - yVal);
294 if (diff > maxResidual) {
297 maxResidualOld = fabs(yValuesToFit[k] / yVal);
302 if (maxResidual > resLimit) {
304 B2INFO(
"FAIL: CellID_i " << CellID_i <<
" maxResidual " << maxResidual <<
" removing entry: " <<
305 NtupleEntries[maxResidualWaveformID] <<
306 " which was waveform number " << maxResidualWaveformID <<
" resLimit was " << resLimit <<
" , resLimit started at " <<
308 B2INFO(
"Old maxResidual of Data/Fit was " << maxResidualOld);
310 B2INFO(
"Iter Time = " << std::time(0) - t << std::endl);
313 std::cout <<
"FAIL: CellID_i " << CellID_i <<
" maxResidual " << maxResidual <<
" removing entry: " <<
314 NtupleEntries[maxResidualWaveformID] <<
315 " which was waveform number " << maxResidualWaveformID <<
" resLimit was " << resLimit <<
" , resLimit started at " <<
316 resLimitOriginal << std::endl;
317 std::cout <<
"wave = [";
318 for (
int k = 0; k < npts; k++) {
319 std::cout << yValuesToFit[k];
320 if (k < (npts - 1)) {
323 std::cout <<
"]" << std::endl;
326 std::cout <<
"fitRes = [";
327 for (
int k = 0; k < npts; k++) {
328 std::cout << TotalFitFunction->Eval(xValuesToFit[k]);
329 if (k < (npts - 1)) {
332 std::cout <<
"]" << std::endl;
337 EntriesToSkip.push_back(NtupleEntries[maxResidualWaveformID]);
349 B2INFO(
"AttemptCounter reach limit: " << AttemptCounter <<
" counterWaveforms: " << counterWaveforms);
353 EntriesToSkip.clear();
359 B2INFO(
"Increasing resLimit to " << resLimit);
366 B2INFO(
"PASS: CellID_i " << CellID_i <<
" maxResidual " << maxResidual <<
" number of waveforms used was " << counterWaveforms <<
367 " resLimit was " << resLimit);
371 limitResidualArray.push_back(resLimit);
378 auto gFitResult =
new TGraph(FitResultX.size(), FitResultX.data(), FitResultY.data());
379 gFitResult->SetName(Form(
"gFitResult_%d",
int(CellID_i)));
382 cellIDArray.push_back(CellID_i);
383 maxResidualArray.push_back(maxResidual);
387 gWaveformToFit->Write();
391 float tempPhotonPar11[11];
392 tempPhotonPar11[0] = ParMin11t[0];
393 for (
unsigned int k = 0; k < 10; k++) tempPhotonPar11[k + 1] = TotalFitFunction->GetParameter((3 * FFsize) + k);
396 FitFunctions[0]->SetParameter(0, 0);
397 FitFunctions[0]->SetParameter(1, 0);
398 FitFunctions[0]->SetParameter(2, 1);
399 for (
int k = 0; k < 10; k++) {
400 FitFunctions[0]->SetParameter(4 + k, tempPhotonPar11[k + 1]);
401 FitFunctions[0]->SetParameter(10 + 4 + k, tempPhotonPar11[k + 1]);
403 FitFunctions[0]->FixParameter(24, ParMin11t[0]);
404 FitFunctions[0]->FixParameter(25, 1);
407 double MaxVal = -1.0;
408 const double cnpts = 2000;
409 for (
int k = 0; k < cnpts; k++) {
411 double yVal = FitFunctions[0]->Eval(xVal);
412 if (yVal > MaxVal) MaxVal = yVal;
414 B2INFO(
"MaxVal " << MaxVal);
415 tempPhotonPar11[0] /= MaxVal;
416 FitFunctions[0]->FixParameter(24, tempPhotonPar11[0]);
419 PhotonParameters->
setTemplateParameters(CellID_i, tempPhotonPar11, tempPhotonPar11, tempPhotonPar11);
422 for (
unsigned int k = 0; k < PhotonWaveformArray.size();
423 k++) PhotonWaveformArray[k] = FitFunctions[0]->Eval(((
double)k) * (1. / 1000.)) ;
427 for (
int w = 0; w < (int)FitFunctions.size(); w++) FitFunctions[w]->Delete();
428 TotalFitFunction->Delete() ;
429 gWaveformToFit->Delete();
435 auto gmaxResidual =
new TGraph(cellIDArray.size(), cellIDArray.data(), maxResidualArray.data());
436 gmaxResidual->SetName(
"gmaxResidual");
437 auto glimitResidualArray =
new TGraph(cellIDArray.size(), cellIDArray.data(), limitResidualArray.data());
438 glimitResidualArray->SetName(
"glimitResidualArray");
439 auto gparLimitFactorArray =
new TGraph(cellIDArray.size(), cellIDArray.data(), parLimitFactorArray.data());
440 gparLimitFactorArray->SetName(
"gparLimitFactorArray");
442 gmaxResidual->Write();
443 glimitResidualArray->Write();
444 gparLimitFactorArray->Write();
449 f_PhotonTemplateOutput->cd();
451 f_PhotonTemplateOutput->Write();
452 f_PhotonTemplateOutput->Close();
453 delete f_PhotonTemplateOutput;
457 B2INFO(
"eclWaveformTemplateCalibrationC2Algorithm: successfully stored " << Form(
"PhotonParameters_CellID%d_CellID%d",
Base class for calibration algorithms.
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
void updateDBObjPtrs(const unsigned int event, const int run, const int experiment)
Updates any DBObjPtrs by calling update(event) for DBStore.
void setDescription(const std::string &description)
Set algorithm description (in constructor)
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
EResult
The result of calibration.
@ c_OK
Finished successfully =0 in Python.
Class for accessing objects in the database.
Abstract base class for different kinds of events.
Struct containing exp number and run number.