GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TGainMatch.cxx
Go to the documentation of this file.
1#include "TGainMatch.h"
2
3#include <algorithm>
4#include <map>
5
6#include "TError.h"
7#include "TRandom2.h"
8#include "Math/Minimizer.h"
9#include "Math/Factory.h"
10#include "Math/Functor.h"
11
13
15 : TCal(copy), fCoarseRange(fDefaultCoarseRange)
16{
17 std::cout << "THIS ---> " << fCoarseRange << std::endl;
18 copy.Copy(*this);
19}
20
21void TGainMatch::Copy(TObject& obj) const
22{
23 static_cast<TGainMatch&>(obj).fCoarseMatch = fCoarseMatch;
24 static_cast<TGainMatch&>(obj).fCoarseRange = fCoarseRange;
25 TCal::Copy(obj);
26}
27
28void TGainMatch::CalculateGain(Double_t cent1, Double_t cent2, Double_t eng1, Double_t eng2)
29{
30 // Put the peaks in order for ease (if the user put them in the wrong order)
31 // Apparantly there is a TGraph Sort method. Might look into this later.
32 if(eng1 > eng2) {
33 std::swap(eng1, eng2);
34 std::swap(cent1, cent2);
35 }
36
37 SetPoint(0, cent1, eng1);
38 SetPoint(1, cent2, eng2);
39
40 auto* gainFit = new TF1("gain", "pol1");
41
42 TFitResultPtr res = Fit(gainFit, "SC0");
43 SetFitFunction(GetFunction("gain")); // Have to do this because I want to delete gainfit
44 fGainCoeffs[0] = res->Parameter(0);
45 fGainCoeffs[1] = res->Parameter(1);
46
47 delete gainFit;
48}
49
50Bool_t TGainMatch::CoarseMatch(TH1* hist, Int_t chanNum, Double_t energy1, Double_t energy2)
51{
52 // This functions is used to perform a rough gain matching on a 60Co
53 // source by default. This makes gain matching over a wide range much easier to do afterwards
54 // This might have to be changed slightly if someone wants a different type of gaussian fitted
55 // for gain matching purposes.
56 if(hist == nullptr) {
57 return false;
58 }
59
60 fHist = hist;
61
62 // I might want to do a clear of the gainmatching parameters at this point.
63
64 // Check to see that the histogram isn't empty
65 if(hist->GetEntries() < 1) {
66 Error("CoarseMatch", "The histogram is empty");
67 return false;
68 }
69
70 // See if the channel exists. There is no point in finding the gains if we don't have anywhere to write it
72 if(chan == nullptr) {
73 if(chanNum != 9999) {
74 Warning("CoarseMatch", "Channel Number %d does not exist in current memory.", chanNum);
75 }
76 } else {
77 // Set the channel number
78 SetChannel(chan);
79 }
80
81 std::vector<Double_t> engVec; // This is the vector of expected energies
82 // We do coarse gain matching on 60Co. So I have hardcoded the energies for now
83 engVec.push_back(energy1);
84 engVec.push_back(energy2);
85
86 // Sort these in case the peak is returned in the wrong order
87 std::sort(engVec.begin(), engVec.end());
88
89 // Change the histogram range so that the search ignores potential low energy noise.
90 Int_t high = hist->GetXaxis()->GetLast();
91 hist->GetXaxis()->SetRangeUser(100, hist->GetXaxis()->GetBinCenter(high));
92 // Use a TSpectrum to find the two largest peaks in the spectrum
93 auto* spec = new TSpectrum; // This might not have to be allocated
94 Int_t nFound = spec->Search(hist, 2, "", 0.50); // This returns peaks in order of their height in the spectrum.
95
96 // return histogram to proper range
97 hist->GetXaxis()->UnZoom();
98 // If we didn't find two peaks, it is likely we gave it garbage
99 if(nFound < 2) {
100 Error("CoarseMatch", "Did not find enough peaks");
101 delete spec;
102 return false;
103 }
104
105 // We want to store the centroids of the found peaks
106 std::vector<Double_t> foundBin;
107 for(int x = 0; x < 2;
108 x++) { // I have hard-coded this to 2 because I'm assuming the rough match peaks will be by far the largest.
109 foundBin.push_back(spec->GetPositionX()[x]);
110 std::cout << "Found peak at bin " << foundBin[x] << std::endl;
111 }
112 std::sort(foundBin.begin(), foundBin.end());
113
114 // Get Bin Width for use later. I am using the first peak found to set the bin width
115 // If you are using antisymmetric binning... god help you.
116 // Double_t binWidth = hist->GetBinWidth(foundBin[0]);
117
118 // Set the number of data points to 2. In the gain matching graph.
119 Set(2);
120
121 // We now want to create a peak for each one we found (2) and fit them.
122 for(int x = 0; x < 2; x++) {
123 TPeak tmpPeak(foundBin[x], foundBin[x] - fCoarseRange, foundBin[x] + fCoarseRange);
124 tmpPeak.SetName(Form("GM_Cent_%lf", foundBin[x])); // Change the name of the TPeak to know it's origin
125 tmpPeak.SetLineColor(static_cast<Color_t>(2 * x + 2));
126 tmpPeak.Fit(hist, "Q+");
127 tmpPeak.ReleaseParameter(3);
128 tmpPeak.ReleaseParameter(4);
129 tmpPeak.SetLineColor(static_cast<Color_t>(2 * x + 3));
130 tmpPeak.Fit(hist, "Q+");
131 SetPoint(x, tmpPeak.GetParameter("centroid"), engVec[x]);
132 }
133
134 auto* gainFit = new TF1("gain", "pol1");
135
136 TFitResultPtr res = Fit(gainFit, "SC0");
137 SetFitFunction(GetFunction("gain")); // Have to do this because I want to delete gainFit.
138 fGainCoeffs[0] = res->Parameter(0);
139 fGainCoeffs[1] = res->Parameter(1);
140
141 delete gainFit;
142 // We have finished gain matching so let the TGainMatch know that it is a coarse gain
143 fCoarseMatch = true;
144 delete spec;
145 return true;
146}
147
148Bool_t TGainMatch::FineMatchFast(TH1* hist1, TPeak* peak1, TH1* hist2, TPeak* peak2, Int_t channelNum)
149{
150 // You need to pass a TPeak with the centroid and range set to the real energy centroid and ranges.
151 // The function uses the coarse gain parameters to find the peak and gain matches the raw spectrum.
152 // This is more useful as it allows a script to find all of the peaks.
153 if((hist1 == nullptr) || (hist2 == nullptr)) {
154 Error("FineMatchFast", "No histogram being pointed to");
155 return false;
156 }
157
158 // Check to see that the histogram isn't empty
159 if(hist1->GetEntries() < 1 || hist2->GetEntries() < 1) {
160 Error("FineMatchFast", "Histogram is empty");
161 return false;
162 }
163
164 // See if the channel exists. There is no point in finding the gains if we don't have anywhere to write it
165 Double_t gain = 0.;
166 Double_t offset = 0.;
167 TChannel* chan = TChannel::GetChannelByNumber(channelNum);
168 if(chan == nullptr) {
169 if(channelNum != 9999) {
170 Warning("FineMatchFast", "Channel Number %d does not exist in current memory.", channelNum);
171 }
172 if(GetFitFunction() != nullptr) {
173 gain = GetParameter(1);
174 offset = GetParameter(0);
175 } else {
176 Error("Fine Match", "There needs to be a coarse gain set to do a fine gain");
177 return false;
178 }
179 } else {
180 SetChannel(chan);
181 // The first thing we need to do is "un-gain correct" the centroids of the TPeaks.
182 // What we are actually doing is applying the recirpocal gain to the energy of the real peak
183 // to figure out where the centroid of the charge is spectrum is roughly going to be
184 // First read in the rough gain coefficients
185 std::vector<Float_t> roughCoeffs = chan->GetENGCoeff();
186 gain = roughCoeffs.at(1);
187 offset = roughCoeffs.at(0);
188 }
189
190 // The reason I'm using TPeak here is we might want to gain match "TPhotopeaks", or "TElectronPeaks", or
191 //"TCrazyNonGaussian" Peak. All we care about is that it has a centroid.
192 if((peak1 == nullptr) || (peak2 == nullptr)) {
193 Error("FineMatchFast", "No TPeak being pointed to");
194 return false;
195 }
196
197 // Set the channel number
198 Set(2);
199
200 // Find the energy of the peak that we want to use
201 std::array<Double_t, 2> energy = {peak1->GetParameter("centroid"), peak2->GetParameter("centroid")};
202 std::cout << peak1->GetParameter("centroid") << " ENERGIES " << energy[1] << std::endl;
203 // Offsets are very small right now so I'm not including them until they become a problem.
204 peak1->SetParameter("centroid", (energy[0] - offset) / gain);
205 peak2->SetParameter("centroid", (energy[1] - offset) / gain);
206 // Change the range for the fit to be in the gain corrected spectrum
207
208 peak1->SetRange((peak1->GetXmin() - offset) / gain, (peak1->GetXmax() - offset) / gain);
209 peak2->SetRange((peak2->GetXmin() - offset) / gain, (peak2->GetXmax() - offset) / gain);
210
211 // The gains won't be perfect, so we need to search for the peak within a range.
212 hist1->GetXaxis()->SetRangeUser(peak1->GetXmin() - 20., peak1->GetXmax() + 20.);
213 TSpectrum spec;
214 Int_t nFound = spec.Search(hist1, 2, "", 0.3);
215
216 for(int x = 0; x < nFound; x++) {
217 std::cout << spec.GetPositionX()[x] << std::endl;
218 }
219 Double_t closestPeak = 0;
220 Double_t closestDiff = 10000;
221 for(int x = 0; x < nFound; x++) {
222 if(fabs(peak1->GetCentroid() - spec.GetPositionX()[x]) < closestDiff) {
223 closestPeak = spec.GetPositionX()[x];
224 closestDiff = fabs(peak1->GetCentroid() - spec.GetPositionX()[x]);
225 }
226 }
227
228 Double_t rangeWidth = (peak1->GetXmax() - peak1->GetXmin()) / 2.;
229 peak1->SetParameter("centroid", closestPeak);
230 peak1->SetRange(peak1->GetCentroid() - rangeWidth, peak1->GetCentroid() + rangeWidth);
231 std::cout << "Centroid Guess " << peak1->GetCentroid() << std::endl;
232 std::cout << "Range Low " << peak1->GetXmin() << " " << peak1->GetXmax() << std::endl;
233
234 closestPeak = 0;
235 closestDiff = 10000;
236 hist2->GetXaxis()->SetRangeUser(peak2->GetXmin() - 20., peak2->GetXmax() + 20.);
237 TSpectrum spec2;
238 nFound = spec2.Search(hist2, 2, "", 0.3); // Search the next histogram
239 for(int x = 0; x < nFound; x++) {
240 std::cout << spec2.GetPositionX()[x] << std::endl;
241 }
242
243 for(int x = 0; x < nFound; x++) {
244 if(fabs(peak2->GetCentroid() - spec2.GetPositionX()[x]) < closestDiff) {
245 closestPeak = spec2.GetPositionX()[x];
246 closestDiff = fabs(peak2->GetCentroid() - spec2.GetPositionX()[x]);
247 }
248 }
249 Double_t rangeWidth2 = (peak2->GetXmax() - peak2->GetXmin()) / 2.;
250 peak2->SetParameter("centroid", closestPeak);
251 peak2->SetRange(peak2->GetCentroid() - rangeWidth2, peak2->GetCentroid() + rangeWidth2);
252
253 std::cout << "Centroid Guess " << peak2->GetCentroid() << std::endl;
254 std::cout << "Range High " << peak2->GetXmin() << " " << peak2->GetXmax() << std::endl;
255
256 hist1->GetXaxis()->UnZoom();
257 hist2->GetXaxis()->UnZoom();
258 peak1->Fit(hist1, "MS+");
259 peak2->Fit(hist2, "MS+");
260
261 hist1->Draw();
262 peak1->Draw("same");
263 peak2->Draw("same");
264
265 std::array<Double_t, 2> centroid = {peak1->GetCentroid(), peak2->GetCentroid()};
266
267 // Put the peaks in order for ease (if the user put them in the wrong order)
268 // Apparantly there is a TGraph Sort method. Might look into this later.
269 if(energy[0] > energy[1]) {
270 std::swap(energy[0], energy[1]);
271 std::swap(centroid[0], centroid[1]);
272 }
273
274 SetPoint(0, centroid[0], energy[0]);
275 SetPoint(1, centroid[1], energy[1]);
276
277 auto* gainFit = new TF1("gain", "pol1");
278
279 TFitResultPtr res = Fit(gainFit, "SC0");
280 SetFitFunction(GetFunction("gain")); // Have to do this because I want to delete gainFit
281 fGainCoeffs[0] = res->Parameter(0);
282 fGainCoeffs[1] = res->Parameter(1);
283
284 delete gainFit;
285
286 fCoarseMatch = false;
287 return true;
288}
289
290Bool_t TGainMatch::FineMatchFast(TH1* hist, TPeak* peak1, TPeak*, Int_t channelNum)
291{
292 // You need to pass a TPeak with the centroid and range set to the real energy centroid and ranges.
293 // The function uses the coarse gain parameters to find the peak and gain matches the raw spectrum.
294 // This is more useful as it allows a script to find all of the peaks.
295 return FineMatchFast(hist, peak1, hist, peak1, channelNum);
296}
297
298Bool_t TGainMatch::FineMatchFast(TH1* hist, Double_t energy1, Double_t energy2, Int_t channelNum)
299{
300 // Performs fine gain matching on one histogram once we have set the rough energy
301 // coefficients.
302 return FineMatchFast(hist, energy1, hist, energy2, channelNum);
303}
304
305Bool_t TGainMatch::FineMatchFast(TH1* hist1, Double_t energy1, TH1* hist2, Double_t energy2, Int_t channelNum)
306{
307 // Performs fine gain matching on two histograms once we have set the rough energy
308 // coefficients. You use this if you have a two different sources giving your full range
309 // of energy. This histograms should be binned the same. NOT IMPLEMENTED
310
311 // using an automatic range of 10 keV for testing purposes.
312 auto* peak1 = new TPeak(energy1, energy1 - 10.0, energy1 + 10.0);
313 auto* peak2 = new TPeak(energy2, energy2 - 10.0, energy2 + 10.0);
314
315 Bool_t result = FineMatchFast(hist1, peak1, hist2, peak2, channelNum);
316
317 delete peak1;
318 delete peak2;
319 return result;
320 return false;
321}
322
324{
325 if(GetChannel() == nullptr) {
326 Error("WriteToChannel", "No Channel Set");
327 return;
328 }
330 std::cout << "Writing to channel " << GetChannel()->GetNumber() << std::endl;
331 std::cout << "p0 = " << GetParameter(0) << " \t p1 = " << GetParameter(1) << std::endl;
332 // Set the energy parameters based on the fitted gains.
333 GetChannel()->AddENGCoefficient(static_cast<float>(GetParameter(0)));
334 GetChannel()->AddENGCoefficient(static_cast<float>(GetParameter(1)));
335}
336
337void TGainMatch::Print(Option_t*) const
338{
339 std::cout << std::endl;
340 std::cout << "GainMatching: " << std::endl;
341 if(fCoarseMatch) {
342 std::cout << "COARSE" << std::endl;
343 } else {
344 std::cout << "FINE" << std::endl;
345 }
346
347 if(fAligned) {
348 std::cout << "Aligned" << std::endl;
349 } else {
350 std::cout << "NOT Aligned" << std::endl;
351 }
352
353 std::cout << "Gain Coefficients: " << fGainCoeffs[0] << "\t" << fGainCoeffs[1] << std::endl;
354 std::cout << "Align Coefficients: " << fAlignCoeffs[0] << "\t" << fAlignCoeffs[1] << std::endl;
355
356 TCal::Print();
357}
358
359Bool_t TGainMatch::CoarseMatchAll(TCalManager* cm, TH2* mat, Double_t, Double_t)
360{
361 // If you supply this function with a matrix of Channel vs. energy it will automatically slice,
362 // and figure out the coarse gains. I might add a channel range option later
363 std::vector<Int_t> badlist;
364 auto* gm = new TGainMatch;
365 if(cm == nullptr) {
366 gm->Error("CoarseMatchAll", "CalManager Pointer is nullptr");
367 return false;
368 }
369 if(mat == nullptr) {
370 gm->Error("CoarseMatchAll", "TH2 Pointer is nullptr");
371 return false;
372 }
373 // Find the range of channels provided
374 Int_t first_chan = mat->GetXaxis()->GetFirst();
375 Int_t last_chan = mat->GetXaxis()->GetLast();
376 // The first thing we need to do is slice the matrix into it's channel vs energy.
377 for(int chan = first_chan; chan <= last_chan; chan++) {
378 gm->Clear();
379 std::cout << std::endl
380 << "Now fitting channel " << chan - 1 << std::endl;
381 auto* h1 = static_cast<TH1D*>(mat->ProjectionY(Form("Channel%d", chan), chan, chan, "o"));
382 std::cout << "BIN WIDTH " << h1->GetXaxis()->GetBinWidth(h1->GetXaxis()->GetFirst() + 1) << std::endl;
383 if(h1->Integral() < 100) {
384 continue;
385 }
386
387 if(!(gm->CoarseMatch(h1, chan - 1))) {
388 badlist.push_back(chan - 1);
389 continue;
390 }
391 gm->SetName(Form("gm_chan_%d", chan - 1));
392 cm->AddToManager(gm, chan - 1);
393 }
394 if(!badlist.empty()) {
395 std::cout << "The following channels did not gain match properly: ";
396 }
397 for(int i : badlist) {
398 std::cout << i << "\t";
399 }
400
401 delete gm;
402
403 return true;
404}
405
406Bool_t TGainMatch::FineMatchFastAll(TCalManager* cm, TH2* mat1, Double_t energy1, TH2* mat2, Double_t energy2)
407{
408
409 // using an automatic range of 10 keV for testing purposes.
410 // We need to include bin width likely?
411 auto* peak1 = new TPeak(energy1, energy1 - 10.0, energy1 + 10.0);
412 auto* peak2 = new TPeak(energy2, energy2 - 10.0, energy2 + 10.0);
413
414 Bool_t result = TGainMatch::FineMatchFastAll(cm, mat1, peak1, mat2, peak2);
415
416 delete peak1;
417 delete peak2;
418 return result;
419}
420
421Bool_t TGainMatch::FineMatchFastAll(TCalManager* cm, TH2* mat, Double_t energy1, Double_t energy2)
422{
423 return TGainMatch::FineMatchFastAll(cm, mat, energy1, mat, energy2);
424}
425
426Bool_t TGainMatch::FineMatchFastAll(TCalManager* cm, TH2* mat1, TPeak* peak1, TH2* mat2, TPeak* peak2)
427{
428 // If you supply this function with a matrix of Channel vs. energy it will automatically slice,
429 // and figure out the fine gains. I might add a channel range option later
430 std::vector<Int_t> badlist;
431 auto* gm = new TGainMatch;
432 if(cm == nullptr) {
433 gm->Error("FineMatchFastAll", "CalManager Pointer is nullptr");
434 return false;
435 }
436 if((mat1 == nullptr) || (mat2 == nullptr)) {
437 gm->Error("FineMatchFastAll", "TH2 Pointer is nullptr");
438 return false;
439 }
440 if((peak1 == nullptr) || (peak2 == nullptr)) {
441 gm->Error("FineMatchFastAll", "No TPeak being pointed to");
442 return false;
443 }
444
445 // Find the range of channels provided, we want to use the largest range provided
446 Int_t first_chan1 = mat1->GetXaxis()->GetFirst();
447 Int_t last_chan1 = mat1->GetXaxis()->GetLast();
448 Int_t first_chan2 = mat2->GetXaxis()->GetFirst();
449 Int_t last_chan2 = mat2->GetXaxis()->GetLast();
450 Int_t first_chan = std::min(first_chan1, first_chan2);
451 Int_t last_chan = std::max(last_chan1, last_chan2);
452 // The first thing we need to do is slice the matrix into it's channel vs energy.
453 for(int chan = first_chan; chan <= last_chan; chan++) {
454 gm->Clear();
455 // Make a copy of the TPeaks so that we can fit each of them
456 // TPeak* copyPeak1 = static_cast<TPeak*>(peak1->Clone());//Clone creates smart pointers so these will be
457 // dealoccated automatically.
458 // TPeak* copyPeak2 = static_cast<TPeak*>(peak2->Clone());
459 auto* copyPeak1 = new TPeak(*peak1);
460 auto* copyPeak2 = new TPeak(*peak2);
461
462 std::cout << std::endl
463 << "Now fitting channel: " << chan - 1 << std::endl;
464 auto* h1 = static_cast<TH1D*>(mat1->ProjectionY(Form("Channel%d_mat1", chan - 1), chan, chan, "o"));
465 auto* h2 = static_cast<TH1D*>(mat2->ProjectionY(Form("Channel%d_mat2", chan - 1), chan, chan, "o"));
466 if(h1->Integral() < 100 || h2->Integral() < 100) {
467 gm->Warning("FineMatchFastAll", "Empty channel = %d", chan - 1);
468 continue;
469 }
470 if(!(gm->FineMatchFast(h1, copyPeak1, h2, copyPeak2, chan - 1))) {
471 badlist.push_back(chan - 1);
472 continue;
473 }
474 cm->AddToManager(gm);
475
476 delete copyPeak1;
477 delete copyPeak2;
478 }
479 if(!badlist.empty()) {
480 std::cout << "The following channels did not gain match properly: ";
481 }
482 for(int i : badlist) {
483 std::cout << i << "\t";
484 }
485
486 delete gm;
487
488 return true;
489}
490
491Bool_t TGainMatch::FineMatchFastAll(TCalManager* cm, TH2* mat, TPeak* peak1, TPeak* peak2)
492{
493 return TGainMatch::FineMatchFastAll(cm, mat, peak1, mat, peak2);
494}
495
496void TGainMatch::Clear(Option_t*)
497{
498 fCoarseMatch = true;
499 fAligned = false;
500 fAlignCoeffs[0] = 0.0;
501 fAlignCoeffs[1] = 1.0;
502 fGainCoeffs[0] = 0.0;
503 fGainCoeffs[1] = 1.0;
504 TCal::Clear();
505}
506
507Bool_t TGainMatch::Align(TH1* test, TH1* hist, Int_t low_range, Int_t high_range)
508{
509 // Minimizes the chi^2 between the bin contents of the test histogram and the bin contents of the histogram to be
510 // matched'
511 if(!((test != nullptr) && (hist != nullptr))) {
512 std::cout << "Unassigned histogram" << std::endl;
513 return false;
514 }
515 fHist = hist; // Need this histogram to be seen by ftotal...Don't know how else to do this right now.
516 // TF1* tmpfunc = new
517 // TF1("tmpfunc",this,&TGainMatch::HistCompare,test->GetMinimumStored()+1,test->GetMaximumStored()-1,3);
518
519 /* //Initial Guesses
520
521 TSpectrum stest,shist;
522 TGraph guess_graph;
523 Int_t nFoundtest = stest.Search(test);
524 Int_t nFoundhist = shist.Search(hist);
525 TF1* shift_guess = new TF1("shift_guess","pol1");
526 if(nFoundtest >1 && nFoundhist > 1) {
527 stest.GetPositionX()[0]/shist.GetPositionX()[0];
528 guess_graph.SetPoint(0,shist.GetPositionX()[0],stest.GetPositionX()[0]);
529 guess_graph.SetPoint(1,shist.GetPositionX()[1],stest.GetPositionX()[1]);
530
531 guess_graph.Print();
532
533 TFitResultPtr res = guess_graph.Fit(shift_guess,"SC0");
534 }
535 */
536
537 TF1* tmpfunc = new TF1("tmpfunc", this, &TGainMatch::HistCompare, low_range, high_range, 3);
538 tmpfunc->SetNpx(10000);
539 tmpfunc->SetParameters(1.0, 1.0, 1.0);
540
541 // hist->Sumw2();
542
543 /* for(int i =0; i<hist->GetXaxis()->GetNbins();i++) {
544 // hist->SetBinError(i,TMath::Sqrt(TMath::Sqrt(hist->GetBinContent(i))));
545 hist->SetBinError(i,1./(1.+hist->GetBinContent(i)));
546 }
547 */
548 // ftot->SetParLimits(0,.3*norm,norm);
549
550 // const char* minName = "Minuit2";
551 // const char* algoName = "Scan";
552 // ROOT::Math::MinimizerOptions::SetDefaultMinimizer("TMinuit2","Simplex");
553 // TVirtualFitter::SetDefaultFitter("Fumili"); //USE COMBINATION!!!!
554 ROOT::Math::MinimizerOptions::SetDefaultMinimizer("Minuit2", "Combination");
555 TVirtualFitter::SetPrecision(1.0e-10);
556 TVirtualFitter::SetMaxIterations(10000);
557 TFitResultPtr res = test->Fit("tmpfunc", "RSILV");
558 fAlignCoeffs[0] = res->Parameter(1);
559 fAlignCoeffs[1] = res->Parameter(2);
560 std::cout << "Chi2: " << res->Chi2() / res->Ndf() << std::endl;
561 if(std::abs(res->Parameter(0) - 1.00) > 20.) {
562 return false;
563 }
564 if(std::abs(res->Parameter(2) - 1.00) > 5.) {
565 return false;
566 }
567 fAligned = true;
568 // delete shift_guess;
569 return true;
570}
571
572Bool_t TGainMatch::AlignAll(TCalManager* cm, TH1* hist, TH2* mat, Int_t low_range, Int_t high_range)
573{
574 std::vector<Int_t> badlist;
575 auto* gm = new TGainMatch;
576 if(cm == nullptr) {
577 gm->Error("AlignAll", "CalManager Pointer is nullptr");
578 return false;
579 }
580 if(mat == nullptr) {
581 gm->Error("AlignAll", "TH2 Pointer is nullptr");
582 return false;
583 }
584 if(hist == nullptr) {
585 gm->Error("AlignAll", "TH1 Pointer is nullptr");
586 return false;
587 }
588 // Find the range of channels provided
589 Int_t first_chan = mat->GetXaxis()->GetFirst();
590 // Int_t last_chan = mat->GetXaxis()->GetLast();
591 // The first thing we need to do is slice the matrix into it's channel vs energy.
592 // for(int chan=first_chan; chan<=last_chan;chan++) {
593 for(int chan = first_chan; chan <= 2; chan++) {
594 gm->Clear();
595 std::cout << std::endl
596 << "Now fitting channel: " << chan << std::endl;
597 auto* h1 = static_cast<TH1D*>(mat->ProjectionY(Form("Channel%d", chan), chan + 1, chan + 1, "o"));
598 std::cout << "BIN WIDTH " << h1->GetXaxis()->GetBinWidth(h1->GetXaxis()->GetFirst() + 1) << std::endl;
599 if(h1->Integral() < 100) {
600 continue;
601 }
602
603 if(!(gm->Align(hist, h1, low_range, high_range))) {
604 badlist.push_back(chan);
605 continue;
606 }
607 cm->AddToManager(gm);
608 }
609 if(!badlist.empty()) {
610 std::cout << "The following channels did not gain match properly: ";
611 }
612 for(int i : badlist) {
613 std::cout << i << "\t";
614 }
615
616 delete gm;
617
618 return true;
619}
620
621Bool_t TGainMatch::FineMatchAll(TCalManager* cm, TH2* charge_mat, TH2* eng_mat, Int_t testchan, Double_t energy1,
622 Double_t energy2, Int_t low_range, Int_t high_range)
623{
624 // Should be run to correct drifting gains. Requires a previously gain matched spectrum, and ideally a roughly energy
625 // calibrated matrix
626
627 auto* gm = new TGainMatch;
628 if(cm == nullptr) {
629 gm->Error("FineMatchAll", "CalManager Pointer is nullptr");
630 return false;
631 }
632 if((charge_mat == nullptr) || (eng_mat == nullptr)) {
633 gm->Error("FineMatchAll", "TH2 Pointer is nullptr");
634 return false;
635 }
636 auto binwidth = static_cast<Int_t>(0.5 + 1. / eng_mat->GetYaxis()->GetBinWidth(100));
637 eng_mat->RebinY(binwidth);
638
639 std::vector<Int_t> badlist;
640 // Find the range of channels provided, we want to use the largest range provided
641 Int_t first_chan = eng_mat->GetXaxis()->GetFirst();
642 Int_t last_chan = eng_mat->GetXaxis()->GetLast();
643 // The first thing we need to do is slice the matrix into it's channel vs energy.
644 auto* testhist = static_cast<TH1D*>(eng_mat->ProjectionY(Form("Test%d_mat", testchan), testchan + 1, testchan + 1, "o"));
645
646 for(int chan = first_chan; chan <= last_chan; chan++) {
647 std::cout << std::endl
648 << "Now fitting channel: " << chan - 1 << std::endl;
649 auto* chargeh = static_cast<TH1D*>(charge_mat->ProjectionY(Form("Charge%d_mat", chan - 1), chan, chan, "o"));
650 auto* engh = static_cast<TH1D*>(eng_mat->ProjectionY(Form("Energy%d_mat", chan - 1), chan, chan, "o"));
651 if(chargeh->Integral() < 100 || chargeh->Integral() < 100) {
652 gm->Warning("FineMatchAll", "Empty channel = %d", chan - 1);
653 continue;
654 }
655 // Needs to align before matching
656 if(!(gm->FineMatch(engh, testhist, chargeh, energy1, energy2, low_range, high_range, chan - 1))) {
657 badlist.push_back(chan - 1);
658 continue;
659 }
660 cm->AddToManager(gm);
661 }
662 if(!badlist.empty()) {
663 std::cout << "The following channels did not gain match properly: ";
664 }
665 for(int i : badlist) {
666 std::cout << i << "\t";
667 }
668
669 delete testhist;
670 delete gm;
671
672 return true;
673}
674
675Bool_t TGainMatch::FineMatch(TH1* energyHist, TH1* testhist, TH1* chargeHist, Double_t energy1, Double_t energy2,
676 Int_t low_range, Int_t high_range, Int_t channelNum)
677{
678 // Should be run to correct drifting gains. Requires a previously gain matched spectrum, and ideally a roughly energy
679 // calibrated spectrum
680 if(!Align(testhist, energyHist, low_range, high_range)) {
681 return false;
682 }
683 TH1* hist2 = chargeHist; // Cheating for easier modification later
684 if((chargeHist == nullptr) || (testhist == nullptr) || (energyHist == nullptr)) {
685 Error("FineMatch", "No histogram being pointed to");
686 return false;
687 }
688
689 // Check to see that the histogram isn't empty
690 if(chargeHist->GetEntries() < 1 || testhist->GetEntries() < 1 || energyHist->GetEntries() < 1) {
691 Error("FineMatchFast", "Histogram is empty");
692 return false;
693 }
694
695 // See if the channel exists. There is no point in finding the gains if we don't have anywhere to write it
696 Double_t gain = 0.;
697 Double_t offset = 0.;
698 TChannel* chan = TChannel::GetChannelByNumber(channelNum);
699 if(chan == nullptr) {
700 if(channelNum != 9999) {
701 Warning("FineMatch", "Channel Number %d does not exist in current memory.", channelNum);
702 }
703 if(GetFitFunction() != nullptr) {
704 gain = GetParameter(1);
705 offset = GetParameter(0);
706 } else {
707 Error("FineMatch", "There needs to be a coarse gain set to do a fine gain");
708 return false;
709 }
710 } else {
711 SetChannel(chan);
712 // The first thing we need to do is "un-gain correct" the centroids of the TPeaks.
713 // What we are actually doing is applying the recirpocal gain to the energy of the real peak
714 // to figure out where the centroid of the charge is spectrum is roughly going to be
715 // First read in the rough gain coefficients
716 std::vector<Float_t> roughCoeffs = chan->GetENGCoeff();
717 gain = roughCoeffs.at(1);
718 offset = roughCoeffs.at(0);
719 }
720
721 auto* peak1 = new TPeak(energy1, energy1 - 15.0, energy1 + 15.0);
722 auto* peak2 = new TPeak(energy2, energy2 - 15.0, energy2 + 15.0);
723
724 // Set the channel number
725 // Graph()->Set(2);
726 Set(2);
727
728 // Find the energy of the peak that we want to use
729 std::array<Double_t, 2> energy = {peak1->GetParameter("centroid"), peak2->GetParameter("centroid")};
730 std::cout << peak1->GetParameter("centroid") << " ENERGIES " << energy[1] << std::endl;
731 // Offsets are very small right now so I'm not including them until they become a problem.
732 // peak1->SetParameter("centroid",((energy[0]-offset)/gain)*fAlignCoeffs[1] + fAlignCoeffs[0]);
733 // peak2->SetParameter("centroid",((energy[0]-offset)/gain)*fAlignCoeffs[1] + fAlignCoeffs[0]);
734 peak1->SetParameter("centroid", (energy[0] * fAlignCoeffs[1] + fAlignCoeffs[0] - offset) / gain);
735 peak2->SetParameter("centroid", (energy[1] * fAlignCoeffs[1] + fAlignCoeffs[0] - offset) / gain);
736 // Change the range for the fit to be in the gain corrected spectrum
737
738 std::cout << " Should be: " << ((energy[0] * fAlignCoeffs[1] + fAlignCoeffs[0] - offset) / gain) << " from (" << energy[0] << " * " << fAlignCoeffs[1] << " + " << fAlignCoeffs[0] << " - " << offset << ")/" << gain << std::endl;
739 std::cout << " Should be: " << ((energy[1] * fAlignCoeffs[1] + fAlignCoeffs[0] - offset) / gain) << " from (" << energy[1] << " * " << fAlignCoeffs[1] << " + " << fAlignCoeffs[0] << " - " << offset << ")/" << gain << std::endl;
740
741 // peak1->SetRange(((peak1->GetXmin()-offset)/gain)*fAlignCoeffs[1] +
742 // fAlignCoeffs[0],((peak1->GetXmax()-offset)/gain)*fAlignCoeffs[1] + fAlignCoeffs[0]);
743 // peak2->SetRange(((peak2->GetXmin()-offset)/gain)*fAlignCoeffs[1] +
744 // fAlignCoeffs[0],((peak2->GetXmax()-offset)/gain)*fAlignCoeffs[1] + fAlignCoeffs[0]);
745 peak1->SetRange(((peak1->GetXmin() * fAlignCoeffs[1] + fAlignCoeffs[0] - offset) / gain),
746 ((peak1->GetXmax() * fAlignCoeffs[1] + fAlignCoeffs[0] - offset) / gain));
747 peak2->SetRange(((peak2->GetXmin() * fAlignCoeffs[1] + fAlignCoeffs[0] - offset) / gain),
748 ((peak2->GetXmax() * fAlignCoeffs[1] + fAlignCoeffs[0] - offset) / gain));
749
750 // The gains won't be perfect, so we need to search for the peak within a range.
751 TSpectrum spec;
752 chargeHist->GetXaxis()->SetRangeUser(peak1->GetXmin() - 20., peak1->GetXmax() + 20.);
753 Int_t nFound = spec.Search(chargeHist, 2, "", 0.3);
754
755 for(int x = 0; x < nFound; x++) {
756 std::cout << spec.GetPositionX()[x] << std::endl;
757 }
758 Double_t closestPeak = 0;
759 Double_t closestDiff = 10000;
760 for(int x = 0; x < nFound; x++) {
761 if(fabs(peak1->GetCentroid() - spec.GetPositionX()[x]) < closestDiff) {
762 closestPeak = spec.GetPositionX()[x];
763 closestDiff = fabs(peak1->GetCentroid() - spec.GetPositionX()[x]);
764 }
765 }
766
767 Double_t rangeWidth = (peak1->GetXmax() - peak1->GetXmin()) / 2.;
768 peak1->SetParameter("centroid", closestPeak);
769 peak1->SetRange(peak1->GetCentroid() - rangeWidth, peak1->GetCentroid() + rangeWidth);
770 std::cout << "Centroid Guess " << peak1->GetCentroid() << std::endl;
771 std::cout << "Range Low " << peak1->GetXmin() << " " << peak1->GetXmax() << std::endl;
772
773 closestDiff = 10000;
774 TSpectrum spec2;
775 hist2->GetXaxis()->SetRangeUser(peak2->GetXmin() - 20., peak2->GetXmax() + 20.);
776 std::cout << "RANGE: " << peak2->GetXmin() - 20. << " " << peak2->GetXmax() + 20. << std::endl;
777 nFound = spec2.Search(hist2, 2, "", 0.3); // Search the next histogram
778 if(nFound == 0) {
779 std::cout << "RANGE: " << peak2->GetXmin() - 40. << " " << peak2->GetXmax() + 40. << std::endl;
780 nFound = spec2.Search(hist2, 2, "", 0.3); // Search the next histogram
781 }
782
783 for(int x = 0; x < nFound; x++) {
784 std::cout << spec2.GetPositionX()[x] << std::endl;
785 }
786
787 Double_t largestPeak = spec2.GetPositionX()[0];
788 for(int x = 0; x < nFound; x++) {
789 if(fabs(peak2->GetCentroid() - spec2.GetPositionX()[x]) < closestDiff) {
790 closestDiff = fabs(peak2->GetCentroid() - spec2.GetPositionX()[x]);
791 }
792 }
793 Double_t rangeWidth2 = (peak2->GetXmax() - peak2->GetXmin()) / 2.;
794 peak2->SetParameter("centroid", largestPeak);
795 peak2->SetRange(peak2->GetCentroid() - rangeWidth2, peak2->GetCentroid() + rangeWidth2);
796
797 std::cout << "Centroid Guess " << peak2->GetCentroid() << std::endl;
798 std::cout << "Range High " << peak2->GetXmin() << " " << peak2->GetXmax() << std::endl;
799
800 chargeHist->GetXaxis()->UnZoom();
801 hist2->GetXaxis()->UnZoom();
802 peak1->InitParams(chargeHist);
803 peak2->InitParams(chargeHist);
804 peak1->SetParameter("sigma", TMath::Sqrt(9.0 + 4. * peak1->GetParameter("centroid") / 1000.) / 2.35);
805 peak2->SetParameter("sigma", TMath::Sqrt(9.0 + 4. * peak2->GetParameter("centroid") / 1000.) / 2.35);
806 peak1->Fit(chargeHist, "MSL+");
807 peak2->Fit(hist2, "MSL+");
808
809 chargeHist->Draw();
810 peak1->Draw("same");
811 peak2->Draw("same");
812
813 std::array<Double_t, 2> centroid = {peak1->GetCentroid(), peak2->GetCentroid()};
814
815 // Put the peaks in order for ease (if the user put them in the wrong order)
816 // Apparantly there is a TGraph Sort method. Might look into this later.
817 if(energy[0] > energy[1]) {
818 std::swap(energy[0], energy[1]);
819 std::swap(centroid[0], centroid[1]);
820 }
821
822 SetPoint(0, centroid[0], energy[0]);
823 SetPoint(1, centroid[1], energy[1]);
824
825 auto* gainFit = new TF1("gain", "pol1");
826
827 TFitResultPtr res = Fit(gainFit, "SC0");
828 SetFitFunction(GetFunction("gain")); // Have to do this because I want to delete gainFit
829 fGainCoeffs[0] = res->Parameter(0);
830 fGainCoeffs[1] = res->Parameter(1);
831
832 delete gainFit;
833
834 fCoarseMatch = false;
835 delete peak1;
836 delete peak2;
837 return true;
838}
839
840Double_t TGainMatch::HistCompare(Double_t* x, Double_t* par) // NOLINT(readability-non-const-parameter)
841{
842 Int_t bin = fHist->GetXaxis()->FindBin(x[0] * par[2] + par[1]);
843 Double_t content = fHist->GetBinContent(bin);
844
845 return par[0] * content;
846}
TH2D * mat
Definition UserFillObj.h:12
TH1D * hist
Definition UserFillObj.h:3
Definition TCal.h:44
void Print(Option_t *opt="") const override
Definition TCal.cxx:145
virtual void SetFitFunction(TF1 *func)
Definition TCal.h:60
virtual Double_t GetParameter(size_t parameter) const
Definition TCal.cxx:99
void Clear(Option_t *opt="") override
Definition TCal.cxx:137
Bool_t SetChannel(TChannel *chan)
Definition TCal.cxx:51
void Copy(TObject &obj) const override
Definition TCal.cxx:37
virtual TF1 * GetFitFunction() const
Definition TCal.h:59
TChannel * GetChannel() const
Definition TCal.cxx:121
Bool_t AddToManager(TCal *cal, UInt_t chanNum, Option_t *opt="")
std::vector< Float_t > GetENGCoeff(size_t range=0) const
Definition TChannel.h:192
void DestroyENGCal()
Definition TChannel.cxx:554
int GetNumber() const
Definition TChannel.h:169
static TChannel * GetChannelByNumber(int temp_num)
Definition TChannel.cxx:482
void AddENGCoefficient(Float_t temp, size_t range=0)
Definition TChannel.h:204
std::array< Double_t, 2 > fGainCoeffs
Definition TGainMatch.h:85
Bool_t FineMatchFast(TH1 *hist1, TPeak *peak1, TH1 *hist2, TPeak *peak2, Int_t channelNum=9999)
TH1 * fHist
Definition TGainMatch.h:83
static Double_t fDefaultCoarseRange
Definition TGainMatch.h:89
std::array< Double_t, 2 > fAlignCoeffs
Definition TGainMatch.h:84
void Clear(Option_t *opt="") override
static Bool_t FineMatchFastAll(TCalManager *cm, TH2 *mat, Double_t energy1, Double_t energy2)
Double_t HistCompare(Double_t *x, Double_t *par)
void WriteToChannel() const override
void Copy(TObject &obj) const override
Bool_t CoarseMatch(TH1 *hist, Int_t chanNum=9999, Double_t energy1=1173.228, Double_t energy2=1332.492)
void CalculateGain(Double_t cent1, Double_t cent2, Double_t eng1, Double_t eng2)
Bool_t fCoarseMatch
Definition TGainMatch.h:81
Bool_t Align(TH1 *test, TH1 *hist, Int_t low_range=100, Int_t high_range=600)
static Bool_t CoarseMatchAll(TCalManager *cm, TH2 *mat, Double_t energy1=1173.228, Double_t energy2=1332.492)
static Bool_t AlignAll(TCalManager *cm, TH1 *hist, TH2 *mat, Int_t low_range=100, Int_t high_range=600)
static Bool_t FineMatchAll(TCalManager *cm, TH2 *charge_mat, TH2 *eng_mat, Int_t testchan, Double_t energy1, Double_t energy2, Int_t low_range=100, Int_t high_range=600)
Bool_t FineMatch(TH1 *energyHist, TH1 *testhist, TH1 *chargeHist, Double_t energy1, Double_t energy2, Int_t low_range=100, Int_t high_range=600, Int_t channelNum=9999)
void Print(Option_t *opt="") const override
Double_t fCoarseRange
Definition TGainMatch.h:86
Bool_t fAligned
Definition TGainMatch.h:82
Definition TPeak.h:28
Bool_t Fit(TH1 *fitHist, Option_t *opt="")
Definition TPeak.cxx:191
Double_t GetCentroid() const
Definition TPeak.h:51