gambit is hosted by Hepforge, IPPP Durham
GAMBIT  v1.5.0-2191-ga4742ac
a Global And Modular Bsm Inference Tool
DarkBit_standalone_WIMP.cpp
Go to the documentation of this file.
1 // GAMBIT: Global and Modular BSM Inference Tool
2 // *********************************************
22 
23 #include <iostream>
24 #include <fstream>
25 
30 
31 #include <boost/multi_array.hpp>
32 
33 using namespace DarkBit::Functown; // Functors wrapping the module's actual module functions
34 using namespace BackendIniBit::Functown; // Functors wrapping the backend initialisation functions
35 
36 QUICK_FUNCTION(DarkBit, TH_ProcessCatalog, OLD_CAPABILITY, TH_ProcessCatalog_WIMP, TH_ProcessCatalog, ())
38 QUICK_FUNCTION(DarkBit, DarkMatterConj_ID, OLD_CAPABILITY, DarkMatterConj_ID_WIMP, std::string, ())
40 
41 
42 void dump_array_to_file(const std::string & filename, const
43  boost::multi_array<double, 2> & a, const std::vector<double> & x, const
44  std::vector<double> & y)
45 {
46  std::fstream file;
47  file.open(filename, std::ios_base::out);
48  file << "0.0 ";
49  for (size_t i = 0; i < x.size(); i++)
50  file << x[i] << " ";
51  file << std::endl;
52  for (size_t j = 0; j < y.size(); j++)
53  {
54  file << y[j] << " ";
55  for (size_t i = 0; i < x.size(); i++)
56  {
57  file << a[i][j] << " ";
58  }
59  file << std::endl;
60  }
61  file.close();
62 }
63 
64 void dumpSpectrum(std::string filename, double mWIMP, double sv, std::vector<double> brList, double mPhi = -1)
65 {
66  DarkMatter_ID_WIMP.reset_and_calculate();
67  DarkMatterConj_ID_WIMP.reset_and_calculate();
68  TH_ProcessCatalog_WIMP.setOption<std::vector<double>>("brList", brList);
69  TH_ProcessCatalog_WIMP.setOption<double>("mWIMP", mWIMP);
70  TH_ProcessCatalog_WIMP.setOption<double>("sv", sv);
71  if (mPhi != -1)
72  TH_ProcessCatalog_WIMP.setOption<double>("mPhi", mPhi);
73  TH_ProcessCatalog_WIMP.reset_and_calculate();
74  DM_process_from_ProcessCatalog.reset_and_calculate();
75  RD_fraction_one.reset_and_calculate();
76  SimYieldTable_DarkSUSY.reset_and_calculate();
77  SimYieldTable_MicrOmegas.reset_and_calculate();
78  GA_missingFinalStates.reset_and_calculate();
79  cascadeMC_FinalStates.reset_and_calculate();
80  cascadeMC_DecayTable.reset_and_calculate();
81  cascadeMC_LoopManager.reset_and_calculate();
82  cascadeMC_gammaSpectra.reset_and_calculate();
83  GA_AnnYield_General.reset_and_calculate();
84  dump_GammaSpectrum.setOption<std::string>("filename", filename);
85  dump_GammaSpectrum.reset_and_calculate();
86 }
87 
88 // ---- Set up basic internal structures for direct & indirect detection ----
89 
90 namespace Gambit
91 {
92  namespace DarkBit
93  {
94 
96  {
97  using namespace Pipes::TH_ProcessCatalog_WIMP;
98  using std::vector;
99  using std::string;
100 
101  // Initialize empty catalog and main annihilation process
102  TH_ProcessCatalog catalog;
103  TH_Process process_ann("WIMP", "WIMP");
104  TH_Process process_dec("phi");
105  TH_Process process_dec1("phi1");
106  TH_Process process_dec2("phi2");
107 
109  // Import particle masses and couplings
111 
112  #define addParticle(Name, Mass, spinX2) \
113  catalog.particleProperties.insert(std::pair<string, TH_ParticleProperty> \
114  (Name , TH_ParticleProperty(Mass, spinX2)));
115 
117  double mWIMP = runOptions->getValue<double>("mWIMP");
119  double sv = runOptions->getValue<double>("sv");
120  double b = 0; // defined as sv(v) = sv(v=0) + b*(sv=0)*v**2
122  auto brList = runOptions->getValue<std::vector<double>>("brList");
124  double mPhi = runOptions->getValueOrDef<double>(59.0, "mPhi");
125 
126  addParticle("gamma", 0.0, 2)
127  addParticle("Z0", 91.2, 2)
128  addParticle("W+", 80.39, 2)
129  addParticle("W-", 80.39, 2)
130  addParticle("e+_3", 1.8, 1)
131  addParticle("e-_3", 1.8, 1)
132  addParticle("e+_1", 0.00051, 1)
133  addParticle("e-_1", 0.00051, 1)
134  addParticle("b", 4.9, 1)
135  addParticle("bbar", 4.9, 1)
136  addParticle("d_3", 4.9, 1)
137  addParticle("dbar_3", 4.9, 1)
138 
139  addParticle("WIMP", mWIMP, 0)
140  addParticle("phi", mPhi, 0)
141  addParticle("phi1", 100., 0)
142  addParticle("phi2", 100., 0)
143  #undef addParticle
144 
145  TH_Channel dec_channel(daFunk::vec<string>("gamma", "gamma"), daFunk::cnst(1.));
146  process_dec.channelList.push_back(dec_channel);
147 
148  TH_Channel dec_channel1(daFunk::vec<string>("e+_3", "e-_3"), daFunk::cnst(1.));
149  process_dec1.channelList.push_back(dec_channel1);
150 
151  TH_Channel dec_channel2(daFunk::vec<string>("d_3", "dbar_3"), daFunk::cnst(1.));
152  process_dec2.channelList.push_back(dec_channel2);
153 
154  process_ann.resonances_thresholds.threshold_energy.push_back(2*mWIMP);
155  auto p1 = daFunk::vec<string>("d_3", "gamma", "gamma", "e-_3", "W-", "e-_1", "phi");
156  auto p2 = daFunk::vec<string>("dbar_3", "Z0", "gamma", "e+_3", "W+", "e+_1", "phi2");
157  {
158  for ( unsigned int i = 0; i < brList.size()-1; i++ )
159  {
160  double mtot_final =
161  catalog.getParticleProperty(p1[i]).mass +
162  catalog.getParticleProperty(p2[i]).mass;
163  if ( mWIMP*2 > mtot_final && brList[i]!= 0.)
164  {
165  // std::cout << p1[i] << " " << p2[i] << " " << brList[i] << std::endl;
166  daFunk::Funk kinematicFunction = (daFunk::one("v")+pow(daFunk::var("v"), 2)*b)*sv*brList[i];
167  TH_Channel new_channel(
168  daFunk::vec<string>(p1[i], p2[i]), kinematicFunction
169  );
170  process_ann.channelList.push_back(new_channel);
171  }
172  else
173  {
175  push_back(mtot_final);
176  }
177  }
178  }
179 
180  if ( brList[7] > 0. )
181  {
182  auto E = daFunk::var("E");
183  // Note:: The below is an arbitrary form of the differential section for demonstration purposes
184  daFunk::Funk kinematicFunction = daFunk::one("v", "E1")/(pow(E-50, 4)+1)*sv*brList[7];
185  // Note: In the three body final states, the gamma yield from AnnYield currently is just the contribution
186  // from the first particle in the list (here the photon):
187  TH_Channel new_channel(daFunk::vec<string>("gamma", "e+_1", "e-_1"), kinematicFunction);
188  process_ann.channelList.push_back(new_channel);
189  }
190 
191  catalog.processList.push_back(process_ann);
192  catalog.processList.push_back(process_dec);
193  catalog.processList.push_back(process_dec1);
194  catalog.processList.push_back(process_dec2);
195 
196  catalog.validate();
197 
198  result = catalog;
199  } // function TH_ProcessCatalog_WIMP
200 
201  // Identifier for DM particle
202  void DarkMatter_ID_WIMP(std::string& result)
203  {
204  result = "WIMP";
205  }
206 
207  // Identifier for conjugate DM particle
208  void DarkMatterConj_ID_WIMP(std::string& result)
209  {
210  result = "WIMP";
211  }
212 
214  {
215  using namespace Pipes::DD_couplings_WIMP;
217  result.gps = runOptions->getValueOrDef<double>(0., "gps");
219  result.gns = runOptions->getValueOrDef<double>(0., "gns");
221  result.gpa = runOptions->getValueOrDef<double>(0., "gpa");
223  result.gna = runOptions->getValueOrDef<double>(0., "gna");
224  //std::cout << "DD_coupling says" << std::endl;
225  //std::cout << result.gps << std::endl;
226  }
227  }
228 }
229 
230 int main(int argc, char* argv[])
231 {
232  std::cout << std::endl;
233  std::cout << "Welcome to the DarkBit Generic WIMP standalone program!" << std::endl;
234  std::cout << std::endl;
235  std::cout << "**************************************************************************************" << std::endl;
236  std::cout << "This standalone example demonstrates how to calculate a range of observables and " << std::endl;
237  std::cout << "likelihoods for a generic WIMP model defined by the WIMP mass and an annihilation (or " << std::endl;
238  std::cout << "scattering) cross section. The model also contains three scalar particles which decay:" << std::endl;
239  std::cout << "phi -> gamma gamma phi_1 -> tau+ tau- phi_2 -> b bbar" << std::endl;
240  std::cout << std::endl;
241  std::cout << "Usage: DarkBit_standalone_WIMP mode" << std::endl;
242  std::cout << std::endl;
243  std::cout << "Mode Options: " << std::endl;
244  std::cout << " 0: Outputs spectrum of gamma rays from WIMP annihilation to b bbar (dPhi_dE0.dat)" << std::endl;
245  std::cout << " 1: Outputs spectrum of gamma rays from WIMP annihilation to gamma Z_0 (dPhi_dE1.dat)" << std::endl;
246  std::cout << " 2: Outputs spectrum of gamma rays from WIMP annihilation to gamma gamma (dPhi_dE2.dat)" << std::endl;
247  std::cout << " 3: Outputs spectrum of gamma rays from WIMP annihilation to tau+ tau- (dPhi_dE3.dat)" << std::endl;
248  std::cout << " 4: Outputs spectrum of gamma rays from WIMP annihilation to W+ W- (dPhi_dE4.dat)" << std::endl;
249  std::cout << " 5: Outputs spectrum of gamma rays from WIMP annihilation to gamma e+ e- " << std::endl;
250  std::cout << " (dPhi_dE5.dat)" << std::endl;
251  std::cout << " 6: Outputs tables of gamma-ray likelihoods and the relic density" << std::endl;
252  std::cout << " in <sigma v> / m_WIMP parameter space." << std::endl;
253  std::cout << " 7: Outputs tables of direct detection likelihoods in sigma / m_WIMP parameter" << std::endl;
254  std::cout << " space." << std::endl;
255  std::cout << " >=10: Outputs spectrum of gamma rays from WIMP annihilation to phi phi_2. The" << std::endl;
256  std::cout << " mode value is m_phi while m_phi_2=100 GeV (dPhi_dE_FCMC_(mode).dat)" << std::endl;
257  std::cout << " N.B. Here dPhi/dE = sigma v / m_chi^2 * dN/dE" << std::endl;
258  std::cout << "**************************************************************************************" << std::endl;
259  std::cout << std::endl;
260 
261  try
262  {
263  if (argc==1)
264  {
265  std::cout << "Please select test mode>=0" << std::endl;
266  exit(1);
267  }
268  int mode = std::stoi((std::string)argv[1]);
269  std::cout << "Starting with mode " << mode << std::endl;
270 
271 
272  // ---- Initialise logging and exceptions ----
273 
274  initialise_standalone_logs("runs/DarkBit_standalone_WIMP/logs/");
275  logger()<<"Running DarkBit standalone example"<<LogTags::info<<EOM;
276  model_warning().set_fatal(true);
277 
278 
279  // ---- Check that required backends are present ----
280 
281  if (not Backends::backendInfo().works["DarkSUSY_generic_wimp6.2.5"]) backend_error().raise(LOCAL_INFO, "DarkSUSY_generic_wimp_6.2.5 is missing!");
282  if (not Backends::backendInfo().works["gamLike1.0.1"]) backend_error().raise(LOCAL_INFO, "gamLike 1.0.1 is missing!");
283  if (not Backends::backendInfo().works["DDCalc2.2.0"]) backend_error().raise(LOCAL_INFO, "DDCalc 2.2.0 is missing!");
284  if (not Backends::backendInfo().works["MicrOmegas_MSSM3.6.9.2"]) backend_error().raise(LOCAL_INFO, "MicrOmegas 3.6.9.2 for MSSM is missing!");
285 
286  // ---- Initialize models ----
287 
288  // Initialize halo model
289  ModelParameters* Halo_primary_parameters = Models::Halo_Einasto::Functown::primary_parameters.getcontentsPtr();
290  Halo_primary_parameters->setValue("vrot", 235.); // Local properties
291  Halo_primary_parameters->setValue("v0", 235.);
292  Halo_primary_parameters->setValue("vesc", 550.);
293  Halo_primary_parameters->setValue("rho0", 0.4);
294  Halo_primary_parameters->setValue("r_sun", 8.5);
295 
296  Halo_primary_parameters->setValue("rs", 20.); // Global properties
297  Halo_primary_parameters->setValue("rhos", 0.08);
298  Halo_primary_parameters->setValue("alpha", 0.17);
299 
300 
301  // --- Resolve halo dependencies ---
302  ExtractLocalMaxwellianHalo.notifyOfModel("Halo_Einasto");
303  ExtractLocalMaxwellianHalo.resolveDependency(&Models::Halo_Einasto::Functown::primary_parameters);
304  ExtractLocalMaxwellianHalo.reset_and_calculate();
305 
306  GalacticHalo_Einasto.notifyOfModel("Halo_Einasto");
307  GalacticHalo_Einasto.resolveDependency(&Models::Halo_Einasto::Functown::primary_parameters);
308  GalacticHalo_Einasto.reset_and_calculate();
309 
310  // ---- Initialize backends ----
311 
312  // Assume for direct and indirect detection likelihoods that dark matter
313  // density is always the measured one (despite relic density results)
314  RD_fraction_one.reset_and_calculate();
315 
316  // Set up DDCalc backend initialization
317  Backends::DDCalc_2_2_0::Functown::DDCalc_CalcRates_simple.setStatus(2);
318  Backends::DDCalc_2_2_0::Functown::DDCalc_Experiment.setStatus(2);
319  Backends::DDCalc_2_2_0::Functown::DDCalc_LogLikelihood.setStatus(2);
320  DDCalc_2_2_0_init.resolveDependency(&ExtractLocalMaxwellianHalo);
321  // Assume for direct and indirect detection likelihoods that dark matter
322  // density is always the measured one (despite relic density results)
323  DDCalc_2_2_0_init.resolveDependency(&RD_fraction_one);
324  DDCalc_2_2_0_init.resolveDependency(&mwimp_generic);
325  DDCalc_2_2_0_init.resolveDependency(&DD_couplings_WIMP);
326 
327  // Initialize gamLike backend
328  gamLike_1_0_1_init.reset_and_calculate();
329 
330  // Initialize DarkSUSY backend
331  DarkSUSY_generic_wimp_6_2_5_init.reset_and_calculate();
332 
333  // Initialize MicrOmegas backend
334  // The below allows us to initialise MicrOmegas_MSSM without a particular MSSM model.
335  MicrOmegas_MSSM_3_6_9_2_init.notifyOfModel("Halo_Einasto");
336  MicrOmegas_MSSM_3_6_9_2_init.reset_and_calculate();
337 
338  // ---- Gamma-ray yields ----
339 
340  // Initialize tabulated gamma-ray yields
341  SimYieldTable_DarkSUSY.resolveBackendReq(&Backends::DarkSUSY_generic_wimp_6_2_5::Functown::dsanyield_sim);
342  SimYieldTable_MicrOmegas.resolveBackendReq(&Backends::MicrOmegas_MSSM_3_6_9_2::Functown::dNdE);
343  SimYieldTable_DarkSUSY.setOption<bool>("allow_yield_extrapolation", true);
344  SimYieldTable_MicrOmegas.setOption<bool>("allow_yield_extrapolation", true);
345 
346  // Select SimYieldTable
347  //auto SimYieldTablePointer = &SimYieldTable_MicrOmegas;
348  auto SimYieldTablePointer = &SimYieldTable_DarkSUSY;
349 
350  // Identify process as annihilation rather than decay
353 
354  // Collect missing final states for simulation in cascade MC
355  GA_missingFinalStates.resolveDependency(&TH_ProcessCatalog_WIMP);
356  GA_missingFinalStates.resolveDependency(SimYieldTablePointer);
357  GA_missingFinalStates.resolveDependency(&DarkMatter_ID_WIMP);
358  GA_missingFinalStates.resolveDependency(&DarkMatterConj_ID_WIMP);
360 
361  // Infer for which type of final states particles MC should be performed
362  cascadeMC_FinalStates.setOption<std::vector<std::string>>("cMC_finalStates", daFunk::vec((std::string)"gamma"));
363 
364  // Collect decay information for cascade MC
365  cascadeMC_DecayTable.resolveDependency(&TH_ProcessCatalog_WIMP);
366  cascadeMC_DecayTable.resolveDependency(SimYieldTablePointer);
367 
368  // Set up MC loop manager for cascade MC
369  cascadeMC_LoopManager.setOption<int>("cMC_maxEvents", 20000);
370  cascadeMC_Histograms.setOption<double>("cMC_endCheckFrequency", 25);
371  cascadeMC_Histograms.setOption<double>("cMC_gammaRelError", .05);
372  cascadeMC_Histograms.setOption<int>("cMC_numSpecSamples", 25);
373  cascadeMC_Histograms.setOption<int>("cMC_NhistBins", 300);
374  cascadeMC_LoopManager.resolveDependency(&GA_missingFinalStates);
375  std::vector<functor*> nested_functions = initVector<functor*>(
377  cascadeMC_LoopManager.setNestedList(nested_functions);
378 
379  // Set up initial state for cascade MC step
380  cascadeMC_InitialState.resolveDependency(&GA_missingFinalStates);
381  cascadeMC_InitialState.resolveLoopManager(&cascadeMC_LoopManager);
382  //cascadeMC_InitialState.reset_and_calculate();
383 
384  // Perform MC step for cascade MC
385  cascadeMC_GenerateChain.resolveDependency(&cascadeMC_InitialState);
386  cascadeMC_GenerateChain.resolveDependency(&cascadeMC_DecayTable);
387  cascadeMC_GenerateChain.resolveLoopManager(&cascadeMC_LoopManager);
388  //cascadeMC_GenerateChain.reset_and_calculate();
389 
390  // Generate histogram for cascade MC
391  cascadeMC_Histograms.resolveDependency(&cascadeMC_InitialState);
392  cascadeMC_Histograms.resolveDependency(&cascadeMC_GenerateChain);
393  cascadeMC_Histograms.resolveDependency(&TH_ProcessCatalog_WIMP);
394  cascadeMC_Histograms.resolveDependency(SimYieldTablePointer);
395  cascadeMC_Histograms.resolveDependency(&cascadeMC_FinalStates);
396  cascadeMC_Histograms.resolveLoopManager(&cascadeMC_LoopManager);
397  //cascadeMC_Histograms.reset_and_calculate();
398 
399  // Check convergence of cascade MC
400  cascadeMC_EventCount.resolveDependency(&cascadeMC_InitialState);
401  cascadeMC_EventCount.resolveLoopManager(&cascadeMC_LoopManager);
402  //cascadeMC_EventCount.reset_and_calculate();
403 
404  // Start cascade MC loop
405 
406  // Infer gamma-ray spectra for recorded MC results
407  cascadeMC_gammaSpectra.resolveDependency(&GA_missingFinalStates);
408  cascadeMC_gammaSpectra.resolveDependency(&cascadeMC_FinalStates);
409  cascadeMC_gammaSpectra.resolveDependency(&cascadeMC_Histograms);
410  cascadeMC_gammaSpectra.resolveDependency(&cascadeMC_EventCount);
411 
412  // Calculate total gamma-ray yield (cascade MC + tabulated results)
413  GA_AnnYield_General.resolveDependency(&TH_ProcessCatalog_WIMP);
414  GA_AnnYield_General.resolveDependency(SimYieldTablePointer);
415  GA_AnnYield_General.resolveDependency(&DarkMatter_ID_WIMP);
416  GA_AnnYield_General.resolveDependency(&DarkMatterConj_ID_WIMP);
417  GA_AnnYield_General.resolveDependency(&cascadeMC_gammaSpectra);
418 
419  dump_GammaSpectrum.resolveDependency(&GA_AnnYield_General);
420 
421  // Resolve Galactic halo requirements for gamLike
422  set_gamLike_GC_halo.resolveDependency(&GalacticHalo_Einasto);
423  set_gamLike_GC_halo.resolveBackendReq(&Backends::gamLike_1_0_1::Functown::set_halo_profile);
424 
425  // Calculate Fermi LAT dwarf likelihood
427  // Assume for direct and indirect detection likelihoods that dark matter
428  // density is always the measured one (despite relic density results)
429  lnL_FermiLATdwarfs_gamLike.resolveDependency(&RD_fraction_one);
430  lnL_FermiLATdwarfs_gamLike.resolveBackendReq(&Backends::gamLike_1_0_1::Functown::lnL);
431 
432  lnL_HESSGC_gamLike.resolveDependency(&GA_AnnYield_General);
433  lnL_HESSGC_gamLike.resolveDependency(&RD_fraction_one);
434  lnL_HESSGC_gamLike.resolveBackendReq(&Backends::gamLike_1_0_1::Functown::lnL);
435 
436  lnL_CTAGC_gamLike.resolveDependency(&GA_AnnYield_General);
437  lnL_CTAGC_gamLike.resolveDependency(&RD_fraction_one);
438  lnL_CTAGC_gamLike.resolveBackendReq(&Backends::gamLike_1_0_1::Functown::lnL);
439 
440  lnL_FermiGC_gamLike.resolveDependency(&GA_AnnYield_General);
441  lnL_FermiGC_gamLike.resolveDependency(&RD_fraction_one);
442  lnL_FermiGC_gamLike.resolveBackendReq(&Backends::gamLike_1_0_1::Functown::lnL);
443 
444 
445  // -- Calculate relic density --
446  // *any* of the models listed by "ALLOW_MODELS" in DarkBit_rollcall.hpp will work here
447  RD_eff_annrate_from_ProcessCatalog.notifyOfModel("ScalarSingletDM_Z2");
448  RD_eff_annrate_from_ProcessCatalog.resolveDependency(&TH_ProcessCatalog_WIMP);
449  RD_eff_annrate_from_ProcessCatalog.resolveDependency(&DarkMatter_ID_WIMP);
450  RD_eff_annrate_from_ProcessCatalog.resolveDependency(&DarkMatterConj_ID_WIMP);
451 
453  RD_eff_annrate_from_ProcessCatalog.resolveDependency(&DarkMatter_ID_WIMP);
454  RD_eff_annrate_from_ProcessCatalog.resolveDependency(&DarkMatterConj_ID_WIMP);
455 
457 
458 
459  RD_oh2_DS_general.resolveDependency(&RD_spectrum_ordered_func);
460  RD_oh2_DS_general.resolveDependency(&RD_eff_annrate_from_ProcessCatalog);
461  RD_oh2_DS_general.resolveBackendReq(&Backends::DarkSUSY_generic_wimp_6_2_5::Functown::rdpars);
462  RD_oh2_DS_general.resolveBackendReq(&Backends::DarkSUSY_generic_wimp_6_2_5::Functown::rdtime);
463  RD_oh2_DS_general.resolveBackendReq(&Backends::DarkSUSY_generic_wimp_6_2_5::Functown::dsrdcom);
464  RD_oh2_DS_general.resolveBackendReq(&Backends::DarkSUSY_generic_wimp_6_2_5::Functown::dsrdstart);
465  RD_oh2_DS_general.resolveBackendReq(&Backends::DarkSUSY_generic_wimp_6_2_5::Functown::dsrdens);
466 
467 
468  // ---- Calculate direct detection constraints ----
469 
470  // Calculate direct detection rates for LZ, PandaX 2017, Xenon 1T and PICO-60
471  LZ_Calc.resolveBackendReq(&Backends::DDCalc_2_2_0::Functown::DDCalc_Experiment);
472  LZ_Calc.resolveBackendReq(&Backends::DDCalc_2_2_0::Functown::DDCalc_CalcRates_simple);
473  PandaX_2017_Calc.resolveBackendReq(&Backends::DDCalc_2_2_0::Functown::DDCalc_Experiment);
474  PandaX_2017_Calc.resolveBackendReq(&Backends::DDCalc_2_2_0::Functown::DDCalc_CalcRates_simple);
475  PICO_60_2017_Calc.resolveBackendReq(&Backends::DDCalc_2_2_0::Functown::DDCalc_Experiment);
476  PICO_60_2017_Calc.resolveBackendReq(&Backends::DDCalc_2_2_0::Functown::DDCalc_CalcRates_simple);
477  XENON1T_2017_Calc.resolveBackendReq(&Backends::DDCalc_2_2_0::Functown::DDCalc_Experiment);
478  XENON1T_2017_Calc.resolveBackendReq(&Backends::DDCalc_2_2_0::Functown::DDCalc_CalcRates_simple);
479 
480  // Calculate direct detection likelihood for LZ, PandaX 2017, Xenon 1T and PICO-60
481  LZ_GetLogLikelihood.resolveDependency(&LZ_Calc);
482  LZ_GetLogLikelihood.resolveBackendReq(&Backends::DDCalc_2_2_0::Functown::DDCalc_Experiment);
483  LZ_GetLogLikelihood.resolveBackendReq(&Backends::DDCalc_2_2_0::Functown::DDCalc_LogLikelihood);
484  PandaX_2017_GetLogLikelihood.resolveDependency(&PandaX_2017_Calc);
485  PandaX_2017_GetLogLikelihood.resolveBackendReq(&Backends::DDCalc_2_2_0::Functown::DDCalc_Experiment);
486  PandaX_2017_GetLogLikelihood.resolveBackendReq(&Backends::DDCalc_2_2_0::Functown::DDCalc_LogLikelihood);
487  XENON1T_2017_GetLogLikelihood.resolveDependency(&XENON1T_2017_Calc);
488  XENON1T_2017_GetLogLikelihood.resolveBackendReq(&Backends::DDCalc_2_2_0::Functown::DDCalc_Experiment);
489  XENON1T_2017_GetLogLikelihood.resolveBackendReq(&Backends::DDCalc_2_2_0::Functown::DDCalc_LogLikelihood);
490  PICO_60_2017_GetLogLikelihood.resolveDependency(&PICO_60_2017_Calc);
491  PICO_60_2017_GetLogLikelihood.resolveBackendReq(&Backends::DDCalc_2_2_0::Functown::DDCalc_Experiment);
492  PICO_60_2017_GetLogLikelihood.resolveBackendReq(&Backends::DDCalc_2_2_0::Functown::DDCalc_LogLikelihood);
493 
494  // Provide bin number in LZ
495  LZ_GetBinSignal.resolveDependency(&LZ_Calc);
496  LZ_GetBinSignal.resolveBackendReq(&Backends::DDCalc_2_2_0::Functown::DDCalc_Experiment);
497  LZ_GetBinSignal.resolveBackendReq(&Backends::DDCalc_2_2_0::Functown::DDCalc_Bins);
498  LZ_GetBinSignal.resolveBackendReq(&Backends::DDCalc_2_2_0::Functown::DDCalc_BinSignal);
499 
500  // Set generic WIMP mass object
501  mwimp_generic.resolveDependency(&TH_ProcessCatalog_WIMP);
502  mwimp_generic.resolveDependency(&DarkMatter_ID_WIMP);
503  sigma_SI_p_simple.resolveDependency(&DD_couplings_WIMP);
504  sigma_SI_p_simple.resolveDependency(&mwimp_generic);
505 
506  // Generate gamma-ray spectra for various final states
507  if ( (mode >= 0) and (mode < 6) )
508  {
509  std::cout << "Producing test spectra." << std::endl;
510  double mass = 100.;
511  double sv = 3e-26;
512  // The array that is being passed to dumpSpectrum give the branching fraction to various final states.
513  // They are (as defined in TH_ProcessCatalog_WIMP):
514  // 0: b bbar
515  // 1: gamma Z_0
516  // 2: gamma gamma
517  // 3: tau+ tau-
518  // 4: W+ W-
519  // 5: e+ e-
520  // 6: phi phi2
521  // 7: gamma e+ e-
522  if (mode==5) dumpSpectrum("dPhi_dE5.dat", mass, sv*0.1, daFunk::vec<double>(0., 0., 0., 0., 0., 0., 0., 1.));
523  if (mode==0) dumpSpectrum("dPhi_dE0.dat", mass, sv, daFunk::vec<double>(1., 0., 0., 0., 0., 0., 0., 0.));
524  if (mode==1) dumpSpectrum("dPhi_dE1.dat", mass, sv, daFunk::vec<double>(0., 1., 0., 0., 0., 0., 0., 0.));
525  if (mode==2) dumpSpectrum("dPhi_dE2.dat", mass, sv, daFunk::vec<double>(0., 0., 1., 0., 0., 0., 0., 0.));
526  if (mode==3) dumpSpectrum("dPhi_dE3.dat", mass, sv, daFunk::vec<double>(0., 0., 0., 1., 0., 0., 0., 0.));
527  if (mode==4) dumpSpectrum("dPhi_dE4.dat", mass, sv, daFunk::vec<double>(0., 0., 0., 0., 1., 0., 0., 0.));
528  }
529 
530  // Generate gamma-ray spectra for various masses
531  if (mode >= 10)
532  {
533  std::cout << "Producing test spectra." << std::endl;
534  double mass = 100.;
535  double sv = 3e-26;
536  std::string filename = "dPhi_dE_FCMC_" + std::to_string(mode) + ".dat";
537  dumpSpectrum(filename, mass, sv, daFunk::vec<double>(0., 0., 0., 0., 0., 0., 1., 0.), mode);
538  }
539 
540  // Systematic parameter maps annihilation
541  if (mode==6)
542  {
543  std::cout << "Producing gamma ray test maps." << std::endl;
544  int mBins = 60;
545  int svBins = 60;
546  double oh2, lnL;
547  std::vector<double> sv_list, m_list;
548 
549  GalacticHalo_Einasto.reset_and_calculate();
550  set_gamLike_GC_halo.reset_and_calculate();
551 
552  boost::multi_array<double, 2>
553  lnL_b_array{boost::extents[mBins][svBins]},
554  lnL_b_array2{boost::extents[mBins][svBins]},
555  lnL_b_array3{boost::extents[mBins][svBins]},
556  lnL_b_array4{boost::extents[mBins][svBins]},
557  lnL_tau_array{boost::extents[mBins][svBins]};
558  boost::multi_array<double, 2> oh2_array{boost::extents[mBins][svBins]};
559 
560  sv_list = daFunk::logspace(-28.0, -22.0, svBins);
561 
562  std::cout << "Calculating gamma-ray likelihood tables for annihilation to b bbar." << std::endl;
563  m_list = daFunk::logspace(log10(5.), 4., mBins);
564  for (size_t i = 0; i < m_list.size(); i++)
565  {
566  for (size_t j = 0; j < sv_list.size(); j++)
567  {
568  TH_ProcessCatalog_WIMP.setOption<double>("mWIMP", m_list[i]);
569  TH_ProcessCatalog_WIMP.setOption<double>("sv", sv_list[j]);
570  //std::cout << "Parameters: " << m_list[i] << " " << sv_list[j] << std::endl;
571 
572  TH_ProcessCatalog_WIMP.setOption<std::vector<double>>("brList", daFunk::vec<double>(1., 0., 0., 0., 0., 0., 0., 0.));
573  DarkMatter_ID_WIMP.reset_and_calculate();
574  DarkMatterConj_ID_WIMP.reset_and_calculate();
575  TH_ProcessCatalog_WIMP.reset_and_calculate();
576  DM_process_from_ProcessCatalog.reset_and_calculate();
577  RD_fraction_one.reset_and_calculate();
578  SimYieldTable_DarkSUSY.reset_and_calculate();
579  SimYieldTable_MicrOmegas.reset_and_calculate();
580  GA_missingFinalStates.reset_and_calculate();
581  cascadeMC_FinalStates.reset_and_calculate();
582  cascadeMC_DecayTable.reset_and_calculate();
583  cascadeMC_LoopManager.reset_and_calculate();
584  cascadeMC_gammaSpectra.reset_and_calculate();
585  GA_AnnYield_General.reset_and_calculate();
586  lnL_FermiLATdwarfs_gamLike.setOption<std::string>("version", "pass8");
587  lnL_FermiLATdwarfs_gamLike.reset_and_calculate();
589  //std::cout << "Fermi dwarf likelihood: " << lnL << std::endl;
590  lnL_b_array[i][j] = lnL;
591  lnL_HESSGC_gamLike.setOption<std::string>("version", "integral_fixedJ");
592  lnL_HESSGC_gamLike.reset_and_calculate();
593  lnL = lnL_HESSGC_gamLike(0);
594  //std::cout << "HESS GC likelihood: " << lnL << std::endl;
595  lnL_b_array2[i][j] = lnL;
596  lnL_CTAGC_gamLike.reset_and_calculate();
597  lnL = lnL_CTAGC_gamLike(0);
598  //std::cout << "CTA GC likelihood: " << lnL << std::endl;
599  lnL_b_array3[i][j] = lnL;
600  lnL_FermiGC_gamLike.setOption<std::string>("version", "fixedJ");
601  lnL_FermiGC_gamLike.reset_and_calculate();
602  lnL = lnL_FermiGC_gamLike(0);
603  lnL_b_array4[i][j] = lnL;
604  //std::cout << "Fermi GC likelihood: " << lnL << std::endl;
605  }
606  }
607 
608  dump_array_to_file("FermiD_b_table.dat", lnL_b_array, m_list, sv_list);
609  dump_array_to_file("HESSGC_b_table.dat", lnL_b_array2, m_list, sv_list);
610  dump_array_to_file("CTAGC_b_table.dat", lnL_b_array3, m_list, sv_list);
611  dump_array_to_file("FermiGC_b_table.dat", lnL_b_array4, m_list, sv_list);
612 
613  std::cout << "Calculating Fermi-LAT dwarf spheroidal likehood table for annihilation to tau+ tau-." << std::endl;
614  m_list = daFunk::logspace(log10(1.9), 4., mBins);
615  for (size_t i = 0; i < m_list.size(); i++)
616  {
617  for (size_t j = 0; j < sv_list.size(); j++)
618  {
619  TH_ProcessCatalog_WIMP.setOption<double>("mWIMP", m_list[i]);
620  TH_ProcessCatalog_WIMP.setOption<double>("sv", sv_list[j]);
621  //std::cout << "Parameters: " << m_list[i] << " " << sv_list[j] << std::endl;
622 
623  TH_ProcessCatalog_WIMP.setOption<std::vector<double>>("brList", daFunk::vec<double>(0., 0., 0., 1., 0., 0., 0., 0.));
624  DarkMatter_ID_WIMP.reset_and_calculate();
625  DarkMatterConj_ID_WIMP.reset_and_calculate();
626  TH_ProcessCatalog_WIMP.reset_and_calculate();
627  DM_process_from_ProcessCatalog.reset_and_calculate();
628  RD_fraction_one.reset_and_calculate();
629  SimYieldTable_DarkSUSY.reset_and_calculate();
630  SimYieldTable_MicrOmegas.reset_and_calculate();
631  GA_missingFinalStates.reset_and_calculate();
632  cascadeMC_FinalStates.reset_and_calculate();
633  cascadeMC_DecayTable.reset_and_calculate();
634  cascadeMC_LoopManager.reset_and_calculate();
635  cascadeMC_gammaSpectra.reset_and_calculate();
636  GA_AnnYield_General.reset_and_calculate();
637  lnL_FermiLATdwarfs_gamLike.reset_and_calculate();
639  //std::cout << "Fermi LAT likelihood: " << lnL << std::endl;
640  lnL_tau_array[i][j] = lnL;
641  }
642  }
643 
644  dump_array_to_file("FermiD_tau_table.dat", lnL_tau_array, m_list, sv_list);
645 
646  std::cout << "Calculating table of Omega h^2 values." << std::endl;
647  m_list = daFunk::logspace(-1.0, 4., mBins);
648  for (size_t i = 0; i < m_list.size(); i++)
649  {
650  for (size_t j = 0; j < sv_list.size(); j++)
651  {
652  TH_ProcessCatalog_WIMP.setOption<double>("mWIMP", m_list[i]);
653  TH_ProcessCatalog_WIMP.setOption<double>("sv", sv_list[j]);
654  //std::cout << "Parameters: " << m_list[i] << " " << sv_list[j] << std::endl;
655 
656  TH_ProcessCatalog_WIMP.setOption<std::vector<double>>("brList", daFunk::vec<double>(0., 0., 0., 0., 0., 1., 0., 0.));
657  DarkMatter_ID_WIMP.reset_and_calculate();
658  DarkMatterConj_ID_WIMP.reset_and_calculate();
659  TH_ProcessCatalog_WIMP.reset_and_calculate();
660  RD_eff_annrate_from_ProcessCatalog.reset_and_calculate();
661  RD_spectrum_from_ProcessCatalog.reset_and_calculate();
662  RD_spectrum_ordered_func.reset_and_calculate();
663  RD_oh2_DS_general.reset_and_calculate();
664  oh2 = RD_oh2_DS_general(0);
665  //std::cout << "Omega h^2 = " << oh2 << std::endl;
666  oh2_array[i][j] = oh2;
667  }
668  }
669 
670  dump_array_to_file("oh2_table.dat", oh2_array, m_list, sv_list);
671  }
672 
673  // Systematic parameter maps scattering
674  if (mode==7)
675  {
676  std::cout << "Producing direct detection test maps." << std::endl;
677  double lnL1, lnL2, lnL3, lnL4;
678  int nbins;
679  double g, reduced_mass;
680  //int mBins = 300;
681  //int sBins = 200;
682  int mBins = 120;
683  int sBins = 80;
684  const double mN = (m_proton + m_neutron) / 2;
685  std::vector<double> m_list = daFunk::logspace(0.0, 4.0, mBins);
686  std::vector<double> s_list;
687  boost::multi_array<double, 2> lnL_array1{boost::extents[mBins][sBins]},
688  lnL_array2{boost::extents[mBins][sBins]}, lnL_array3{boost::extents[mBins][sBins]},
689  lnL_array4{boost::extents[mBins][sBins]};
690  TH_ProcessCatalog_WIMP.setOption<double>("sv", 0.);
691  TH_ProcessCatalog_WIMP.setOption<std::vector<double>>("brList", daFunk::vec<double>(1., 0., 0., 0., 0., 0., 0., 0.));
692 
693  s_list = daFunk::logspace(-47., -40., sBins);
694  // Calculate array of sigma_SI and lnL values for LZ, PandaX, XENON1T and PICO-60
695  // assuming gps=gns
696 
697  std::cout << "Calculating tables of SI likelihoods." << std::endl;
698  for (size_t i = 0; i < m_list.size(); i++)
699  {
700  for (size_t j = 0; j < s_list.size(); j++)
701  {
702  // Re-initialize DDCalc with LZ/Xenon/PandaX halo parameters
703  Halo_primary_parameters->setValue("rho0", 0.3);
704  Halo_primary_parameters->setValue("vrot", 232.7); // v_Earth = 245 km/s
705  Halo_primary_parameters->setValue("v0", 220.);
706  Halo_primary_parameters->setValue("vesc", 544.);
707  ExtractLocalMaxwellianHalo.reset_and_calculate();
708 
709  TH_ProcessCatalog_WIMP.setOption<double>("mWIMP", m_list[i]);
710  //std::cout << "Parameters: " << m_list[i] << " " << s_list[j] << std::endl;
711  reduced_mass = (m_list[i] * mN) / (mN + m_list[i]);
712  g = sqrt(s_list[j]*pi/gev2cm2) / (reduced_mass);
713  DarkMatter_ID_WIMP.reset_and_calculate();
714  DarkMatterConj_ID_WIMP.reset_and_calculate();
715  TH_ProcessCatalog_WIMP.reset_and_calculate();
716  // Assume for direct and indirect detection likelihoods that dark matter
717  // density is always the measured one (despite relic density results)
718  RD_fraction_one.reset_and_calculate();
719  DD_couplings_WIMP.setOption<double>("gps", g);
720  DD_couplings_WIMP.setOption<double>("gns", g);
721  DD_couplings_WIMP.setOption<double>("gpa", 0.);
722  DD_couplings_WIMP.setOption<double>("gna", 0.);
723  DD_couplings_WIMP.reset_and_calculate();
724  mwimp_generic.reset_and_calculate();
725 
726  DDCalc_2_2_0_init.reset_and_calculate();
727  LZ_Calc.reset_and_calculate();
728  LZ_GetLogLikelihood.reset_and_calculate();
729 
730  XENON1T_2017_Calc.reset_and_calculate();
731  XENON1T_2017_GetLogLikelihood.reset_and_calculate();
732  PandaX_2017_Calc.reset_and_calculate();
733  PandaX_2017_GetLogLikelihood.reset_and_calculate();
734 
735  lnL1 = LZ_GetLogLikelihood(0);
736  lnL2 = PandaX_2017_GetLogLikelihood(0);
737  lnL3 = XENON1T_2017_GetLogLikelihood(0);
738 
739  // Set LocalHalo Model parameters to PICO-60 values
740  Halo_primary_parameters->setValue("rho0", 0.3);
741  Halo_primary_parameters->setValue("vrot", 220.);
742  Halo_primary_parameters->setValue("v0", 220.);
743  Halo_primary_parameters->setValue("vesc", 544.);
744  ExtractLocalMaxwellianHalo.reset_and_calculate();
745 
746  DDCalc_2_2_0_init.reset_and_calculate();
747  PICO_60_2017_Calc.reset_and_calculate();
748  PICO_60_2017_GetLogLikelihood.reset_and_calculate();
749  lnL4 = PICO_60_2017_GetLogLikelihood(0);
750 
751  //std::cout << "LZ SI lnL = " << lnL1 << std::endl;
752  //std::cout << "PandaX_2017 SI lnL = " << lnL2 << std::endl;
753  //std::cout << "XENON1T_2017 SI lnL = " << lnL3 << std::endl;
754  //std::cout << "PICO_60_2017 SI lnL = " << lnL4 << std::endl;
755 
756  DDCalc_2_2_0_init.reset_and_calculate();
757  LZ_Calc.reset_and_calculate();
758  std::vector<double> events;
759  LZ_GetBinSignal.reset_and_calculate();
760  events = LZ_GetBinSignal(0);
761  nbins = events.size();
762  std::cout << "Number of LZ bins: " << nbins << std::endl;
763  std::cout << "Predicted signal: ";
764  for (int ibin=0;ibin<=nbins-1;ibin++) {
765  std::cout << events[ibin] << " ";
766  }
767  std::cout << std::endl;
768 
769  lnL_array1[i][j] = lnL1;
770  lnL_array2[i][j] = lnL2;
771  lnL_array3[i][j] = lnL3;
772  lnL_array4[i][j] = lnL4;
773  }
774  }
775 
776  dump_array_to_file("LZ_SI_table.dat", lnL_array1, m_list, s_list);
777  dump_array_to_file("PandaX_2017_SI_table.dat", lnL_array2, m_list, s_list);
778  dump_array_to_file("XENON1T_2017_SI_table.dat", lnL_array3, m_list, s_list);
779  dump_array_to_file("PICO_60_2017_SI_table.dat", lnL_array4, m_list, s_list);
780 
781  s_list = daFunk::logspace(-42., -35., sBins);
782  // Calculate array of sigma_SI and lnL values for LZ, PandaX, XENON1T and PICO-60
783  // assuming gna=0 (proton-only)
784 
785  std::cout << "Calculating tables of SD likelihoods." << std::endl;
786  for (size_t i = 0; i < m_list.size(); i++)
787  {
788  for (size_t j = 0; j < s_list.size(); j++)
789  {
790  // Re-initialize DDCalc with LZ/Xenon/PandaX halo parameters
791  Halo_primary_parameters->setValue("rho0", 0.3);
792  Halo_primary_parameters->setValue("vrot", 232.7); // v_Earth = 245 km/s
793  Halo_primary_parameters->setValue("v0", 220.);
794  Halo_primary_parameters->setValue("vesc", 544.);
795  ExtractLocalMaxwellianHalo.reset_and_calculate();
796 
797  TH_ProcessCatalog_WIMP.setOption<double>("mWIMP", m_list[i]);
798  //std::cout << "Parameters: " << m_list[i] << " " << s_list[j] << std::endl;
799  reduced_mass = (m_list[i] * m_proton) / (m_proton + m_list[i]);
800  g = sqrt(s_list[j]*pi/(3*gev2cm2)) / (reduced_mass);
801  DarkMatter_ID_WIMP.reset_and_calculate();
802  DarkMatterConj_ID_WIMP.reset_and_calculate();
803  TH_ProcessCatalog_WIMP.reset_and_calculate();
804  RD_fraction_one.reset_and_calculate();
805  DD_couplings_WIMP.setOption<double>("gps", 0.);
806  DD_couplings_WIMP.setOption<double>("gns", 0.);
807  DD_couplings_WIMP.setOption<double>("gpa", g);
808  DD_couplings_WIMP.setOption<double>("gna", 0.);
809  DD_couplings_WIMP.reset_and_calculate();
810  mwimp_generic.reset_and_calculate();
811 
812  DDCalc_2_2_0_init.reset_and_calculate();
813  LZ_Calc.reset_and_calculate();
814  LZ_GetLogLikelihood.reset_and_calculate();
815  XENON1T_2017_Calc.reset_and_calculate();
816  XENON1T_2017_GetLogLikelihood.reset_and_calculate();
817  PandaX_2017_Calc.reset_and_calculate();
818  PandaX_2017_GetLogLikelihood.reset_and_calculate();
819  lnL1 = LZ_GetLogLikelihood(0);
820  lnL2 = PandaX_2017_GetLogLikelihood(0);
821  lnL3 = XENON1T_2017_GetLogLikelihood(0);
822 
823  // Set LocalHalo Model parameters to PICO-60 values
824  Halo_primary_parameters->setValue("rho0", 0.3);
825  Halo_primary_parameters->setValue("vrot", 220.);
826  Halo_primary_parameters->setValue("v0", 220.);
827  Halo_primary_parameters->setValue("vesc", 544.);
828  ExtractLocalMaxwellianHalo.reset_and_calculate();
829 
830  DDCalc_2_2_0_init.reset_and_calculate();
831  PICO_60_2017_Calc.reset_and_calculate();
832  PICO_60_2017_GetLogLikelihood.reset_and_calculate();
833  lnL4 = PICO_60_2017_GetLogLikelihood(0);
834 
835  //std::cout << "LZ SD lnL = " << lnL1 << std::endl;
836  //std::cout << "PandaX_2017 SD lnL = " << lnL2 << std::endl;
837  //std::cout << "XENON1T_2017 SD lnL = " << lnL3 << std::endl;
838  //std::cout << "PICO_60_2017 SD lnL = " << lnL4 << std::endl;
839 
840  lnL_array1[i][j] = lnL1;
841  lnL_array2[i][j] = lnL2;
842  lnL_array3[i][j] = lnL3;
843  lnL_array4[i][j] = lnL4;
844  }
845  }
846 
847  dump_array_to_file("LZ_SD_table.dat", lnL_array1, m_list, s_list);
848  dump_array_to_file("PandaX_2017_SD_table.dat", lnL_array2, m_list, s_list);
849  dump_array_to_file("XENON1T_2017_SD_table.dat", lnL_array3, m_list, s_list);
850  dump_array_to_file("PICO_60_2017_SD_table.dat", lnL_array4, m_list, s_list);
851 
852  // Reset halo parameters to DarkBit defaults.
853  Halo_primary_parameters->setValue("rho0", 0.4);
854  Halo_primary_parameters->setValue("vrot", 235.);
855  Halo_primary_parameters->setValue("v0", 235.);
856  Halo_primary_parameters->setValue("vesc", 550.);
857  ExtractLocalMaxwellianHalo.reset_and_calculate();
858  GalacticHalo_Einasto.reset_and_calculate();
859 
860  }
861  }
862 
863  catch (std::exception& e)
864  {
865  std::cout << "DarkBit_standalone_WIMP has exited with fatal exception: " << e.what() << std::endl;
866  }
867 
868  return 0;
869 }
void lnL_FermiGC_gamLike(double &result)
Fermi LAT galactic center likelihoods, using gamLike backend.
std::vector< TH_Channel > channelList
List of channels.
void cascadeMC_FinalStates(std::vector< std::string > &list)
Function for retrieving list of final states for cascade decays.
Definition: Cascades.cpp:41
void GA_missingFinalStates(std::vector< std::string > &result)
Identification of final states that are not yet tabulated.
Definition: GamYields.cpp:61
void dumpSpectrum(std::string filename, double mWIMP, double sv, std::vector< double > brList, double mPhi=-1)
TH_ParticleProperty getParticleProperty(str) const
Retrieve properties of a given particle involved in one or more processes in this catalog...
void validate()
Validate kinematics and entries.
const double m_proton
void TH_ProcessCatalog_WIMP(TH_ProcessCatalog &result)
void sigma_SI_p_simple(double &result)
Simple calculator of the spin-independent WIMP-proton cross-section.
Includes everything needed to use a GAMBIT module as a standalone analysis code.
const double m_neutron
void RD_fraction_one(double &result)
#define LOCAL_INFO
Definition: local_info.hpp:34
void SimYieldTable_DarkSUSY(SimYieldTable &result)
SimYieldTable based on DarkSUSY6 tabulated results.
Definition: GamYields.cpp:640
STL namespace.
void cascadeMC_InitialState(std::string &pID)
Function selecting initial state for decay chain.
Definition: Cascades.cpp:140
void dump_GammaSpectrum(double &result)
Helper function to dump gamma-ray spectra.
void SimYieldTable_MicrOmegas(SimYieldTable &result)
SimYieldTable based on MicrOmegas tabulated results.
Definition: GamYields.cpp:773
Helpers for making simple spectra from SLHAea objects.
Funk var(std::string arg)
Definition: daFunk.hpp:921
START_MODEL b
Definition: demo.hpp:270
warning & model_warning()
Model warnings.
void RD_spectrum_from_ProcessCatalog(RD_spectrum_type &result)
Collects information about resonances and threshold energies directly from the ProcessCatalog [NB: th...
General small utility functions.
error & backend_error()
Backend errors.
double mass
Particle mass (GeV)
void cascadeMC_EventCount(std::map< std::string, int > &counts)
Event counter for cascade decays.
Definition: Cascades.cpp:175
void setValue(std::string const &inkey, double const &value)
Set single parameter value.
int main(int argc, char *argv[])
void cascadeMC_GenerateChain(DarkBit::DecayChain::ChainContainer &chain)
Function for generating decay chains.
Definition: Cascades.cpp:198
const double pi
void GA_AnnYield_General(daFunk::Funk &result)
General routine to derive annihilation yield.
Definition: GamYields.cpp:420
TH_resonances_thresholds resonances_thresholds
List of resonances and thresholds.
void initialise_standalone_logs(str)
Logger setup standalone utility function.
void DarkMatterConj_ID_WIMP(std::string &result)
void lnL_FermiLATdwarfs_gamLike(double &result)
Fermi LAT dwarf likelihoods, using gamLike backend.
A container holding all annihilation and decay initial states relevant for DarkBit.
std::vector< double > logspace(double x0, double x1, unsigned int n)
Definition: daFunk.hpp:186
const Logging::endofmessage EOM
Explicit const instance of the end of message struct in Gambit namespace.
Definition: logger.hpp:100
void set_gamLike_GC_halo(bool &result)
Fermi LAT dwarf likelihoods, based on arXiv:1108.2914.
QUICK_FUNCTION(DarkBit, TH_ProcessCatalog, OLD_CAPABILITY, TH_ProcessCatalog_WIMP, TH_ProcessCatalog,()) QUICK_FUNCTION(DarkBit
EXPORT_SYMBOLS Logging::LogMaster & logger()
Function to retrieve a reference to the Gambit global log object.
Definition: logger.cpp:95
void dump_array_to_file(const std::string &filename, const boost::multi_array< double, 2 > &a, const std::vector< double > &x, const std::vector< double > &y)
All data on a single annihilation or decay channel, e.g.
const double gev2cm2
void ExtractLocalMaxwellianHalo(LocalMaxwellianHalo &result)
Module function providing local density and velocity dispersion parameters.
Definition: DarkBit.cpp:177
Rollcall header for module DarkBit.
void lnL_HESSGC_gamLike(double &result)
Funk cnst(double x, Args... argss)
Definition: daFunk.hpp:606
double pow(const double &a)
Outputs a^i.
void mwimp_generic(double &result)
Retrieve the DM mass in GeV for generic models (GeV)
Definition: DarkBit.cpp:60
void cascadeMC_LoopManager()
Loop manager for cascade decays.
Definition: Cascades.cpp:80
void cascadeMC_Histograms(std::map< std::string, std::map< std::string, SimpleHist > > &result)
Function responsible for histogramming, and evaluating end conditions for event loop.
Definition: Cascades.cpp:363
void cascadeMC_DecayTable(DarkBit::DecayChain::DecayTable &table)
Function setting up the decay table used in decay chains.
Definition: Cascades.cpp:58
A container for a single process.
void DM_process_from_ProcessCatalog(std::string &result)
Information about the nature of the DM process in question (i.e.
Definition: DarkBit.cpp:114
std::vector< T > vec(std::vector< T > vector)
Definition: daFunk.hpp:142
void RD_oh2_DS_general(double &result)
General routine for calculation of relic density, using DarkSUSY 6+ Boltzmann solver.
shared_ptr< FunkBase > Funk
Definition: daFunk.hpp:113
std::vector< TH_Process > processList
Vector of all processes in this catalog.
void GalacticHalo_Einasto(GalacticHaloProperties &result)
Module function to generate GalacticHaloProperties for Einasto profile.
Definition: DarkBit.cpp:165
void lnL_CTAGC_gamLike(double &result)
#define addParticle(Name, Mass, spinX2)
void RD_spectrum_ordered_func(RD_spectrum_type &result)
Order RD_spectrum object and derive coannihilation thresholds.
TODO: see if we can use this one:
Definition: Analysis.hpp:33
double E(const double x, const double y)
Funk one(Args... argss)
Definition: daFunk.hpp:602