source.cc
1 #include <map>
2 #include <string>
3 #include <iostream>
4 
5 #include "TMath.h"
6 #include "TObjString.h"
7 #include "TTimeStamp.h"
8 #include "TH1.h"
9 #include "TFile.h"
10 #include "TTree.h"
11 #include "EnvironmentVariable.h"
12 #include "Math/WrappedTF1.h"
13 #include "Math/GaussIntegrator.h"
14 
15 #include "blazars/fava.h"
16 #include "source.hh"
17 #include "icemc_random.h"
18 
19 // SourceModel
20 
21 
22 SourceModel::SourceModel(const char * n)
23  : name(n)
24 {
25  weight_Emin = 0;
26  weight_Emax = 0;
27  average_flux = -1;
28  average_nonzero_flux= -1;
29 }
30 
31 SourceModel::~SourceModel()
32 {
33  for (unsigned i = 0; i < sources.size(); i++) delete sources[i];
34 }
35 
36 
37 // convert redshift to MPC
38 static double luminosity_distance(double z)
39 {
40  static TF1 * f_comoving = 0;
41  if (!f_comoving)
42  {
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);
45  }
46 
47  double d_M = f_comoving->Integral(0,z);
48  return (1+z) * d_M;
49 }
50 
51 
52 // TXS blazar (from icecube)
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);
59 
60 
61 // M77
62 const double m77_index = 3.16;
63 
64 
65 
66 // SNe
67 const double gamma_index = 2;
68 
69 double SourceModel::Restriction::fromString(const char * the_time)
70 {
71  printf("PARSING \"%s\"\n",the_time);
72  // Start times
73  if(strcasestr(the_time, "A4_MIN"))
74  {
75  return 1480707643;
76  }
77  else if(strcasestr(the_time, "A3_MIN"))
78  {
79  return 1418938406;
80  }
81  else if(strcasestr(the_time,"A4_MAX"))
82  {
83  return 1483004563;
84  }
85  else if(strcasestr(the_time, "A3_MAX"))
86  {
87  return 1420777814;
88  }
89  else if(strcasestr(the_time,"MAX"))
90  {
91  return 2147483646;
92  }
93 
94  return atof(the_time);
95 }
96 
97 
98 
100 {
101 
102  if (!key || !strcasecmp(key,"NONE")) return NULL;
103 
104  SourceModel * m = new SourceModel(key);
105 
106  TString str(key);
107  TObjArray * tokens = str.Tokenize("+");
108 
109 
110  for (int i =0; i < tokens->GetEntries(); i++)
111  {
112  TObjString * tok = (TObjString*) tokens->At(i);
113 
114  TString stripped( tok->String().Strip( TString::kBoth));
115 
116  // The flux from TXS
117  if (stripped == "TXS")
118  {
119 
120  double txs_ra = 77.3582; // decimal degrees
121  double txs_dec = 5.69314; // decimal degrees
122 
123  m->addSource(new Source("TXS 0506+056", txs_ra/15., txs_dec,
124  new ConstantExponentialSourceFlux(txs_index,txs_flux,txs_E0) //from IceCube paper, assume it's constant over flight for now
125  )) ;
126  }
127 
128  else if (stripped == "M77")
129  {
130  m->addSource(new Source("M77", 2 + 42./60 + 40.771/3600, -47.84/3600, new ConstantExponentialSourceFlux(m77_index, txs_flux, txs_E0))); //flux and E0 don't matter
131  }
132 
133  // Supernovae
134  else if(stripped=="SN")
135  {
136  TFile *file = new TFile("./supernovae/data/supernovaeTNS.root");
137  TTree *tree = (TTree*) file->Get("tnsTree");
138 
139  // Tree vars
140  float ra, dec;
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;
145 
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);
152 
153  int two_weeks = 24*3600*14;
154  //int two_days = 24*3600*2;
155  for (unsigned int i = 0; i < tree->GetEntries(); i++)
156  {
157  tree->GetEntry(i);
158 
160  // Time cut for finding sources within a certain time period
161  // NOTE we are working with _neutrinos_ instead of photons. Neutrinos reach Earth before photons because they escape the star before the shockwaves of the supernova explosion.
162  // assume photons arrive, at max, 2 weeks after the neutrinos
163  // discoveryUnixTime is the discovery time of the supernova using photons, not neutrinos, thus:
164  if( (discoveryUnixTime < r.whichStartTime) || (discoveryUnixTime > r.whichEndTime + two_weeks) ){continue;}
165  // Declination cut (ANITA won't see neutrinos from these sources)
166  if( abs(dec)>r.maxDec){continue;}
167  // Search for specific subtype
168  if(strcasecmp(r.whichSources,"All"))
169  {
170  // Account for the chosen one
171  if(strcasestr(r.whichSources, objName->c_str()))
172  {
173  std::cout << "Using the specified source: " << objName << std::endl;
174  }
175  else{continue;}
176  }
177  // Only look at a single subtype (REMEMBER this searches for a complete string
178  // but SNII can take the formats of: SN IIa ,SN IIb, SN IIa-91bg-like, etc)
179  if(strcasecmp(r.whichSubtype,"All"))
180  {
181  if(strcasestr(r.whichSubtype,objSubtype->c_str()))
182  {
183  //std::cout << "Using the specified subtype: " << whichSubtype << std::endl;
184  }
185  else{continue;}
186  }
187 
188  m->addSource(new Source(objName->c_str(), ra/=15., dec,
189  new TimeWindowedExponentialSourceFlux(discoveryUnixTime - two_weeks, discoveryUnixTime, 2,txs_flux,txs_E0) // Just use these params as the first step
190  )) ;
191  }
192 
193  }
194 
195  else if (stripped.BeginsWith("GRB"))
196  {
197  TFile *file = new TFile("./GRBs/data/GRBIceCube.root");
198  TTree *tree = (TTree*) file->Get("iceCubeTree");
199 
200  // Tree vars
201  float RA, dec;
202  std::string *objName = new std::string;
203  int unixTriggerTime = 0;
204  float t90, fluence, redshift;
205 
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);
213 
214  for (unsigned int i = 0; i < tree->GetEntries(); i++)
215  {
216  tree->GetEntry(i);
217 
218  if (fluence < 0) fluence = 1e-6;
219 
221  // Time cut for finding sources within a certain time period
222  if( (unixTriggerTime < r.whichStartTime) || (unixTriggerTime > r.whichEndTime) ){continue;}
223  // Declination cut (ANITA won't see neutrinos from these sources)
224  if( abs(dec)>r.maxDec){continue;}
225 
226  printf("%d\n", unixTriggerTime);
227 
228  // Search for specific subtype
229  if(strcasecmp(r.whichSources,"All"))
230  {
231  // Account for the chosen one
232  if(strcasestr(r.whichSources,objName->c_str()))
233  {
234  std::cout << "Using the specified source: " << objName << std::endl;
235  }
236  else{continue;}
237  }
238 
239  // Calculate the explicit maximum neutrino energy via fireball model
240  bool modelCalculation = true;
241  double energyCutoff = 1e10;
242  if(modelCalculation==true)
243  {
244 
245  // Cosmology
246  double omegaRad = 9.236e-5;
247  double omegaM = 0.315;
248  double omegaLambda = 0.679;
249  double H0 = 67.4 * pow(10,3); // We want this in ms^-1/Mpc, not the standard kms^-1/Mpc
250  double c = 299792458; // ms^-1, as we we calculate c/H0
251  double EGammaIso = 0.;
252  double dL = 0.; // we will calculate in Mpc
253  double MpcToCm = 3.08e+24;
254 
255  // GRB model vars
256  double r = 0.;
257 
258  // Model assumptions
259  double lorentzFactor = 300.;
260  double eB = 0.1;
261  double eE = 0.1;
263  double nEx = 1.; // cm^-3
264  double EKinIso = 0.;
265  double shockLorentzFactor = 0.;
266  double energyCGamma = 0.;
267  double fPiAG = 0.;
268  double energyBreakNeutrinoAG = 0.;
269  double shockRadius = 0.;
270  double massProton = 0.0015; // erg/c^2
271  double magFieldShock = 0.;
272  double maxShockNuE = 0.;
273 
274  // GRB spectra calculation starts
275  // Some T90s are not recorded, set to 20s.
276  if(t90==-999){t90 = 20;}
277  // Some redshifts are not calculated, set to 2.
278  if(redshift==-999){redshift = 2;}
279  // Just set the fluence to this for now
280  if(fluence==-999){fluence=1e-6;}
281 
282  // Cosmology stuff
283  // Assume flat Universe
284  // Curvature is k = 0, omegaK = 0
285  TF1 Ez("E(z)", "1/(pow([0] * pow((1+x),3) + [1] + [2] * pow((1+x),4),0.5))", 0, redshift);
286  //TF1 Ez("E(z)", "1/(pow([0] * pow((1+x),3) + [1],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;
292  ig.SetFunction(wEz);
293  ig.SetRelTolerance(0.001);
294  double integral = ig.Integral(0, redshift);
295 
296  // Now calc lum dist
297  dL = (1+redshift) * c/H0 * integral; // Mpc
298 
299  // Calc gamma-ray bolometric energy
300  EGammaIso = 4 * TMath::Pi() * dL*dL * fluence * 1/(1+redshift); // Mpc^2 * erg cm^-2
301  EGammaIso*=pow(MpcToCm,2); // erg
302 
304  // Afterglow specific calcs
306 
307  // Calculate total jet kinetic energy and bulk Lorentz factor of GRB afterglow
308  EKinIso = EGammaIso/eE; //erg
309 
310  // Radii, depths
311  shockRadius = pow(((3*EKinIso)/(4*TMath::Pi() * nEx * massProton * c*c * lorentzFactor*lorentzFactor)),(1./3.)); // in cm
312  shockRadius*=100; // in meters, to compare with internal rad
313  if(shockRadius < r)
314  {
315  std::cout << "Afterglow shock radius < burst radius. Something might be wrong." << std::endl;
316  }
317 
318  shockRadius/=100;
319  magFieldShock = pow(8 * TMath::Pi() * eB * nEx * massProton,0.5);
320 
321  shockLorentzFactor = 195 * pow(((EKinIso)/(pow(10,54))),0.125) * pow(t90/10,-0.375) * pow(nEx/1.,-0.125);
322 
323  // Peak sync gamma energy radiated by electrons in B field:
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.));
325 
326  // Calculate proton-to-pion conversion factor
327  fPiAG = 0.2;
328  //fPiAG = 0.2 * pow(((EKinIso)/(pow(10,54))),(33./24.)) * pow(t90/10,(-9./8.)) * pow(nEx/1,(9./8.));
329 
330  // convert to GeV
331  energyCGamma/=1e+9;
332  // break energies
333  energyBreakNeutrinoAG = 0.015 * shockLorentzFactor * 1/(1+redshift) * pow(energyCGamma,-1); // GeV
334  // max shock neutrino energy
335  maxShockNuE = fPiAG/(4*(1+redshift)) * eE * shockLorentzFactor * magFieldShock * shockRadius;
336 
337  // normal case
338  if(maxShockNuE > energyBreakNeutrinoAG)
339  {
340  // All good
341  }
342  else if(maxShockNuE < 1e9) // account for minimum possible energy
343  {
344  maxShockNuE = 1e9;
345  }
346  else
347  {
348  maxShockNuE = energyBreakNeutrinoAG+energyBreakNeutrinoAG*0.01;
349  }
350 
351  energyCutoff = maxShockNuE;
352 
353  }
354 
355  if (strcasestr(stripped.Data(),"AFTERGLOW") )
356  {
357  m->addSource(new Source(objName->c_str(), RA/=15., dec,
358  new TimeWindowedExponentialSourceFlux(unixTriggerTime,
359  unixTriggerTime + 3600,1.5,
360  fluence,1e9, energyCutoff)
361  ));
362 
363  }
364  else if (strcasestr(stripped.Data(),"PROMPT") )
365  {
366  m->addSource(new Source(objName->c_str(), RA/=15., dec,
367  //fireball model for E > 1e18 eV, I think
368  new TimeWindowedExponentialSourceFlux(unixTriggerTime,
369  unixTriggerTime + t90,4,
370  fluence,1e9)
371  )) ;
372 
373  }
374 
375  // conservative case... 5 minutes before, one hour after, use afterglow flux.
376  else
377  {
378  m->addSource(new Source(objName->c_str(), RA/=15., dec,
379  //fireball model for E > 1e18 eV, I think
380  new TimeWindowedExponentialSourceFlux(unixTriggerTime-300,
381  unixTriggerTime + 3600,1.5,
382  fluence,1e9, energyCutoff)
383  )) ;
384 
385  }
386 
387  }
388  }
389 
390  else if (stripped.BeginsWith("FAVA"))
391  {
392  // Let's load all the blazars from FAVA that occurred during a flight
393 
394  bool flux_weighted = strcasestr(stripped.Data(),"FLUX");
395  bool z_weighted = strcasestr(stripped.Data(),"REDSHIFT");
396 
397  const char * dir = EnvironmentVariable::ICEMC_SRC_DIR() ?: ".";
398 
399  TFile ffava(Form("%s/blazars/fava.root", dir));
400  TTree * fava_tree = (TTree*) ffava.Get("fava");
401  FAVAEntry *fava = 0;
402 
403  fava_tree->SetBranchAddress("fava",&fava);
404 
405  for (int i = 0; i < fava_tree->GetEntries(); i++)
406  {
407  fava_tree->GetEntry(i);
408 
409  // time cut
410  if( (fava->unix_tmax < r.whichStartTime) || (fava->unix_tmin > r.whichEndTime) ){continue;}
411 
412  const char * fava_name = fava->association.GetString().Data();
413 
414  // Skip flares that are not associated with a source
415  if(strcasestr(fava_name,"None"))
416  {
417  continue;
418  }
419 
420  // If we didn't specific all sources
421 
422  if(strcasecmp(r.whichSources,"All"))
423  {
424  // Account for the chosen one
425  if(strcasestr(r.whichSources,fava_name))
426  {
427  std::cout << "Using the specified source: " << fava_name << std::endl;
428  }
429  else{continue;}
430  }
431 
432 
433  //say no to |dec| >30
434  //if (fabs(fava->dec) > r.maxDec) continue;
435 
436  //say no to low HE flux
437 
438  //if (fava->he_sigma < 4 && fava->he_flux < 0) continue;
439  //printf("passed flux cut \n");
440 
441  //only blazars
442  if (fava->source_class.GetString() == "bcu" || fava->source_class.GetString() == "fsrq" || fava->source_class.GetString() == "bll")
443  {
444 
445  m->addSource(new Source(fava_name, (fava->ra)/15., fava->dec,
446  new TimeWindowedExponentialSourceFlux( fava->unix_tmin, fava->unix_tmax, txs_index,
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))
448  );
449 
450  }
451  }
452  }
453  }
454 
455  delete tokens;
456 
457  std::cout << "----------------------" << std::endl;
458 
459  if (!m->getNSources())
460  {
461  fprintf(stderr,"WARNING: no sources added for key %s\n", key);
462  delete m;
463  m = NULL;
464  }
465  else
466  {
467  std::cout << "Sources added: " << m->getNSources() << std::endl;
468  }
469 
470  return m;
471 }
472 
473 TH1 * SourceModel::estimateFlux(double tmin, double tmax, double Emin, double Emax, int nbins, int N)
474 {
475  double from = log10(Emin);
476  double to = log10(Emax);
477 
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++)
483  {
484  double t = rng->Uniform(tmin,tmax);
485  for (unsigned j = 0; j < sources.size(); j++)
486  {
487  //figure out the flux between each energy bin
488  for (int E = 1; E<= nbins; E++)
489  {
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));
493  }
494  }
495  }
496 
497  spectrum->Scale(1./N);
498  return spectrum;
499 }
500 
501 
502 void SourceModel::computeFluxTimeChanges(std::vector<double> * changes) const
503 {
504 
505  changes->clear();
506 
507  for (unsigned i = 0; i < sources.size(); i++)
508  {
509  sources[i]->getFlux()->getFluxTimeChanges(changes);
510  }
511 
512  //sort and remove dupes
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));
516 }
517 
518 int SourceModel::getDirectionAndEnergy( Vector * nudir, double t, double & nuE, double minE, double maxE)
519 {
520  if (!sources.size()) return -1;
521 
522  bool fixedEnergy = minE==maxE;
523  if (fixedEnergy) nuE = minE;
524 
525  //pick a random source, weighted by the flux
526  double total_flux = 0 ;
527  std::vector<double> fluxes(sources.size());
528 
529  for (unsigned i = 0; i < sources.size(); i++)
530  {
531  double f = fixedEnergy ? sources[i]->getFlux()->getFlux(nuE,t) : sources[i]->getFlux()->getFluxBetween(minE,maxE,t);
532  total_flux += f;
533  fluxes[i] = total_flux;
534  }
535 
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();
539 
540 
541  // printf("random: %g total_flux%g, index:%u \n",random, total_flux, index);
542  if (total_flux == 0)
543  {
544  nuE = minE; // do something...
545  return -1;
546  }
547  // This is the chosen source, so we should retain some info about it
548  // Then we can easily access info about *each* simulated neutrino's origin
549  const Source * which = sources[index];
550  //std::cout << objName << std::endl;
551  //std::cout << RA << std::endl;
552  //std::cout << dec << std::endl;
553 
554  if (nudir) *nudir = which->getDirection(t);
555 
556  if (!fixedEnergy) nuE = which->getFlux()->pickEnergy(minE,maxE,t,rng);
557 
558  return index;
559 }
560 
561 
562 Source::Source(const char * nm, double ra, double dc, SourceFlux * f)
563 : name(nm), flux(f) , RA(ra), dec(dc*TMath::DegToRad())
564 {
565 }
566 
567 Vector Source::getDirection(double t) const
568 {
569 
570 #if ROOT_VERSION_CODE < ROOT_VERSION(6,0,0)
571  std::cerr << "ROOT 5 doesn't have AsLMST. Returning nonsense." << std::endl;
572  return Vector(1,0,0);
573 #else
574 
575  time_t secs = t;
576  int nsecs = 1e9 *(t-secs);
577  TTimeStamp ts(secs, nsecs);
578  double lst = ts.AsLMST(-90);
579  double h = lst - RA;
580  h *= (TMath::Pi()/12);
581  return Vector(cos(h) * cos(dec), sin(h) * cos(dec), sin(dec));
582 #endif
583 }
584 
585 
586 ConstantExponentialSourceFlux::ConstantExponentialSourceFlux(double e, double norm, double normE)
587  : gamma(e) , f("f", "[0] * x^-[1]",1e9,1e12)
588 {
589  //figure out what A is
590  A = norm * TMath::Power(normE,gamma);
591  f.SetParameters(A,gamma);
592 }
593 
594 double ConstantExponentialSourceFlux::pickEnergy(double Emin, double Emax, double t, TRandom * rng) const
595 {
596  (void) t;
597  TRandom * old = gRandom;
598  gRandom = rng;
599  if (Emin == Emax) return Emin;
600 
601  //I'm too bad at math to figure out the proper quantile function here. It's probably trivial, but this works and isn't THAT slow.
602  if (f.GetXmin() != Emin || f.GetXmax() !=Emax) f.SetRange(Emin,Emax);
603 
604  double E = f.GetRandom();
605  gRandom = old;
606  return E;
607 }
608 
609 
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)
612 {
613  //figure out what A is
614  A = norm * TMath::Power(normE,gamma);
615  f.SetParameters(A,gamma);
616 }
617 
618 double TimeWindowedExponentialSourceFlux::pickEnergy(double Emin, double Emax, double t, TRandom * rng) const
619 {
620  if (t < t0 || t > t1) return 0;
621  TRandom * old = gRandom;
622  gRandom = rng;
623  if (cutoff && Emax > cutoff) Emax = cutoff;
624 
625  //I'm too bad at math to figure out the proper quantile function here. It's probably trivial, but this works and isn't THAT slow.
626  if (f.GetXmin() != Emin || f.GetXmax() !=Emax) f.SetRange(Emin,Emax);
627 
628  double E = f.GetRandom();
629  gRandom = old;
630  return E;
631 }
632 
633 
634 void SourceModel::setUpWeights(double t0, double t1, double minE, double maxE, int N)
635 {
636 
637  average_flux = 0;
638  average_nonzero_flux = 0;
639 
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());
645  int Nnonzero = 0;
646 
647  weight_Emin = minE;
648  weight_Emax = maxE;
649 
650  TRandom * rng = getRNG(RNG_SOURCE);
651  for (int i =0; i < N; i++)
652  {
653  double t = rng->Uniform(t0,t1);
654  double sumFlux = 0;
655  for (unsigned j = 0; j < sources.size(); j++)
656  {
657  double dFlux = minE == maxE ? sources[j]->getFlux()->getFlux(minE, t)/N : sources[j]->getFlux()->getFluxBetween(weight_Emin,weight_Emax,t);
658  sumFlux += dFlux;
659 
660  average_flux += dFlux/N;
661  per_source_average_flux[j]+=dFlux/N;
662  average_nonzero_flux += dFlux;
663 
664  if (dFlux!=0)
665  {
666  Nnonzero_per_source[j]++;
667  per_source_average_nonzero_flux[j]+= dFlux;
668  }
669  }
670  if (sumFlux != 0)
671  {
672  Nnonzero++;
673  }
674  }
675 
676  for (unsigned j = 0; j < sources.size(); j++)
677  {
678  if (per_source_average_nonzero_flux[j]) per_source_average_nonzero_flux[j] /= Nnonzero_per_source[j];
679  }
680 
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);
683 }
684 
685 double SourceModel::getPerSourceTimeWeight(double t, int i, bool use_average_nonzero_flux) const
686 {
687  double weight = weight_Emin == weight_Emax ? sources[i]->getFlux()->getFlux(weight_Emin,t)
688  : sources[i]->getFlux()->getFluxBetween(weight_Emin,weight_Emax,t);
689 
690  return weight / (use_average_nonzero_flux ? per_source_average_nonzero_flux[i] : per_source_average_flux[i]);
691 }
692 
693 double SourceModel::getTimeWeight(double t, bool use_average_nonzero_flux) const
694 {
695  static int nnag = 0;
696  if (average_flux < 0)
697  {
698  if (nnag++ < 100) fprintf(stderr,"WARNING: setUpWeights() hasn't been called yet. Will return crazy value!\n");
699  }
700  double weight = 0;
701  for (unsigned i = 0; i < sources.size(); i++)
702  {
703  weight += weight_Emin == weight_Emax ? sources[i]->getFlux()->getFlux(weight_Emin,t) : sources[i]->getFlux()->getFluxBetween(weight_Emin,weight_Emax,t);
704  }
705  return weight / (use_average_nonzero_flux ? average_nonzero_flux : average_flux);
706 }
707 
int getDirectionAndEnergy(Vector *nudir, double t, double &nuE, double minE=1e9, double maxE=1e12)
Definition: source.cc:518
const char * ICEMC_SRC_DIR()
void computeFluxTimeChanges(std::vector< double > *changes) const
Definition: source.cc:502
static SourceModel * getSourceModel(const char *key, Restriction r=Restriction())
Definition: source.cc:99
void addSource(Source *source)
Definition: source.hh:65
This class represents a three-vector. Operators are overloaded to provide for the familiar operations...
Definition: vector.hh:27
void setUpWeights(double t0, double t1, double minE=1e9, double maxE=1e12, int N=1e6)
Definition: source.cc:634