REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
Loading...
Searching...
No Matches
TRestDataSetGainMap.cxx
1/*************************************************************************
2 * This file is part of the REST software framework. *
3 * *
4 * Copyright (C) 2016 GIFNA/TREX (University of Zaragoza) *
5 * For more information see https://gifna.unizar.es/trex *
6 * *
7 * REST is free software: you can redistribute it and/or modify *
8 * it under the terms of the GNU General Public License as published by *
9 * the Free Software Foundation, either version 3 of the License, or *
10 * (at your option) any later version. *
11 * *
12 * REST is distributed in the hope that it will be useful, *
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15 * GNU General Public License for more details. *
16 * *
17 * You should have a copy of the GNU General Public License along with *
18 * REST in $REST_PATH/LICENSE. *
19 * If not, see https://www.gnu.org/licenses/. *
20 * For the list of contributors see $REST_PATH/CREDITS. *
21 *************************************************************************/
22
139
140#include "TRestDataSetGainMap.h"
141
142ClassImp(TRestDataSetGainMap);
147
162TRestDataSetGainMap::TRestDataSetGainMap(const char* configFilename, std::string name)
163 : TRestMetadata(configFilename) {
166}
167
172
179 this->Initialize();
181
182 // Load Module from RML
183 TiXmlElement* moduleDefinition = GetElement("module");
184 while (moduleDefinition != nullptr) {
185 fModulesCal.push_back(Module(*this));
186 fModulesCal.back().LoadConfigFromTiXmlElement(moduleDefinition);
187
188 moduleDefinition = GetNextElement(moduleDefinition);
189 }
190
191 fCut = (TRestCut*)InstantiateChildMetadata("TRestCut");
192
194}
195
200 for (auto& mod : fModulesCal) {
201 RESTInfo << "Generating gain map of plane " << mod.GetPlaneId() << " module " << mod.GetModuleId()
202 << RESTendl;
203 mod.GenerateGainMap();
205 mod.DrawSpectrum();
206 mod.DrawFullSpectrum();
207 mod.DrawGainMap();
208 }
209 }
210}
211
224void TRestDataSetGainMap::CalibrateDataSet(const std::string& dataSetFileName, std::string outputFileName,
225 std::vector<std::string> excludeColumns) {
226 if (fModulesCal.empty()) {
227 RESTError << "TRestDataSetGainMap::CalibrateDataSet: No modules defined." << RESTendl;
228 return;
229 }
230
231 TRestDataSet dataSet;
232 dataSet.EnableMultiThreading(true);
233 dataSet.Import(dataSetFileName);
234 auto dataFrame = dataSet.GetDataFrame();
235
236 // Define a new column with the identifier (pmID) of the module for each row (event)
237 std::string pmIDname = (std::string)GetName() + "_pmID";
238 std::string modCut = fModulesCal[0].GetModuleDefinitionCut();
239 if (modCut.empty()) modCut = "1"; // if no cut is defined, use "1" (all events)
240 int pmID = fModulesCal[0].GetPlaneId() * 10 + fModulesCal[0].GetModuleId();
241
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");
245 else
246 dataFrame = dataFrame.Redefine(pmIDname, modCut + " ? " + std::to_string(pmID) + " : -1");
247
248 for (size_t n = 1; n < fModulesCal.size(); n++) {
249 modCut = fModulesCal[n].GetModuleDefinitionCut();
250 if (modCut.empty()) modCut = "1"; // if no cut is defined, use "1" (all events)
251 pmID = fModulesCal[n].GetPlaneId() * 10 + fModulesCal[n].GetModuleId();
252 dataFrame = dataFrame.Redefine(pmIDname, (modCut + " ? " + std::to_string(pmID) + " : " + pmIDname));
253 }
254
255 // Define a new column with the calibrated observable
256 auto calibrate = [this](double val, double x, double y, int pmID) {
257 for (auto& m : fModulesCal) {
258 if (pmID == m.GetPlaneId() * 10 + m.GetModuleId())
259 return m.GetSlope(x, y) * val + m.GetIntercept(x, y);
260 }
261 // RESTError << "TRestDataSetGainMap::CalibrateDataSet: Module not found for pmID " << pmID <<
262 // RESTendl;
263 return std::numeric_limits<double>::quiet_NaN();
264 };
265 std::string calibObsName = (std::string)GetName() + "_";
266 calibObsName += GetObservable().erase(0, GetObservable().find("_") + 1); // remove the "rawAna_" part
267 dataFrame = dataFrame.Define(calibObsName, calibrate,
269
270 // Define a new column with the calibrated observable for the whole module calibration
271 auto calibrateFullSpc = [this](double val, int pmID) {
272 for (auto& m : fModulesCal) {
273 if (pmID == m.GetPlaneId() * 10 + m.GetModuleId())
274 return m.GetSlopeFullSpc() * val + m.GetInterceptFullSpc();
275 }
276 // RESTError << "TRestDataSetGainMap::CalibrateDataSet: Module not found for pmID " << pmID <<
277 // RESTendl;
278 return std::numeric_limits<double>::quiet_NaN();
279 };
280 std::string calibObsNameFullSpc = (std::string)GetName() + "_";
281 calibObsNameFullSpc +=
282 GetObservable().erase(0, GetObservable().find("_") + 1); // remove the "rawAna_" part
283 calibObsNameFullSpc += "_NoSegmentation";
284 dataFrame = dataFrame.Define(calibObsNameFullSpc, calibrateFullSpc, {fObservable, pmIDname});
285
286 dataSet.SetDataFrame(dataFrame);
287
288 // Format the output file name and export the dataSet
289 if (outputFileName.empty()) outputFileName = dataSetFileName;
290 if (outputFileName == dataSetFileName) { // TRestDataSet cannot be overwritten
291 std::string gmName = GetName();
292 outputFileName = outputFileName.substr(0, outputFileName.find_last_of(".")); // remove extension
293 outputFileName += "_" + gmName + "." + TRestTools::GetFileNameExtension(dataSetFileName);
294 }
295
296 // Export dataset. Exclude columns if requested.
297 std::set<std::string> excludeCol; // vector with the explicit column names to be excluded
298 auto columns = dataSet.GetDataFrame().GetColumnNames();
299 // Get the columns to be excluded from the list of columns. It accepts wildcards "*" and "?"
300 for (auto& eC : excludeColumns) {
301 if (eC.find("*") != std::string::npos || eC.find("?") != std::string::npos) {
302 for (auto& c : columns)
303 if (MatchString(c, eC)) excludeCol.insert(c);
304 } else if (std::find(columns.begin(), columns.end(), eC) != columns.end())
305 excludeCol.insert(eC);
306 }
307 // Remove the calibObsName, calibObsNameFullSpc and pmIDname from the list of columns to be excluded
308 excludeCol.erase(calibObsName);
309 excludeCol.erase(calibObsNameFullSpc);
310 excludeCol.erase(pmIDname);
311
312 RESTDebug << "Excluding columns: ";
313 for (auto& c : excludeCol) RESTDebug << c << ", ";
314 RESTDebug << RESTendl;
315
316 dataSet.Export(outputFileName, std::vector<std::string>(excludeCol.begin(), excludeCol.end()));
317
318 // Add this TRestDataSetGainMap metadata to the output file
319 TFile* f = TFile::Open(outputFileName.c_str(), "UPDATE");
320 this->Write();
321 f->Close();
322 delete f;
323}
324
330 if (index < fModulesCal.size()) return &fModulesCal[index];
331
332 RESTError << "No ModuleCalibration with index " << index;
333 if (fModulesCal.empty())
334 RESTError << ". There are no modules defined." << RESTendl;
335 else
336 RESTError << ". Max index is " << fModulesCal.size() - 1 << RESTendl;
337 return nullptr;
338}
339
344TRestDataSetGainMap::Module* TRestDataSetGainMap::GetModule(const int planeID, const int moduleID) {
345 for (auto& i : fModulesCal) {
346 if (i.GetPlaneId() == planeID && i.GetModuleId() == moduleID) return &i;
347 }
348 RESTError << "No ModuleCalibration with planeID " << planeID << " and moduleID " << moduleID << RESTendl;
349 return nullptr;
350}
351
357double TRestDataSetGainMap::GetSlopeParameter(const int planeID, const int moduleID, const double x,
358 const double y) {
359 Module* moduleCal = GetModule(planeID, moduleID);
360 if (moduleCal == nullptr) return 0; // return numeric_limits<double>::quiet_NaN()
361 return moduleCal->GetSlope(x, y);
362}
363
368double TRestDataSetGainMap::GetSlopeParameterFullSpc(const int planeID, const int moduleID) {
369 Module* moduleCal = GetModule(planeID, moduleID);
370 if (moduleCal == nullptr) return 0; // return numeric_limits<double>::quiet_NaN()
371 return moduleCal->GetSlopeFullSpc();
372}
373
379double TRestDataSetGainMap::GetInterceptParameter(const int planeID, const int moduleID, const double x,
380 const double y) {
381 Module* moduleCal = GetModule(planeID, moduleID);
382 if (moduleCal == nullptr) return 0; // return numeric_limits<double>::quiet_NaN()
383 return moduleCal->GetIntercept(x, y);
384}
385
390double TRestDataSetGainMap::GetInterceptParameterFullSpc(const int planeID, const int moduleID) {
391 Module* moduleCal = GetModule(planeID, moduleID);
392 if (moduleCal == nullptr) return 0; // return numeric_limits<double>::quiet_NaN()
393 return moduleCal->GetInterceptFullSpc();
394}
395
400std::set<int> TRestDataSetGainMap::GetPlaneIDs() const {
401 std::set<int> planeIDs;
402 for (const auto& mc : fModulesCal) planeIDs.insert(mc.GetPlaneId());
403 return planeIDs;
404}
405
410std::set<int> TRestDataSetGainMap::GetModuleIDs(const int planeId) const {
411 std::set<int> moduleIDs;
412 for (const auto& mc : fModulesCal)
413 if (mc.GetPlaneId() == planeId) moduleIDs.insert(mc.GetModuleId());
414 return moduleIDs;
415}
416
420std::map<int, std::set<int>> TRestDataSetGainMap::GetModuleIDs() const {
421 std::map<int, std::set<int>> moduleIds;
422 for (const int planeId : GetPlaneIDs())
423 moduleIds.insert(std::pair<int, std::set<int>>(planeId, GetModuleIDs(planeId)));
424 return moduleIds;
425}
426
427TRestDataSetGainMap& TRestDataSetGainMap::operator=(TRestDataSetGainMap& src) {
428 SetName(src.GetName());
429 fOutputFileName = src.GetOutputFileName();
430 fObservable = src.GetObservable();
431 fSpatialObservableX = src.GetSpatialObservableX();
432 fSpatialObservableY = src.GetSpatialObservableY();
433 fCut = src.GetCut();
434 fModulesCal.clear();
435 for (auto pID : src.GetPlaneIDs())
436 for (auto mID : src.GetModuleIDs(pID)) fModulesCal.push_back(*src.GetModule(pID, mID));
437 return *this;
438}
439
445 for (auto& i : fModulesCal) {
446 if (i.GetPlaneId() == moduleCal.GetPlaneId() && i.GetModuleId() == moduleCal.GetModuleId()) {
447 i = moduleCal;
448 return;
449 }
450 }
451 fModulesCal.push_back(moduleCal);
452}
453
458void TRestDataSetGainMap::Import(const std::string& fileName) {
459 if (fileName.empty()) {
460 RESTError << "No input calibration file defined" << RESTendl;
461 return;
462 }
463
464 if (TRestTools::isRootFile(fileName)) {
465 RESTInfo << "Opening " << fileName << RESTendl;
466 TFile* f = TFile::Open(fileName.c_str(), "READ");
467 if (f == nullptr) {
468 RESTError << "Cannot open calibration file " << fileName << RESTendl;
469 return;
470 }
471
472 TRestDataSetGainMap* cal = nullptr;
473 if (f != nullptr) {
474 TIter nextkey(f->GetListOfKeys());
475 TKey* key;
476 while ((key = (TKey*)nextkey())) {
477 std::string kName = key->GetClassName();
478 if (REST_Reflection::GetClassQuick(kName.c_str()) != nullptr &&
479 REST_Reflection::GetClassQuick(kName.c_str())->InheritsFrom("TRestDataSetGainMap")) {
480 cal = f->Get<TRestDataSetGainMap>(key->GetName());
481 *this = *cal;
482 }
483 }
484 }
485 } else
486 RESTError << "File extension not supported for " << fileName << RESTendl;
488}
489
496void TRestDataSetGainMap::Export(const std::string& fileName) {
497 if (!fileName.empty()) fOutputFileName = fileName;
498 if (fOutputFileName.empty()) {
499 RESTError << "No output file defined" << RESTendl;
500 return;
501 }
502
504 TFile* f = TFile::Open(fOutputFileName.c_str(), "UPDATE");
505 this->Write(GetName());
506 f->Close();
507 delete f;
508 RESTInfo << "Calibration saved to " << fOutputFileName << RESTendl;
509 } else
510 RESTError << "File extension for " << fOutputFileName << "is not supported." << RESTendl;
511}
512
518 RESTMetadata << " Calibration dataset: " << fCalibFileName << RESTendl;
519 if (fCut) {
520 RESTMetadata << " Cuts applied: ";
521 /* print only cutStrings and paramCut because
522 TRestDataSet::MakeCut() uses fCut->GetCutStrings() and fCut->GetParamCut() */
523 for (const auto& cut : fCut->GetCutStrings()) RESTMetadata << cut << ", " << RESTendl;
524 for (const auto& cut : fCut->GetParamCut())
525 RESTMetadata << cut.first << " " << cut.second << ", " << RESTendl;
526 }
527 RESTMetadata << " Output file: " << fOutputFileName << RESTendl;
528 RESTMetadata << RESTendl;
529 RESTMetadata << " Number of planes: " << GetNumberOfPlanes() << RESTendl;
530 RESTMetadata << " Number of modules: " << GetNumberOfModules() << RESTendl;
531 RESTMetadata << " Calibration observable: " << fObservable << RESTendl;
532 RESTMetadata << " Spatial observable X: " << fSpatialObservableX << RESTendl;
533 RESTMetadata << " Spatial observable Y: " << fSpatialObservableY << RESTendl;
534 RESTMetadata << "-----------------------------------------------" << RESTendl;
535 for (auto& i : fModulesCal) i.Print();
536 RESTMetadata << "***********************************************" << RESTendl;
537 RESTMetadata << RESTendl;
538}
539
547std::pair<int, int> TRestDataSetGainMap::Module::GetIndexMatrix(const double x, const double y) const {
548 int index_x = -1, index_y = -1;
549
550 if (fSplitX.upper_bound(x) != fSplitX.end()) {
551 index_x = std::distance(fSplitX.begin(), fSplitX.upper_bound(x)) - 1;
552 if (index_x < 0) {
553 RESTWarning << "index_x < 0 for x = " << x << " and fSplitX[0]=" << *fSplitX.begin()
554 << p->RESTendl;
555 index_x = 0;
556 }
557 } else {
558 RESTWarning << "x is out of split for x = " << x << p->RESTendl;
559 index_x = fSplitX.size() - 2;
560 }
561
562 if (fSplitY.upper_bound(y) != fSplitY.end()) {
563 index_y = std::distance(fSplitY.begin(), fSplitY.upper_bound(y)) - 1;
564 if (index_y < 0) {
565 RESTWarning << "index_y < 0 for y = " << y << " and fSplitY[0]=" << *fSplitY.begin()
566 << p->RESTendl;
567 index_y = 0;
568 }
569 } else {
570 RESTWarning << "y is out of split for y = " << y << p->RESTendl;
571 index_y = fSplitY.size() - 2;
572 }
573
574 return std::make_pair(index_x, index_y);
575}
583double TRestDataSetGainMap::Module::GetSlope(const double x, const double y) const {
584 auto [index_x, index_y] = GetIndexMatrix(x, y);
585 if (fSlope.empty()) {
586 RESTError << "Calibration slope matrix is empty. Returning 0" << p->RESTendl;
587 return 0;
588 }
589
590 if (index_x > (int)fSlope.size() || index_y > (int)fSlope.at(0).size()) {
591 RESTError << "Index out of range. Returning 0" << p->RESTendl;
592 return 0;
593 }
594
595 return fSlope[index_x][index_y];
596}
597
605double TRestDataSetGainMap::Module::GetIntercept(const double x, const double y) const {
606 auto [index_x, index_y] = GetIndexMatrix(x, y);
607 if (fIntercept.empty()) {
608 RESTError << "Calibration constant matrix is empty. Returning 0" << p->RESTendl;
609 return 0;
610 }
611
612 if (index_x > (int)fIntercept.size() || index_y > (int)fIntercept.at(0).size()) {
613 RESTError << "Index out of range. Returning 0" << p->RESTendl;
614 return 0;
615 }
616
617 return fIntercept[index_x][index_y];
618}
619
624 SetSplitX();
625 SetSplitY();
626}
634 if (fNumberOfSegmentsX < 1) {
635 RESTError << "SetSplitX: fNumberOfSegmentsX must be >=1." << p->RESTendl;
636 return;
637 }
638 std::set<double> split;
639 for (int i = 0; i <= fNumberOfSegmentsX; i++) { // <= so that the last segment is included
640 double x =
641 fReadoutRange.X() + ((fReadoutRange.Y() - fReadoutRange.X()) / (float)fNumberOfSegmentsX) * i;
642 split.insert(x);
643 }
644 SetSplitX(split);
645}
646
647void TRestDataSetGainMap::Module::SetSplitX(const std::set<double>& splitX) {
648 if (splitX.size() < 2) {
649 RESTError << "SetSplitX: split size must be >=2 (start and end of range must be included)."
650 << p->RESTendl;
651 return;
652 }
653 if (!fSlope.empty())
654 RESTWarning << "SetSplitX: changing split but current gain map and calibration paremeters correspond "
655 "to previous splitting. Use GenerateGainMap() to update them."
656 << p->RESTendl;
657 fSplitX = splitX;
658 fNumberOfSegmentsX = fSplitX.size() - 1;
659}
667 if (fNumberOfSegmentsY < 1) {
668 RESTError << "SetSplitY: fNumberOfSegmentsY must be >=1." << p->RESTendl;
669 return;
670 }
671 std::set<double> split;
672 for (int i = 0; i <= fNumberOfSegmentsY; i++) { // <= so that the last segment is included
673 double y =
674 fReadoutRange.X() + ((fReadoutRange.Y() - fReadoutRange.X()) / (float)fNumberOfSegmentsY) * i;
675 split.insert(y);
676 }
677 SetSplitY(split);
678}
679
680void TRestDataSetGainMap::Module::SetSplitY(const std::set<double>& splitY) {
681 if (splitY.size() < 2) {
682 RESTError << "SetSplitY: split size must be >=2 (start and end of range must be included)."
683 << p->RESTendl;
684 return;
685 }
686 if (!fSlope.empty())
687 RESTWarning << "SetSplitY: changing split but current gain map and calibration paremeters correspond "
688 "to previous splitting. Use GenerateGainMap() to update them."
689 << p->RESTendl;
690 fSplitY = splitY;
691 fNumberOfSegmentsY = fSplitY.size() - 1;
692}
693
711 //--- Initial checks and settings ---
712 std::string dsFileName = fDataSetFileName;
713 if (dsFileName.empty()) dsFileName = p->GetCalibrationFileName();
714 if (dsFileName.empty()) {
715 RESTError << "No calibration file defined" << p->RESTendl;
716 return;
717 }
718
719 TRestDataSet dataSet;
720 dataSet.EnableMultiThreading(true);
721 if (TRestTools::isDataSet(dsFileName)) {
722 dataSet.Import(dsFileName);
723 fDataSetFileName = dsFileName;
724 } else {
725 RESTWarning << dsFileName << " is not a dataset. Generating a temporal one..." << p->RESTendl;
726 // get all the observables needed for the gain map
727 std::vector<std::string> obsList;
728
729 obsList.push_back(p->GetObservable());
730 obsList.push_back(p->GetSpatialObservableX());
731 obsList.push_back(p->GetSpatialObservableY());
732
733 // look for observables (characterized by having a _ in the name) in the definition cut
734 auto modDefCutObs = TRestTools::GetObservablesInString(fDefinitionCut, true);
735 obsList.insert(obsList.end(), modDefCutObs.begin(), modDefCutObs.end());
736
737 // look for observables in the cut
738 for (const auto& cut : p->GetCut()->GetCutStrings()) {
739 auto cutObs = TRestTools::GetObservablesInString(cut, true);
740 obsList.insert(obsList.end(), cutObs.begin(), cutObs.end());
741 }
742 for (const auto& [variable, condition] : p->GetCut()->GetParamCut()) {
743 auto cutObs = TRestTools::GetObservablesInString(variable, true);
744 obsList.insert(obsList.end(), cutObs.begin(), cutObs.end());
745 // not sure if any obs can be in the condition. Just in case...
746 cutObs = TRestTools::GetObservablesInString(condition, true);
747 obsList.insert(obsList.end(), cutObs.begin(), cutObs.end());
748 }
749
750 // remove duplicates
751 std::sort(obsList.begin(), obsList.end());
752 obsList.erase(std::unique(obsList.begin(), obsList.end()), obsList.end());
753
754 // generate the dataset with the needed observables
755 dataSet.SetFilePattern(dsFileName);
756 dataSet.SetObservablesList(obsList);
757 dataSet.GenerateDataSet();
758 fDataSetFileName = dsFileName;
759 }
760
761 dataSet.SetDataFrame(dataSet.MakeCut(p->GetCut()));
762
763 if (fSplitX.empty()) SetSplitX();
764 if (fSplitY.empty()) SetSplitY();
765
766 //--- Get the calibration range if not provided (default is 0,0) ---
767 if (fCalibRange.X() >= fCalibRange.Y()) {
768 // Get spectrum for this file
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()));
773 // double xMin = hpunt->GetXaxis()->GetXmin();
774 double xMax = hpunt->GetXaxis()->GetXmax();
775
776 // Reduce the range to avoid the possible empty (nCounts<1%) end part of the spectrum
777 double fraction = 1, nAtEndSpc = 0, step = 0.66;
778 while (nAtEndSpc * 1. / hpunt->Integral() < 0.01 && fraction > 0.001) {
779 fraction *= step;
780 nAtEndSpc = hpunt->Integral(hpunt->FindFixBin(hpunt->GetXaxis()->GetXmax() * fraction),
781 hpunt->FindFixBin(hpunt->GetXaxis()->GetXmax()));
782 }
783 xMax = hpunt->GetXaxis()->GetXmax() * fraction /
784 step; // previous step is the last one that nAtEndSpc<1%
785 // Set the calibration range if needed
786 // fCalibRange.SetX(xMin);
787 fCalibRange.SetY(xMax);
788 hpunt.reset(); // delete hpunt;
789 RESTDebug << "Calibration range (auto)set to (" << fCalibRange.X() << "," << fCalibRange.Y() << ")"
790 << p->RESTendl;
791 }
792
793 // --- Definition of histogram whole module ---
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());
797
798 // build the spectrum for the whole module
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());
805
806 //--- Definition of histogram matrix ---
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(),
812 fCalibRange.Y()); // h[column][row] equivalent to h[x][y]
813 }
814 }
815
816 // build the spectrum for each segment
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++) {
821 // Get the segment limits from the splits
822 auto xLower = *itX;
823 auto xUpper = *std::next(itX);
824 auto yLower = *itY;
825 auto yUpper = *std::next(itY);
826
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;
837 auto histo = dataSet.GetDataFrame()
838 .Filter(segment_cut)
839 .Histo1D({"temp", "", h[i][j]->GetNbinsX(), h[i][j]->GetXaxis()->GetXmin(),
840 h[i][j]->GetXaxis()->GetXmax()},
841 GetObservable());
842 std::unique_ptr<TH1F> hpunt = std::unique_ptr<TH1F>(static_cast<TH1F*>(histo->Clone()));
843 h[i][j]->Add(hpunt.get());
844 hpunt.reset(); // delete hpunt;
845 itY++;
846 }
847 itX++;
848 }
849
850 //--- Fit every peak energy for every segment ---
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;
862 }
863 fSlope = calParSlope;
864 fIntercept = calParIntercept;
865 fSegSpectra = h;
866
867 //--- Fit every peak energy for the whole module ---
868 delete fFullLinearFit;
869 fFullLinearFit = new TGraph();
870 auto [intercept, slope] = FitPeaks(fFullSpectrum, fFullLinearFit);
871 fFullSlope = slope;
872 fFullIntercept = intercept;
873}
874
875std::pair<double, double> TRestDataSetGainMap::Module::FitPeaks(TH1F* hSeg, TGraph* gr) {
876 if (!hSeg) {
877 RESTError << "No histogram for fitting" << p->RESTendl;
878 return std::make_pair(0, 0);
879 }
880 if (hSeg->Integral() == 0) {
881 RESTError << "Empty spectrum " << hSeg->GetName() << p->RESTendl;
882 return std::make_pair(0, 0);
883 }
884 std::shared_ptr<TGraph> graph = std::shared_ptr<TGraph>(new TGraph());
885 RESTExtreme << "Fitting peaks for " << hSeg->GetName() << p->RESTendl;
886
887 // Search for peaks --> peakPos
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>());
893 const double ratio =
894 peakPos.size() == 0 ? 1 : peakPos.front() / fEnergyPeaks.front(); // to estimate peak position
895
896 // Initialize graph for linear fit
897 graph->SetName("grFit");
898 graph->SetTitle((";" + GetObservable() + ";energy").c_str());
899
900 // Fit every energy peak
901 int c = 0;
902 double mu = 0;
903 for (const auto& energy : fEnergyPeaks) {
904 RESTExtreme << "\t fitting energy " << DoubleToString(energy, "%g") << p->RESTendl;
905 // estimation of the peak position is between start and end
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()) { // if range is defined use it
910 start = fRangePeaks.at(c).X();
911 end = fRangePeaks.at(c).Y();
912 }
913
914 do {
915 if (fAutoRangePeaks) {
916 if (peakPos.size() > 0) {
917 // Find the peak position that is between start and end
918 pos = peakPos.at(0);
919 while (!(start < pos && pos < end)) {
920 // if none of the peak position is
921 // between start and end, use the greater one.
922 if (pos == peakPos.back()) {
923 pos = peakPos.at(0);
924 break;
925 }
926 pos = *std::next(std::find(peakPos.begin(), peakPos.end(),
927 pos)); // get next peak position
928 }
929 peakPos.erase(std::find(peakPos.begin(), peakPos.end(),
930 pos)); // remove this peak position from the list
931 // final estimation of the peak range (idem fitting window) with this peak
932 // position pos
933 start = pos * 0.8;
934 end = pos * 1.2;
935 const double relDist = peakPos.size() > 0 ? (pos - peakPos.front()) / pos : 999;
936 if (relDist < 0.2) { // if the next peak is too close reduce the window width
937 start = pos * (1 - relDist / 2);
938 end = pos * (1 + relDist / 2);
939 }
940 }
941 }
942
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("
946 << DoubleToString(start, "%.3g") << ", " << DoubleToString(end, "%.3g") << ")"
947 << p->RESTendl;
948
949 if (hSeg->GetFunction(name.c_str())) // remove previous fit
950 hSeg->GetListOfFunctions()->Remove(hSeg->GetFunction(name.c_str()));
951
952 hSeg->Fit(g, "R+Q0"); // use 0 to not draw the fit but save it
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)); // avoid small peaks on main peak tail
957 graph->SetPoint(c++, mu, energy);
958 }
959 s.reset(); // delete s;
960
961 if (fZeroPoint) graph->SetPoint(c++, 0, 0);
962 while (graph->GetN() < 2) { // minimun 2 points needed for linear fit
963 graph->SetPoint(c++, 0, 0);
964 SetZeroPoint(true);
965 RESTDebug << "Not enough points for linear fit. Adding and setting zero point to true" << p->RESTendl;
966 }
967
968 // Linear fit
969 std::unique_ptr<TF1> linearFit;
970 linearFit = std::unique_ptr<TF1>(new TF1("linearFit", "pol1"));
971 graph->Fit("linearFit", "SQ"); // Q for quiet mode
972
973 if (gr) *gr = *(TGraph*)graph->Clone(); // if nullptr is passed, do not copy the graph
974 return std::make_pair(linearFit->GetParameter(0), linearFit->GetParameter(1));
975}
976
984void TRestDataSetGainMap::Module::Refit(const TVector2& position, const double energyPeak,
985 const TVector2& range) {
986 auto [index_x, index_y] = GetIndexMatrix(position.X(), position.Y());
987 int peakNumber = -1;
988 for (size_t i = 0; i < fEnergyPeaks.size(); i++)
989 if (fEnergyPeaks.at(i) == energyPeak) {
990 peakNumber = i;
991 break;
992 }
993 if (peakNumber == -1) {
994 RESTError << "Energy " << energyPeak << " not found in the list of energy peaks" << p->RESTendl;
995 return;
996 }
997 Refit((size_t)index_x, (size_t)index_y, (size_t)peakNumber, range);
998}
999
1009void TRestDataSetGainMap::Module::Refit(const size_t x, const size_t y, const size_t peakNumber,
1010 const TVector2& range) {
1011 if (fSegSpectra.empty()) {
1012 RESTError << "No gain map found. Use GenerateGainMap() first." << p->RESTendl;
1013 return;
1014 }
1015 if (x >= fSegSpectra.size() || y >= fSegSpectra.at(0).size()) {
1016 RESTError << "Segment with index (" << x << ", " << y << ") not found" << p->RESTendl;
1017 return;
1018 }
1019 if (peakNumber >= fEnergyPeaks.size()) {
1020 RESTError << "Peak with index " << peakNumber << " not found" << p->RESTendl;
1021 return;
1022 }
1023
1024 // Refit the desired peak
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())) // clear previous fits for this peakNumber
1029 h->GetListOfFunctions()->Remove(h->GetFunction(name.c_str()));
1030 h->Fit(g, "R+Q0"); // use 0 to not draw the fit but save it
1031
1032 // Change the point of the graph
1033 UpdateCalibrationFits(x, y);
1034}
1035
1043void TRestDataSetGainMap::Module::RefitFullSpc(const double energyPeak, const TVector2& range) {
1044 int peakNumber = -1;
1045 for (size_t i = 0; i < fEnergyPeaks.size(); i++)
1046 if (fEnergyPeaks.at(i) == energyPeak) {
1047 peakNumber = i;
1048 break;
1049 }
1050 if (peakNumber == -1) {
1051 RESTError << "Energy " << energyPeak << " not found in the list of energy peaks" << p->RESTendl;
1052 return;
1053 }
1054 RefitFullSpc((size_t)peakNumber, range);
1055}
1056
1064void TRestDataSetGainMap::Module::RefitFullSpc(const size_t peakNumber, const TVector2& range) {
1065 if (!fFullSpectrum) {
1066 RESTError << "No gain map found. Use GenerateGainMap() first." << p->RESTendl;
1067 return;
1068 }
1069 if (peakNumber >= fEnergyPeaks.size()) {
1070 RESTError << "Peak with index " << peakNumber << " not found" << p->RESTendl;
1071 return;
1072 }
1073
1074 // Refit the desired peak
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())) // clear previous fits for this peakNumber
1078 fFullSpectrum->GetListOfFunctions()->Remove(fFullSpectrum->GetFunction(name.c_str()));
1079 fFullSpectrum->Fit(g, "R+Q0"); // use 0 to not draw the fit but save it
1080
1081 // Change the point of the graph
1082 UpdateCalibrationFitsFullSpc();
1083}
1084
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;
1096 return;
1097 }
1098 if (x >= fSegSpectra.size() || y >= fSegSpectra.at(0).size()) {
1099 RESTError << "Segment with index (" << x << ", " << y << ") not found" << p->RESTendl;
1100 return;
1101 }
1102
1103 TH1F* h = fSegSpectra.at(x).at(y);
1104 TGraph* gr = fSegLinearFit.at(x).at(y);
1105
1106 auto [intercept, slope] = UpdateCalibrationFits(h, gr);
1107 fSlope[x][y] = slope;
1108 fIntercept[x][y] = intercept;
1109}
1110
1111std::pair<double, double> TRestDataSetGainMap::Module::UpdateCalibrationFits(TH1F* h, TGraph* gr) {
1112 if (!h) {
1113 RESTError << "No histogram for updating fits" << p->RESTendl;
1114 return std::make_pair(0, 0);
1115 }
1116 if (!gr) {
1117 RESTError << "No graph for updating fits" << p->RESTendl;
1118 return std::make_pair(0, 0);
1119 }
1120 if (h->Integral() == 0) {
1121 RESTError << "Empty spectrum " << h->GetName() << p->RESTendl;
1122 return std::make_pair(0, 0);
1123 }
1124
1125 // Clear the points of the graph
1126 for (size_t i = 0; i < fEnergyPeaks.size(); i++) gr->RemovePoint(i);
1127 // Add the new points to the graph
1128 int c = 0;
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());
1132 if (!g) {
1133 RESTWarning << "No fit ( " << fitName << " ) found for energy peak " << fEnergyPeaks[i]
1134 << " in histogram " << h->GetName() << p->RESTendl;
1135 continue;
1136 }
1137 gr->SetPoint(c++, g->GetParameter(1), fEnergyPeaks[i]);
1138 }
1139
1140 // Add zero points if needed (if there are less than 2 points)
1141 while (gr->GetN() < 2) {
1142 gr->SetPoint(c++, 0, 0);
1143 }
1144
1145 // Refit the calibration curve
1146 TF1* lf = nullptr;
1147 if (gr->GetFunction("linearFit"))
1148 lf = gr->GetFunction("linearFit");
1149 else
1150 lf = new TF1("linearFit", "pol1");
1151 gr->Fit(lf, "SQ"); // Q for quiet mode
1152
1153 return std::make_pair(lf->GetParameter(0), lf->GetParameter(1));
1154}
1155
1162 auto [intercept, slope] = UpdateCalibrationFits(fFullSpectrum, fFullLinearFit);
1163 fFullSlope = slope;
1164 fFullIntercept = intercept;
1165}
1166
1182 if (module == nullptr) {
1183 RESTError << "TRestDataSetGainMap::Module::LoadConfigFromTiXmlElement: module is nullptr"
1184 << p->RESTendl;
1185 return;
1186 }
1187
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));
1192
1193 el = !module->Attribute("moduleDefinitionCut") ? "Not defined" : module->Attribute("moduleDefinitionCut");
1194 if (!(el.empty() || el == "Not defined")) this->SetModuleDefinitionCut(el);
1195
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));
1202
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));
1207
1208 el = !module->Attribute("dataSetFileName") ? "Not defined" : module->Attribute("dataSetFileName");
1209 if (!(el.empty() || el == "Not defined")) this->SetDataSetFileName(el);
1210
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");
1215
1216 // Get peaks energy and range
1217 TiXmlElement* peakDefinition = (TiXmlElement*)module->FirstChildElement("peak");
1218 while (peakDefinition != nullptr) {
1219 double energy = 0;
1220 TVector2 range = TVector2(0, 0);
1221
1222 std::string ell =
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;
1226 exit(1);
1227 }
1228 energy = StringToDouble(ell);
1229
1230 ell = !peakDefinition->Attribute("range") ? "Not defined" : peakDefinition->Attribute("range");
1231 if (!(ell.empty() || ell == "Not defined")) range = StringTo2DVector(ell);
1232
1233 this->AddPeak(energy, range);
1234 peakDefinition = (TiXmlElement*)peakDefinition->NextSiblingElement();
1235 }
1236}
1237
1238void TRestDataSetGainMap::Module::DrawSpectrum(const TVector2& position, bool drawFits, int color,
1239 TCanvas* c) {
1240 std::pair<size_t, size_t> index = GetIndexMatrix(position.X(), position.Y());
1241 DrawSpectrum(index.first, index.second, drawFits, color, c);
1242}
1243
1244void TRestDataSetGainMap::Module::DrawSpectrum(const int index_x, const int index_y, bool drawFits, int color,
1245 TCanvas* c) {
1246 if (fSegSpectra.size() == 0) {
1247 RESTError << "Spectra matrix is empty." << p->RESTendl;
1248 return;
1249 }
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;
1253 return;
1254 }
1255 if (!fSegSpectra[index_x][index_y]) {
1256 RESTError << "No Spectrum for segment (" << index_x << ", " << index_y << ")." << p->RESTendl;
1257 return;
1258 }
1259
1260 if (!c) {
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());
1264 }
1265
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);
1270 std::string tH = "Spectrum x=[" + DoubleToString(xLower, "%g") + ", " + DoubleToString(xUpper, "%g") +
1271 ") y=[" + DoubleToString(yLower, "%g") + ", " + DoubleToString(yUpper, "%g") + ");" +
1272 GetObservable() + ";counts";
1273 fSegSpectra[index_x][index_y]->SetTitle(tH.c_str());
1274
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");
1278
1279 if (drawFits)
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());
1282 if (!fit)
1283 RESTWarning << "Fit for energy peak " << fEnergyPeaks[c] << " not found." << p->RESTendl;
1284 if (!fit) continue;
1285 fit->SetLineColor(c + 2 != colorT ? c + 2 : ++colorT); /* does not work with kRed, kBlue, etc.
1286 as they are not defined with the same number as the first 10 basic colors. See
1287 https://root.cern.ch/doc/master/classTColor.html#C01 and
1288 https://root.cern.ch/doc/master/classTColor.html#C02 */
1289 fit->Draw("same");
1290 }
1291}
1292
1312void TRestDataSetGainMap::Module::DrawSpectrum(const bool drawFits, const int color, TCanvas* c) {
1313 if (fSegSpectra.size() == 0) {
1314 RESTError << "Spectra matrix is empty." << p->RESTendl;
1315 return;
1316 }
1317 if (!c) {
1318 std::string t = "spectrum_" + std::to_string(fPlaneId) + "_" + std::to_string(fModuleId);
1319 c = new TCanvas(t.c_str(), t.c_str());
1320 }
1321
1322 size_t nPads = 0;
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;
1328 return;
1329 } else if (nPads == 0)
1330 c->Divide(fSegSpectra.size(), fSegSpectra.at(0).size());
1331
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;
1335 c->cd(pad);
1336 DrawSpectrum(i, j, drawFits, color, c);
1337 }
1338 }
1339}
1340void TRestDataSetGainMap::Module::DrawFullSpectrum(const bool drawFits, const int color, TCanvas* c) {
1341 if (!fFullSpectrum) {
1342 RESTError << "Spectrum is empty." << p->RESTendl;
1343 return;
1344 }
1345
1346 if (!c) {
1347 std::string t = "fullSpc_" + std::to_string(fPlaneId) + "_" + std::to_string(fModuleId);
1348 c = new TCanvas(t.c_str(), t.c_str());
1349 }
1350 c->cd();
1351
1352 fFullSpectrum->SetTitle(("Full spectrum;" + GetObservable() + ";counts").c_str());
1353
1354 if (color > 0) fFullSpectrum->SetLineColor(color);
1355 size_t colorT = fFullSpectrum->GetLineColor();
1356 fFullSpectrum->Draw("same");
1357
1358 if (drawFits)
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;
1362 if (!fit) continue;
1363 fit->SetLineColor(c + 2 != colorT ? c + 2 : ++colorT); /* does not work with kRed, kBlue, etc.
1364 as they are not defined with the same number as the first 10 basic colors. See
1365 https://root.cern.ch/doc/master/classTColor.html#C01 and
1366 https://root.cern.ch/doc/master/classTColor.html#C02 */
1367 fit->Draw("same");
1368 }
1369}
1370
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);
1374}
1375
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;
1379 return;
1380 }
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;
1384 return;
1385 }
1386 if (!fSegLinearFit[index_x][index_y]) {
1387 RESTError << "No linear fit for segment (" << index_x << ", " << index_y << ")." << p->RESTendl;
1388 return;
1389 }
1390
1391 if (!c) {
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());
1395 }
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);
1400 std::string tH = "Linear Fit x=[" + DoubleToString(xLower, "%g") + ", " + DoubleToString(xUpper, "%g") +
1401 ") y=[" + DoubleToString(yLower, "%g") + ", " + DoubleToString(yUpper, "%g") + ");" +
1402 GetObservable() + ";energy";
1403 fSegLinearFit[index_x][index_y]->SetTitle(tH.c_str());
1404 fSegLinearFit[index_x][index_y]->Draw("AL*");
1405}
1406
1407void TRestDataSetGainMap::Module::DrawLinearFit(TCanvas* c) {
1408 if (fSegLinearFit.size() == 0) {
1409 RESTError << "Spectra matrix is empty." << p->RESTendl;
1410 return;
1411 }
1412 if (!c) {
1413 std::string t = "linearFits_" + std::to_string(fPlaneId) + "_" + std::to_string(fModuleId);
1414 c = new TCanvas(t.c_str(), t.c_str());
1415 }
1416
1417 size_t nPads = 0;
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;
1423 return;
1424 } else if (nPads == 0)
1425 c->Divide(fSegLinearFit.size(), fSegLinearFit.at(0).size());
1426
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;
1430 c->cd(pad);
1431 DrawLinearFit(i, j, c);
1432 }
1433 }
1434}
1435
1444void TRestDataSetGainMap::Module::DrawGainMap(const int peakNumber, bool fullModuleAsRef) {
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;
1448 return;
1449 }
1450 if (fSegLinearFit.size() == 0) {
1451 RESTError << "Linear fit matrix is empty." << p->RESTendl;
1452 return;
1453 }
1454 if (!fFullLinearFit) {
1455 RESTError << "Full linear fit is empty." << p->RESTendl;
1456 return;
1457 }
1458
1459 double peakEnergy = fEnergyPeaks[peakNumber];
1460 std::string title = "Gain map for energy " + DoubleToString(peakEnergy, "%g") + ";" +
1461 GetSpatialObservableX() + ";" + GetSpatialObservableY(); // + " keV";
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());
1465 gainMap->cd();
1466 TH2F* hGainMap = new TH2F(("h" + t).c_str(), title.c_str(), fNumberOfSegmentsX, fReadoutRange.X(),
1467 fReadoutRange.Y(), fNumberOfSegmentsY, fReadoutRange.X(), fReadoutRange.Y());
1468
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);
1474 }
1475
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++) {
1480 auto xLower = *itX;
1481 auto xUpper = *std::next(itX);
1482 auto yLower = *itY;
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);
1489 itY++;
1490 }
1491 itX++;
1492 }
1493 hGainMap->SetStats(0);
1494 hGainMap->Draw("colz");
1495 hGainMap->SetBarOffset(0.2);
1496 hGainMap->Draw("TEXT SAME");
1497}
1498
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;
1509
1510 RESTMetadata << " Calibration dataset: " << fDataSetFileName << p->RESTendl;
1511 RESTMetadata << p->RESTendl;
1512
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: ";
1520 anyRange = true;
1521 break;
1522 }
1523 if (anyRange)
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() << " )"
1529 << p->RESTendl;
1530 RESTMetadata << " Number of bins: " << fNBins << p->RESTendl;
1531 RESTMetadata << p->RESTendl;
1532
1533 RESTMetadata << " Number of segments X: " << fNumberOfSegmentsX << p->RESTendl;
1534 RESTMetadata << " Number of segments Y: " << fNumberOfSegmentsY << p->RESTendl;
1535 // RESTMetadata << " Draw: " << (fDrawVar ? "true" : "false") << p->RESTendl;
1536 RESTMetadata << " Readout range (" << fReadoutRange.X() << ", " << fReadoutRange.Y() << " )"
1537 << p->RESTendl;
1538 RESTMetadata << "SplitX: ";
1539 for (auto& i : fSplitX) {
1540 RESTMetadata << " " << i;
1541 }
1542 RESTMetadata << p->RESTendl;
1543 RESTMetadata << "SplitY: ";
1544 for (auto& i : fSplitY) {
1545 RESTMetadata << " " << i;
1546 }
1547 RESTMetadata << p->RESTendl;
1548 RESTMetadata << p->RESTendl;
1549
1550 RESTMetadata << " Slope: " << p->RESTendl;
1551 size_t maxSize = 0;
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") << " ";
1558 else
1559 RESTMetadata << " ";
1560 }
1561 RESTMetadata << p->RESTendl;
1562 }
1563 RESTMetadata << " Intercept: " << p->RESTendl;
1564 maxSize = 0;
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") << " ";
1571 else
1572 RESTMetadata << " ";
1573 }
1574 RESTMetadata << p->RESTendl;
1575 }
1576 RESTMetadata << p->RESTendl;
1577 RESTMetadata << " Full slope: " << DoubleToString(fFullSlope, "%.3e") << p->RESTendl;
1578 RESTMetadata << " Full intercept: " << DoubleToString(fFullIntercept, "%+.3e") << p->RESTendl;
1579
1580 RESTMetadata << "-----------------------------------------------" << p->RESTendl;
1581}
A class to help on cuts definitions. To be used with TRestAnalysisTree.
Definition TRestCut.h:31
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...
A base class for any REST metadata class.
virtual void PrintMetadata()
Implemented it in the derived metadata class to print out specific metadata information.
endl_t RESTendl
Termination flag object for TRestStringOutput.
TiXmlElement * GetElement(std::string eleDeclare, TiXmlElement *e=nullptr)
Get an xml element from a given parent element, according to its declaration.
Int_t LoadConfigFromFile(const std::string &configFilename, const std::string &sectionName="")
Give the file name, find out the corresponding section. Then call the main starter.
TRestMetadata * InstantiateChildMetadata(int index, std::string pattern="")
This method will retrieve a new TRestMetadata instance of a child element of the present TRestMetadat...
virtual void InitFromConfigFile()
To make settings from rml file. This method must be implemented in the derived class.
TRestStringOutput::REST_Verbose_Level GetVerboseLevel()
returns the verboselevel in type of REST_Verbose_Level enumerator
void SetSectionName(std::string sName)
set the section name, clear the section content
std::string fConfigFileName
Full name of the rml file.
virtual Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0)
overwriting the write() method with fStore considered
TiXmlElement * GetNextElement(TiXmlElement *e)
Get the next sibling xml element of this element, with same eleDeclare.
@ REST_Info
+show most of the information for each steps
@ REST_Debug
+show the defined debug messages
static std::string GetFileNameExtension(const std::string &fullname)
Gets the file extension as the substring found after the latest ".".
static bool isDataSet(const std::string &filename)
It checks if the file contains a dataset object.
static std::vector< std::string > GetObservablesInString(const std::string &observablesStr, bool removeDuplicates=true)
Returns a vector with the observables names found in the input string. The observables names must con...
static bool isRootFile(const std::string &filename)
Returns true if the filename has *.root* extension.
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.