REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestDetectorPositionMappingProcess.cxx
1
16
17#include "TRestDetectorPositionMappingProcess.h"
18
19#include <TLegend.h>
20#include <TPaveText.h>
21#include <TRandom.h>
22
23using namespace std;
24
26
27TRestDetectorPositionMappingProcess::TRestDetectorPositionMappingProcess() { Initialize(); }
28
29TRestDetectorPositionMappingProcess::~TRestDetectorPositionMappingProcess() {}
30
32 SetSectionName(this->ClassName());
33 SetLibraryVersion(LIBRARY_VERSION);
34
35 fHitsEvent = nullptr;
36
37 fReadout = nullptr;
38}
39
41 fReadout = GetMetadata<TRestDetectorReadout>();
42 fCalib = GetMetadata<TRestDetectorGainMap>();
43 fGas = GetMetadata<TRestDetectorGas>();
44 if (fReadout == nullptr) {
45 if (fCreateGainMap) {
46 RESTError << "You must set a TRestDetectorReadout metadata object to create gain map!"
47 << RESTendl;
48 abort();
49 }
50 } else {
51 double xmin, xmax, ymin, ymax;
52 fReadout->GetReadoutPlane(0)->GetBoundaries(xmin, xmax, ymin, ymax);
53 fAreaGainMap = new TH2F("areaGainMap", "areaGainMap", fNBinsX, xmin, xmax, fNBinsY, ymin, ymax);
54 fAreaCounts = new TH2D("areaCounts", "areaCounts", fNBinsX, xmin, xmax, fNBinsY, ymin, ymax);
55 fAreaThrIntegralSum =
56 new TH2D("areaThrIntegralSum", "areaThrIntegralSum", fNBinsX, xmin, xmax, fNBinsY, ymin, ymax);
57 }
58
59 if (fApplyGainCorrection) {
60 if (fCalib == nullptr || fCalib->f2DGainMapping == nullptr) {
61 RESTError << "You must set a TRestDetectorGainMap metadata object to apply gain correction!"
62 << RESTendl;
63 abort();
64 }
65 }
66}
67
69 fHitsEvent = (TRestDetectorHitsEvent*)inputEvent;
70
71 Double_t newEnergy = 0;
72
73 SetObservableValue("AreaCorrEnergy", newEnergy);
74
75 double x = fHitsEvent->GetMeanPositionX();
76 double y = fHitsEvent->GetMeanPositionY();
77 double z = fHitsEvent->GetMeanPositionZ();
78
79 double e = fHitsEvent->GetTotalEnergy();
80 double n = fHitsEvent->GetNumberOfHits();
81
82 // cout << x << " " << y << " " << e << " " << n << endl;
83
84 if (TMath::IsNaN(x) || TMath::IsNaN(y)) {
85 return fHitsEvent;
86 }
87
88 if (fCreateGainMap) {
89 if (e > fEnergyCutRange.X() && e < fEnergyCutRange.Y() && n > fNHitsCutRange.X() &&
90 n < fNHitsCutRange.Y()) {
91 // Calculate mean energy for different small areas as the area gain
92 int bin = fAreaThrIntegralSum->FindBin(x, y);
93 fAreaThrIntegralSum->AddBinContent(bin, e);
94 fAreaCounts->AddBinContent(bin, 1);
95 }
96 }
97
98 else if (fApplyGainCorrection) {
99 double correction = GetCorrection2(x, y);
100 correction *= GetCorrection3(x, y, z);
101
102 newEnergy = e * correction;
103 }
104
105 SetObservableValue("AreaCorrEnergy", newEnergy);
106
107 return fHitsEvent;
108}
109
110double TRestDetectorPositionMappingProcess::GetCorrection2(double x, double y) {
111 int bin = fCalib->f2DGainMapping->FindBin(x, y);
112 return fCalib->f2DGainMapping->GetBinContent(bin);
113}
114
115double TRestDetectorPositionMappingProcess::GetCorrection3(double x, double y, double z) {
116 double result = 1;
117 if (fCalib->f3DGainMapping != nullptr && fCalib->f3DGainMapping->GetEntries() > 0) {
118 int bin = fCalib->f3DGainMapping->FindBin(x, y, z);
119 result *= fCalib->f2DGainMapping->GetBinContent(bin);
120 }
121 if (fGas != nullptr && fGas->GetElectronLifeTime() != 0) {
122 double dt = z / fGas->GetDriftVelocity();
123 result *= exp(dt / fGas->GetElectronLifeTime());
124 }
125 return result;
126}
127
129 if (fCreateGainMap) {
130 // Calculate the mean of each bin's spectrum
131 double sum = 0;
132 double n = 0;
133 for (int i = 1; i <= fAreaGainMap->GetNbinsX(); i++) {
134 for (int j = 1; j <= fAreaGainMap->GetNbinsY(); j++) {
135 if (fAreaCounts->GetBinContent(i, j) > 100) {
136 double meanthrintegral =
137 fAreaThrIntegralSum->GetBinContent(i, j) / fAreaCounts->GetBinContent(i, j);
138 fAreaGainMap->SetBinContent(i, j, meanthrintegral);
139 sum += meanthrintegral;
140 n++;
141 }
142 }
143 }
144
145 // the mean value of all the valued bins
146 double meanmean = sum / n;
147
148 // normalize and fill the result
149 for (int i = 1; i <= fAreaGainMap->GetNbinsX(); i++) {
150 for (int j = 1; j <= fAreaGainMap->GetNbinsY(); j++) {
151 if (fAreaGainMap->GetBinContent(i, j) == 0) {
152 fAreaGainMap->SetBinContent(i, j, 1);
153 } else {
154 fAreaGainMap->SetBinContent(i, j, meanmean / fAreaGainMap->GetBinContent(i, j));
155 }
156 }
157 }
158
159 if (fCalib == nullptr) {
160 fCalib = new TRestDetectorGainMap();
161 }
162 fCalib->f2DGainMapping = fAreaGainMap;
163 fCalib->SetName("PositionCalibration");
164
165 TRestRun* r = new TRestRun();
166 r->SetOutputFileName(fMappingSave);
167 r->AddMetadata(fCalib);
168 r->AddMetadata(fReadout);
169 r->FormOutputFile();
170 if (fAreaGainMap != nullptr) fAreaGainMap->Write();
171 }
172}
173
174// setting amplification:
175// <parameter name="modulesAmp" value = "2-1:5-1.2:6-0.8:8-0.9" />
176// setting readout modules to draw:
177// <parameter name="modulesHist" value="2:5:6:8"/>
179 string mode = GetParameter("mode", "apply");
180 if (mode == "create") {
181 fCreateGainMap = true;
182 fApplyGainCorrection = false;
183 fSingleThreadOnly = true;
184 } else if (mode == "apply") {
185 fCreateGainMap = false;
186 fApplyGainCorrection = true;
187 fSingleThreadOnly = false;
188 } else {
189 RESTError << "illegal mode definition! supported: create, apply" << RESTendl;
190 abort();
191 }
192
193 fEnergyCutRange = StringTo2DVector(GetParameter("energyRange", "(0,1e9)"));
194 fNHitsCutRange = StringTo2DVector(GetParameter("nHitsRange", "(4,14)"));
195 fMappingSave = GetParameter("save", "calib.root");
196
197 fNBinsX = StringToDouble(GetParameter("nBinsX", "100"));
198 fNBinsY = StringToDouble(GetParameter("nBinsY", "100"));
199 fNBinsZ = StringToDouble(GetParameter("nBinsZ", "100"));
200}
void InitProcess() override
To be executed at the beginning of the run (outside event loop)
void Initialize() override
Making default settings.
void EndProcess() override
To be executed at the end of the run (outside event loop)
void InitFromConfigFile() override
To make settings from rml file. This method must be implemented in the derived class.
TRestEvent * ProcessEvent(TRestEvent *inputEvent) override
Process one event.
void GetBoundaries(double &xmin, double &xmax, double &ymin, double &ymax)
Finds the xy boundaries of the readout plane delimited by the modules.
TRestDetectorReadoutPlane * GetReadoutPlane(int p)
Returns a pointer to the readout plane by index.
void SetObservableValue(const std::string &name, const T &value)
Set observable value for AnalysisTree.
A base class for any REST event.
Definition: TRestEvent.h:38
endl_t RESTendl
Termination flag object for TRestStringOutput.
void SetLibraryVersion(TString version)
Set the library version of this metadata class.
void SetSectionName(std::string sName)
set the section name, clear the section content
std::string GetParameter(std::string parName, TiXmlElement *e, TString defaultValue=PARAMETER_NOT_FOUND_STR)
Returns the value for the parameter named parName in the given section.
Data provider and manager in REST.
Definition: TRestRun.h:18
TFile * FormOutputFile()
Create a new TFile as REST output file. Writing metadata objects into it.
Definition: TRestRun.cxx:1036
void AddMetadata(TRestMetadata *meta)
add metadata object to the metadata list
Definition: TRestRun.h:112
Double_t StringToDouble(std::string in)
Gets a double from a string.
TVector2 StringTo2DVector(std::string in)
Gets a 2D-vector from a string.