OscProb
Utils.h
Go to the documentation of this file.
1
2#include "colormod.h"
3
4#include "../tutorial/SetNiceStyle.C"
5
6#include "TMath.h"
7#include "TFile.h"
8
9#include "PMNS_OQS.h"
10#include "PMNS_LIV.h"
11#include "PMNS_Deco.h"
12#include "PMNS_Iter.h"
13#include "PMNS_NUNM.h"
14#include "PMNS_SNSI.h"
15#include "PMNS_Decay.h"
16#include "PMNS_Sterile.h"
17#include "PMNS_Avg.h"
18
19//.............................................................................
21
22 // Set to NuFIT 5.2 values (NO w/ SK)
23 p->SetDm(2, 7.41e-5);
24 p->SetDm(3, 2.507e-3);
25 p->SetAngle(1,2, asin(sqrt(0.303)));
26 p->SetAngle(1,3, asin(sqrt(0.02225)));
27 p->SetAngle(2,3, asin(sqrt(0.451)));
28 p->SetDelta(1,3, 232 * TMath::DegToRad());
29
30}
31
32//.............................................................................
33OscProb::PMNS_Fast* GetFast(bool is_nominal){
34
37 return p;
38
39}
40
41//.............................................................................
42OscProb::PMNS_Iter* GetIter(bool is_nominal){
43
46 return p;
47
48}
49
50//.............................................................................
51OscProb::PMNS_Deco* GetDeco(bool is_nominal){
52
55 if(!is_nominal){
56 p->SetGamma(2, 1e-23);
57 p->SetGamma(3, 1e-22);
58 }
59
60 return p;
61
62}
63
64//.............................................................................
65OscProb::PMNS_OQS* GetOQS(bool is_nominal){
66
69 if(!is_nominal){
70 p->SetDecoElement(3, sqrt(2e-23));
71 p->SetDecoElement(8, sqrt(4e-23));
72 p->SetDecoAngle(3,8, acos(0.5));
73 }
74
75 return p;
76
77}
78
79//.............................................................................
81
84 if(!is_nominal){
85 p->SetDm(4, 0.1);
86 p->SetAngle(1,4, 0.1);
87 p->SetAngle(2,4, 0.1);
88 p->SetAngle(3,4, 0.1);
89 }
90
91 return p;
92
93}
94
95//.............................................................................
96OscProb::PMNS_Decay* GetDecay(bool is_nominal){
97
100 if(!is_nominal){
101 p->SetAlpha3(1e-4);
102 }
103
104 return p;
105
106}
107
108//.............................................................................
109OscProb::PMNS_NSI* GetNSI(bool is_nominal){
110
113 if(!is_nominal){
114 p->SetEps(0,0, 0.1, 0);
115 p->SetEps(0,1, 0.2, 0);
116 p->SetEps(0,2, 0.3, 0);
117 p->SetEps(1,1, 0.4, 0);
118 p->SetEps(1,2, 0.5, 0);
119 p->SetEps(2,2, 0.6, 0);
120 }
121
122 return p;
123
124}
125
126//.............................................................................
127OscProb::PMNS_SNSI* GetSNSI(bool is_nominal){
128
131 if(!is_nominal){
132 p->SetEps(0,0, 0.1, 0);
133 p->SetEps(0,1, 0.2, 0);
134 p->SetEps(0,2, 0.3, 0);
135 p->SetEps(1,1, 0.4, 0);
136 p->SetEps(1,2, 0.5, 0);
137 p->SetEps(2,2, 0.6, 0);
138 }
139
140 return p;
141
142}
143
144//.............................................................................
145OscProb::PMNS_LIV* GetLIV(bool is_nominal){
146
149 if(!is_nominal){
150 p->SetaT(0,0, 0.1e-22, 0);
151 p->SetaT(0,1, 0.2e-22, 0);
152 p->SetaT(0,2, 0.3e-22, 0);
153 p->SetaT(1,1, 0.4e-22, 0);
154 p->SetaT(1,2, 0.5e-22, 0);
155 p->SetaT(2,2, 0.6e-22, 0);
156 p->SetcT(0,0, 0.1e-22, 0);
157 p->SetcT(0,1, 0.2e-22, 0);
158 p->SetcT(0,2, 0.3e-22, 0);
159 p->SetcT(1,1, 0.4e-22, 0);
160 p->SetcT(1,2, 0.5e-22, 0);
161 p->SetcT(2,2, 0.6e-22, 0);
162 }
163
164 return p;
165
166}
167
168//.............................................................................
169OscProb::PMNS_NUNM* GetNUNM(bool is_nominal){
170
173 if(!is_nominal){
174 p->SetAlpha(0,0, 0.05, 0);
175 p->SetAlpha(1,0, 0.06, 0);
176 p->SetAlpha(2,0, 0.07, 0);
177 p->SetAlpha(1,1, 0.08, 0);
178 p->SetAlpha(2,1, 0.09, 0);
179 p->SetAlpha(2,2, 0.1, 0);
180 }
181
182 return p;
183
184}
185
186//.............................................................................
187OscProb::PMNS_Avg* GetAvg(bool is_nominal){
188
191 return p;
192
193}
194
195//.............................................................................
196OscProb::PMNS_Base* GetModel(string model, bool is_nominal = false){
197
198 if(model == "Iter") return GetIter(is_nominal);
199 if(model == "Deco") return GetDeco(is_nominal);
200 if(model == "Sterile") return GetSterile(is_nominal);
201 if(model == "Decay") return GetDecay(is_nominal);
202 if(model == "NSI") return GetNSI(is_nominal);
203 if(model == "LIV") return GetLIV(is_nominal);
204 if(model == "SNSI") return GetSNSI(is_nominal);
205 if(model == "NUNM") return GetNUNM(is_nominal);
206 if(model == "OQS") return GetOQS(is_nominal);
207 if(model == "Avg") return GetAvg(is_nominal);
208
209 return GetFast(is_nominal);
210
211}
212
213//.............................................................................
214vector<string> GetListOfModels(){
215
216 return {"Fast", "Iter", "Sterile", "NSI",
217 "Deco", "Decay", "LIV", "SNSI",
218 "NUNM", "OQS", "Avg"};
219
220}
221
222//.............................................................................
224
225 p->SetPath(1000, 2);
226 p->AddPath(1000, 4);
227 p->AddPath(1000, 2);
228
229}
230
231//.............................................................................
232void SaveTestFile(OscProb::PMNS_Base* p, TString filename){
233
234 SetTestPath(p);
235
236 int nbins = 100;
237 vector<double> xbins = GetLogAxis(nbins, 0.1, 10);
238 TH1D* h = 0;
239
240 TFile* f = new TFile("data/"+filename, "recreate");
241
242 for(int flvi=0; flvi<3; flvi++){
243 for(int flvf=0; flvf<3; flvf++){
244 for(int isnb=0; isnb<2; isnb++){
245 p->SetIsNuBar(isnb);
246 TString hname = TString::Format("h%d%d%d",flvi,flvf,isnb);
247 h = new TH1D(hname, "", nbins, &xbins[0]);
248 for(int i=1; i<=nbins; i++){
249 double energy = h->GetBinCenter(i);
250 double dE = h->GetBinWidth(i);
251 h->SetBinContent(i, p->AvgProb(flvi, flvf, energy, dE));
252 }
253 h->Write();
254 delete h;
255 }}}
256
257 f->Close();
258 cout << "Saved new test file: data/" + filename << endl;
259
260}
261
262//.............................................................................
263int CheckProb(OscProb::PMNS_Base* p, TString filename){
264
265 SetTestPath(p);
266
267 TFile* f = TFile::Open("data/"+filename, "read");
268
269 if(!f){
270 printf((Color::FAILED + " data/%s not found\n").c_str(), filename.Data());
271 return 1;
272 }
273
274 int ntests = 0;
275 int fails = 0;
276
277 TCanvas* c1 = 0;
278 TH1D* h0 = 0;
279 TH1D* h = 0;
280
281 for(int flvi=0; flvi<3; flvi++){
282 for(int flvf=0; flvf<3; flvf++){
283 for(int isnb=0; isnb<2; isnb++){
284 p->SetIsNuBar(isnb);
285 TString hname = TString::Format("h%d%d%d",flvi,flvf,isnb);
286 h0 = (TH1D*)f->Get(hname);
287 h = (TH1D*)h0->Clone();
288 bool plot = false;
289 for(int i=1; i<=h0->GetNbinsX(); i++){
290 double energy = h->GetBinCenter(i);
291 double dE = h->GetBinWidth(i);
292 double p0 = h0->GetBinContent(i);
293 double p1 = p->AvgProb(flvi, flvf, energy, dE);
294 ntests++;
295 if(abs(p0-p1)>1e-12){
296 plot = true;
297 fails++;
298 }
299 h->SetBinContent(i, p1);
300 }
301 if(plot){
302 c1 = new TCanvas();
303 c1->Divide(1,2);
304 c1->cd(1);
305 SetHist(h0, kBlue);
306 SetHist(h, kRed);
307 TString nu_lab = isnb ? "#bar{#nu}" : "#nu";
308 TString flv_lab[3] = {"e","#mu","#tau"};
309 TString ylab = "P(" + nu_lab + "_{" + flv_lab[flvi] +
310 "}#rightarrow" + nu_lab + "_{" + flv_lab[flvf] + "})";
311 h0->SetTitle(";Energy [GeV];"+ylab+";");
312 h->SetLineStyle(7);
313 double ymax = max(h->GetMaximum(),
314 h0->GetMaximum());
315 double ymin = min(h->GetMinimum(),
316 h0->GetMinimum());
317 h0->GetYaxis()->SetRangeUser(ymin, ymax);
318 h0->DrawCopy("hist");
319 h->DrawCopy("hist same");
320 SetTH1Margin();
321 gPad->SetLogx();
322 c1->cd(2);
323 h->SetTitle(";Energy [GeV];#Delta"+ylab+";");
324 h->Add(h0,-1);
325 h->SetLineStyle(1);
326 h->DrawCopy("hist");
327 SetTH1Margin();
328 gPad->SetLogx();
329 MiscText(0.55,0.8,0.1,filename,kGray,1);
330 c1->DrawClone();
331 TString pngfile = filename;
332 pngfile.ReplaceAll(".root",".png");
333 c1->SaveAs("plots/Failed_"+hname+"_"+pngfile);
334 delete c1;
335 }
336 if(h0) delete h0;
337 if(h) delete h;
338 }}}
339
340 if(fails>0){
341 printf((Color::FAILED + " Found %d differences in %d tests (%.3g%%) in %s\n").c_str(),
342 fails, ntests, 100.*fails/ntests, filename.Data());
343 }
344 else {
345 cout << Color::PASSED << " No differences found in " << filename << endl;
346 }
347
348 return fails;
349
350}
OscProb::PMNS_NUNM * GetNUNM(bool is_nominal)
Definition: Utils.h:169
OscProb::PMNS_LIV * GetLIV(bool is_nominal)
Definition: Utils.h:145
OscProb::PMNS_Fast * GetFast(bool is_nominal)
Definition: Utils.h:33
OscProb::PMNS_Deco * GetDeco(bool is_nominal)
Definition: Utils.h:51
OscProb::PMNS_Sterile * GetSterile(bool is_nominal)
Definition: Utils.h:80
vector< string > GetListOfModels()
Definition: Utils.h:214
void SetNominalPars(OscProb::PMNS_Base *p)
Definition: Utils.h:20
OscProb::PMNS_Iter * GetIter(bool is_nominal)
Definition: Utils.h:42
OscProb::PMNS_NSI * GetNSI(bool is_nominal)
Definition: Utils.h:109
OscProb::PMNS_Base * GetModel(string model, bool is_nominal=false)
Definition: Utils.h:196
void SetTestPath(OscProb::PMNS_Base *p)
Definition: Utils.h:223
OscProb::PMNS_Decay * GetDecay(bool is_nominal)
Definition: Utils.h:96
void SaveTestFile(OscProb::PMNS_Base *p, TString filename)
Definition: Utils.h:232
OscProb::PMNS_SNSI * GetSNSI(bool is_nominal)
Definition: Utils.h:127
OscProb::PMNS_OQS * GetOQS(bool is_nominal)
Definition: Utils.h:65
int CheckProb(OscProb::PMNS_Base *p, TString filename)
Definition: Utils.h:263
OscProb::PMNS_Avg * GetAvg(bool is_nominal)
Definition: Utils.h:187
Implementation of oscillations of neutrinos in matter in a three-neutrino framework with a first orde...
Definition: PMNS_Avg.h:31
Base class implementing general functions for computing neutrino oscillations.
Definition: PMNS_Base.h:26
virtual void SetDm(int j, double dm)
Set the mass-splitting dm_j1 in eV^2.
Definition: PMNS_Base.cxx:674
virtual void SetDelta(int i, int j, double delta)
Set the CP phase delta_ij.
Definition: PMNS_Base.cxx:602
virtual void SetIsNuBar(bool isNuBar)
Set the anti-neutrino flag.
Definition: PMNS_Base.cxx:243
virtual double AvgProb(vectorC nu_in, int flvf, double E, double dE=0)
Compute the average probability over a bin of energy.
Definition: PMNS_Base.cxx:1568
virtual void AddPath(NuPath p)
Add a path to the sequence.
Definition: PMNS_Base.cxx:307
virtual void SetPath(NuPath p)
Set a single path.
Definition: PMNS_Base.cxx:330
virtual void SetAngle(int i, int j, double th)
Set the mixing angle theta_ij.
Definition: PMNS_Base.cxx:539
Implementation of neutrino decay in a three-neutrino framework.
Definition: PMNS_Decay.h:30
virtual void SetAlpha3(double alpha3)
Definition: PMNS_Decay.cxx:62
Implementation of oscillations of neutrinos in matter in a three-neutrino framework with decoherence.
Definition: PMNS_Deco.h:27
virtual void SetGamma(int j, double val)
Set any given decoherence parameter.
Definition: PMNS_Deco.cxx:49
Implementation of oscillations of neutrinos in matter in a three-neutrino framework.
Definition: PMNS_Fast.h:40
Implementation of oscillations of neutrinos in matter in a three-neutrino framework.
Definition: PMNS_Iter.h:35
Implements oscillations with LIV as modelled by SME.
Definition: PMNS_LIV.h:28
virtual void SetaT(int flvi, int flvj, int dim, double val, double phase)
Definition: PMNS_LIV.cxx:67
virtual void SetcT(int flvi, int flvj, int dim, double val, double phase)
Definition: PMNS_LIV.cxx:119
Implementation of oscillations of neutrinos in matter in a three-neutrino framework with NSI.
Definition: PMNS_NSI.h:28
virtual void SetEps(int flvi, int flvj, double val, double phase)
Set any given NSI parameter.
Definition: PMNS_NSI.cxx:86
Implementation of oscillations of neutrinos in matter in a three-neutrino framework with Non unitary ...
Definition: PMNS_NUNM.h:30
virtual void SetAlpha(int i, int j, double val, double phase)
Set any given NUNM parameter.
Definition: PMNS_NUNM.cxx:93
Implements neutrino oscillations using an open quantum system approach.
Definition: PMNS_OQS.h:28
virtual void SetDecoAngle(int i, int j, double th)
Set mixing angle between two decoherence parameters a_i, a_j.
Definition: PMNS_OQS.cxx:75
virtual void SetDecoElement(int i, double val)
Set value of the a_i decoherence element in Gell-Mann basis.
Definition: PMNS_OQS.cxx:58
Implementation of oscillations of neutrinos in matter in a three-neutrino framework with scalar NSI.
Definition: PMNS_SNSI.h:26
Implementation of oscillations of neutrinos in matter in a N-neutrino framework.
Definition: PMNS_Sterile.h:34
static const string PASSED
Definition: colormod.h:23
static const string FAILED
Definition: colormod.h:22