MaCh3 DUNE 1.0.0
Reference Guide
Loading...
Searching...
No Matches
splinesDUNE.cpp
Go to the documentation of this file.
1#include "splinesDUNE.h"
3
4#include "TKey.h"
5#include "TROOT.h"
6
7splinesDUNE::splinesDUNE(covarianceXsec* xsec_cov) : splineFDBase(xsec_cov) {
8 MACH3LOG_INFO("Created splinesDUNE object");
9}
10
12 MACH3LOG_INFO("Deleting splineSKBase object");
13}
14
15//****************************************
16void splinesDUNE::FillSampleArray(std::string SampleName, std::vector<std::string> OscChanFileNames)
17// Performs two jobs
18// 1. Fills indexing/each sample
19// 2. Creates the big spline vector
20//****************************************
21{
22
23 int iSample = getSampleIndex(SampleName);
24
25 enum ModeState {
26 kInConfig = 0,
27 kInConfigAndInFile,
28 kInFile
29 };
30
31 std::vector<std::map<std::string,ModeState>> ModeStatus;
32
33 //Create map of Mode Status
34 for (unsigned iSyst = 0; iSyst < SplineFileParPrefixNames[iSample].size(); iSyst++)
35 {
36 auto modes = SplineModeVecs[iSample][iSyst];
37 ModeStatus.emplace_back();
38
39 for (auto const & mode : modes)
40 {
41 //Add Modes from config
42 ModeStatus.back()[MaCh3mode_ToDUNEString((MaCh3_Mode)mode).c_str()] = kInConfig;
43 }
44 }
45
46 int nOscChannels = nOscChans[iSample];
47
48 for (int iOscChan = 0; iOscChan < nOscChannels; iOscChan++)
49 {
50 std::cout << "Processing:" << OscChanFileNames[iOscChan] << std::endl;
51
52 TFile *File = new TFile(OscChanFileNames[iOscChan].c_str(), "READ");
53 if (!File || File->IsZombie())
54 {
55 std::cerr << "File " << OscChanFileNames[iOscChan] << " not found" << std::endl;
56 throw;
57 }
58
59 //This is the MC specific part of the code
60 //i.e. we always assume that the splines are just store in single TDirectory and they're all in there as single objects
61 TIter Next(File->GetListOfKeys());
62 TKey *Key;
63
64 std::set<std::string> unique_spline_names;
65 int nb_splines = 0;
66
67 while ((Key = (TKey *)Next()))
68 {
69 TClass *Class = gROOT->GetClass(Key->GetClassName());
70
71 //Skip the TGraphs also in the spline files
72 if (!Class->InheritsFrom("TSpline3")){continue;}
73 const char* keyName = Key->GetName();
74
75 char* SplineName = new char[strlen(keyName) + 1];
76 strcpy(SplineName, keyName);
77
78 nb_splines += 1;
79 if(unique_spline_names.count(std::string(SplineName)) > 0){
80 if (std::string(SplineName).find("unknown") == std::string::npos){
81 //std::cout << "Repeated entry for spline named: " << std::string(SplineName) << std::endl;
82 continue;
83 }
84 }
85 unique_spline_names.insert(std::string(SplineName));
86
87 char *Syst;
88 char *Mode;
89 int Var1Bin;
90 int Var2Bin;
91 int Var3Bin;
92
93 char *Token = strtok(SplineName, "_");
94 Token = strtok(NULL, "_");
95 Syst = Token;
96
97 int SystNum = -1;
98 for (unsigned iSyst = 0; iSyst < SplineFileParPrefixNames[iSample].size(); iSyst++) {
99 if (strcmp(Syst, SplineFileParPrefixNames[iSample][iSyst].c_str()) == 0) {
100 SystNum = iSyst;
101 break;
102 }
103 }
104
105 // If the syst doesn't match any of the spline names then skip it
106 if (SystNum == -1){
107 continue;
108 }
109
110 int ModeNum = -1;
111 Mode = strtok(NULL, "_");
112 for (unsigned int iMode = 0; iMode < SplineModeVecs[iSample][SystNum].size(); iMode++) {
113 if (strcmp(Mode, MaCh3mode_ToDUNEString((MaCh3_Mode)SplineModeVecs[iSample][SystNum][iMode]).c_str()) == 0) {
114 ModeNum = iMode;
115 break;
116 }
117 }
118
119 //Check if mode has been registered already
120 if(ModeStatus[SystNum].count(Mode))
121 {
122 //Chech if mode has been found in config
123 if(ModeStatus[SystNum][Mode]!=kInFile)
124 {
125 ModeStatus[SystNum][Mode]=kInConfigAndInFile;
126 }
127 else
128 {
129 continue; //Skip if mode has been found in file but not config
130 }
131 }
132 else
133 {
134 ModeStatus[SystNum][Mode]=kInFile; //If mode hasn't been registered then skip because it's not in the config
135 continue;
136 }
137
138 TSpline3 *Obj = (TSpline3 *)Key->ReadObj();
139 TSpline3_red *Spline = new TSpline3_red(Obj);
140 delete Obj;
141
142 Token = strtok(NULL, "_"); // DB Needed to remove sp from spline name
143
144 Var1Bin = atoi(strtok(NULL, "_"));
145 Var2Bin = atoi(strtok(NULL, "_"));
146
147 char *Var3Bin_Char = strtok(NULL, "_");
148 if (Var3Bin_Char == NULL)
149 {
150 Var3Bin = 0;
151 }
152 else
153 {
154 Var3Bin = atoi(Var3Bin_Char);
155 }
156
157 if (isValidSplineIndex(SampleName, iOscChan, SystNum, ModeNum, Var1Bin, Var2Bin, Var3Bin))
158 {
159 // loop over all the spline knots and check their value
160 // if the value is 1 then set the flat bool to false
161 int nKnots = Spline->GetNp();
162 bool isFlat = true;
163
164 for (int iKnot = 0; iKnot < nKnots; iKnot++)
165 {
166 double x = -999;
167 double y = -999;
168 Spline->GetKnot(iKnot, x, y);
169
170 if (x == -999 || y == -999)
171 {
172 std::cerr << "Something has gone wrong... knot position is at -999" << std::endl;
173 std::cerr << "This error brought you by the folks at : "<<__FILE__<<" : "<<__LINE__<<std::endl;
174 throw;
175 }
176
177 double Eval = Spline->Eval(x);
178 if (Eval < 0.99999 || Eval > 1.00001)
179 {
180 isFlat = false;
181 break;
182 }
183 }
184
185 //Rather than keeping a mega vector of splines then converting, this should just keep everything nice in memory!
186 indexvec[iSample][iOscChan][SystNum][ModeNum][Var1Bin][Var2Bin][Var3Bin]=MonolithIndex;
187
188 coeffindexvec.push_back(CoeffIndex);
189 // Should save memory rather saving [x_i_0 ,... x_i_maxknots] for every spline!
190 if (isFlat)
191 {
192 splinevec_Monolith.push_back(NULL);
193 delete Spline;
194 }
195 else{
196 splinevec_Monolith.push_back(Spline);
197 int np=Spline->GetNp();
198 uniquecoeffindices.push_back(MonolithIndex); //So we can get the unique coefficients and skip flat splines later on!
199 CoeffIndex+=np;
200 }
201
202 MonolithIndex+=1;
203 }
204 }//End of loop over all TKeys in file
205 //ETA - I have no idea why but this breaks in ROOT 6.24 :/
206 std::cout << "Got " << nb_splines << " total splines with " << unique_spline_names.size() << " unique names." << std::endl;
207 delete File;
208 } //End of oscillation channel loop
209
210 //Find all modes which have been found in the spline file but have not been specified in the config
211 for (unsigned iSyst = 0; iSyst < SplineFileParPrefixNames[iSample].size(); iSyst++)
212 {
213 std::vector<std::string> MissedModes;
214 for (auto const & ModeUsage : ModeStatus[iSyst])
215 {
216 if(ModeUsage.second==kInFile)
217 {
218 MissedModes.push_back(ModeUsage.first);
219 }
220 }
221
222 if(MissedModes.size()!=0)
223 {
224 MACH3LOG_INFO("Parameter {} has splines for {} modes which have not been read in!", SplineFileParPrefixNames[iSample][iSyst].c_str(), MissedModes.size());
225 std::cout << "Modes: ";
226 for (auto const & Mode : MissedModes)
227 {
228 std::cout << Mode << " ";
229 }
230 std::cout << std::endl;
231 }
232 }
233
234 return;
235}
236
237//****************************************
238std::vector< std::vector<int> > splinesDUNE::GetEventSplines(std::string SampleName, int iOscChan, int EventMode, double Var1Val, double Var2Val, double Var3Val)
239//****************************************
240{
241 std::vector<std::vector<int>> ReturnVec;
242 int SampleIndex = -1;
243 for (unsigned int iSample = 0; iSample < SampleNames.size(); iSample++) {
244 if (SampleName == SampleNames[iSample]) {
245 SampleIndex = iSample;
246 }
247 }
248
249 if (SampleIndex == -1)
250 {
251 MACH3LOG_ERROR("Sample not found: {}", SampleName);
252 throw MaCh3Exception(__FILE__, __LINE__);
253 }
254
255 int nSplineSysts = (int)indexvec[SampleIndex][iOscChan].size();
256 //ETA- this is already a MaCh3 mode
257 //int Mode = MaCh3Mode_to_SplineMode(MaCh3_Mode(EventMode));
258 int Mode = MaCh3_Mode(EventMode);
259
260 int Var1Bin = SplineBinning[SampleIndex][iOscChan][0]->FindBin(Var1Val)-1;
261 if (Var1Bin < 0 || Var1Bin >= SplineBinning[SampleIndex][iOscChan][0]->GetNbins()){
262 //Explicitly push back with an empty vector
263 ReturnVec.push_back(std::vector<int>());
264 return ReturnVec;
265 }
266
267 int Var2Bin = SplineBinning[SampleIndex][iOscChan][1]->FindBin(Var2Val)-1;
268 if (Var2Bin < 0 || Var2Bin >= SplineBinning[SampleIndex][iOscChan][1]->GetNbins()){
269 //Explicitly push back with an empty vector
270 ReturnVec.push_back(std::vector<int>());
271 return ReturnVec;
272 }
273
274 int Var3Bin = SplineBinning[SampleIndex][iOscChan][2]->FindBin(Var3Val)-1;
275
276 if (Var3Bin < 0 || Var3Bin >= SplineBinning[SampleIndex][iOscChan][2]->GetNbins()){
277 //Explicitly push back with an empty vector
278 ReturnVec.push_back(std::vector<int>());
279 return ReturnVec;
280 }
281
282 for(int iSyst=0; iSyst<nSplineSysts; iSyst++){
283 std::vector<int> spline_modes = SplineModeVecs[SampleIndex][iSyst];
284 int nSampleModes = (int)spline_modes.size();
285
286 //ETA - look here at the length of spline_modes and what you're actually comparing against
287 for(int iMode = 0; iMode<nSampleModes ; iMode++){
288 if(Mode == spline_modes[iMode]){
289 std::vector<int> event_vec(7);
290 event_vec[0]=SampleIndex;
291 event_vec[1]=iOscChan;
292 event_vec[2]=iSyst;
293 event_vec[3]=iMode;
294 event_vec[4]=Var1Bin;
295 event_vec[5]=Var2Bin;
296 event_vec[6]=Var3Bin;
297 int splineID=indexvec[SampleIndex][iOscChan][iSyst][iMode][Var1Bin][Var2Bin][Var3Bin];
298 if(!isflatarray[splineID]){
299 ReturnVec.push_back(event_vec);
300 }
301 }
302 }
303 }
304
305 return ReturnVec;
306}
307
308// Converts MaCh3 modes to unique modes that splines apply to
309// (for example if CCRES and CCCoherent are treated as one spline mode)
310std::vector< std::vector<int> > splinesDUNE::StripDuplicatedModes(std::vector< std::vector<int> > InputVector) {
311
312 //ETA - this is of size nPars from the xsec model
313 int InputVectorSize = InputVector.size();
314 std::vector< std::vector<int> > ReturnVec(InputVectorSize);
315
316 //ETA - loop over all systematics
317 for (int iVec=0;iVec<InputVectorSize;iVec++) {
318 std::vector<int> TmpVec;
319
320 //Loop over the modes that we've listed in xsec cov
321 for (unsigned int iMode = 0 ; iMode < InputVector[iVec].size() ; iMode++){
322 //Convert the MaCh3 mode to spline mode (could use MaCh3Mode_SplineMode_Map here)
323 int iSplineMode = MaCh3Mode_to_SplineMode(InputVector[iVec][iMode]);
324 bool IncludeMode = true;
325
326 //Now check to see if we've already included this unique spline mode
327 for(unsigned int entry_i = 0 ; entry_i < TmpVec.size() ; entry_i++){
328 //Check to see if we've already included this unique SK spline mode
329 if(iSplineMode == TmpVec[entry_i]){
330 IncludeMode = false;
331 }
332 }
333
334 if(IncludeMode){
335 //Push back with the spline mode a systematic applies to
336 TmpVec.push_back(iSplineMode);
337 }
338
339 }//end of loop over modes for a syst
340
341 ReturnVec[iVec] = TmpVec;
342 }//end of loop over systs
343
344 return ReturnVec;
345}
std::string MaCh3mode_ToDUNEString(MaCh3_Mode i)
int MaCh3Mode_to_SplineMode(int iMode)
MaCh3_Mode
splinesDUNE(covarianceXsec *xsec_cov)
Constructor.
virtual std::vector< std::vector< int > > StripDuplicatedModes(std::vector< std::vector< int > > InputVector) override
Find sand strips dubplicate modes to ensure everything corresponds to MaCh3 modes.
virtual std::vector< std::vector< int > > GetEventSplines(std::string SampleName, int iOscChan, int EventMode, double Var1Val, double Var2Val, double Var3Val) override
Getter method for each spline.
virtual ~splinesDUNE()
Destructor.
virtual void FillSampleArray(std::string SampleName, std::vector< std::string > OscChanFileNames) override
Fills indexing for each sample and generates a large spline vector.