gambit is hosted by Hepforge, IPPP Durham
GAMBIT  v1.5.0-2191-ga4742ac
a Global And Modular Bsm Inference Tool
FlavBit_standalone_example.cpp File Reference
Include dependency graph for FlavBit_standalone_example.cpp:

Go to the source code of this file.

Namespaces

 Gambit
 TODO: see if we can use this one:
 
 Gambit::FlavBit
 

Functions

std::string infile ("FlavBit/data/example.slha")
 
 QUICK_FUNCTION (FlavBit, unimproved_MSSM_spectrum, NEW_CAPABILITY, createSpectrum, Spectrum,(MSSM30atQ, MSSM30atMGUT)) QUICK_FUNCTION(FlavBit
 
Spectrum QUICK_FUNCTION (FlavBit, Z_decay_rates, NEW_CAPABILITY, GammaZ, DecayTable::Entry) QUICK_FUNCTION(FlavBit
 
void Gambit::FlavBit::createSpectrum (Spectrum &outSpec)
 Make an unimproved GAMBIT spectrum object from an SLHA file. More...
 
void Gambit::FlavBit::relabelSpectrum (Spectrum &outSpec)
 Relabel it as a complete spectrum. More...
 
void Gambit::FlavBit::GammaW (DecayTable::Entry &result)
 W decays – only need the total width for SuperIso. More...
 
void Gambit::FlavBit::GammaZ (DecayTable::Entry &result)
 Z decays – only need the total width for SuperIso. More...
 
int main (int argc, char **argv)
 

Variables

 MSSM_spectrum
 
 NEW_CAPABILITY
 
 relabelSpectrum
 
 Spectrum
 
 MSSM30atQ
 
 MSSM30atMGUT
 
 unimproved_MSSM_spectrum
 
Spectrum W_plus_decay_rates
 
Spectrum GammaW
 

Function Documentation

◆ infile()

std::string infile ( "FlavBit/data/example.slha"  )

Referenced by Gambit::FlavBit::createSpectrum(), Gambit::ColliderBit::getHepMCEvent(), and main().

Here is the caller graph for this function:

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 72 of file FlavBit_standalone_example.cpp.

References Gambit::FlavBit::b2ll_likelihood(), Gambit::FlavBit::b2ll_measurements(), Gambit::FlavBit::b2sgamma_likelihood(), Gambit::FlavBit::b2sll_likelihood(), Gambit::FlavBit::b2sll_measurements(), Gambit::FlavBit::createSpectrum(), Gambit::FlavBit::deltaMB_likelihood(), Gambit::Utils::ensure_path_exists(), Gambit::EOM, Gambit::FlavBit::FH_DeltaMs(), Gambit::FlavBit::FH_FlavourObs(), GammaW, Gambit::FlavBit::GammaZ(), infile(), Gambit::LogTags::info, Gambit::Logging::LogMaster::initialise(), Gambit::logger(), MSSM_spectrum, relabelSpectrum, Gambit::FlavBit::SI_BDmunu(), Gambit::FlavBit::SI_BDstarmunu(), Gambit::FlavBit::SI_BDtaunu(), Gambit::FlavBit::SI_Bmumu(), Gambit::FlavBit::SI_bsgamma(), Gambit::FlavBit::SI_Bsmumu_untag(), Gambit::FlavBit::SI_Btaunu(), Gambit::FlavBit::SI_Dmunu(), Gambit::FlavBit::SI_Dsmunu(), Gambit::FlavBit::SI_Dstaunu(), Gambit::FlavBit::SI_fill(), Gambit::FlavBit::SI_RD(), Gambit::FlavBit::SI_RDstar(), Gambit::FlavBit::SL_likelihood(), and Gambit::FlavBit::SL_measurements().

73 {
74 
75  cout << "Starting FlavBit_standalone" << endl;
76 
77  try
78  {
79 
80  // Get the SLHA filename from the command line, if it has been given.
81  if (argc >= 2) infile = argv[1];
82 
83  // Make a logging object
84  std::map<std::string, std::string> loggerinfo;
85 
86  // Define where the logs will end up
87  std::string prefix("runs/FlavBit_standalone/logs/");
88 
89  // Ensure that the above directory exists
91 
92  // Add entries to the loggerinfo map
93  loggerinfo["Core, Error"] = prefix+"core_errors.log";
94  loggerinfo["Default"] = prefix+"default.log";
95  loggerinfo["Warning"] = prefix+"warnings.log";
96  loggerinfo["FlavBit, Info"] = prefix+"FlavBit_info.log";
97 
98  // Initialise global LogMaster object
99  logger().initialise(loggerinfo);
100 
101  logger()<<"Running FlavBit standalone example"<<LogTags::info<<EOM;
102 
103  // Notify all module functions that care of the model being scanned.
104  createSpectrum.notifyOfModel("MSSM30atQ");
105  relabelSpectrum.notifyOfModel("MSSM30atQ");
106  SI_fill.notifyOfModel("MSSM30atQ");
107 
108  // Arrange the spectrum chain
109  relabelSpectrum.resolveDependency(&createSpectrum);
110 
111  // Set up the deltaMB_LL likelihood
112  // Have to resolve dependencies by hand
113  // deltaMB_likelihood depends on:
114  // - deltaMs
115  deltaMB_likelihood.resolveDependency(&FH_DeltaMs);
116 
117  // FH_deltaMs depends on:
118  // - FH_FlavourObs
119  FH_DeltaMs.resolveDependency(&FH_FlavourObs);
120 
121  // FH_FlavourObs has only one backend requirement:
122  // - FHFlavour
123  FH_FlavourObs.resolveBackendReq(&Backends::FeynHiggs_2_11_3::Functown::FHFlavour);
124 
125  // The FeynHiggs initialisation function depends on:
126  // - unimproved_MSSM_spectrum
127  FeynHiggs_2_11_3_init.resolveDependency(&createSpectrum);
128 
129  // Set up the b2sll_LL likelihood
130  // Have to resolve dependencies by hand
131  // b2sll_likelihood depends on:
132  // - b2sll_M
133  b2sll_likelihood.resolveDependency(&b2sll_measurements);
134 
135  //SI_fill needs on:
136  // - BEreq Init_param
137  // - BEreq slha_adjust
138  // - BEreq mb_1S
139  // - MSSM_spectrum (or SM_spectrum)
140  // - W_plus_decay_rates
141  // - Z_decay_rates
142  SI_fill.resolveDependency(&relabelSpectrum);
143  SI_fill.resolveDependency(&GammaZ);
144  SI_fill.resolveDependency(&GammaW);
145  SI_fill.resolveBackendReq(&Backends::SuperIso_3_6::Functown::Init_param);
146  SI_fill.resolveBackendReq(&Backends::SuperIso_3_6::Functown::slha_adjust);
147  SI_fill.resolveBackendReq(&Backends::SuperIso_3_6::Functown::mb_1S);
148 
149  // b2sll_measurements depends on:
150  // - BKstarmumu_11_25
151  // - BKstarmumu_25_40
152  // - BKstarmumu_40_60
153  // - BKstarmumu_60_80
154  // - BKstarmumu_15_17
155  // - BKstarmumu_17_19
156  b2sll_measurements.resolveDependency(&SI_BKstarmumu_11_25);
157  b2sll_measurements.resolveDependency(&SI_BKstarmumu_25_40);
158  b2sll_measurements.resolveDependency(&SI_BKstarmumu_40_60);
159  b2sll_measurements.resolveDependency(&SI_BKstarmumu_60_80);
160  b2sll_measurements.resolveDependency(&SI_BKstarmumu_15_17);
161  b2sll_measurements.resolveDependency(&SI_BKstarmumu_17_19);
162 
163  // Now resolve dependencies of the BKstar mu mu measurements
164  // Each depends on:
165  // - SI_fill
166  // Each has a BE requirement on:
167  // - BKstarmumu_CONV (from SuperIso)
168  SI_BKstarmumu_11_25.resolveDependency(&SI_fill);
169  SI_BKstarmumu_11_25.resolveBackendReq(&Backends::SuperIso_3_6::Functown::BKstarmumu_CONV);
170  SI_BKstarmumu_25_40.resolveDependency(&SI_fill);
171  SI_BKstarmumu_25_40.resolveBackendReq(&Backends::SuperIso_3_6::Functown::BKstarmumu_CONV);
172  SI_BKstarmumu_40_60.resolveDependency(&SI_fill);
173  SI_BKstarmumu_40_60.resolveBackendReq(&Backends::SuperIso_3_6::Functown::BKstarmumu_CONV);
174  SI_BKstarmumu_60_80.resolveDependency(&SI_fill);
175  SI_BKstarmumu_60_80.resolveBackendReq(&Backends::SuperIso_3_6::Functown::BKstarmumu_CONV);
176  SI_BKstarmumu_15_17.resolveDependency(&SI_fill);
177  SI_BKstarmumu_15_17.resolveBackendReq(&Backends::SuperIso_3_6::Functown::BKstarmumu_CONV);
178  SI_BKstarmumu_17_19.resolveDependency(&SI_fill);
179  SI_BKstarmumu_17_19.resolveBackendReq(&Backends::SuperIso_3_6::Functown::BKstarmumu_CONV);
180 
181  // Now do the b2ll likelihood
182  b2ll_likelihood.resolveDependency(&b2ll_measurements);
183 
184  // b2ll_measurements depends on:
185  // - Bsmumu_untag
186  // - Bmumu
187  b2ll_measurements.resolveDependency(&SI_Bsmumu_untag);
188  b2ll_measurements.resolveDependency(&SI_Bmumu);
189 
190  // Resolve dependencies of SI_Bsmumu_untag
191  // These are:
192  // - SI_fill
193  // Plus BE reqs:
194  // - Bsll_untag_CONV
195  SI_Bsmumu_untag.resolveDependency(&SI_fill);
196  SI_Bsmumu_untag.resolveBackendReq(&Backends::SuperIso_3_6::Functown::Bsll_untag_CONV);
197 
198  // Resolve dependencies of SI_Bmumu
199  // These are:
200  // - SI_fill
201  // Plus BE reqs:
202  // - Bmumu
203  // - CW_calculator
204  // - C_calculator_base1
205  // - CQ_calculator
206  SI_Bmumu.resolveDependency(&SI_fill);
207  SI_Bmumu.resolveBackendReq(&Backends::SuperIso_3_6::Functown::Bll_CONV);
208 
209  // Now do the semi-leptonic likelihood SL_LL
210  // This depends on:
211  // - SL_M
212  SL_likelihood.resolveDependency(&SL_measurements);
213 
214  // Resolve dependencies of SL_measurements
215  // which are:
216  // - RD
217  // - RDstar
218  // - BDmunu
219  // - BDstarmunu
220  // - Btaunu
221  // - Dstaunu
222  // - Dsmunu
223  // - Dmunu
224  SL_measurements.resolveDependency(&SI_RD);
225  SL_measurements.resolveDependency(&SI_RDstar);
226  SL_measurements.resolveDependency(&SI_BDmunu);
227  SL_measurements.resolveDependency(&SI_BDstarmunu);
228  SL_measurements.resolveDependency(&SI_Btaunu);
229  SL_measurements.resolveDependency(&SI_Dsmunu);
230  SL_measurements.resolveDependency(&SI_Dstaunu);
231  SL_measurements.resolveDependency(&SI_Dmunu);
232 
233  // Resolve all of the individual dependencies and backend reqs
234  // These are:
235  // - SI_fill
236  // BE Req: BDtaunu, etc
237  SI_Btaunu.resolveDependency(&SI_fill);
238  SI_Btaunu.resolveBackendReq(&Backends::SuperIso_3_6::Functown::Btaunu);
239  SI_BDtaunu.resolveDependency(&SI_fill);
240  SI_BDtaunu.resolveBackendReq(&Backends::SuperIso_3_6::Functown::BRBDlnu);
241  SI_RD.resolveDependency(&SI_fill);
242  SI_RD.resolveBackendReq(&Backends::SuperIso_3_6::Functown::BDtaunu_BDenu);
243  SI_RDstar.resolveDependency(&SI_fill);
244  SI_RDstar.resolveBackendReq(&Backends::SuperIso_3_6::Functown::BDstartaunu_BDstarenu);
245  SI_Dstaunu.resolveDependency(&SI_fill);
246  SI_Dstaunu.resolveBackendReq(&Backends::SuperIso_3_6::Functown::Dstaunu);
247  SI_Dsmunu.resolveDependency(&SI_fill);
248  SI_Dsmunu.resolveBackendReq(&Backends::SuperIso_3_6::Functown::Dsmunu);
249  SI_Dmunu.resolveDependency(&SI_fill);
250  SI_Dmunu.resolveBackendReq(&Backends::SuperIso_3_6::Functown::Dmunu);
251 
252  // Now resolve dependencies for the b->sgamma likelihoods
253  b2sgamma_likelihood.resolveDependency(&SI_bsgamma);
254 
255  // Resolve dependencies and backend requirements of SI_bsgamma:
256  // - SI_fill
257  // Plus BE reqs:
258  // - bsgamma_CONV
259  SI_bsgamma.resolveDependency(&SI_fill);
260  SI_bsgamma.resolveBackendReq(&Backends::SuperIso_3_6::Functown::bsgamma_CONV);
261 
262  // Double-check which backend requirements have been filled with what
263  std::cout << std::endl << "My function SI_fill has had its backend requirement on Init_param filled by:" << std::endl;
264  std::cout << FlavBit::Pipes::SI_fill::BEreq::Init_param.origin() << "::";
265  std::cout << FlavBit::Pipes::SI_fill::BEreq::Init_param.name() << std::endl;
266  std::cout << std::endl << "My function SI_fill has had its backend requirement on slha_adjust filled by:" << std::endl;
267  std::cout << FlavBit::Pipes::SI_fill::BEreq::slha_adjust.origin() << "::";
268  std::cout << FlavBit::Pipes::SI_fill::BEreq::slha_adjust.name() << std::endl;
269 
270  // Double-check which backend requirements have been filled with what
271  std::cout << std::endl << "My function SI_Bmumu has had its backend requirement on Bll_CONV filled by:" << std::endl;
272  std::cout << FlavBit::Pipes::SI_Bmumu::BEreq::Bll_CONV.origin() << "::";
273  std::cout << FlavBit::Pipes::SI_Bmumu::BEreq::Bll_CONV.name() << std::endl;
274 
275  // Double-check which dependencies have been filled with whatever (not every combination is shown)
276  std::cout << std::endl << "My function SI_fill has had its dependency on MSSM_spectrum filled by:" << endl;
277  std::cout << FlavBit::Pipes::SI_fill::Dep::MSSM_spectrum.origin() << "::";
278  std::cout << FlavBit::Pipes::SI_fill::Dep::MSSM_spectrum.name() << std::endl;
279 
280  std::cout << std::endl << "My function b2sll_likelihood has had its dependency on b2sll_M filled by:" << endl;
281  std::cout << FlavBit::Pipes::b2sll_likelihood::Dep::b2sll_M.origin() << "::";
282  std::cout << FlavBit::Pipes::b2sll_likelihood::Dep::b2sll_M.name() << std::endl;
283 
284  std::cout << std::endl << "My function b2sll_measurements has had its dependency on BKstarmumu_11_25 filled by:" << endl;
285  std::cout << FlavBit::Pipes::b2sll_measurements::Dep::BKstarmumu_11_25.origin() << "::";
286  std::cout << FlavBit::Pipes::b2sll_measurements::Dep::BKstarmumu_11_25.name() << std::endl;
287  std::cout << std::endl << "My function b2sll_measurements has had its dependency on BKstarmumu_25_40 filled by:" << endl;
288  std::cout << FlavBit::Pipes::b2sll_measurements::Dep::BKstarmumu_25_40.origin() << "::";
289  std::cout << FlavBit::Pipes::b2sll_measurements::Dep::BKstarmumu_25_40.name() << std::endl;
290  std::cout << std::endl << "My function b2sll_measurements has had its dependency on BKstarmumu_40_60 filled by:" << endl;
291  std::cout << FlavBit::Pipes::b2sll_measurements::Dep::BKstarmumu_40_60.origin() << "::";
292  std::cout << FlavBit::Pipes::b2sll_measurements::Dep::BKstarmumu_40_60.name() << std::endl;
293  std::cout << std::endl << "My function b2sll_measurements has had its dependency on BKstarmumu_60_80 filled by:" << endl;
294  std::cout << FlavBit::Pipes::b2sll_measurements::Dep::BKstarmumu_60_80.origin() << "::";
295  std::cout << FlavBit::Pipes::b2sll_measurements::Dep::BKstarmumu_60_80.name() << std::endl;
296  std::cout << std::endl << "My function b2sll_measurements has had its dependency on BKstarmumu_15_17 filled by:" << endl;
297  std::cout << FlavBit::Pipes::b2sll_measurements::Dep::BKstarmumu_15_17.origin() << "::";
298  std::cout << FlavBit::Pipes::b2sll_measurements::Dep::BKstarmumu_15_17.name() << std::endl;
299  std::cout << std::endl << "My function b2sll_measurements has had its dependency on BKstarmumu_17_19 filled by:" << endl;
300  std::cout << FlavBit::Pipes::b2sll_measurements::Dep::BKstarmumu_17_19.origin() << "::";
301  std::cout << FlavBit::Pipes::b2sll_measurements::Dep::BKstarmumu_17_19.name() << std::endl;
302 
303  // Now we start the actual calculations (which we would loop over if doing a scan).
304  double loglike;
305  std::cout << std::endl;
306 
307  // Initialise the spectra
308  createSpectrum.reset_and_calculate();
309  relabelSpectrum.reset_and_calculate();
310 
311  // Initialise the backends
312  SuperIso_3_6_init.reset_and_calculate();
313  FeynHiggs_2_11_3_init.reset_and_calculate();
314 
315  // Now call the module functions in an appropriate order
316  SI_fill.reset_and_calculate();
317 
318  // Calculate the B meson mass asymmetry likelihood
319  FH_FlavourObs.reset_and_calculate();
320  FH_DeltaMs.reset_and_calculate();
321  deltaMB_likelihood.reset_and_calculate();
322  loglike = deltaMB_likelihood(0);
323  std::cout << std::endl << "B meson mass asymmetry log-likelihood: " << loglike << std::endl;
324 
325  // Calculate the B -> ll likelihood
326  SI_Bsmumu_untag.reset_and_calculate();
327  SI_Bmumu.reset_and_calculate();
328  b2ll_measurements.reset_and_calculate();
329  b2ll_likelihood.reset_and_calculate();
330  loglike = b2ll_likelihood(0);
331  std::cout << "Fully leptonic B decays (B->ll) joint log-likelihood: " << loglike << std::endl;
332 
333  // Calculate the B -> Xs ll likelihood
334  SI_BKstarmumu_11_25.reset_and_calculate();
335  SI_BKstarmumu_25_40.reset_and_calculate();
336  SI_BKstarmumu_40_60.reset_and_calculate();
337  SI_BKstarmumu_60_80.reset_and_calculate();
338  SI_BKstarmumu_15_17.reset_and_calculate();
339  SI_BKstarmumu_17_19.reset_and_calculate();
340  b2sll_measurements.reset_and_calculate();
341  b2sll_likelihood.reset_and_calculate();
342  loglike = b2sll_likelihood(0);
343  std::cout << "Leptonic penguin B decays (B->X_s ll) joint log-likelihood: " << loglike << std::endl;
344 
345  // Calculate the semi-leptonic (SL) likelihood
346  SI_Btaunu.reset_and_calculate();
347  SI_BDtaunu.reset_and_calculate();
348  SI_RD.reset_and_calculate();
349  SI_RDstar.reset_and_calculate();
350  SI_Dstaunu.reset_and_calculate();
351  SI_Dsmunu.reset_and_calculate();
352  SI_Dmunu.reset_and_calculate();
353  SL_measurements.reset_and_calculate();
354  SL_likelihood.reset_and_calculate();
355  loglike = SL_likelihood(0);
356  std::cout << "Semi-leptonic B decays (B->D l nu) joint log-likelihood: " << loglike << std::endl;
357 
358  // Calculate the B -> X_s gamma likelihood
359  SI_bsgamma.reset_and_calculate();
360  b2sgamma_likelihood.reset_and_calculate();
361  loglike = b2sgamma_likelihood(0);
362  std::cout << "B->X_s gamma log-likelihood: " << loglike << std::endl;
363 
364  std::cout << endl;
365 
366  }
367 
368  catch (std::exception& e)
369  {
370  std::cout << "FlavBit_standalone example has exited with fatal exception: " << e.what() << std::endl;
371  }
372 
373 }
void deltaMB_likelihood(double &result)
Likelihood for Delta Ms.
Definition: FlavBit.cpp:1271
void b2ll_likelihood(double &result)
Likelihood for rare purely leptonic B decays.
Definition: FlavBit.cpp:1413
void SI_BDstarmunu(double &result)
Br B -> D* mu nu.
Definition: FlavBit.cpp:683
void SI_RDstar(double &result)
B->D* tau nu / B-> D* e nu decays.
Definition: FlavBit.cpp:719
void SI_Bsmumu_untag(double &result)
Br Bs->mumu decays for the untagged case (CP-averaged)
Definition: FlavBit.cpp:516
void createSpectrum(Spectrum &outSpec)
Make an unimproved GAMBIT spectrum object from an SLHA file.
void SI_BDtaunu(double &result)
Br B -> D tau nu.
Definition: FlavBit.cpp:617
void b2sgamma_likelihood(double &result)
Likelihood for b->s gamma.
Definition: FlavBit.cpp:1308
void SI_RD(double &result)
B-> D tau nu / B-> D e nu decays.
Definition: FlavBit.cpp:705
void SI_Dmunu(double &result)
Br D -> mu nu.
Definition: FlavBit.cpp:603
void SI_fill(parameters &result)
Fill SuperIso model info structure.
Definition: FlavBit.cpp:86
void FH_DeltaMs(double &result)
Definition: FlavBit.cpp:1109
void b2ll_measurements(predictions_measurements_covariances &pmc)
Measurements for rare purely leptonic B decays.
Definition: FlavBit.cpp:1347
void SL_likelihood(double &result)
Likelihood for tree-level leptonic and semileptonic B decays.
Definition: FlavBit.cpp:1550
void SI_bsgamma(double &result)
Br b-> s gamma decays.
Definition: FlavBit.cpp:501
void b2sll_measurements(predictions_measurements_covariances &pmc)
Measurements for electroweak penguin decays.
Definition: FlavBit.cpp:1117
void b2sll_likelihood(double &result)
Likelihood for electroweak penguin decays.
Definition: FlavBit.cpp:1231
void GammaZ(DecayTable::Entry &result)
Z decays – only need the total width for SuperIso.
const Logging::endofmessage EOM
Explicit const instance of the end of message struct in Gambit namespace.
Definition: logger.hpp:100
std::string infile("FlavBit/data/example.slha")
EXPORT_SYMBOLS Logging::LogMaster & logger()
Function to retrieve a reference to the Gambit global log object.
Definition: logger.cpp:95
EXPORT_SYMBOLS const str & ensure_path_exists(const str &)
Ensure that a path exists (and then return the path, for chaining purposes)
void SI_Dstaunu(double &result)
Br B->D_s tau nu.
Definition: FlavBit.cpp:575
void initialise(std::vector< std::pair< std::set< std::string >, std::string >> &)
Function to construct loggers according to blueprint.
Definition: logmaster.cpp:227
void SI_BDmunu(double &result)
Br B -> D mu nu.
Definition: FlavBit.cpp:639
void SL_measurements(predictions_measurements_covariances &pmc)
Measurements for tree-level leptonic and semileptonic B decays.
Definition: FlavBit.cpp:1453
void SI_Bmumu(double &result)
Br B0->mumu decays.
Definition: FlavBit.cpp:546
Spectrum GammaW
void FH_FlavourObs(fh_FlavourObs &result)
Flavour observables from FeynHiggs: B_s mass asymmetry, Br B_s -> mu mu, Br B -> X_s gamma...
Definition: FlavBit.cpp:1068
void SI_Btaunu(double &result)
Br B->tau nu_tau decays.
Definition: FlavBit.cpp:561
void SI_Dsmunu(double &result)
Br B->D_s mu nu.
Definition: FlavBit.cpp:589
Here is the call graph for this function:

◆ QUICK_FUNCTION() [1/2]

QUICK_FUNCTION ( FlavBit  ,
unimproved_MSSM_spectrum  ,
NEW_CAPABILITY  ,
createSpectrum  ,
Spectrum  ,
(MSSM30atQ, MSSM30atMGUT  
)

◆ QUICK_FUNCTION() [2/2]

Spectrum QUICK_FUNCTION ( FlavBit  ,
Z_decay_rates  ,
NEW_CAPABILITY  ,
GammaZ  ,
DecayTable::Entry   
)

Variable Documentation

◆ GammaW

Spectrum GammaW

Definition at line 38 of file FlavBit_standalone_example.cpp.

Referenced by main().

◆ MSSM30atMGUT

MSSM30atMGUT

Definition at line 36 of file FlavBit_standalone_example.cpp.

◆ MSSM30atQ

MSSM30atQ

Definition at line 36 of file FlavBit_standalone_example.cpp.

◆ MSSM_spectrum

MSSM_spectrum

Definition at line 36 of file FlavBit_standalone_example.cpp.

Referenced by main().

◆ NEW_CAPABILITY

Spectrum NEW_CAPABILITY

Definition at line 36 of file FlavBit_standalone_example.cpp.

◆ relabelSpectrum

relabelSpectrum

Definition at line 36 of file FlavBit_standalone_example.cpp.

Referenced by main().

◆ Spectrum

Definition at line 36 of file FlavBit_standalone_example.cpp.

◆ unimproved_MSSM_spectrum

◆ W_plus_decay_rates