140#include "TRestDataSetGainMap.h"
183 TiXmlElement* moduleDefinition =
GetElement(
"module");
184 while (moduleDefinition !=
nullptr) {
186 fModulesCal.back().LoadConfigFromTiXmlElement(moduleDefinition);
201 RESTInfo <<
"Generating gain map of plane " << mod.GetPlaneId() <<
" module " << mod.GetModuleId()
203 mod.GenerateGainMap();
206 mod.DrawFullSpectrum();
225 std::vector<std::string> excludeColumns) {
227 RESTError <<
"TRestDataSetGainMap::CalibrateDataSet: No modules defined." <<
RESTendl;
232 dataSet.EnableMultiThreading(
true);
233 dataSet.
Import(dataSetFileName);
237 std::string pmIDname = (std::string)GetName() +
"_pmID";
238 std::string modCut =
fModulesCal[0].GetModuleDefinitionCut();
239 if (modCut.empty()) modCut =
"1";
242 auto columnList = dataFrame.GetColumnNames();
243 if (std::find(columnList.begin(), columnList.end(), pmIDname) == columnList.end())
244 dataFrame = dataFrame.Define(pmIDname, modCut +
" ? " + std::to_string(pmID) +
" : -1");
246 dataFrame = dataFrame.Redefine(pmIDname, modCut +
" ? " + std::to_string(pmID) +
" : -1");
250 if (modCut.empty()) modCut =
"1";
252 dataFrame = dataFrame.Redefine(pmIDname, (modCut +
" ? " + std::to_string(pmID) +
" : " + pmIDname));
256 auto calibrate = [
this](
double val,
double x,
double y,
int pmID) {
258 if (pmID == m.GetPlaneId() * 10 + m.GetModuleId())
259 return m.GetSlope(x, y) * val + m.GetIntercept(x, y);
263 return std::numeric_limits<double>::quiet_NaN();
265 std::string calibObsName = (std::string)GetName() +
"_";
266 calibObsName += GetObservable().erase(0, GetObservable().find(
"_") + 1);
267 dataFrame = dataFrame.Define(calibObsName, calibrate,
271 auto calibrateFullSpc = [
this](
double val,
int pmID) {
273 if (pmID == m.GetPlaneId() * 10 + m.GetModuleId())
274 return m.GetSlopeFullSpc() * val + m.GetInterceptFullSpc();
278 return std::numeric_limits<double>::quiet_NaN();
280 std::string calibObsNameFullSpc = (std::string)GetName() +
"_";
281 calibObsNameFullSpc +=
282 GetObservable().erase(0, GetObservable().find(
"_") + 1);
283 calibObsNameFullSpc +=
"_NoSegmentation";
284 dataFrame = dataFrame.Define(calibObsNameFullSpc, calibrateFullSpc, {
fObservable, pmIDname});
286 dataSet.SetDataFrame(dataFrame);
289 if (outputFileName.empty()) outputFileName = dataSetFileName;
290 if (outputFileName == dataSetFileName) {
291 std::string gmName = GetName();
292 outputFileName = outputFileName.substr(0, outputFileName.find_last_of(
"."));
297 std::set<std::string> excludeCol;
300 for (
auto& eC : excludeColumns) {
301 if (eC.find(
"*") != std::string::npos || eC.find(
"?") != std::string::npos) {
302 for (
auto& c : columns)
304 }
else if (std::find(columns.begin(), columns.end(), eC) != columns.end())
305 excludeCol.insert(eC);
308 excludeCol.erase(calibObsName);
309 excludeCol.erase(calibObsNameFullSpc);
310 excludeCol.erase(pmIDname);
312 RESTDebug <<
"Excluding columns: ";
313 for (
auto& c : excludeCol) RESTDebug << c <<
", ";
316 dataSet.
Export(outputFileName, std::vector<std::string>(excludeCol.begin(), excludeCol.end()));
319 TFile* f = TFile::Open(outputFileName.c_str(),
"UPDATE");
332 RESTError <<
"No ModuleCalibration with index " << index;
334 RESTError <<
". There are no modules defined." <<
RESTendl;
346 if (i.GetPlaneId() == planeID && i.GetModuleId() == moduleID)
return &i;
348 RESTError <<
"No ModuleCalibration with planeID " << planeID <<
" and moduleID " << moduleID <<
RESTendl;
360 if (moduleCal ==
nullptr)
return 0;
370 if (moduleCal ==
nullptr)
return 0;
371 return moduleCal->GetSlopeFullSpc();
382 if (moduleCal ==
nullptr)
return 0;
392 if (moduleCal ==
nullptr)
return 0;
393 return moduleCal->GetInterceptFullSpc();
401 std::set<int> planeIDs;
402 for (
const auto& mc :
fModulesCal) planeIDs.insert(mc.GetPlaneId());
411 std::set<int> moduleIDs;
413 if (mc.GetPlaneId() == planeId) moduleIDs.insert(mc.GetModuleId());
421 std::map<int, std::set<int>> moduleIds;
423 moduleIds.insert(std::pair<
int, std::set<int>>(planeId,
GetModuleIDs(planeId)));
428 SetName(src.GetName());
446 if (i.GetPlaneId() == moduleCal.GetPlaneId() && i.GetModuleId() == moduleCal.GetModuleId()) {
459 if (fileName.empty()) {
460 RESTError <<
"No input calibration file defined" <<
RESTendl;
465 RESTInfo <<
"Opening " << fileName <<
RESTendl;
466 TFile* f = TFile::Open(fileName.c_str(),
"READ");
468 RESTError <<
"Cannot open calibration file " << fileName <<
RESTendl;
474 TIter nextkey(f->GetListOfKeys());
476 while ((key = (TKey*)nextkey())) {
477 std::string kName = key->GetClassName();
486 RESTError <<
"File extension not supported for " << fileName <<
RESTendl;
499 RESTError <<
"No output file defined" <<
RESTendl;
505 this->
Write(GetName());
520 RESTMetadata <<
" Cuts applied: ";
523 for (
const auto& cut :
fCut->GetCutStrings()) RESTMetadata << cut <<
", " <<
RESTendl;
524 for (
const auto& cut :
fCut->GetParamCut())
525 RESTMetadata << cut.first <<
" " << cut.second <<
", " <<
RESTendl;
529 RESTMetadata <<
" Number of planes: " << GetNumberOfPlanes() <<
RESTendl;
530 RESTMetadata <<
" Number of modules: " << GetNumberOfModules() <<
RESTendl;
534 RESTMetadata <<
"-----------------------------------------------" <<
RESTendl;
536 RESTMetadata <<
"***********************************************" <<
RESTendl;
548 int index_x = -1, index_y = -1;
551 index_x = std::distance(
fSplitX.begin(),
fSplitX.upper_bound(x)) - 1;
553 RESTWarning <<
"index_x < 0 for x = " << x <<
" and fSplitX[0]=" << *
fSplitX.begin()
558 RESTWarning <<
"x is out of split for x = " << x <<
p->
RESTendl;
563 index_y = std::distance(
fSplitY.begin(),
fSplitY.upper_bound(y)) - 1;
565 RESTWarning <<
"index_y < 0 for y = " << y <<
" and fSplitY[0]=" << *
fSplitY.begin()
570 RESTWarning <<
"y is out of split for y = " << y <<
p->
RESTendl;
574 return std::make_pair(index_x, index_y);
584 auto [index_x, index_y] = GetIndexMatrix(x, y);
585 if (fSlope.empty()) {
586 RESTError <<
"Calibration slope matrix is empty. Returning 0" << p->RESTendl;
590 if (index_x > (
int)fSlope.size() || index_y > (
int)fSlope.at(0).size()) {
591 RESTError <<
"Index out of range. Returning 0" << p->RESTendl;
595 return fSlope[index_x][index_y];
606 auto [index_x, index_y] = GetIndexMatrix(x, y);
607 if (fIntercept.empty()) {
608 RESTError <<
"Calibration constant matrix is empty. Returning 0" << p->RESTendl;
612 if (index_x > (
int)fIntercept.size() || index_y > (
int)fIntercept.at(0).size()) {
613 RESTError <<
"Index out of range. Returning 0" << p->RESTendl;
617 return fIntercept[index_x][index_y];
634 if (fNumberOfSegmentsX < 1) {
635 RESTError <<
"SetSplitX: fNumberOfSegmentsX must be >=1." << p->RESTendl;
638 std::set<double> split;
639 for (
int i = 0; i <= fNumberOfSegmentsX; i++) {
641 fReadoutRange.X() + ((fReadoutRange.Y() - fReadoutRange.X()) / (
float)fNumberOfSegmentsX) * i;
648 if (splitX.size() < 2) {
649 RESTError <<
"SetSplitX: split size must be >=2 (start and end of range must be included)."
654 RESTWarning <<
"SetSplitX: changing split but current gain map and calibration paremeters correspond "
655 "to previous splitting. Use GenerateGainMap() to update them."
658 fNumberOfSegmentsX = fSplitX.size() - 1;
667 if (fNumberOfSegmentsY < 1) {
668 RESTError <<
"SetSplitY: fNumberOfSegmentsY must be >=1." << p->RESTendl;
671 std::set<double> split;
672 for (
int i = 0; i <= fNumberOfSegmentsY; i++) {
674 fReadoutRange.X() + ((fReadoutRange.Y() - fReadoutRange.X()) / (
float)fNumberOfSegmentsY) * i;
681 if (splitY.size() < 2) {
682 RESTError <<
"SetSplitY: split size must be >=2 (start and end of range must be included)."
687 RESTWarning <<
"SetSplitY: changing split but current gain map and calibration paremeters correspond "
688 "to previous splitting. Use GenerateGainMap() to update them."
691 fNumberOfSegmentsY = fSplitY.size() - 1;
712 std::string dsFileName = fDataSetFileName;
713 if (dsFileName.empty()) dsFileName = p->GetCalibrationFileName();
714 if (dsFileName.empty()) {
715 RESTError <<
"No calibration file defined" << p->RESTendl;
720 dataSet.EnableMultiThreading(
true);
722 dataSet.
Import(dsFileName);
723 fDataSetFileName = dsFileName;
725 RESTWarning << dsFileName <<
" is not a dataset. Generating a temporal one..." << p->RESTendl;
727 std::vector<std::string> obsList;
729 obsList.push_back(p->GetObservable());
730 obsList.push_back(p->GetSpatialObservableX());
731 obsList.push_back(p->GetSpatialObservableY());
735 obsList.insert(obsList.end(), modDefCutObs.begin(), modDefCutObs.end());
738 for (
const auto& cut : p->GetCut()->GetCutStrings()) {
740 obsList.insert(obsList.end(), cutObs.begin(), cutObs.end());
742 for (
const auto& [variable, condition] : p->GetCut()->GetParamCut()) {
744 obsList.insert(obsList.end(), cutObs.begin(), cutObs.end());
747 obsList.insert(obsList.end(), cutObs.begin(), cutObs.end());
751 std::sort(obsList.begin(), obsList.end());
752 obsList.erase(std::unique(obsList.begin(), obsList.end()), obsList.end());
755 dataSet.SetFilePattern(dsFileName);
756 dataSet.SetObservablesList(obsList);
758 fDataSetFileName = dsFileName;
761 dataSet.SetDataFrame(dataSet.
MakeCut(p->GetCut()));
763 if (fSplitX.empty()) SetSplitX();
764 if (fSplitY.empty()) SetSplitY();
767 if (fCalibRange.X() >= fCalibRange.Y()) {
769 std::string cut = fDefinitionCut;
770 if (cut.empty()) cut =
"1";
771 auto histo = dataSet.
GetDataFrame().Filter(cut).Histo1D({
"temp",
"", fNBins, 0, 0}, GetObservable());
772 std::unique_ptr<TH1F> hpunt = std::unique_ptr<TH1F>(
static_cast<TH1F*
>(histo->Clone()));
774 double xMax = hpunt->GetXaxis()->GetXmax();
777 double fraction = 1, nAtEndSpc = 0, step = 0.66;
778 while (nAtEndSpc * 1. / hpunt->Integral() < 0.01 && fraction > 0.001) {
780 nAtEndSpc = hpunt->Integral(hpunt->FindFixBin(hpunt->GetXaxis()->GetXmax() * fraction),
781 hpunt->FindFixBin(hpunt->GetXaxis()->GetXmax()));
783 xMax = hpunt->GetXaxis()->GetXmax() * fraction /
787 fCalibRange.SetY(xMax);
789 RESTDebug <<
"Calibration range (auto)set to (" << fCalibRange.X() <<
"," << fCalibRange.Y() <<
")"
794 std::string hModuleName =
"hSpc_" + std::to_string(fPlaneId) +
"_" + std::to_string(fModuleId);
795 delete fFullSpectrum;
796 fFullSpectrum =
new TH1F(hModuleName.c_str(),
"", fNBins, fCalibRange.X(), fCalibRange.Y());
799 std::string cut = fDefinitionCut;
800 if (cut.empty()) cut =
"1";
801 auto histoMod = dataSet.
GetDataFrame().Filter(cut).Histo1D(
802 {
"tempMod",
"", fNBins, fCalibRange.X(), fCalibRange.Y()}, GetObservable());
803 std::unique_ptr<TH1F> hpuntMod = std::unique_ptr<TH1F>(
static_cast<TH1F*
>(histoMod->Clone()));
804 fFullSpectrum->Add(hpuntMod.get());
807 std::vector<std::vector<TH1F*>> h(fNumberOfSegmentsX, std::vector<TH1F*>(fNumberOfSegmentsY,
nullptr));
808 for (
size_t i = 0; i < h.size(); i++) {
809 for (
size_t j = 0; j < h.at(0).size(); j++) {
810 std::string name = hModuleName +
"_" + std::to_string(i) +
"_" + std::to_string(j);
811 h[i][j] =
new TH1F(name.c_str(),
"", fNBins, fCalibRange.X(),
817 auto itX = fSplitX.begin();
818 for (
size_t i = 0; i < h.size(); i++) {
819 auto itY = fSplitY.begin();
820 for (
size_t j = 0; j < h.at(0).size(); j++) {
823 auto xUpper = *std::next(itX);
825 auto yUpper = *std::next(itY);
827 std::string segment_cut =
"";
828 if (!GetSpatialObservableX().empty())
829 segment_cut += GetSpatialObservableX() +
">=" + std::to_string(xLower) +
"&&" +
830 GetSpatialObservableX() +
"<" + std::to_string(xUpper);
831 if (!GetSpatialObservableY().empty())
832 segment_cut +=
"&&" + GetSpatialObservableY() +
">=" + std::to_string(yLower) +
"&&" +
833 GetSpatialObservableY() +
"<" + std::to_string(yUpper);
834 if (!fDefinitionCut.empty()) segment_cut +=
"&&" + fDefinitionCut;
835 if (segment_cut.empty()) segment_cut =
"1";
836 RESTExtreme <<
"Segment[" << i <<
"][" << j <<
"] cut: " << segment_cut << p->RESTendl;
839 .Histo1D({
"temp",
"", h[i][j]->GetNbinsX(), h[i][j]->GetXaxis()->GetXmin(),
840 h[i][j]->GetXaxis()->GetXmax()},
842 std::unique_ptr<TH1F> hpunt = std::unique_ptr<TH1F>(
static_cast<TH1F*
>(histo->Clone()));
843 h[i][j]->Add(hpunt.get());
851 std::vector<std::vector<double>> calParSlope(fNumberOfSegmentsX,
852 std::vector<double>(fNumberOfSegmentsY, 0));
853 std::vector<std::vector<double>> calParIntercept(fNumberOfSegmentsX,
854 std::vector<double>(fNumberOfSegmentsY, 0));
855 fSegLinearFit = std::vector(h.size(), std::vector<TGraph*>(h.at(0).size(),
nullptr));
856 for (
size_t i = 0; i < h.size(); i++)
857 for (
int j = 0; j < (int)h.at(0).size(); j++) {
858 fSegLinearFit[i][j] =
new TGraph();
859 auto [intercept, slope] = FitPeaks(h[i][j], fSegLinearFit[i][j]);
860 calParSlope[i][j] = slope;
861 calParIntercept[i][j] = intercept;
863 fSlope = calParSlope;
864 fIntercept = calParIntercept;
868 delete fFullLinearFit;
869 fFullLinearFit =
new TGraph();
870 auto [intercept, slope] = FitPeaks(fFullSpectrum, fFullLinearFit);
872 fFullIntercept = intercept;
875std::pair<double, double> TRestDataSetGainMap::Module::FitPeaks(TH1F* hSeg, TGraph* gr) {
877 RESTError <<
"No histogram for fitting" << p->RESTendl;
878 return std::make_pair(0, 0);
880 if (hSeg->Integral() == 0) {
881 RESTError <<
"Empty spectrum " << hSeg->GetName() << p->RESTendl;
882 return std::make_pair(0, 0);
884 std::shared_ptr<TGraph> graph = std::shared_ptr<TGraph>(
new TGraph());
885 RESTExtreme <<
"Fitting peaks for " << hSeg->GetName() << p->RESTendl;
888 std::unique_ptr<TSpectrum> s(
new TSpectrum(2 * fEnergyPeaks.size() + 1));
889 std::vector<double> peakPos;
890 s->Search(hSeg, 2,
"goff", 0.1);
891 for (
int k = 0; k < s->GetNPeaks(); k++) peakPos.push_back(s->GetPositionX()[k]);
892 std::sort(peakPos.begin(), peakPos.end(), std::greater<double>());
894 peakPos.size() == 0 ? 1 : peakPos.front() / fEnergyPeaks.front();
897 graph->SetName(
"grFit");
898 graph->SetTitle((
";" + GetObservable() +
";energy").c_str());
903 for (
const auto& energy : fEnergyPeaks) {
904 RESTExtreme <<
"\t fitting energy " <<
DoubleToString(energy,
"%g") << p->RESTendl;
906 double pos = energy * ratio;
907 double start = pos * 0.8;
908 double end = pos * 1.2;
909 if (fRangePeaks.at(c).X() < fRangePeaks.at(c).Y()) {
910 start = fRangePeaks.at(c).X();
911 end = fRangePeaks.at(c).Y();
915 if (fAutoRangePeaks) {
916 if (peakPos.size() > 0) {
919 while (!(start < pos && pos < end)) {
922 if (pos == peakPos.back()) {
926 pos = *std::next(std::find(peakPos.begin(), peakPos.end(),
929 peakPos.erase(std::find(peakPos.begin(), peakPos.end(),
935 const double relDist = peakPos.size() > 0 ? (pos - peakPos.front()) / pos : 999;
937 start = pos * (1 - relDist / 2);
938 end = pos * (1 + relDist / 2);
943 std::string name =
"g" + std::to_string(c);
944 TF1* g =
new TF1(name.c_str(),
"gaus", start, end);
945 RESTExtreme <<
"\t\tat " <<
DoubleToString(pos,
"%.3g") <<
". Range("
949 if (hSeg->GetFunction(name.c_str()))
950 hSeg->GetListOfFunctions()->Remove(hSeg->GetFunction(name.c_str()));
952 hSeg->Fit(g,
"R+Q0");
953 mu = g->GetParameter(1);
954 RESTExtreme <<
"\t\tgaus mean " <<
DoubleToString(mu,
"%g") << p->RESTendl;
955 }
while (fAutoRangePeaks && peakPos.size() > 0 &&
956 !(start < mu && mu < end));
957 graph->SetPoint(c++, mu, energy);
961 if (fZeroPoint) graph->SetPoint(c++, 0, 0);
962 while (graph->GetN() < 2) {
963 graph->SetPoint(c++, 0, 0);
965 RESTDebug <<
"Not enough points for linear fit. Adding and setting zero point to true" << p->RESTendl;
969 std::unique_ptr<TF1> linearFit;
970 linearFit = std::unique_ptr<TF1>(
new TF1(
"linearFit",
"pol1"));
971 graph->Fit(
"linearFit",
"SQ");
973 if (gr) *gr = *(TGraph*)graph->Clone();
974 return std::make_pair(linearFit->GetParameter(0), linearFit->GetParameter(1));
985 const TVector2& range) {
986 auto [index_x, index_y] = GetIndexMatrix(position.X(), position.Y());
988 for (
size_t i = 0; i < fEnergyPeaks.size(); i++)
989 if (fEnergyPeaks.at(i) == energyPeak) {
993 if (peakNumber == -1) {
994 RESTError <<
"Energy " << energyPeak <<
" not found in the list of energy peaks" << p->RESTendl;
997 Refit((
size_t)index_x, (
size_t)index_y, (
size_t)peakNumber, range);
1010 const TVector2& range) {
1011 if (fSegSpectra.empty()) {
1012 RESTError <<
"No gain map found. Use GenerateGainMap() first." << p->RESTendl;
1015 if (x >= fSegSpectra.size() || y >= fSegSpectra.at(0).size()) {
1016 RESTError <<
"Segment with index (" << x <<
", " << y <<
") not found" << p->RESTendl;
1019 if (peakNumber >= fEnergyPeaks.size()) {
1020 RESTError <<
"Peak with index " << peakNumber <<
" not found" << p->RESTendl;
1025 std::string name =
"g" + std::to_string(peakNumber);
1026 TF1* g =
new TF1(name.c_str(),
"gaus", range.X(), range.Y());
1027 TH1F* h = fSegSpectra.at(x).at(y);
1028 while (h->GetFunction(name.c_str()))
1029 h->GetListOfFunctions()->Remove(h->GetFunction(name.c_str()));
1033 UpdateCalibrationFits(x, y);
1044 int peakNumber = -1;
1045 for (
size_t i = 0; i < fEnergyPeaks.size(); i++)
1046 if (fEnergyPeaks.at(i) == energyPeak) {
1050 if (peakNumber == -1) {
1051 RESTError <<
"Energy " << energyPeak <<
" not found in the list of energy peaks" << p->RESTendl;
1054 RefitFullSpc((
size_t)peakNumber, range);
1065 if (!fFullSpectrum) {
1066 RESTError <<
"No gain map found. Use GenerateGainMap() first." << p->RESTendl;
1069 if (peakNumber >= fEnergyPeaks.size()) {
1070 RESTError <<
"Peak with index " << peakNumber <<
" not found" << p->RESTendl;
1075 std::string name =
"g" + std::to_string(peakNumber);
1076 TF1* g =
new TF1(name.c_str(),
"gaus", range.X(), range.Y());
1077 while (fFullSpectrum->GetFunction(name.c_str()))
1078 fFullSpectrum->GetListOfFunctions()->Remove(fFullSpectrum->GetFunction(name.c_str()));
1079 fFullSpectrum->Fit(g,
"R+Q0");
1082 UpdateCalibrationFitsFullSpc();
1093void TRestDataSetGainMap::Module::UpdateCalibrationFits(
const size_t x,
const size_t y) {
1094 if (fSegSpectra.empty()) {
1095 RESTError <<
"No gain map found. Use GenerateGainMap() first." << p->RESTendl;
1098 if (x >= fSegSpectra.size() || y >= fSegSpectra.at(0).size()) {
1099 RESTError <<
"Segment with index (" << x <<
", " << y <<
") not found" << p->RESTendl;
1103 TH1F* h = fSegSpectra.at(x).at(y);
1104 TGraph* gr = fSegLinearFit.at(x).at(y);
1106 auto [intercept, slope] = UpdateCalibrationFits(h, gr);
1107 fSlope[x][y] = slope;
1108 fIntercept[x][y] = intercept;
1111std::pair<double, double> TRestDataSetGainMap::Module::UpdateCalibrationFits(TH1F* h, TGraph* gr) {
1113 RESTError <<
"No histogram for updating fits" << p->RESTendl;
1114 return std::make_pair(0, 0);
1117 RESTError <<
"No graph for updating fits" << p->RESTendl;
1118 return std::make_pair(0, 0);
1120 if (h->Integral() == 0) {
1121 RESTError <<
"Empty spectrum " << h->GetName() << p->RESTendl;
1122 return std::make_pair(0, 0);
1126 for (
size_t i = 0; i < fEnergyPeaks.size(); i++) gr->RemovePoint(i);
1129 for (
size_t i = 0; i < fEnergyPeaks.size(); i++) {
1130 std::string fitName = (std::string)
"g" + std::to_string(i);
1131 TF1* g = h->GetFunction(fitName.c_str());
1133 RESTWarning <<
"No fit ( " << fitName <<
" ) found for energy peak " << fEnergyPeaks[i]
1134 <<
" in histogram " << h->GetName() << p->RESTendl;
1137 gr->SetPoint(c++, g->GetParameter(1), fEnergyPeaks[i]);
1141 while (gr->GetN() < 2) {
1142 gr->SetPoint(c++, 0, 0);
1147 if (gr->GetFunction(
"linearFit"))
1148 lf = gr->GetFunction(
"linearFit");
1150 lf =
new TF1(
"linearFit",
"pol1");
1153 return std::make_pair(lf->GetParameter(0), lf->GetParameter(1));
1162 auto [intercept, slope] = UpdateCalibrationFits(fFullSpectrum, fFullLinearFit);
1164 fFullIntercept = intercept;
1182 if (module ==
nullptr) {
1183 RESTError <<
"TRestDataSetGainMap::Module::LoadConfigFromTiXmlElement: module is nullptr"
1188 std::string el = !
module->Attribute("planeId") ? "Not defined" : module->Attribute("planeId");
1189 if (!(el.empty() || el ==
"Not defined")) this->SetPlaneId(
StringToInteger(el));
1190 el = !
module->Attribute("moduleId") ? "Not defined" : module->Attribute("moduleId");
1191 if (!(el.empty() || el ==
"Not defined")) this->SetModuleId(
StringToInteger(el));
1193 el = !
module->Attribute("moduleDefinitionCut") ? "Not defined" : module->Attribute("moduleDefinitionCut");
1194 if (!(el.empty() || el ==
"Not defined")) this->SetModuleDefinitionCut(el);
1196 el = !
module->Attribute("numberOfSegmentsX") ? "Not defined" : module->Attribute("numberOfSegmentsX");
1197 if (!(el.empty() || el ==
"Not defined")) this->SetNumberOfSegmentsX(
StringToInteger(el));
1198 el = !
module->Attribute("numberOfSegmentsY") ? "Not defined" : module->Attribute("numberOfSegmentsY");
1199 if (!(el.empty() || el ==
"Not defined")) this->SetNumberOfSegmentsY(
StringToInteger(el));
1200 el = !
module->Attribute("readoutRange") ? "Not defined" : module->Attribute("readoutRange");
1201 if (!(el.empty() || el ==
"Not defined")) this->SetReadoutRange(
StringTo2DVector(el));
1203 el = !
module->Attribute("calibRange") ? "Not defined" : module->Attribute("calibRange");
1204 if (!(el.empty() || el ==
"Not defined")) this->SetCalibrationRange(
StringTo2DVector(el));
1205 el = !
module->Attribute("nBins") ? "Not defined" : module->Attribute("nBins");
1206 if (!(el.empty() || el ==
"Not defined")) this->SetNBins(
StringToInteger(el));
1208 el = !
module->Attribute("dataSetFileName") ? "Not defined" : module->Attribute("dataSetFileName");
1209 if (!(el.empty() || el ==
"Not defined")) this->SetDataSetFileName(el);
1211 el = !
module->Attribute("zeroPoint") ? "Not defined" : module->Attribute("zeroPoint");
1212 if (!(el.empty() || el ==
"Not defined")) this->SetZeroPoint(
ToLower(el) ==
"true");
1213 el = !
module->Attribute("autoRangePeaks") ? "Not defined" : module->Attribute("autoRangePeaks");
1214 if (!(el.empty() || el ==
"Not defined")) this->SetAutoRangePeaks(
ToLower(el) ==
"true");
1217 TiXmlElement* peakDefinition = (TiXmlElement*)module->FirstChildElement(
"peak");
1218 while (peakDefinition !=
nullptr) {
1220 TVector2 range = TVector2(0, 0);
1223 !peakDefinition->Attribute(
"energy") ?
"Not defined" : peakDefinition->Attribute(
"energy");
1224 if (ell.empty() || ell ==
"Not defined") {
1225 RESTError <<
"< peak variable key does not contain energy!" << p->RESTendl;
1230 ell = !peakDefinition->Attribute(
"range") ?
"Not defined" : peakDefinition->Attribute(
"range");
1233 this->AddPeak(energy, range);
1234 peakDefinition = (TiXmlElement*)peakDefinition->NextSiblingElement();
1240 std::pair<size_t, size_t> index = GetIndexMatrix(position.X(), position.Y());
1241 DrawSpectrum(index.first, index.second, drawFits, color, c);
1246 if (fSegSpectra.size() == 0) {
1247 RESTError <<
"Spectra matrix is empty." << p->RESTendl;
1250 if (index_x < 0 || index_y < 0 || index_x >= (
int)fSegSpectra.size() ||
1251 index_y >= (
int)fSegSpectra.at(index_x).size()) {
1252 RESTError <<
"Index out of range." << p->RESTendl;
1255 if (!fSegSpectra[index_x][index_y]) {
1256 RESTError <<
"No Spectrum for segment (" << index_x <<
", " << index_y <<
")." << p->RESTendl;
1261 std::string t =
"spectrum_" + std::to_string(fPlaneId) +
"_" + std::to_string(fModuleId) +
"_" +
1262 std::to_string(index_x) +
"_" + std::to_string(index_y);
1263 c =
new TCanvas(t.c_str(), t.c_str());
1266 auto xLower = *std::next(fSplitX.begin(), index_x);
1267 auto xUpper = *std::next(fSplitX.begin(), index_x + 1);
1268 auto yLower = *std::next(fSplitY.begin(), index_y);
1269 auto yUpper = *std::next(fSplitY.begin(), index_y + 1);
1272 GetObservable() +
";counts";
1273 fSegSpectra[index_x][index_y]->SetTitle(tH.c_str());
1275 if (color > 0) fSegSpectra[index_x][index_y]->SetLineColor(color);
1276 size_t colorT = fSegSpectra[index_x][index_y]->GetLineColor();
1277 fSegSpectra[index_x][index_y]->Draw(
"same");
1280 for (
size_t c = 0; c < fEnergyPeaks.size(); c++) {
1281 auto fit = fSegSpectra[index_x][index_y]->GetFunction((
"g" + std::to_string(c)).c_str());
1283 RESTWarning <<
"Fit for energy peak " << fEnergyPeaks[c] <<
" not found." << p->RESTendl;
1285 fit->SetLineColor(c + 2 != colorT ? c + 2 : ++colorT);
1313 if (fSegSpectra.size() == 0) {
1314 RESTError <<
"Spectra matrix is empty." << p->RESTendl;
1318 std::string t =
"spectrum_" + std::to_string(fPlaneId) +
"_" + std::to_string(fModuleId);
1319 c =
new TCanvas(t.c_str(), t.c_str());
1323 for (
const auto&
object : *c->GetListOfPrimitives())
1324 if (object->InheritsFrom(TVirtualPad::Class())) ++nPads;
1325 if (nPads != 0 && nPads != fSegSpectra.size() * fSegSpectra.at(0).size()) {
1326 RESTError <<
"Canvas " << c->GetName() <<
" has " << nPads <<
" pads, but "
1327 << fSegSpectra.size() * fSegSpectra.at(0).size() <<
" are needed." << p->RESTendl;
1329 }
else if (nPads == 0)
1330 c->Divide(fSegSpectra.size(), fSegSpectra.at(0).size());
1332 for (
size_t i = 0; i < fSegSpectra.size(); i++) {
1333 for (
size_t j = 0; j < fSegSpectra[i].size(); j++) {
1334 int pad = fSegSpectra.size() * (fSegSpectra[i].size() - 1) + 1 + i - fSegSpectra.size() * j;
1336 DrawSpectrum(i, j, drawFits, color, c);
1340void TRestDataSetGainMap::Module::DrawFullSpectrum(
const bool drawFits,
const int color, TCanvas* c) {
1341 if (!fFullSpectrum) {
1342 RESTError <<
"Spectrum is empty." << p->RESTendl;
1347 std::string t =
"fullSpc_" + std::to_string(fPlaneId) +
"_" + std::to_string(fModuleId);
1348 c =
new TCanvas(t.c_str(), t.c_str());
1352 fFullSpectrum->SetTitle((
"Full spectrum;" + GetObservable() +
";counts").c_str());
1354 if (color > 0) fFullSpectrum->SetLineColor(color);
1355 size_t colorT = fFullSpectrum->GetLineColor();
1356 fFullSpectrum->Draw(
"same");
1359 for (
size_t c = 0; c < fEnergyPeaks.size(); c++) {
1360 auto fit = fFullSpectrum->GetFunction((
"g" + std::to_string(c)).c_str());
1361 if (!fit) RESTWarning <<
"Fit for energy peak" << fEnergyPeaks[c] <<
" not found." << p->RESTendl;
1363 fit->SetLineColor(c + 2 != colorT ? c + 2 : ++colorT);
1371void TRestDataSetGainMap::Module::DrawLinearFit(
const TVector2& position, TCanvas* c) {
1372 std::pair<size_t, size_t> index = GetIndexMatrix(position.X(), position.Y());
1373 DrawLinearFit(index.first, index.second, c);
1376void TRestDataSetGainMap::Module::DrawLinearFit(
const int index_x,
const int index_y, TCanvas* c) {
1377 if (fSegLinearFit.size() == 0) {
1378 RESTError <<
"Spectra matrix is empty." << p->RESTendl;
1381 if (index_x < 0 || index_y < 0 || index_x >= (
int)fSegLinearFit.size() ||
1382 index_y >= (
int)fSegLinearFit.at(index_x).size()) {
1383 RESTError <<
"Index out of range." << p->RESTendl;
1386 if (!fSegLinearFit[index_x][index_y]) {
1387 RESTError <<
"No linear fit for segment (" << index_x <<
", " << index_y <<
")." << p->RESTendl;
1392 std::string t =
"linearFit_" + std::to_string(fPlaneId) +
"_" + std::to_string(fModuleId) +
"_" +
1393 std::to_string(index_x) +
"_" + std::to_string(index_y);
1394 c =
new TCanvas(t.c_str(), t.c_str());
1396 auto xLower = *std::next(fSplitX.begin(), index_x);
1397 auto xUpper = *std::next(fSplitX.begin(), index_x + 1);
1398 auto yLower = *std::next(fSplitY.begin(), index_y);
1399 auto yUpper = *std::next(fSplitY.begin(), index_y + 1);
1402 GetObservable() +
";energy";
1403 fSegLinearFit[index_x][index_y]->SetTitle(tH.c_str());
1404 fSegLinearFit[index_x][index_y]->Draw(
"AL*");
1407void TRestDataSetGainMap::Module::DrawLinearFit(TCanvas* c) {
1408 if (fSegLinearFit.size() == 0) {
1409 RESTError <<
"Spectra matrix is empty." << p->RESTendl;
1413 std::string t =
"linearFits_" + std::to_string(fPlaneId) +
"_" + std::to_string(fModuleId);
1414 c =
new TCanvas(t.c_str(), t.c_str());
1418 for (
const auto&
object : *c->GetListOfPrimitives())
1419 if (object->InheritsFrom(TVirtualPad::Class())) ++nPads;
1420 if (nPads != 0 && nPads != fSegLinearFit.size() * fSegLinearFit.at(0).size()) {
1421 RESTError <<
"Canvas " << c->GetName() <<
" has " << nPads <<
" pads, but "
1422 << fSegLinearFit.size() * fSegLinearFit.at(0).size() <<
" are needed." << p->RESTendl;
1424 }
else if (nPads == 0)
1425 c->Divide(fSegLinearFit.size(), fSegLinearFit.at(0).size());
1427 for (
size_t i = 0; i < fSegLinearFit.size(); i++) {
1428 for (
size_t j = 0; j < fSegLinearFit[i].size(); j++) {
1429 int pad = fSegLinearFit.size() * (fSegLinearFit[i].size() - 1) + 1 + i - fSegLinearFit.size() * j;
1431 DrawLinearFit(i, j, c);
1445 if (peakNumber < 0 || peakNumber >= (
int)fEnergyPeaks.size()) {
1446 RESTError <<
"Peak number out of range (peakNumber should be between 0 and "
1447 << fEnergyPeaks.size() - 1 <<
" )" << p->RESTendl;
1450 if (fSegLinearFit.size() == 0) {
1451 RESTError <<
"Linear fit matrix is empty." << p->RESTendl;
1454 if (!fFullLinearFit) {
1455 RESTError <<
"Full linear fit is empty." << p->RESTendl;
1459 double peakEnergy = fEnergyPeaks[peakNumber];
1460 std::string title =
"Gain map for energy " +
DoubleToString(peakEnergy,
"%g") +
";" +
1461 GetSpatialObservableX() +
";" + GetSpatialObservableY();
1462 std::string t =
"gainMap" + std::to_string(peakNumber) +
"_" + std::to_string(fPlaneId) +
"_" +
1463 std::to_string(fModuleId);
1464 TCanvas* gainMap =
new TCanvas(t.c_str(), t.c_str());
1466 TH2F* hGainMap =
new TH2F((
"h" + t).c_str(), title.c_str(), fNumberOfSegmentsX, fReadoutRange.X(),
1467 fReadoutRange.Y(), fNumberOfSegmentsY, fReadoutRange.X(), fReadoutRange.Y());
1469 double peakPosRef = fFullLinearFit->GetPointX(peakNumber);
1470 if (!fullModuleAsRef) {
1471 int index_x = fNumberOfSegmentsX > 0 ? (fNumberOfSegmentsX - 1) / 2 : 0;
1472 int index_y = fNumberOfSegmentsY > 0 ? (fNumberOfSegmentsY - 1) / 2 : 0;
1473 peakPosRef = fSegLinearFit[index_x][index_y]->GetPointX(peakNumber);
1476 auto itX = fSplitX.begin();
1477 for (
size_t i = 0; i < fSegLinearFit.size(); i++) {
1478 auto itY = fSplitY.begin();
1479 for (
size_t j = 0; j < fSegLinearFit.at(0).size(); j++) {
1481 auto xUpper = *std::next(itX);
1483 auto yUpper = *std::next(itY);
1484 float xMean = (xUpper + xLower) / 2.;
1485 float yMean = (yUpper + yLower) / 2.;
1486 auto [index_x, index_y] = GetIndexMatrix(xMean, yMean);
1487 if (!fSegLinearFit[index_x][index_y])
continue;
1488 hGainMap->Fill(xMean, yMean, fSegLinearFit[index_x][index_y]->GetPointX(peakNumber) / peakPosRef);
1493 hGainMap->SetStats(0);
1494 hGainMap->Draw(
"colz");
1495 hGainMap->SetBarOffset(0.2);
1496 hGainMap->Draw(
"TEXT SAME");
1504 RESTMetadata <<
"-----------------------------------------------" << p->RESTendl;
1505 RESTMetadata <<
" Plane ID: " << fPlaneId << p->RESTendl;
1506 RESTMetadata <<
" Module ID: " << fModuleId << p->RESTendl;
1507 RESTMetadata <<
" Definition cut: " << fDefinitionCut << p->RESTendl;
1508 RESTMetadata << p->RESTendl;
1510 RESTMetadata <<
" Calibration dataset: " << fDataSetFileName << p->RESTendl;
1511 RESTMetadata << p->RESTendl;
1513 RESTMetadata <<
" Energy peaks: ";
1514 for (
const auto& peak : fEnergyPeaks) RESTMetadata << peak <<
", ";
1515 RESTMetadata << p->RESTendl;
1516 bool anyRange =
false;
1517 for (
const auto& r : fRangePeaks)
1518 if (r.X() < r.Y()) {
1519 RESTMetadata <<
" Range peaks: ";
1524 for (
const auto& r : fRangePeaks) RESTMetadata <<
"(" << r.X() <<
", " << r.Y() <<
") ";
1525 if (anyRange) RESTMetadata << p->RESTendl;
1526 RESTMetadata <<
" Auto range peaks: " << (fAutoRangePeaks ?
"true" :
"false") << p->RESTendl;
1527 RESTMetadata <<
" Zero point: " << (fZeroPoint ?
"true" :
"false") << p->RESTendl;
1528 RESTMetadata <<
" Calibration range: (" << fCalibRange.X() <<
", " << fCalibRange.Y() <<
" )"
1530 RESTMetadata <<
" Number of bins: " << fNBins << p->RESTendl;
1531 RESTMetadata << p->RESTendl;
1533 RESTMetadata <<
" Number of segments X: " << fNumberOfSegmentsX << p->RESTendl;
1534 RESTMetadata <<
" Number of segments Y: " << fNumberOfSegmentsY << p->RESTendl;
1536 RESTMetadata <<
" Readout range (" << fReadoutRange.X() <<
", " << fReadoutRange.Y() <<
" )"
1538 RESTMetadata <<
"SplitX: ";
1539 for (
auto& i : fSplitX) {
1540 RESTMetadata <<
" " << i;
1542 RESTMetadata << p->RESTendl;
1543 RESTMetadata <<
"SplitY: ";
1544 for (
auto& i : fSplitY) {
1545 RESTMetadata <<
" " << i;
1547 RESTMetadata << p->RESTendl;
1548 RESTMetadata << p->RESTendl;
1550 RESTMetadata <<
" Slope: " << p->RESTendl;
1552 for (
auto& x : fSlope)
1553 if (maxSize < x.size()) maxSize = x.size();
1554 for (
size_t j = 0; j < maxSize; j++) {
1555 for (
size_t k = 0; k < fSlope.size(); k++) {
1556 if (j < fSlope[k].size())
1557 RESTMetadata <<
DoubleToString(fSlope[k][fSlope[k].size() - 1 - j],
"%.3e") <<
" ";
1559 RESTMetadata <<
" ";
1561 RESTMetadata << p->RESTendl;
1563 RESTMetadata <<
" Intercept: " << p->RESTendl;
1565 for (
auto& x : fIntercept)
1566 if (maxSize < x.size()) maxSize = x.size();
1567 for (
size_t j = 0; j < maxSize; j++) {
1568 for (
size_t k = 0; k < fIntercept.size(); k++) {
1569 if (j < fIntercept[k].size())
1570 RESTMetadata <<
DoubleToString(fIntercept[k][fIntercept[k].size() - 1 - j],
"%+.3e") <<
" ";
1572 RESTMetadata <<
" ";
1574 RESTMetadata << p->RESTendl;
1576 RESTMetadata << p->RESTendl;
1577 RESTMetadata <<
" Full slope: " <<
DoubleToString(fFullSlope,
"%.3e") << p->RESTendl;
1578 RESTMetadata <<
" Full intercept: " <<
DoubleToString(fFullIntercept,
"%+.3e") << p->RESTendl;
1580 RESTMetadata <<
"-----------------------------------------------" << p->RESTendl;
A class to help on cuts definitions. To be used with TRestAnalysisTree.
void SetSplitY()
Function to set the class members for segmentation of the detector plane along the Y axis.
void SetSplits()
Function to set the class members for segmentation of the detector plane along the X and Y axis.
void SetSplitX()
Function to set the class members for segmentation of the detector plane along the X axis.
void DrawSpectrum(const bool drawFits=true, const int color=-1, TCanvas *c=nullptr)
Function to draw the spectrum for each segment of the module on the same canvas. The canvas is divide...
void LoadConfigFromTiXmlElement(const TiXmlElement *module)
Function to read the parameters from the RML element (TiXmlElement) and set those class members.
void GenerateGainMap()
Function that calculates the calibration parameters for each segment defined at fSplitX and fSplitY a...
double GetSlope(const double x, const double y) const
Function to get the calibration parameter slope for a given x and y position on the detector plane.
std::set< double > fSplitX
Split points in the x direction.
void Print() const
Prints on screen the information about the members of Module.
void Refit(const TVector2 &position, const double energy, const TVector2 &range)
Function to fit again manually a peak for a given segment of the module.
std::set< double > fSplitY
Split points in the y direction.
void RefitFullSpc(const double energy, const TVector2 &range)
Function to fit again manually a peak for the whole module spectrum. The calibration curve is updated...
const TRestDataSetGainMap * p
Pointer to the parent class.
void UpdateCalibrationFitsFullSpc()
Function to update the calibration curve for the whole module. The calibration curve is cleared and t...
void DrawGainMap(const int peakNumber=0, const bool fullModuleAsRef=true)
Function to draw the relative gain map for a given energy peak of the module.
double GetIntercept(const double x, const double y) const
Function to get the calibration parameter intercept for a given x and y position on the detector plan...
std::pair< int, int > GetIndexMatrix(const double x, const double y) const
Function to get the index of the matrix of calibration parameters for a given x and y position on the...
Metadata class to calculate,store and apply the gain corrected calibration of a group of detectors.
double GetInterceptParameter(const int planeID, const int moduleID, const double x, const double y)
Function to get the intercept parameter of the module with planeID and moduleID at physical position ...
double GetSlopeParameter(const int planeID, const int moduleID, const double x, const double y)
Function to get the slope parameter of the module with planeID and moduleID at physical position (x,...
std::string fObservable
Observable that will be used to calculate the gain map.
double GetInterceptParameterFullSpc(const int planeID, const int moduleID)
Function to get the intercept parameter of the whole module with planeID and moduleID.
void CalibrateDataSet(const std::string &dataSetFileName, std::string outputFileName="", std::vector< std::string > excludeColumns={})
Function to calibrate a dataset with this gain map.
~TRestDataSetGainMap()
Default destructor.
void GenerateGainMap()
Function to calculate the calibration parameters of all modules.
TRestDataSetGainMap()
Default constructor.
void Initialize() override
Making default settings.
std::vector< Module > fModulesCal
List of modules.
void SetModule(const Module &moduleCal)
Function to set a module calibration. If the module calibration already exists (same planeId and modu...
void InitFromConfigFile() override
Initialization of TRestDataSetGainMap members through a RML file.
double GetSlopeParameterFullSpc(const int planeID, const int moduleID)
Function to get the slope parameter of the whole module with planeID and moduleID.
void Export(const std::string &fileName="")
Function to export the calibration to the file fileName.
void PrintMetadata() override
Prints on screen the information about the metadata members.
std::map< int, std::set< int > > GetModuleIDs() const
Function to get the map of the module IDs for each plane ID.
std::string fSpatialObservableY
Observable that will be used to segmentize the gain map in the y direction.
std::string fCalibFileName
Name of the file that contains the calibration data.
TRestCut * fCut
Cut to be applied to the calibration data.
std::string fOutputFileName
Name of the file where the gain map was (or will be) exported.
std::set< int > GetPlaneIDs() const
Function to get a list (set) of the plane IDs.
std::string fSpatialObservableX
Observable that will be used to segmentize the gain map in the x direction.
Module * GetModule(const size_t index=0)
Function to retrieve the module calibration by index. Default is 0.
void Import(const std::string &fileName)
Function to import the calibration parameters from the root file fileName.
It allows to group a number of runs that satisfy given metadata conditions.
void Import(const std::string &fileName)
This function imports metadata from a root file it import metadata info from the previous dataSet whi...
ROOT::RDF::RNode GetDataFrame() const
Gives access to the RDataFrame.
ROOT::RDF::RNode MakeCut(const TRestCut *cut)
This function applies a TRestCut to the dataframe and returns a dataframe with the applied cuts....
void GenerateDataSet()
This function generates the data frame with the filelist and column names (or observables) that have ...
void Export(const std::string &filename, std::vector< std::string > excludeColumns={})
It will generate an output file with the dataset compilation. Only the selected branches and the file...
@ REST_Info
+show most of the information for each steps
@ REST_Debug
+show the defined debug messages
TClass * GetClassQuick()
Get the type of a "class" object, returning the wrapped type identifier "TClass".
Bool_t MatchString(std::string str, std::string matcher)
This method matches a string with certain matcher. Returns true if matched. Supports wildcard charact...
Double_t StringToDouble(std::string in)
Gets a double from a string.
Int_t StringToInteger(std::string in)
Gets an integer from a string.
std::string DoubleToString(Double_t d, std::string format="%8.6e")
Gets a string from a double.
TVector2 StringTo2DVector(std::string in)
Gets a 2D-vector from a string.
std::string ToLower(std::string in)
Convert string to its lower case. Alternative of TString::ToLower.