Belle II Software development
eclWaveformTemplateCalibrationC4Algorithm.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9/* Own header. */
10#include <ecl/calibration/eclWaveformTemplateCalibrationC4Algorithm.h>
11
12/* ECL headers. */
13#include <ecl/dataobjects/ECLElementNumbers.h>
14#include <ecl/dbobjects/ECLDigitWaveformParameters.h>
15
16/* Basf2 headers. */
17#include <framework/database/DBObjPtr.h>
18
19/* ROOT headers. */
20#include <TFile.h>
21#include <TGraph.h>
22
23using namespace std;
24using namespace Belle2;
25using namespace ECL;
26using namespace Calibration;
27
30 CalibrationAlgorithm("DummyCollector")
31{
33 "Collects results from C3 to produce final payload, which contains new waveform templates"
34 );
35}
36
38{
39 B2INFO("eclWaveformTemplateCalibrationC4Algorithm m_firstCellID, m_lastCellID " << m_firstCellID << " " << m_lastCellID);
40
41 //save to db
42 ECLDigitWaveformParameters* PhotonHadronDiodeParameters = new ECLDigitWaveformParameters();
43
44 std::vector<float> cellIDs;
45 std::vector<float> photonNorms;
46 std::vector<float> hadronNorms;
47 std::vector<float> diodeNorms;
48
49 //int experimentNumber = -999;
50
51 for (int i = 0; i < m_numBatches; i++) {
52
53 int first = (i * m_batchsize) + 1;
54 int last = ((i + 1) * m_batchsize);
56
57 B2INFO("eclWaveformTemplateCalibrationC4Algorithm i " << i << " " << first << " " << last);
58 B2INFO("eclWaveformTemplateCalibrationC4Algorithm " << Form("PhotonParameters_CellID%d_CellID%d", first, last));
59 B2INFO("eclWaveformTemplateCalibrationC4Algorithm " << Form("HadronDiodeParameters_CellID%d_CellID%d", first, last));
60
61 DBObjPtr<ECLDigitWaveformParameters> tempexistingPhotonWaveformParameters(Form("PhotonParameters_CellID%d_CellID%d", first,
62 last));
63 DBObjPtr<ECLDigitWaveformParameters> tempexistingHadronDiodeWaveformParameters(Form("HadronDiodeParameters_CellID%d_CellID%d",
64 first, last));
65
66 auto runs = getRunList();
67 // Take the first run
68 ExpRun chosenRun = runs.front();
69 B2INFO("merging using the ExpRun (" << chosenRun.second << "," << chosenRun.first << ")");
70 // After here your DBObjPtrs are correct
71 updateDBObjPtrs(1, chosenRun.second, chosenRun.first);
72
73 for (int j = first; j <= last; j++) {
74 B2INFO("Check Norm Parms CellID " << j);
75 B2INFO("P " << j << " " << tempexistingPhotonWaveformParameters->getPhotonParameters(j)[0]);
76 B2INFO("H " << j << " " << tempexistingHadronDiodeWaveformParameters->getHadronParameters(j)[0]);
77 B2INFO("D " << j << " " << tempexistingHadronDiodeWaveformParameters->getDiodeParameters(j)[0]);
78
79 cellIDs.push_back(j);
80 photonNorms.push_back(tempexistingPhotonWaveformParameters->getPhotonParameters(j)[0]);
81 hadronNorms.push_back(tempexistingPhotonWaveformParameters->getHadronParameters(j)[0]);
82 diodeNorms.push_back(tempexistingPhotonWaveformParameters->getDiodeParameters(j)[0]);
83
84 float tempPhotonWaveformParameters[11];
85 float tempHadronWaveformParameters[11];
86 float tempDiodeWaveformParameters[11];
87
88 for (int k = 0; k < 11; k++) {
89 tempPhotonWaveformParameters[k] = tempexistingPhotonWaveformParameters->getPhotonParameters(j)[k];
90 tempHadronWaveformParameters[k] = tempexistingHadronDiodeWaveformParameters->getHadronParameters(j)[k];
91 tempDiodeWaveformParameters[k] = tempexistingHadronDiodeWaveformParameters->getDiodeParameters(j)[k];
92 }
93 PhotonHadronDiodeParameters->setTemplateParameters(j, tempPhotonWaveformParameters, tempHadronWaveformParameters,
94 tempDiodeWaveformParameters);
95 }
96 }
97
98 auto gphotonNorms = new TGraph(cellIDs.size(), cellIDs.data(), photonNorms.data());
99 gphotonNorms->SetName("gphotonNorms");
100 auto ghadronNorms = new TGraph(cellIDs.size(), cellIDs.data(), hadronNorms.data());
101 ghadronNorms->SetName("ghadronNorms");
102 auto gdiodeNorms = new TGraph(cellIDs.size(), cellIDs.data(), diodeNorms.data());
103 gdiodeNorms->SetName("gdiodeNorms");
104
105 TString fName = m_outputName;
106 TFile* histfile = new TFile(fName, "recreate");
107 histfile->cd();
108 gphotonNorms->Write();
109 ghadronNorms->Write();
110 gdiodeNorms->Write();
111 histfile->Close();
112 delete histfile;
113
114
115 B2INFO("eclWaveformTemplateCalibrationC4Algorithm: Successful, now writing DB PAyload");
116 saveCalibration(PhotonHadronDiodeParameters, "ECLDigitWaveformParameters");
117
118 return c_OK;
119}
120
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.
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
DB object to store photon, hadron and diode shape parameters.
void setTemplateParameters(int cellID, const float photonInput[11], const float hadronInput[11], const float diodeInput[11])
Set photon, hadron and diode template parameters for crystal.
const int c_NCrystals
Number of crystals.
Abstract base class for different kinds of events.
STL namespace.
Struct containing exp number and run number.
Definition: Splitter.h:51