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
18//.............................................................................
20
21 // Set to NuFIT 5.2 values (NO w/ SK)
22 p->SetDm(2, 7.41e-5);
23 p->SetDm(3, 2.507e-3);
24 p->SetAngle(1,2, asin(sqrt(0.303)));
25 p->SetAngle(1,3, asin(sqrt(0.02225)));
26 p->SetAngle(2,3, asin(sqrt(0.451)));
27 p->SetDelta(1,3, 232 * TMath::DegToRad());
28
29}
30
31//.............................................................................
32OscProb::PMNS_Fast* GetFast(bool is_nominal){
33
36 return p;
37
38}
39
40//.............................................................................
41OscProb::PMNS_Iter* GetIter(bool is_nominal){
42
45 return p;
46
47}
48
49//.............................................................................
50OscProb::PMNS_Deco* GetDeco(bool is_nominal){
51
54 if(!is_nominal){
55 p->SetGamma(2, 1e-23);
56 p->SetGamma(3, 1e-22);
57 }
58
59 return p;
60
61}
62
63//.............................................................................
64OscProb::PMNS_OQS* GetOQS(bool is_nominal){
65
68 if(!is_nominal){
69 p->SetDecoElement(3, sqrt(2e-23));
70 p->SetDecoElement(8, sqrt(4e-23));
71 p->SetDecoAngle(3,8, acos(0.5));
72 }
73
74 return p;
75
76}
77
78//.............................................................................
80
83 if(!is_nominal){
84 p->SetDm(4, 0.1);
85 p->SetAngle(1,4, 0.1);
86 p->SetAngle(2,4, 0.1);
87 p->SetAngle(3,4, 0.1);
88 }
89
90 return p;
91
92}
93
94//.............................................................................
95OscProb::PMNS_Decay* GetDecay(bool is_nominal){
96
99 if(!is_nominal){
100 p->SetAlpha3(1e-4);
101 }
102
103 return p;
104
105}
106
107//.............................................................................
108OscProb::PMNS_NSI* GetNSI(bool is_nominal){
109
112 if(!is_nominal){
113 p->SetEps(0,0, 0.1, 0);
114 p->SetEps(0,1, 0.2, 0);
115 p->SetEps(0,2, 0.3, 0);
116 p->SetEps(1,1, 0.4, 0);
117 p->SetEps(1,2, 0.5, 0);
118 p->SetEps(2,2, 0.6, 0);
119 }
120
121 return p;
122
123}
124
125//.............................................................................
126OscProb::PMNS_SNSI* GetSNSI(bool is_nominal){
127
130 if(!is_nominal){
131 p->SetEps(0,0, 0.1, 0);
132 p->SetEps(0,1, 0.2, 0);
133 p->SetEps(0,2, 0.3, 0);
134 p->SetEps(1,1, 0.4, 0);
135 p->SetEps(1,2, 0.5, 0);
136 p->SetEps(2,2, 0.6, 0);
137 }
138
139 return p;
140
141}
142
143//.............................................................................
144OscProb::PMNS_LIV* GetLIV(bool is_nominal){
145
148 if(!is_nominal){
149 p->SetaT(0,0, 0.1e-22, 0);
150 p->SetaT(0,1, 0.2e-22, 0);
151 p->SetaT(0,2, 0.3e-22, 0);
152 p->SetaT(1,1, 0.4e-22, 0);
153 p->SetaT(1,2, 0.5e-22, 0);
154 p->SetaT(2,2, 0.6e-22, 0);
155 p->SetcT(0,0, 0.1e-22, 0);
156 p->SetcT(0,1, 0.2e-22, 0);
157 p->SetcT(0,2, 0.3e-22, 0);
158 p->SetcT(1,1, 0.4e-22, 0);
159 p->SetcT(1,2, 0.5e-22, 0);
160 p->SetcT(2,2, 0.6e-22, 0);
161 }
162
163 return p;
164
165}
166
167//.............................................................................
168OscProb::PMNS_NUNM* GetNUNM(bool is_nominal){
169
172 if(!is_nominal){
173 p->SetAlpha(0,0, 0.05, 0);
174 p->SetAlpha(1,0, 0.06, 0);
175 p->SetAlpha(2,0, 0.07, 0);
176 p->SetAlpha(1,1, 0.08, 0);
177 p->SetAlpha(2,1, 0.09, 0);
178 p->SetAlpha(2,2, 0.1, 0);
179 }
180
181 return p;
182
183}
184
185//.............................................................................
186OscProb::PMNS_Base* GetModel(string model, bool is_nominal = false){
187
188 if(model == "Iter") return GetIter(is_nominal);
189 if(model == "Deco") return GetDeco(is_nominal);
190 if(model == "Sterile") return GetSterile(is_nominal);
191 if(model == "Decay") return GetDecay(is_nominal);
192 if(model == "NSI") return GetNSI(is_nominal);
193 if(model == "LIV") return GetLIV(is_nominal);
194 if(model == "SNSI") return GetSNSI(is_nominal);
195 if(model == "NUNM") return GetNUNM(is_nominal);
196 if(model == "OQS") return GetOQS(is_nominal);
197
198 return GetFast(is_nominal);
199
200}
201
202//.............................................................................
203vector<string> GetListOfModels(){
204
205 return {"Fast", "Iter", "Sterile", "NSI",
206 "Deco", "Decay", "LIV", "SNSI",
207 "NUNM", "OQS"};
208
209}
210
211//.............................................................................
213
214 p->SetPath(1000, 2);
215 p->AddPath(1000, 4);
216 p->AddPath(1000, 2);
217
218}
219
220//.............................................................................
221void SaveTestFile(OscProb::PMNS_Base* p, TString filename){
222
223 SetTestPath(p);
224
225 int nbins = 100;
226 vector<double> xbins = GetLogAxis(nbins, 0.1, 10);
227 TH1D* h = 0;
228
229 TFile* f = new TFile("data/"+filename, "recreate");
230
231 for(int flvi=0; flvi<3; flvi++){
232 for(int flvf=0; flvf<3; flvf++){
233 for(int isnb=0; isnb<2; isnb++){
234 p->SetIsNuBar(isnb);
235 TString hname = TString::Format("h%d%d%d",flvi,flvf,isnb);
236 h = new TH1D(hname, "", nbins, &xbins[0]);
237 for(int i=1; i<=nbins; i++){
238 double energy = h->GetBinCenter(i);
239 double dE = h->GetBinWidth(i);
240 h->SetBinContent(i, p->AvgProb(flvi, flvf, energy, dE));
241 }
242 h->Write();
243 delete h;
244 }}}
245
246 f->Close();
247 cout << "Saved new test file: data/" + filename << endl;
248
249}
250
251//.............................................................................
252int CheckProb(OscProb::PMNS_Base* p, TString filename){
253
254 SetTestPath(p);
255
256 TFile* f = TFile::Open("data/"+filename, "read");
257
258 if(!f){
259 printf((Color::FAILED + " data/%s not found\n").c_str(), filename.Data());
260 return 1;
261 }
262
263 int ntests = 0;
264 int fails = 0;
265
266 TCanvas* c1 = 0;
267 TH1D* h0 = 0;
268 TH1D* h = 0;
269
270 for(int flvi=0; flvi<3; flvi++){
271 for(int flvf=0; flvf<3; flvf++){
272 for(int isnb=0; isnb<2; isnb++){
273 p->SetIsNuBar(isnb);
274 TString hname = TString::Format("h%d%d%d",flvi,flvf,isnb);
275 h0 = (TH1D*)f->Get(hname);
276 h = (TH1D*)h0->Clone();
277 bool plot = false;
278 for(int i=1; i<=h0->GetNbinsX(); i++){
279 double energy = h->GetBinCenter(i);
280 double dE = h->GetBinWidth(i);
281 double p0 = h0->GetBinContent(i);
282 double p1 = p->AvgProb(flvi, flvf, energy, dE);
283 ntests++;
284 if(abs(p0-p1)>1e-12){
285 plot = true;
286 fails++;
287 }
288 h->SetBinContent(i, p1);
289 }
290 if(plot){
291 c1 = new TCanvas();
292 c1->Divide(1,2);
293 c1->cd(1);
294 SetHist(h0, kBlue);
295 SetHist(h, kRed);
296 TString nu_lab = isnb ? "#bar{#nu}" : "#nu";
297 TString flv_lab[3] = {"e","#mu","#tau"};
298 TString ylab = "P(" + nu_lab + "_{" + flv_lab[flvi] +
299 "}#rightarrow" + nu_lab + "_{" + flv_lab[flvf] + "})";
300 h0->SetTitle(";Energy [GeV];"+ylab+";");
301 h->SetLineStyle(7);
302 double ymax = max(h->GetMaximum(),
303 h0->GetMaximum());
304 double ymin = min(h->GetMinimum(),
305 h0->GetMinimum());
306 h0->GetYaxis()->SetRangeUser(ymin, ymax);
307 h0->DrawCopy("hist");
308 h->DrawCopy("hist same");
309 SetTH1Margin();
310 gPad->SetLogx();
311 c1->cd(2);
312 h->SetTitle(";Energy [GeV];#Delta"+ylab+";");
313 h->Add(h0,-1);
314 h->SetLineStyle(1);
315 h->DrawCopy("hist");
316 SetTH1Margin();
317 gPad->SetLogx();
318 MiscText(0.55,0.8,0.1,filename,kGray,1);
319 c1->DrawClone();
320 TString pngfile = filename;
321 pngfile.ReplaceAll(".root",".png");
322 c1->SaveAs("plots/Failed_"+hname+"_"+pngfile);
323 delete c1;
324 }
325 if(h0) delete h0;
326 if(h) delete h;
327 }}}
328
329 if(fails>0){
330 printf((Color::FAILED + " Found %d differences in %d tests (%.3g%%) in %s\n").c_str(),
331 fails, ntests, 100.*fails/ntests, filename.Data());
332 }
333 else {
334 cout << Color::PASSED << " No differences found in " << filename << endl;
335 }
336
337 return fails;
338
339}
OscProb::PMNS_NUNM * GetNUNM(bool is_nominal)
Definition: Utils.h:168
OscProb::PMNS_LIV * GetLIV(bool is_nominal)
Definition: Utils.h:144
OscProb::PMNS_Fast * GetFast(bool is_nominal)
Definition: Utils.h:32
OscProb::PMNS_Deco * GetDeco(bool is_nominal)
Definition: Utils.h:50
OscProb::PMNS_Sterile * GetSterile(bool is_nominal)
Definition: Utils.h:79
vector< string > GetListOfModels()
Definition: Utils.h:203
void SetNominalPars(OscProb::PMNS_Base *p)
Definition: Utils.h:19
OscProb::PMNS_Iter * GetIter(bool is_nominal)
Definition: Utils.h:41
OscProb::PMNS_NSI * GetNSI(bool is_nominal)
Definition: Utils.h:108
OscProb::PMNS_Base * GetModel(string model, bool is_nominal=false)
Definition: Utils.h:186
void SetTestPath(OscProb::PMNS_Base *p)
Definition: Utils.h:212
OscProb::PMNS_Decay * GetDecay(bool is_nominal)
Definition: Utils.h:95
void SaveTestFile(OscProb::PMNS_Base *p, TString filename)
Definition: Utils.h:221
OscProb::PMNS_SNSI * GetSNSI(bool is_nominal)
Definition: Utils.h:126
OscProb::PMNS_OQS * GetOQS(bool is_nominal)
Definition: Utils.h:64
int CheckProb(OscProb::PMNS_Base *p, TString filename)
Definition: Utils.h:252
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