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);
225 B2INFO(
"eclWaveformTemplateCalibrationC2Algorithm: warning total entries for cell ID " << CellID_i <<
" is only: " <<
234 auto gWaveformToFit =
new TGraph(xValuesToFit.size(), xValuesToFit.data(), yValuesToFit.data());
235 gWaveformToFit->SetName(Form(
"gWaveformToFit_%d",
int(CellID_i)));
240 FitFunctions.clear();
241 for (
int i = 0; i < counterWaveforms; i++) {
242 FitFunctions.push_back(
new TF1(Form(
"Shp_%d", i), Belle2::ECL::WaveFuncTwoComponent, 0, 30.5, 26));
243 FitFunctions[i]->SetNpx(10000);
244 FitFunctions[i]->FixParameter(3, 0);
245 for (
int k = 0; k < 10; k++) {
246 FitFunctions[i]->SetParameter(4 + k, ParMin11t[k + 1]);
247 FitFunctions[i]->FixParameter(10 + 4 + k, ParMin11t[k + 1]);
249 FitFunctions[i]->FixParameter(24, ParMin11t[0]);
250 FitFunctions[i]->FixParameter(25, 1);
254 TF1* TotalFitFunction =
new TF1(
"TotalFitFunction", fitf, 0, counterWaveforms *
m_NumberofADCPoints,
255 (3 * FitFunctions.size()) + 10);
258 int FFsize = FitFunctions.size();
259 for (
int i = 0; i < FFsize; i++) {
260 TotalFitFunction->SetParameter(i, guessTime[i]);
261 TotalFitFunction->SetParameter(FFsize + i, guessBaseline[i]);
262 TotalFitFunction->SetParameter((2 * FFsize) + i, guessAmp[i]);
263 for (
int k = 0; k < 10; k++) {
264 TotalFitFunction->SetParameter((3 * FFsize) + k, ParMin11t[k + 1]);
266 TotalFitFunction->SetParLimits((3 * FFsize) + k, ParMin11t[k + 1] -
m_ParamLimitFactor * fabs(ParMin11t[k + 1]),
269 TotalFitFunction->ReleaseParameter((3 * FFsize) + k);
275 gWaveformToFit->Fit(
"TotalFitFunction",
"Q M W N 0 R",
"", 0, counterWaveforms *
m_NumberofADCPoints);
278 std::vector<int> FitResultY;
279 std::vector<int> FitResultX;
280 int maxResidualWaveformID = 0;
284 double npts = xValuesToFit.size();
285 double maxResidualOld = 0.0;
286 for (
int k = 0; k < npts; k++) {
287 double xVal = xValuesToFit[k];
288 double yVal = TotalFitFunction->Eval(xVal);
289 FitResultX.push_back(xVal);
290 FitResultY.push_back(yVal);
291 double diff = fabs(yValuesToFit[k] - yVal);
292 if (diff > maxResidual) {
295 maxResidualOld = fabs(yValuesToFit[k] / yVal);
300 if (maxResidual > resLimit) {
302 B2INFO(
"FAIL: CellID_i " << CellID_i <<
" maxResidual " << maxResidual <<
" removing entry: " <<
303 NtupleEntries[maxResidualWaveformID] <<
304 " which was waveform number " << maxResidualWaveformID <<
" resLimit was " << resLimit <<
" , resLimit started at " <<
306 B2INFO(
"Old maxResidual of Data/Fit was " << maxResidualOld);
308 B2INFO(
"Iter Time = " << std::time(0) - t << std::endl);
311 std::cout <<
"FAIL: CellID_i " << CellID_i <<
" maxResidual " << maxResidual <<
" removing entry: " <<
312 NtupleEntries[maxResidualWaveformID] <<
313 " which was waveform number " << maxResidualWaveformID <<
" resLimit was " << resLimit <<
" , resLimit started at " <<
314 resLimitOriginal << std::endl;
315 std::cout <<
"wave = [";
316 for (
int k = 0; k < npts; k++) {
317 std::cout << yValuesToFit[k];
318 if (k < (npts - 1)) {
321 std::cout <<
"]" << std::endl;
324 std::cout <<
"fitRes = [";
325 for (
int k = 0; k < npts; k++) {
326 std::cout << TotalFitFunction->Eval(xValuesToFit[k]);
327 if (k < (npts - 1)) {
330 std::cout <<
"]" << std::endl;
335 EntriesToSkip.push_back(NtupleEntries[maxResidualWaveformID]);
347 B2INFO(
"AttemptCounter reach limit: " << AttemptCounter <<
" counterWaveforms: " << counterWaveforms);
351 EntriesToSkip.clear();
357 B2INFO(
"Increasing resLimit to " << resLimit);
364 B2INFO(
"PASS: CellID_i " << CellID_i <<
" maxResidual " << maxResidual <<
" number of waveforms used was " << counterWaveforms <<
365 " resLimit was " << resLimit);
369 limitResidualArray.push_back(resLimit);
376 auto gFitResult =
new TGraph(FitResultX.size(), FitResultX.data(), FitResultY.data());
377 gFitResult->SetName(Form(
"gFitResult_%d",
int(CellID_i)));
380 cellIDArray.push_back(CellID_i);
381 maxResidualArray.push_back(maxResidual);
385 gWaveformToFit->Write();
389 float tempPhotonPar11[11];
390 tempPhotonPar11[0] = ParMin11t[0];
391 for (
unsigned int k = 0; k < 10; k++) tempPhotonPar11[k + 1] = TotalFitFunction->GetParameter((3 * FFsize) + k);
394 FitFunctions[0]->SetParameter(0, 0);
395 FitFunctions[0]->SetParameter(1, 0);
396 FitFunctions[0]->SetParameter(2, 1);
397 for (
int k = 0; k < 10; k++) {
398 FitFunctions[0]->SetParameter(4 + k, tempPhotonPar11[k + 1]);
399 FitFunctions[0]->SetParameter(10 + 4 + k, tempPhotonPar11[k + 1]);
401 FitFunctions[0]->FixParameter(24, ParMin11t[0]);
402 FitFunctions[0]->FixParameter(25, 1);
405 double MaxVal = -1.0;
406 const double cnpts = 2000;
407 for (
int k = 0; k < cnpts; k++) {
409 double yVal = FitFunctions[0]->Eval(xVal);
410 if (yVal > MaxVal) MaxVal = yVal;
412 B2INFO(
"MaxVal " << MaxVal);
413 tempPhotonPar11[0] /= MaxVal;
414 FitFunctions[0]->FixParameter(24, tempPhotonPar11[0]);
417 PhotonParameters->
setTemplateParameters(CellID_i, tempPhotonPar11, tempPhotonPar11, tempPhotonPar11);
420 for (
unsigned int k = 0; k < PhotonWaveformArray.size();
421 k++) PhotonWaveformArray[k] = FitFunctions[0]->Eval(((
double)k) * (1. / 1000.)) ;
425 for (
int w = 0; w < (int)FitFunctions.size(); w++) FitFunctions[w]->Delete();
426 TotalFitFunction->Delete() ;
427 gWaveformToFit->Delete();
433 auto gmaxResidual =
new TGraph(cellIDArray.size(), cellIDArray.data(), maxResidualArray.data());
434 gmaxResidual->SetName(
"gmaxResidual");
435 auto glimitResidualArray =
new TGraph(cellIDArray.size(), cellIDArray.data(), limitResidualArray.data());
436 glimitResidualArray->SetName(
"glimitResidualArray");
437 auto gparLimitFactorArray =
new TGraph(cellIDArray.size(), cellIDArray.data(), parLimitFactorArray.data());
438 gparLimitFactorArray->SetName(
"gparLimitFactorArray");
440 gmaxResidual->Write();
441 glimitResidualArray->Write();
442 gparLimitFactorArray->Write();
447 f_PhotonTemplateOutput->cd();
449 f_PhotonTemplateOutput->Write();
450 f_PhotonTemplateOutput->Close();
451 delete f_PhotonTemplateOutput;
455 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 successfuly =0 in Python.
@ c_NotEnoughData
Needs more data =2 in Python.
Class for accessing objects in the database.
Abstract base class for different kinds of events.
Struct containing exp number and run number.