6 #include "TObjString.h" 7 #include "TTimeStamp.h" 11 #include "EnvironmentVariable.h" 12 #include "Math/WrappedTF1.h" 13 #include "Math/GaussIntegrator.h" 15 #include "blazars/fava.h" 17 #include "icemc_random.h" 22 SourceModel::SourceModel(
const char * n)
28 average_nonzero_flux= -1;
31 SourceModel::~SourceModel()
33 for (
unsigned i = 0; i < sources.size(); i++)
delete sources[i];
38 static double luminosity_distance(
double z)
40 static TF1 * f_comoving = 0;
43 f_comoving =
new TF1(
"comoving_distance",
"[0] / sqrt([1]*(1+x)^4 + [2] * (1+x)^3 + [3] * (1+x)^2 + [4])", 0,10);
44 f_comoving->SetParameters(67.74, 5e-5, 0.3089, 0, 0.6911);
47 double d_M = f_comoving->Integral(0,z);
53 const double txs_index = 2;
54 const double txs_norm = 2.2e-8;
55 const double txs_E0 = 100e3;
56 const double txs_flux = 3.8;
57 const double txs_z = 0.3365;
58 const double txs_D = luminosity_distance(txs_z);
62 const double m77_index = 3.16;
67 const double gamma_index = 2;
69 double SourceModel::Restriction::fromString(
const char * the_time)
71 printf(
"PARSING \"%s\"\n",the_time);
73 if(strcasestr(the_time,
"A4_MIN"))
77 else if(strcasestr(the_time,
"A3_MIN"))
81 else if(strcasestr(the_time,
"A4_MAX"))
85 else if(strcasestr(the_time,
"A3_MAX"))
89 else if(strcasestr(the_time,
"MAX"))
94 return atof(the_time);
102 if (!key || !strcasecmp(key,
"NONE"))
return NULL;
107 TObjArray * tokens = str.Tokenize(
"+");
110 for (
int i =0; i < tokens->GetEntries(); i++)
112 TObjString * tok = (TObjString*) tokens->At(i);
114 TString stripped( tok->String().Strip( TString::kBoth));
117 if (stripped ==
"TXS")
120 double txs_ra = 77.3582;
121 double txs_dec = 5.69314;
128 else if (stripped ==
"M77")
134 else if(stripped==
"SN")
136 TFile *file =
new TFile(
"./supernovae/data/supernovaeTNS.root");
137 TTree *tree = (TTree*) file->Get(
"tnsTree");
141 std::string *objName =
new std::string;
142 std::string *fullObjType =
new std::string;
143 std::string *objSubtype =
new std::string;
144 int discoveryUnixTime = 0;
146 tree->SetBranchAddress(
"ra",&ra);
147 tree->SetBranchAddress(
"dec",&dec);
148 tree->SetBranchAddress(
"name",&objName);
149 tree->SetBranchAddress(
"fullObjType",&fullObjType);
150 tree->SetBranchAddress(
"objSubtype",&objSubtype);
151 tree->SetBranchAddress(
"discoveryUnixTime",&discoveryUnixTime);
153 int two_weeks = 24*3600*14;
155 for (
unsigned int i = 0; i < tree->GetEntries(); i++)
164 if( (discoveryUnixTime < r.whichStartTime) || (discoveryUnixTime > r.whichEndTime + two_weeks) ){
continue;}
166 if( abs(dec)>r.maxDec){
continue;}
168 if(strcasecmp(r.whichSources,
"All"))
171 if(strcasestr(r.whichSources, objName->c_str()))
173 std::cout <<
"Using the specified source: " << objName << std::endl;
179 if(strcasecmp(r.whichSubtype,
"All"))
181 if(strcasestr(r.whichSubtype,objSubtype->c_str()))
195 else if (stripped.BeginsWith(
"GRB"))
197 TFile *file =
new TFile(
"./GRBs/data/GRBIceCube.root");
198 TTree *tree = (TTree*) file->Get(
"iceCubeTree");
202 std::string *objName =
new std::string;
203 int unixTriggerTime = 0;
204 float t90, fluence, redshift;
206 tree->SetBranchAddress(
"RA",&RA);
207 tree->SetBranchAddress(
"dec",&dec);
208 tree->SetBranchAddress(
"name",&objName);
209 tree->SetBranchAddress(
"unixTriggerTime",&unixTriggerTime);
210 tree->SetBranchAddress(
"T90",&t90);
211 tree->SetBranchAddress(
"fluence",&fluence);
212 tree->SetBranchAddress(
"redshift",&redshift);
214 for (
unsigned int i = 0; i < tree->GetEntries(); i++)
218 if (fluence < 0) fluence = 1e-6;
222 if( (unixTriggerTime < r.whichStartTime) || (unixTriggerTime > r.whichEndTime) ){
continue;}
224 if( abs(dec)>r.maxDec){
continue;}
226 printf(
"%d\n", unixTriggerTime);
229 if(strcasecmp(r.whichSources,
"All"))
232 if(strcasestr(r.whichSources,objName->c_str()))
234 std::cout <<
"Using the specified source: " << objName << std::endl;
240 bool modelCalculation =
true;
241 double energyCutoff = 1e10;
242 if(modelCalculation==
true)
246 double omegaRad = 9.236e-5;
247 double omegaM = 0.315;
248 double omegaLambda = 0.679;
249 double H0 = 67.4 * pow(10,3);
250 double c = 299792458;
251 double EGammaIso = 0.;
253 double MpcToCm = 3.08e+24;
259 double lorentzFactor = 300.;
265 double shockLorentzFactor = 0.;
266 double energyCGamma = 0.;
268 double energyBreakNeutrinoAG = 0.;
269 double shockRadius = 0.;
270 double massProton = 0.0015;
271 double magFieldShock = 0.;
272 double maxShockNuE = 0.;
276 if(t90==-999){t90 = 20;}
278 if(redshift==-999){redshift = 2;}
280 if(fluence==-999){fluence=1e-6;}
285 TF1 Ez(
"E(z)",
"1/(pow([0] * pow((1+x),3) + [1] + [2] * pow((1+x),4),0.5))", 0, redshift);
287 Ez.SetParameter(0,omegaM); Ez.SetParName(0,
"omegaM");
288 Ez.SetParameter(1,omegaLambda); Ez.SetParName(1,
"omegaLambda");
289 Ez.SetParameter(2,omegaRad); Ez.SetParName(2,
"omegaRad");
290 ROOT::Math::WrappedTF1 wEz(Ez);
291 ROOT::Math::GaussIntegrator ig;
293 ig.SetRelTolerance(0.001);
294 double integral = ig.Integral(0, redshift);
297 dL = (1+redshift) * c/H0 * integral;
300 EGammaIso = 4 * TMath::Pi() * dL*dL * fluence * 1/(1+redshift);
301 EGammaIso*=pow(MpcToCm,2);
308 EKinIso = EGammaIso/eE;
311 shockRadius = pow(((3*EKinIso)/(4*TMath::Pi() * nEx * massProton * c*c * lorentzFactor*lorentzFactor)),(1./3.));
315 std::cout <<
"Afterglow shock radius < burst radius. Something might be wrong." << std::endl;
319 magFieldShock = pow(8 * TMath::Pi() * eB * nEx * massProton,0.5);
321 shockLorentzFactor = 195 * pow(((EKinIso)/(pow(10,54))),0.125) * pow(t90/10,-0.375) * pow(nEx/1.,-0.125);
324 energyCGamma = 0.4 * pow(((EKinIso)/(pow(10,54))),(-25./24.)) * pow(((lorentzFactor)/(300)),(4./3.)) * pow(((t90)/(10)),(9./8.)) * pow(((nEx)/(1)),(-11./24.));
333 energyBreakNeutrinoAG = 0.015 * shockLorentzFactor * 1/(1+redshift) * pow(energyCGamma,-1);
335 maxShockNuE = fPiAG/(4*(1+redshift)) * eE * shockLorentzFactor * magFieldShock * shockRadius;
338 if(maxShockNuE > energyBreakNeutrinoAG)
342 else if(maxShockNuE < 1e9)
348 maxShockNuE = energyBreakNeutrinoAG+energyBreakNeutrinoAG*0.01;
351 energyCutoff = maxShockNuE;
355 if (strcasestr(stripped.Data(),
"AFTERGLOW") )
359 unixTriggerTime + 3600,1.5,
360 fluence,1e9, energyCutoff)
364 else if (strcasestr(stripped.Data(),
"PROMPT") )
369 unixTriggerTime + t90,4,
381 unixTriggerTime + 3600,1.5,
382 fluence,1e9, energyCutoff)
390 else if (stripped.BeginsWith(
"FAVA"))
394 bool flux_weighted = strcasestr(stripped.Data(),
"FLUX");
395 bool z_weighted = strcasestr(stripped.Data(),
"REDSHIFT");
399 TFile ffava(Form(
"%s/blazars/fava.root", dir));
400 TTree * fava_tree = (TTree*) ffava.Get(
"fava");
403 fava_tree->SetBranchAddress(
"fava",&fava);
405 for (
int i = 0; i < fava_tree->GetEntries(); i++)
407 fava_tree->GetEntry(i);
410 if( (fava->unix_tmax < r.whichStartTime) || (fava->unix_tmin > r.whichEndTime) ){
continue;}
412 const char * fava_name = fava->association.GetString().Data();
415 if(strcasestr(fava_name,
"None"))
422 if(strcasecmp(r.whichSources,
"All"))
425 if(strcasestr(r.whichSources,fava_name))
427 std::cout <<
"Using the specified source: " << fava_name << std::endl;
442 if (fava->source_class.GetString() ==
"bcu" || fava->source_class.GetString() ==
"fsrq" || fava->source_class.GetString() ==
"bll")
447 txs_norm * (flux_weighted ? fava->he_flux / txs_flux : z_weighted ? txs_flux * pow(txs_D/luminosity_distance(fava->z),2): 1), txs_E0))
457 std::cout <<
"----------------------" << std::endl;
459 if (!m->getNSources())
461 fprintf(stderr,
"WARNING: no sources added for key %s\n", key);
467 std::cout <<
"Sources added: " << m->getNSources() << std::endl;
473 TH1 * SourceModel::estimateFlux(
double tmin,
double tmax,
double Emin,
double Emax,
int nbins,
int N)
475 double from = log10(Emin);
476 double to = log10(Emax);
478 TH1D * spectrum =
new TH1D(
"spec", name, nbins, from,
to);
479 spectrum->SetDirectory(0);
480 spectrum->GetXaxis()->SetTitle(
"log10 (E GeV)");
481 TRandom * rng = getRNG(RNG_SOURCE);
482 for (
int i = 0; i < N; i++)
484 double t = rng->Uniform(tmin,tmax);
485 for (
unsigned j = 0; j < sources.size(); j++)
488 for (
int E = 1; E<= nbins; E++)
490 double l = TMath::Power(10,spectrum->GetBinLowEdge(E));
491 double h = TMath::Power(10,spectrum->GetBinLowEdge(E) + spectrum->GetBinWidth(E));
492 spectrum->Fill(spectrum->GetBinCenter(E), sources[j]->getFlux()->getFluxBetween(l,h,t));
497 spectrum->Scale(1./N);
507 for (
unsigned i = 0; i < sources.size(); i++)
509 sources[i]->getFlux()->getFluxTimeChanges(changes);
513 std::sort(changes->begin(), changes->end());
514 std::vector<double>::iterator it = std::unique(changes->begin(),changes->end());
515 changes->resize(std::distance(changes->begin(),it));
520 if (!sources.size())
return -1;
522 bool fixedEnergy = minE==maxE;
523 if (fixedEnergy) nuE = minE;
526 double total_flux = 0 ;
527 std::vector<double> fluxes(sources.size());
529 for (
unsigned i = 0; i < sources.size(); i++)
531 double f = fixedEnergy ? sources[i]->getFlux()->getFlux(nuE,t) : sources[i]->getFlux()->getFluxBetween(minE,maxE,t);
533 fluxes[i] = total_flux;
536 TRandom * rng = getRNG(RNG_SOURCE);
537 double random = rng->Uniform(0, total_flux);
538 unsigned index = std::upper_bound(fluxes.begin(), fluxes.end(), random) - fluxes.begin();
549 const Source * which = sources[index];
554 if (nudir) *nudir = which->getDirection(t);
556 if (!fixedEnergy) nuE = which->getFlux()->pickEnergy(minE,maxE,t,rng);
562 Source::Source(
const char * nm,
double ra,
double dc,
SourceFlux * f)
563 : name(nm), flux(f) , RA(ra), dec(dc*TMath::DegToRad())
567 Vector Source::getDirection(
double t)
const 570 #if ROOT_VERSION_CODE < ROOT_VERSION(6,0,0) 571 std::cerr <<
"ROOT 5 doesn't have AsLMST. Returning nonsense." << std::endl;
576 int nsecs = 1e9 *(t-secs);
577 TTimeStamp ts(secs, nsecs);
578 double lst = ts.AsLMST(-90);
580 h *= (TMath::Pi()/12);
581 return Vector(cos(h) * cos(dec), sin(h) * cos(dec), sin(dec));
586 ConstantExponentialSourceFlux::ConstantExponentialSourceFlux(
double e,
double norm,
double normE)
587 : gamma(e) , f(
"f",
"[0] * x^-[1]",1e9,1e12)
590 A = norm * TMath::Power(normE,gamma);
591 f.SetParameters(A,gamma);
594 double ConstantExponentialSourceFlux::pickEnergy(
double Emin,
double Emax,
double t, TRandom * rng)
const 597 TRandom * old = gRandom;
599 if (Emin == Emax)
return Emin;
602 if (f.GetXmin() != Emin || f.GetXmax() !=Emax) f.SetRange(Emin,Emax);
604 double E = f.GetRandom();
610 TimeWindowedExponentialSourceFlux::TimeWindowedExponentialSourceFlux(
double t0,
double t1,
double e,
double norm,
double normE,
double cutoff)
611 : gamma(e) , f(
"f",
"[0] * x^-[1]",1e9,cutoff ? cutoff : 1e12) , t0(t0), t1(t1) , cutoff(cutoff)
614 A = norm * TMath::Power(normE,gamma);
615 f.SetParameters(A,gamma);
618 double TimeWindowedExponentialSourceFlux::pickEnergy(
double Emin,
double Emax,
double t, TRandom * rng)
const 620 if (t < t0 || t > t1)
return 0;
621 TRandom * old = gRandom;
623 if (cutoff && Emax > cutoff) Emax = cutoff;
626 if (f.GetXmin() != Emin || f.GetXmax() !=Emax) f.SetRange(Emin,Emax);
628 double E = f.GetRandom();
638 average_nonzero_flux = 0;
640 per_source_average_flux.clear();
641 per_source_average_nonzero_flux.clear();
642 per_source_average_flux.resize(sources.size());
643 per_source_average_nonzero_flux.resize(sources.size());
644 std::vector<int> Nnonzero_per_source(sources.size());
650 TRandom * rng = getRNG(RNG_SOURCE);
651 for (
int i =0; i < N; i++)
653 double t = rng->Uniform(t0,t1);
655 for (
unsigned j = 0; j < sources.size(); j++)
657 double dFlux = minE == maxE ? sources[j]->getFlux()->getFlux(minE, t)/N : sources[j]->getFlux()->getFluxBetween(weight_Emin,weight_Emax,t);
660 average_flux += dFlux/N;
661 per_source_average_flux[j]+=dFlux/N;
662 average_nonzero_flux += dFlux;
666 Nnonzero_per_source[j]++;
667 per_source_average_nonzero_flux[j]+= dFlux;
676 for (
unsigned j = 0; j < sources.size(); j++)
678 if (per_source_average_nonzero_flux[j]) per_source_average_nonzero_flux[j] /= Nnonzero_per_source[j];
681 if (average_nonzero_flux) average_nonzero_flux /= Nnonzero;
682 printf(
"Set up weights: average flux is %g (average non-zero flux: %g)\n", average_flux, average_nonzero_flux);
685 double SourceModel::getPerSourceTimeWeight(
double t,
int i,
bool use_average_nonzero_flux)
const 687 double weight = weight_Emin == weight_Emax ? sources[i]->getFlux()->getFlux(weight_Emin,t)
688 : sources[i]->getFlux()->getFluxBetween(weight_Emin,weight_Emax,t);
690 return weight / (use_average_nonzero_flux ? per_source_average_nonzero_flux[i] : per_source_average_flux[i]);
693 double SourceModel::getTimeWeight(
double t,
bool use_average_nonzero_flux)
const 696 if (average_flux < 0)
698 if (nnag++ < 100) fprintf(stderr,
"WARNING: setUpWeights() hasn't been called yet. Will return crazy value!\n");
701 for (
unsigned i = 0; i < sources.size(); i++)
703 weight += weight_Emin == weight_Emax ? sources[i]->getFlux()->getFlux(weight_Emin,t) : sources[i]->getFlux()->getFluxBetween(weight_Emin,weight_Emax,t);
705 return weight / (use_average_nonzero_flux ? average_nonzero_flux : average_flux);
int getDirectionAndEnergy(Vector *nudir, double t, double &nuE, double minE=1e9, double maxE=1e12)
const char * ICEMC_SRC_DIR()
void computeFluxTimeChanges(std::vector< double > *changes) const
static SourceModel * getSourceModel(const char *key, Restriction r=Restriction())
void addSource(Source *source)
This class represents a three-vector. Operators are overloaded to provide for the familiar operations...
void setUpWeights(double t0, double t1, double minE=1e9, double maxE=1e12, int N=1e6)