GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TAngularCorrelation.h
Go to the documentation of this file.
1#ifndef TANGULARCORRELATION_H
2#define TANGULARCORRELATION_H
3
4/** \addtogroup Fitting Fitting & Analysis
5 * @{
6 */
7
8#include <iostream>
9#include <fstream>
10#include <string>
11#include <vector>
12#include "TH1.h"
13#include "TH2.h"
14#include "THnSparse.h"
15#include "TGraphAsymmErrors.h"
16#include "TPeak.h"
17#include "TCanvas.h"
18
19class TAngularCorrelation : public TObject {
20private:
21 TH1D* fIndexCorrelation; /// 1D plot of counts vs. angular index
22
23 // for diagnostics and re-fitting
24 TH1D* fChi2{nullptr}; /// 1D plot of chi^2 vs. angular index
25 TH1D* fCentroid{nullptr}; /// 1D plot of centroid vs. angular index
26 TH1D* fFWHM{nullptr}; /// 1D plot of FWHM vs. angular index
27 std::map<Int_t, TPeak*> fPeaks; /// array of TPeaks used to create fIndexCorrelations
28 std::map<Int_t, TH1D*> f1DSlices; /// array of 1D histograms used to create fIndexCorrelations
29
30 // mapping information
31 std::map<Int_t, std::map<Int_t, Int_t>>
32 fIndexMap; /// 2D square array correlating array number pairs with angular index
33 Int_t fNumIndices{0}; /// number of angular indices
34 Int_t fIndexMapSize; /// size of fIndexMap
35 std::vector<Double_t> fAngleMap; /// array correlating angular index with opening angle
36 std::vector<Int_t> fWeights; /// array correlating angular index with weight (number of detector pairs at that index)
37
38 // folding and grouping information
39 std::vector<Int_t> fGroups; /// array correlating angular index with group assignment
40 std::vector<Double_t> fGroupAngles; /// array correlating group assignment with their average angles
41 Bool_t fFolded; /// switch to indicate a folded correlation
42 Bool_t fGrouped; /// switch to indicate a grouped correlation
43
44 // modified indices information
45 std::vector<Int_t> fModifiedIndices; // array correlating angular index with modified index
46 std::vector<Int_t> fModifiedWeights; // array correlating modified index with weights
47 std::vector<Double_t> fModifiedAngles; // array correlating modified index with angles
48
49public:
53 TAngularCorrelation& operator=(const TAngularCorrelation&) = default;
54 TAngularCorrelation& operator=(TAngularCorrelation&&) noexcept = default;
55 ~TAngularCorrelation() override;
56
57 //----------------- getters -----------------
59
60 // diagnostics and re-fitting
61 TH1D* GetChi2Hst() { return fChi2; }
62 TH1D* GetCentroidHst() { return fCentroid; }
63 TH1D* GetFWHMHst() { return fFWHM; }
64 TPeak* GetPeak(Int_t index);
65 TH1D* Get1DSlice(Int_t index) { return f1DSlices[index]; }
66
67 // information
68 Int_t GetAngularIndex(Int_t arraynum1, Int_t arraynum2); // returns the angular index for a pair of detectors
69 // TODO: move the next function to implementation file and check if in range
70 Double_t GetAngleFromIndex(Int_t index)
71 {
72 return fAngleMap[index];
73 } // returns the opening angle for a specific angular index
74 // TODO: move the next function to implementation file and check if in range
75 Double_t GetWeightFromIndex(Int_t index)
76 {
77 return fWeights[index];
78 } // returns the weight for a specific angular index
79 Int_t GetGroupFromIndex(Int_t index)
80 {
81 return fGroups[index];
82 } // returns the assigned group for a specific angular index
83 Double_t GetGroupAngleFromIndex(Int_t gindex)
84 {
85 return fGroupAngles[gindex];
86 } // returns the angle for each group index
87 Int_t GetModifiedIndex(Int_t index)
88 {
89 return fModifiedIndices[index];
90 } // returns the modified index from the angular index
91 Int_t GetModifiedWeight(Int_t modindex)
92 {
93 return fModifiedWeights[modindex];
94 } // returns in the weight from the modified index
95 Double_t GetModifiedAngleFromIndex(Int_t modindex)
96 {
97 return fModifiedAngles[modindex];
98 } // returns the angle from the modified index
99 Int_t GetNumGroups() const; // returns the number of groups assigned
100 Int_t GetNumModIndices() const; // returns the number of modified indices
101 Int_t GetWeightsSize() const { return fWeights.size(); } // returns in the size of the fWeights array
102
103 //----------------- setters -----------------
104 void SetIndexCorrelation(TH1D* hst) { fIndexCorrelation = hst; }
105
106 // diagnostics and re-fitting
107 // TODO: move the next function to implementation file and update fIndexCorrelation
108 void SetPeak(Int_t index, TPeak* peak) { fPeaks[index] = peak; }
109 void Set1DSlice(Int_t index, TH1D* slice) { f1DSlices[index] = slice; }
110
111 //----------------- functions that do most of the work
112 TH2D* Create2DSlice(THnSparse* hst, Double_t min, Double_t max, Bool_t fold, Bool_t group);
113 TH2D* Create2DSlice(TObjArray* hstarray, Double_t min, Double_t max, Bool_t fold, Bool_t group);
114 TH2D* Modify2DSlice(TH2* hst, Bool_t fold, Bool_t group);
115 TH1D* IntegralSlices(TH2* hst, Double_t min, Double_t max);
116 TH1D* FitSlices(TH2* hst, TPeak* peak, Bool_t visualization);
117 TH1D* DivideByWeights(TH1* hst, Bool_t fold, Bool_t group);
118 TGraphAsymmErrors* CreateGraphFromHst(TH1* hst, Bool_t fold, Bool_t group);
119 TGraphAsymmErrors* CreateGraphFromHst() { return CreateGraphFromHst(fIndexCorrelation, kFALSE, kFALSE); }
120
121 //----------------- functions for diagnostics, re-fitting, modification of ACs
122 void UpdatePeak(Int_t index, TPeak* peak);
123 void ScaleSingleIndex(TH1* hst, Int_t index, Double_t factor);
125 void UpdateDiagnostics();
126 void DisplayDiagnostics(TCanvas* c_diag);
127
128 //----------------- functions for checking and printing the mapping
129 void PrintIndexMap(); // print the map
130 void PrintAngleMap(); // print the map
131 void PrintWeights(); // print the map
132 void PrintGroupIndexMap(); // print the map
133 void PrintGroupAngleMap(); // print the group angle map
134 void PrintModifiedIndexMap(); // prints a map between angular and modified indices
135 void PrintModifiedAngleMap(); // prints a map of angles for the modified indices
136 void PrintModifiedWeights(); // prints a map of modified weights
137 void PrintModifiedConditions() const; // prints the current folding and grouping conditions
138 Bool_t CheckGroups(std::vector<Int_t>& group) const;
139 Bool_t CheckGroupAngles(std::vector<Double_t>& groupangles) const;
140 Bool_t CheckMaps(Bool_t fold, Bool_t group); // checks to make sure fIndexMap, fAngleMap, and fWeights are consistent
141 Bool_t CheckModifiedHistogram(TH1* hst) const; // checks to make sure histogram is consistent with current settings
142
143 //----------------- functions for generating the mapping
144 // original maps
145 static std::vector<Double_t> GenerateAngleMap(std::vector<Int_t>& arraynumbers, std::vector<Int_t>& distances);
146 static std::map<Int_t, std::map<Int_t, Int_t>> GenerateIndexMap(std::vector<Int_t>& arraynumbers,
147 std::vector<Int_t>& distances,
148 std::vector<Double_t>& anglemap);
149 static std::vector<Int_t> GenerateWeights(
150 std::vector<Int_t>& arraynumbers, std::vector<Int_t>& distances,
151 std::map<Int_t, std::map<Int_t, Int_t>>& indexmap); // with input of array number array (crystals that were
152 // present in data collection), generates the weights for each
153 // angular index (no input generates weights for 16 detectors)
154 Int_t GenerateMaps(std::vector<Int_t>& arraynumbers, std::vector<Int_t>& distances);
155 Int_t GenerateMaps(Int_t detectors, Int_t distance);
156
157 // modified maps
158 Int_t AssignGroupMaps(std::vector<Int_t>& group, std::vector<Double_t>& groupangles);
159 Int_t GenerateGroupMaps(std::vector<Int_t>& arraynumbers, std::vector<Int_t>& distances, std::vector<Int_t>& group,
160 std::vector<Double_t>& groupangles);
161 Int_t GenerateModifiedMaps(Bool_t fold, Bool_t group);
162 std::vector<Int_t> GenerateModifiedIndices(Bool_t fold, Bool_t group);
163 std::vector<Double_t> GenerateModifiedAngles(Bool_t fold, Bool_t group);
164 std::vector<Int_t> GenerateModifiedWeights(std::vector<Int_t>& modindices, std::vector<Int_t>& weights);
165 static std::vector<Int_t> GenerateFoldedIndices(std::vector<Double_t>& folds, std::vector<Double_t>& anglemap);
166 static std::vector<Double_t> GenerateFoldedAngles(std::vector<Double_t>& anglemap);
167 void ClearModifiedMaps();
168
169 /// \cond CLASSIMP
170 ClassDefOverride(TAngularCorrelation, 1) // NOLINT(readability-else-after-return)
171 /// \endcond
172};
173/*! @} */
174#endif
std::vector< Double_t > fGroupAngles
array correlating angular index with group assignment
TGraphAsymmErrors * CreateGraphFromHst()
void Set1DSlice(Int_t index, TH1D *slice)
TAngularCorrelation(TAngularCorrelation &&) noexcept=default
Double_t GetWeightFromIndex(Int_t index)
Int_t GetAngularIndex(Int_t arraynum1, Int_t arraynum2)
TPeak * GetPeak(Int_t index)
TH1D * DivideByWeights(TH1 *hst, Bool_t fold, Bool_t group)
Bool_t CheckModifiedHistogram(TH1 *hst) const
std::vector< Double_t > fAngleMap
size of fIndexMap
Int_t GetModifiedWeight(Int_t modindex)
Int_t AssignGroupMaps(std::vector< Int_t > &group, std::vector< Double_t > &groupangles)
std::vector< Int_t > fGroups
array correlating angular index with weight (number of detector pairs at that index)
static std::vector< Int_t > GenerateFoldedIndices(std::vector< Double_t > &folds, std::vector< Double_t > &anglemap)
Bool_t CheckGroups(std::vector< Int_t > &group) const
std::vector< Double_t > GenerateModifiedAngles(Bool_t fold, Bool_t group)
Int_t GetGroupFromIndex(Int_t index)
TH2D * Create2DSlice(THnSparse *hst, Double_t min, Double_t max, Bool_t fold, Bool_t group)
std::vector< Int_t > GenerateModifiedWeights(std::vector< Int_t > &modindices, std::vector< Int_t > &weights)
Bool_t fGrouped
switch to indicate a folded correlation
void PrintModifiedIndexMap()
Prints map of angular index vs. Folded Index.
Int_t fIndexMapSize
number of angular indices
TH1D * fFWHM
1D plot of centroid vs. angular index
TH1D * IntegralSlices(TH2 *hst, Double_t min, Double_t max)
Bool_t CheckGroupAngles(std::vector< Double_t > &groupangles) const
static std::vector< Double_t > GenerateAngleMap(std::vector< Int_t > &arraynumbers, std::vector< Int_t > &distances)
std::vector< Int_t > fModifiedWeights
TH1D * FitSlices(TH2 *hst, TPeak *peak, Bool_t visualization)
static std::vector< Int_t > GenerateWeights(std::vector< Int_t > &arraynumbers, std::vector< Int_t > &distances, std::map< Int_t, std::map< Int_t, Int_t > > &indexmap)
Bool_t fFolded
array correlating group assignment with their average angles
void ScaleSingleIndex(TH1 *hst, Int_t index, Double_t factor)
Int_t GenerateGroupMaps(std::vector< Int_t > &arraynumbers, std::vector< Int_t > &distances, std::vector< Int_t > &group, std::vector< Double_t > &groupangles)
static std::map< Int_t, std::map< Int_t, Int_t > > GenerateIndexMap(std::vector< Int_t > &arraynumbers, std::vector< Int_t > &distances, std::vector< Double_t > &anglemap)
std::vector< Int_t > GenerateModifiedIndices(Bool_t fold, Bool_t group)
void SetPeak(Int_t index, TPeak *peak)
void UpdatePeak(Int_t index, TPeak *peak)
TAngularCorrelation(const TAngularCorrelation &)=default
TH1D * fChi2
1D plot of counts vs. angular index
Double_t GetGroupAngleFromIndex(Int_t gindex)
Int_t fNumIndices
2D square array correlating array number pairs with angular index
Int_t GenerateMaps(std::vector< Int_t > &arraynumbers, std::vector< Int_t > &distances)
TH1D * Get1DSlice(Int_t index)
Double_t GetModifiedAngleFromIndex(Int_t modindex)
std::map< Int_t, TH1D * > f1DSlices
array of TPeaks used to create fIndexCorrelations
static std::vector< Double_t > GenerateFoldedAngles(std::vector< Double_t > &anglemap)
std::vector< Int_t > fWeights
array correlating angular index with opening angle
std::map< Int_t, TPeak * > fPeaks
1D plot of FWHM vs. angular index
void DisplayDiagnostics(TCanvas *c_diag)
Bool_t CheckMaps(Bool_t fold, Bool_t group)
std::vector< Int_t > fModifiedIndices
switch to indicate a grouped correlation
Int_t GetModifiedIndex(Int_t index)
std::vector< Double_t > fModifiedAngles
Int_t GenerateModifiedMaps(Bool_t fold, Bool_t group)
std::map< Int_t, std::map< Int_t, Int_t > > fIndexMap
array of 1D histograms used to create fIndexCorrelations
TH1D * fCentroid
1D plot of chi^2 vs. angular index
Double_t GetAngleFromIndex(Int_t index)
TH2D * Modify2DSlice(TH2 *hst, Bool_t fold, Bool_t group)
void SetIndexCorrelation(TH1D *hst)
Definition TPeak.h:23