GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
GriffinCTFix.cxx
Go to the documentation of this file.
1#include <iostream>
2
3#include "TFile.h"
4#include "TH2.h"
5#include "TF1.h"
6#include "TCutG.h"
7#include "TGraphErrors.h"
8#include "TString.h"
9#include "TObjString.h"
10#include "TObjArray.h"
11#include "TProfile.h"
12
13#include "TChannel.h"
14#include "TGriffin.h"
15#include "TUserSettings.h"
16
17double CrossTalkFit(double* x, double* par) // NOLINT(readability-non-const-parameter)
18{
19 // This function is the linear fit function, but uses the CT coefficients as parameters instead of slope and
20 // intercept
21 double k0 = par[0] / (1. - par[0]);
22 double k1 = par[1] / (1. - par[1]);
23 double energy = par[2];
24
25 double slope = -(1. + k0) / (1. + k1);
26 double intercept = energy * (1. + k0 / (1. + k1));
27
28 return x[0] * slope + intercept;
29}
30
31double* CrossTalkFix(int det, double energy, TFile* inputFile)
32{
33 // The outfile is implicit since it was the last file that was open.
34 static double largestCorrection = 0.0;
35
36 // Clear all of the CT calibrations for this detector.
37 for(int i = 0; i < 4; ++i) {
38 TChannel* chan = TChannel::FindChannelByName(Form("GRG%02d%sN00A", det, TGriffin::GetColorFromNumber(i)));
39 if(chan == nullptr) {
40 std::cout << DRED << "Couldn't find a channel for "
41 << Form("GRG%02d%sN00A", det, TGriffin::GetColorFromNumber(i)) << RESET_COLOR << std::endl;
42 continue;
43 }
44 chan->DestroyCTCal();
45 }
46
47 static int largestDet = -1;
48 static int largestCrystal1 = -1;
49 static int largestCrystal2 = -1;
50
51 std::string namebase = Form("det_%d", det);
52
53 // This range seems to be working fairly well since no shift should be larger than say 6 or 7 keV
54 double lowCut = energy - 15;
55 double highCut = energy + 15;
56
57 // create diagonal cut
58 std::array<double, 5> xpts = {lowCut, 0, 0, highCut, lowCut};
59 std::array<double, 5> ypts = {0, lowCut, highCut, 0, 0};
60 TCutG cut("cut", 5, xpts.data(), ypts.data());
61
62 auto* d = new double[16]; // matrix of coefficients
63 auto* eD = new double[16]; // matrix of errors
64
65 for(int crystal1 = 0; crystal1 < 4; crystal1++) {
66 for(int crystal2 = crystal1 + 1; crystal2 < 4; crystal2++) {
67 // Load all of the addback matrices in and put them into a vector of TH2*
68 std::string name = Form("%s_%d_%d", namebase.c_str(), crystal1, crystal2);
69 TH2* mat = dynamic_cast<TH2*>(inputFile->Get(name.c_str()));
70 if(mat == nullptr) {
71 std::cout << "can not find: " << name << std::endl;
72 return nullptr;
73 }
74 std::cout << mat->GetName() << std::endl;
75 int xbins = mat->GetNbinsX();
76 int ybins = mat->GetNbinsY();
77
78 TH2* cmat = dynamic_cast<TH2*>(mat->Clone(Form("%s_clone", mat->GetName())));
79 cmat->Reset();
80
81 // I make a graph out of the "addback line" because I don't like the way TProfile handles the empty bins
82 auto* fitGraph = new TGraphErrors;
83 fitGraph->SetNameTitle(Form("%s_graph", mat->GetName()), "Graph");
84
85 // This loop turns the addback plot and TCut into the TGraphErrors
86 for(int i = 1; i <= xbins; i++) {
87 bool insideYet = false;
88 for(int j = 1; j <= ybins; j++) {
89 double xc = mat->GetXaxis()->GetBinCenter(i);
90 double yc = mat->GetYaxis()->GetBinCenter(j);
91 if(cut.IsInside(xc, yc) != 0) {
92 if(!insideYet) {
93 insideYet = true;
94 }
95 cmat->Fill(xc, yc, mat->GetBinContent(i, j));
96 }
97 }
98 cmat->GetXaxis()->SetRange(i, i);
99 // This makes sure that there are at least 4 counts in the "y bin". I'd prefer this to be higher,
100 // but that requires more 60Co statistics. The reason I do this is because RMS and SD of the mean really only
101 // works for us if we have enough counts that the mean is actually a good representation of the true value.
102 // This is something TProfile does not do for us, and seems to skew the result a bit.
103 if(cmat->Integral() > 6) {
104 fitGraph->SetPoint(fitGraph->GetN(), cmat->GetYaxis()->GetBinCenter(i), cmat->GetMean(2));
105 fitGraph->SetPointError(fitGraph->GetN() - 1, cmat->GetXaxis()->GetBinWidth(i), cmat->GetMeanError(2));
106 }
107 }
108 // unzoom the x-axis
109 cmat->GetXaxis()->SetRange(1, cmat->GetXaxis()->GetNbins());
110 cmat->Write();
111
112 // This fits the TGraph
113 auto* fpx = new TF1(Form("pxfit_%i_%i_%i", det, crystal1, crystal2), CrossTalkFit, 6, 1167, 3);
114 fpx->SetParameter(0, 0.0001);
115 fpx->SetParameter(1, 0.0001);
116 fpx->SetParameter(2, energy);
117 fpx->FixParameter(2, energy);
118 fitGraph->Fit(fpx);
119 fitGraph->Write();
120 delete fitGraph;
121
122 // create and write profile
123 TProfile* px = cmat->ProfileX();
124 px->Write();
125 // Make a residuals plot
126 TH1* residualPlot = new TH1D(Form("%s_resid", fpx->GetName()), "residuals", 2000, 0, 2000);
127 for(int i = 0; i < residualPlot->GetNbinsX(); ++i) {
128 if(px->GetBinContent(i) != 0) {
129 residualPlot->SetBinContent(i, px->GetBinContent(i) - fpx->Eval(residualPlot->GetBinCenter(i)));
130 }
131 residualPlot->SetBinError(i, px->GetBinError(i));
132 }
133 residualPlot->Write();
134 delete residualPlot;
135 delete cmat;
136
137 // for some reason from here on out crystal1 and crystal2 are used in reverse?
138 // was this maybe because before we used xind and yind instead which were from the tokenized name using entry N-1 (for xind) and N-2 (for yind)?
139 // temporary swapped the two - except for the largest stuff?
140 std::cout << "=====================" << std::endl;
141 std::cout << mat->GetName() << std::endl;
142 std::cout << "d" << crystal1 << crystal2 << " at zero " << (fpx->Eval(energy)) / energy << std::endl;
143 std::cout << "=====================" << std::endl;
144
145 // Fill the parameter matrix with the parameters from the fit.
146 d[crystal2 * 4 + crystal1] = fpx->GetParameter(0);
147 d[crystal1 * 4 + crystal2] = fpx->GetParameter(1);
148 eD[crystal2 * 4 + crystal1] = fpx->GetParError(0);
149 eD[crystal1 * 4 + crystal2] = fpx->GetParError(1);
150
151 // Keep track of the largest correction and output that to screen,
152 // This helps identify problem channels, or mistakes
153 if(fpx->GetParameter(0) > largestCorrection) {
154 largestCorrection = fpx->GetParameter(0);
155 largestDet = det;
156 largestCrystal2 = crystal2;
157 largestCrystal1 = crystal1;
158 }
159 if(fpx->GetParameter(1) > largestCorrection) {
160 largestCorrection = fpx->GetParameter(1);
161 largestDet = det;
162 largestCrystal2 = crystal2;
163 largestCrystal1 = crystal1;
164 }
165 } //for(int crystal2 = crystal1 + 1; crystal2 < 4; crystal2++)
166 } // for(int crystal1 = 0; crystal1 < 4; crystal1++)
167
168 std::cout << " -------------------- " << std::endl;
169 std::cout << " -------------------- " << std::endl;
170 std::cout << std::endl;
171
172 // Set the diagonal elements to 0 since you can't cross-talk yourself
173 // store current precision and set to fixed 10 precision
174 std::streamsize prec = std::cout.precision();
175 std::cout.precision(10);
176 std::cout.setf(std::ios::fixed, std::ios::floatfield);
177 for(int i = 0; i < 4; i++) {
178 for(int j = 0; j < 4; j++) {
179 if(i == j) {
180 d[i * 4 + j] = 0.0000;
181 eD[i * 4 + j] = 0.0000;
182 }
183 // output a matrix to screen
184 std::cout << d[j * 4 + i] << "\t";
185
186 // Time to find the proper channels and build the corrections xind/i is row number
187 TChannel* chan = TChannel::FindChannelByName(Form("GRG%02d%sN00A", det, TGriffin::GetColorFromNumber(j)));
188 if(chan == nullptr) {
189 std::cout << DRED << "Couldn't find a channel for "
190 << Form("GRG%02d%sN00A", det, TGriffin::GetColorFromNumber(j)) << RESET_COLOR << std::endl;
191 continue;
192 }
193 // Writes the coefficient to the found channel above
194 chan->AddCTCoefficient(d[i * 4 + j]);
195 }
196 std::cout << std::endl;
197 }
198 // set precision back to old value and remove fixed flag
199 std::cout.precision(prec);
200 std::cout.unsetf(std::ios::floatfield);
201
202 std::cout << std::endl;
203 std::cout << "Largest correction = " << largestCorrection << " Shift = " << largestCorrection * energy << std::endl;
204 std::cout << "Largest combo, det = " << largestDet << " Crystals = " << largestCrystal1 << ", " << largestCrystal2 << std::endl;
205 std::cout << " -------------------- " << std::endl;
206 std::cout << " -------------------- " << std::endl;
207 return d;
208}
209
210void FixAll(TFile* inputFile, double energy)
211{
212 // This function loops over 16 clovers (always does) and by default uses the 1332 keV gamma ray in 60Co
213 for(int d = 1; d <= 16; d++) {
214 std::cout << "Starting CrossTalkFix for detector " << d << ", using energy " << energy << std::endl;
215 CrossTalkFix(d, energy, inputFile);
216 }
217}
218
219#ifndef __CINT__
220int main(int argc, char** argv)
221{
222 // Do basic file checks
223 if(argc != 3 && argc != 4) {
224 std::cout << "Usage: " << argv[0] << " <matrix file> <cal file> <user settings (optional)>" << std::endl;
225 return 1;
226 }
227
228 // We need a cal file to find the channels to write the corrections to
229 if(TChannel::ReadCalFile(argv[2]) < 0) {
230 std::cout << "Aborting" << std::endl;
231 exit(1);
232 }
233
234 auto* inputFile = new TFile(argv[1]);
235 if(inputFile == nullptr || !inputFile->IsOpen()) {
236 std::cout << "Failed to open file '" << argv[1] << "'!" << std::endl;
237 return 1;
238 }
239
240 auto* userSettings = new TUserSettings();
241 if(argc == 4) {
242 userSettings->ReadSettings(argv[3]);
243 }
244
245 // Create the output file
246 std::string outputFileName = argv[1];
247 auto lastSlash = outputFileName.find_last_of('/');
248 if(lastSlash != std::string::npos) {
249 outputFileName = std::string("ct_") + outputFileName.substr(lastSlash + 1);
250 } else {
251 outputFileName = std::string("ct_") + std::string(argv[1]);
252 }
253
254 auto* outputFile = new TFile(outputFileName.c_str(), "recreate");
255
256 FixAll(inputFile, userSettings->GetDouble("CrossTalkEnergy", 1332.));
257
258 // This function writes a corrections cal_file which can be loaded in with your normal cal file.
259 TChannel::WriteCTCorrections("ct_correction.cal");
260
261 outputFile->Write();
262 outputFile->Close();
263}
264#endif
#define DRED
Definition Globals.h:18
#define RESET_COLOR
Definition Globals.h:5
int main(int argc, char **argv)
double * CrossTalkFix(int det, double energy, TFile *inputFile)
void FixAll(TFile *inputFile, double energy)
double CrossTalkFit(double *x, double *par)
TH2D * mat
Definition UserFillObj.h:12
void AddCTCoefficient(double temp)
Definition TChannel.h:214
static TChannel * FindChannelByName(const char *ccName)
Definition TChannel.cxx:494
static Int_t ReadCalFile(std::ifstream &infile)
static void WriteCTCorrections(const std::string &outfilename="")
void DestroyCTCal()
Definition TChannel.cxx:587
static const char * GetColorFromNumber(int number)
Definition TGriffin.cxx:664