/bxdecay0

C++ port of the Decay0/GENBB fortran Monte Carlo code for the generation of standard decay or double beta decay processes for various radioactive nuclides of interest

Primary LanguageFortranGNU General Public License v3.0GPL-3.0

BxDecay0 - C++ port of the legacy Decay0 FORTRAN library

Authors: François Mauger
Vladimir Tretyak
Emma Mauger
Date: 2020-01-16
Copyright: Copyright (C) 2017-2020 the BxCppDev group

The BxDecay0 C++ library provides a set of classes and functions for the random generation of simulated nuclear decays. It consists in a C++ port of the original Decay0 Fortran program (also known as GENBB in some other context) written and maintained by Vladimir Tretyak (KINR). Decay0 was created to address the Monte Carlo generation of nuclear decays in the context of double beta decay and dark matter experimental research.

BxDecay0 aims to be used by any C++ simulation program that needs to generate the kinematics of primary particles emitted from specific nuclear decays (example: input for a Geant4 based Monte Carlo simulation program). It can be easily interfaced with the HepMC event/particle model or any other client application.

The first version of the Decay0 program was written in the early 90's by Vladimir Tretyak and collaborators, using the Fortran 77 programming language and the CERNLIB library. It has been maintained and improved up to now (2020).

From 1992 to 2010, the Decay0/GENBB code was embedded in the GEANT3 based simulation software of the NEMO2 and then NEMO3 double beta decay experiments.

From 2005 to 2011, the NEMO3 collaboration has initiated a R&D program to design a new generation double beta decay experiment: SuperNEMO. The choice was made to use C++ as the main programming language. In this context, François Mauger has created a C++ version of GENBB. The first version was a C++ wrapper of the original Fortran code, using binding mechanism based on plain static C structures.

Later, a pure C++ version was released without any dependency on Fortran and CERNLIB. This code was integrated in 2011 as the genbb_help module of the Bayeux C++ library (version < 4).

This release of the BxDecay0 C++ library is extracted from the Bayeux genbb_help module with some changes from the 2017-03-01 release of the Decay0 program by Vladimir Tretyak. It has then been updated to the 2018-12-05 release.

BxDecay0 is a standalone library with very few dependencies (mostly the GSL library for numerical integration and a few special math functions). External random engines can be used through a simple wrapping functor interface, particularly the ones provided by the C++ standard random library can be used by default. However the user is free to provide its own uniform deviates random generator (based on GSL, ROOT or whatever).

Versions:

  • Prerelease 1.0.0 : mixed port from Decay0 2017-03-01 and embedded decay0 C++ code in Bayeux 3 based on Decay0 2013.
  • First release 1.0.0 : updated from Decay0 2018-12-05

The core of the BxDecay0 code does not follow a fully object-oriented approach. In order to ensure the easy synchronization of its low level code with the original Decay0 code, BxDecay0 mimics the layout of the Fortran code (including massive usage of GOTO statements!). BxDecay0 provides a large collection of plain generator functions for about 100 radioactive nuclei split in two categories: double beta decay and background/calibration decay. When a Decay0 fix or improvement is published in the original Fortran code by its author (V.Tretyak), it is thus rather easy to adequately change/adapt the C++ code in the relevant section of BxDecay0.

Hopefully, BxDecay0 gets rid of the original common block based data model in Decay0 which has strong limitations in terms of usability in a modern OOP context (static data structures). The BxDecay0 API introduces its own OOP data model through the bxdecay0::event and bxdecay0::particle classes (see the ex01 example). It is thus easy to use such classes through any C++ client program and/or to interface with some high level event generator library (i.e. HepMC3). See the ex02 example.

More, BxDecay0 provides the bxdecay0::decay0_generator class which wraps low-level functions with a simple OOP interface.

Finally, it is also possible to use the low level C++ functions ported from the original Fortran code. However it is not recommended and should be reserved to experts and developpers of the library.

BxDecay0 is developped on a Ubuntu Linux (18.04 LTS) and should work on any Unix/BSD flavor with a recent C++ compiler with c++11 support (i.e. GNU g++ >= 4.9) including macOS.

The following lines give some hints to prepare your system before the installation of BxDecay0. Some instructions may vary depending on your own system.

  1. Install GNU C++ compiler:

    $ sudo apt-get install g++
  2. Install CMake:

    $ sudo apt-get install cmake
  3. Install the GNU scientific library (development package):

    $ sudo apt-get install libgsl-dev
    $ gsl-config --version
    2.4

Clone the Git development repository on your filesystem:

$ cd /tmp
$ git clone https://github.com/BxCppDev/bxdecay0.git bxdecay0.git
$ cd bxdecay0.git

Or download the archive associated to a released version :

$ cd /tmp
$ wget https://github.com/BxCppDev/bxdecay0/archive/1.0.0.tar.gz
$ tar xvzf bxdecay0-1.0.0.tar.gz
$ cd bxdecay0-1.0.0

Here we use a temporary build directory and choose to install BxDecay0 in our home directory:

$ mkdir /tmp/_build.d
$ cd /tmp/_build.d
$ cmake -DCMAKE_INSTALL_PREFIX=${HOME}/bxdecay0 /tmp/bxdecay0.git

or:

$ cmake -DCMAKE_INSTALL_PREFIX=${HOME}/bxdecay0 /tmp/bxdecay0-1.0.0

From the build directory:

$ make -j4
$ make test
$ make install

Add the following line in your shell startup script (i.e. ~/.bashrc):

$ export PATH=${HOME}/bxdecay0/bin:$PATH

The bxdecay0-query script will be usable from your projects:

$ which bxdecay0-query
  • The bxdecay0-query utility script allows you to fetch informations about your installation of the BxDecay0 library.

    $ bxdecay0-query --help
    $ bxdecay0-query --prefix
    $ bxdecay0-query --version
    $ bxdecay0-query --cmakedir
  • CMake configuration scripts are provided:

    • BxDecay0Config.cmake,
    • BxDecay0ConfigVersion.cmake.

    The find_package(BxDecay0 1.0 CONFIG) CMake command can be given the following variable to locate BxDecay0 on your system from a client project which uses the CMake build system:

    $ cmake -DBxDecay0_DIR="$(bxdecay0-query --cmakedir)" ...

BxDecay0 comes with CMake and pkg-config support. The BxDecay0 installation directory contains dedicated scripts usable by client applications.

The following program is extracted from the BxDecay0's ex00 example. It randomly generates 10 simulated events corresponding to the neutrinoless double beta decay (DBD) process of 100 Mo to the ground state of 100 Ru. The resulting events are printed in the terminal in a very simple format. It is of course possible to adapt this program and use the OOP interface of the bxdecay0::event class in order to extract physical quantities of interest (particles' type and momentum...).

#include <iostream>
#include <bxdecay0/std_random.h>        // Wrapper for the standard random module's PRNG
#include <bxdecay0/event.h>             // Decay event data model
#include <bxdecay0/decay0_generator.h>  // Decay0 generator with OOP interface

int main()
{
  // Declare a PRNG:
  unsigned int seed = 314159;                 // Random seed
  std::default_random_engine generator(seed); // Standard PRNG
  bxdecay0::std_random prng(generator);       // PRNG wrapper

  // Declare a Decay0 generator:
  bxdecay0::decay0_generator decay0;

  // Configure the Decay0 generator:
  decay0.set_decay_category(bxdecay0::decay0_generator::DECAY_CATEGORY_DBD); // Double-beta decay process
  decay0.set_decay_isotope("Mo100");              // Emitter nucleus
  decay0.set_decay_dbd_level(0);                  // Ground state of the daughter nucleus
  decay0.set_decay_dbd_mode(bxdecay0::DBDMODE_1); // Neutrinoless DBD (mass mechanism)
  // or :
  // decay0.set_decay_dbd_mode_by_label("0nubb_mn");

  // Initialize the Decay0 generator.
  // We need to pass some PRNG to pre-compute some quantities from energy distributions:
  decay0.initialize(prng);

  // Shoot some decay events:
  std::size_t nevents = 10;
  for (std::size_t ievent = 0; ievent < nevents; ievent++) {
    bxdecay0::event gendecay;     // Declare an empty decay event
    decay0.shoot(prng, gendecay); // Randomize the decay event
    gendecay.store(std::cout);    // Basic ASCII output
  }

  decay0.reset(); // Terminate the generator
  return 0;
}
  • ex00 : Minimal program for the generation of Mo100 neutrinoless double beta decay events (mass mechanism) with plain ASCII output,
  • ex01 : Generation of Mo100 two neutrino double beta decay events with plain ASCII output,
  • ex02 : Generation of Mo100 two neutrino double beta decay events with HepMC3 formatted ASCII output,
  • ex03 : Generation of Co60 decay events with plain ASCII output,
  • ex04 : Use of the plumbing bxdecay0::genbbsub function (expert/developper only).

BxDecay0 is released under the GNU GENERAL PUBLIC LICENSE, version 3. See the LICENSE.txt file.

  • Vladimir Tretyak (KINR, Kiev Institute for Nuclear Research, Lepton Physics Department, Ukraine) is the original author and maintener of the Fortran Decay0 package,
  • François Mauger (LPC Caen, Laboratoire de Physique Corpusculaire de Caen, Université de Caen Normandie, France) is the author and maintener of the original C++ port of Decay0 within Bayeux and the BxDecay0 library,
  • Emma Mauger (formerly Université de Caen Normandie) has done a large part of the extraction and port to C++ of the standalone BxDecay0 from the original Bayeux genbb library module.
  • Vladimir Tretyak, DECAY0 event generator for initial kinematics of particles in alpha, beta and double beta decays, talk given at Laboratori Nazionali del Gran Sasso, 17 March 2015 :
  • O.A.Ponkratenko, V.I.Tretyak, Yu.G.Zdesenko, Event Generator DECAY4 for Simulating Double-Beta Processes and Decays of Radioactive Nuclei, Phys. At. Nucl. 63 (2000) 1282 (nucl-ex/0104018)

From the dbd_isotopes.lis resource file:

  • Ca40
  • Ca46
  • Ca48
  • Ni58
  • Zn64
  • Zn70
  • Ge76
  • Se74
  • Se82
  • Sr84
  • Zr94
  • Zr96
  • Mo92
  • Mo100
  • Ru96
  • Ru104
  • Cd106
  • Cd108
  • Cd114
  • Cd116
  • Sn112
  • Sn122
  • Sn124
  • Te120
  • Te128
  • Te130
  • Xe136
  • Ce136
  • Ce138
  • Ce142
  • Nd148
  • Nd150
  • Dy156
  • Dy158
  • W180
  • W186
  • Os184
  • Os192
  • Pt190
  • Pt198
  • Bi214 (for Bi214+At214)
  • Pb214 (for Pb214+Po214)
  • Po218 (for Po218+Rn218+Po214)
  • Rn222 (for Rn222+Ra222+Rn218+Po214)
  • Sm144
  • Sm154
  • Er162
  • Er164
  • Er170
  • Yb168
  • Yb176
  • Ca40 -> Ar40 :

    1. 0+ (gs) {0 MeV}
  • Ca46 -> Ti46 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.889 MeV}
  • Ca48 -> Ti48 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.984 MeV}
    3. 2+ (2) {2.421 MeV}
  • Ni58 -> Fe58 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.811 MeV}
    3. 2+ (2) {1.675 MeV}
  • Zn64 -> Ni64 :

    1. 0+ (gs) {0 MeV}
  • Zn70 -> Ge70 :

    1. 0+ (gs) {0 MeV}
  • Ge76 -> Se76 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.559 MeV}
    3. 0+ (1) {1.122 MeV}
    4. 2+ (2) {1.216 MeV}
  • Se74 -> Ge74 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.596 MeV}
    3. 2+ (2) {1.204 MeV}
  • Se82 -> Kr82 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.776 MeV}
    3. 2+ (2) {1.475 MeV}
  • Sr84 -> Kr84 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.882 MeV}
  • Zr94 -> Mo94 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.871 MeV}
  • Zr96 -> Mo96 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.778 MeV}
    3. 0+ (1) {1.148 MeV}
    4. 2+ (2) {1.498 MeV}
    5. 2+ (3) {1.626 MeV}
    6. 2+ (4) {2.096 MeV}
    7. 2+ (5) {2.426 MeV}
    8. 0+ (2) {2.623 MeV}
    9. 2+ (6) {2.700 MeV}
    10. 2+?(7) {2.713 MeV}
  • Mo92 -> Zr92 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.934 MeV}
    3. 0+ (1) {1.383 MeV}
  • Mo100 -> Ru100 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.540 MeV}
    3. 0+ (1) {1.130 MeV}
    4. 2+ (2) {1.362 MeV}
    5. 0+ (2) {1.741 MeV}
  • Ru96 -> Mo96 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.778 MeV}
    3. 0+ (1) {1.148 MeV}
    4. 2+ (2) {1.498 MeV}
    5. 2+ (3) {1.626 MeV}
    6. 2+ (4) {2.096 MeV}
    7. 2+ (5) {2.426 MeV}
    8. 0+ (2) {2.623 MeV}
    9. 2+ (6) {2.700 MeV}
    10. 2+?(7) {2.713 MeV}
  • Ru104 -> Pd104 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.556 MeV}
  • Cd106 -> Pd106 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.512 MeV}
    3. 2+ (2) {1.128 MeV}
    4. 0+ (1) {1.134 MeV}
    5. 2+ (3) {1.562 MeV}
    6. 0+ (2) {1.706 MeV}
  • Cd108 -> Pd108 :

    1. 0+ (gs) {0 MeV}
  • Cd114 -> Sn114 :

    1. 0+ (gs) {0 MeV}
  • Cd116 -> Sn116 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {1.294 MeV}
    3. 0+ (1) {1.757 MeV}
    4. 0+ (2) {2.027 MeV}
    5. 2+ (2) {2.112 MeV}
    6. 2+ (3) {2.225 MeV}
  • Sn112 -> Cd112 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.618 MeV}
    3. 0+ (1) {1.224 MeV}
    4. 2+ (2) {1.312 MeV}
    5. 0+ (2) {1.433 MeV}
    6. 2+ (3) {1.469 MeV}
    7. 0+ (3) {1.871 MeV}
  • Sn122 -> Te122 :

    1. 0+ (gs) {0 MeV}
  • Sn124 -> Te124 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.603 MeV}
    3. 2+ (2) {1.326 MeV}
    4. 0+ (1) {1.657 MeV}
    5. 0+ (2) {1.883 MeV}
    6. 2+ (3) {2.039 MeV}
    7. 2+ (4) {2.092 MeV}
    8. 0+ (3) {2.153 MeV}
    9. 2+ (5) {2.182 MeV}
  • Te120 -> Sn120 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {1.171 MeV}
  • Te128 -> Xe128 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.443 MeV}
  • Te130 -> Xe130 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.536 MeV}
    3. 2+ (2) {1.122 MeV}
    4. 0+ (1) {1.794 MeV}
  • Xe136 -> Ba136 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.819 MeV}
    3. 2+ (2) {1.551 MeV}
    4. 0+ (1) {1.579 MeV}
    5. 2+ (3) (2.080 MeV}
    6. 2+ (4) {2.129 MeV}
    7. 0+ (2) {2.141 MeV}
    8. 2+ (5) {2.223 MeV}
    9. 0+ (3) {2.315 MeV}
    10. 2+ (6) {2.400 MeV}
  • Ce136 -> Ba136 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.819 MeV}
    3. 2+ (2) {1.551 MeV}
    4. 0+ (1) {1.579 MeV}
    5. 2+ (3) (2.080 MeV}
    6. 2+ (4) {2.129 MeV}
    7. 0+ (2) {2.141 MeV}
    8. 2+ (5) {2.223 MeV}
    9. 0+ (3) {2.315 MeV}
    10. 2+ (6) {2.400 MeV}
  • Ce138 -> Ba138 :

    1. 0+ (gs) {0 MeV}
  • Ce142 -> Nd142 :

    1. 0+ (gs) {0 MeV}
  • Nd148 -> Sm148 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.550 MeV}
    3. 0+ (1) {1.424 MeV}
    4. 2+ (2) {1.454 MeV}
    5. 2+ (3) {1.664 MeV}
    6. 0+ (2) {1.921 MeV}
  • Nd150 -> Sm150 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.334 MeV}
    3. 0+ (1) {0.740 MeV}
    4. 2+ (2) {1.046 MeV}
    5. 2+ (3) {1.194 MeV}
    6. 0+ (2) {1.256 MeV}
  • Sm144 -> Nd144 (new : Decay0 2018-12-05) :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.697 MeV}
    3. 2+ (2) {1.561 MeV}
  • Sm154 -> Gd144 (new : Decay0 2018-12-05) :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.123 MeV}
    3. 0+ (1) {0.681 MeV}
    4. 2+ (2) {0.815 MeV}
    5. 2+ (3) {0.996 MeV}
    6. 0+ (2) {1.182 MeV}
  • Dy156 -> Gd156 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.089 MeV}
    3. 0+ (1) {1.050 MeV}
    4. 2+ (2) {1.129 MeV}
    5. 2+ (3) {1.154 MeV}
    6. 0+ (2) {1.168 MeV}
    7. 2+ (4) {1.258 MeV}
    8. 0+ (3) {1.715 MeV}
    9. 2+ (5) {1.771 MeV}
    10. 2+ (6) {1.828 MeV}
    11. 0+ (4) {1.851 MeV}
    12. 2+ (7) {1.915 MeV}
    13. 1- {1.946 MeV}
    14. 0- {1.952 MeV}
    15. 0+ (5) {1.989 MeV}
    16. 2+ (8) {2.004 MeV}
  • Dy158 -> Gd158 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.080 MeV}
    3. 4+ (1) {0.261 MeV}
  • Er162 -> Dy162 (new : Decay0 2018-12-05) :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.081 MeV}
    3. 2+ (2) {0.888 MeV}
    4. 0+ (1) {1.400 MeV}
    5. 2+ (3) {1.453 MeV}
    6. 0+ (2) {1.666 MeV}
    7. 2+ (4) {1.728 MeV}
    8. 2+ (5) {1.783 MeV}
  • Er164 -> Dy164 (new : Decay0 2018-12-05) :

    1. 0+ (gs) {0 MeV};
  • Er170 -> Yb170 (new : Decay0 2018-12-05) :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.084 MeV}
  • Yb168 -> Er168 (new : Decay0 2018-12-05) :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.080 MeV}
    3. 2+ (2) {0.821 MeV}
    4. 0+ (1) {1.217 MeV}
    5. 2+ (3) {1.276 MeV}
  • Yb176 -> Hf176 (new : Decay0 2018-12-05) :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.088 MeV}
  • W180 -> Hf180 :

    1. 0+ (gs) {0 MeV}
  • W186 -> Os186 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.137 MeV}
  • Os184 -> W184 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.111 MeV}
    3. 2+ (2) {0.903 MeV}
    4. 0+ (1) {1.002 MeV}
    5. 2+ (3) {1.121 MeV}
    6. 0+ (2) {1.322 MeV}
    7. 2+ (4) {1.386 MeV}
    8. 2+ (5) {1.431 MeV}
  • Os192 -> Pt192 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.317 MeV}
  • Pt190 -> Os190 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.187 MeV}
    3. 2+ (2) {0.558 MeV}
    4. 0+ (1) {0.912 MeV}
    5. 2+ (3) {1.115 MeV}
    6. 0+ (2) {1.382 MeV}
  • Pt198 -> Hg198 :

    1. 0+ (gs) {0 MeV}
    2. 2+ (1) {0.412 MeV}
  • Bi214 -> At214 :

    1. 1- (gs) {0 MeV}
  • Pb214 -> Po214 :

    1. 0+ (gs) {0 MeV}
  • Po218 -> Rn218 :

    1. 0+ (gs) {0 MeV}
  • Rn222 -> Ra222 :

    1. 0+ (gs) {0 MeV}

From the bxdecay0::bb_utils.h C++ header and the dbd_modes.lis resource file:

BxDecay0 mode Identification label Decay0 mode Description
DBDMODE_1 0nubb_mn 1 0nubb(mn) 0+ -> 0+ {2n} with neutrino mass
DBDMODE_2 0nubb_rhc_lambda_0 2 0nubb(rhc-lambda) 0+ -> 0+ {2n} with RHC lambda
DBDMODE_3 0nubb_rhc_lambda_02 3 0nubb(rhc-lambda) 0+ -> 0+, 2+ {N*} with RHC lambda
DBDMODE_4 2nubb 4 2nubb 0+ -> 0+ {2n}
DBDMODE_5 0nubbM1 5 0nubbM1 0+ -> 0+ {2n} (Majoron, SI=1)
DBDMODE_6 0nubbM3 7 0nubbM3 0+ -> 0+ {2n} (Majoron, SI=3)
DBDMODE_7 0nubb_rhc_lambda_2 9 0nubb(rhc-lambda) 0+ -> 2+ {2n} with RHC lambda
DBDMODE_8 2nubb_2 10 2nubb 0+ -> 2+ {2n}, {N*}
DBDMODE_9 0nuKb+ 11 0nuKb+ 0+ -> 0+, 2+
DBDMODE_10 2nuKb+ 12 2nuKb+ 0+ -> 0+, 2+
DBDMODE_11 0nu2K 13 0nu2K 0+ -> 0+, 2+
DBDMODE_12 2nu2K 14 2nu2K 0+ -> 0+, 2+
DBDMODE_13 0nubbM7 8 0nubbM7 0+ -> 0+ {2n} (Majoron, SI=7)
DBDMODE_14 0nubbM2 6 0nubbM2 0+ -> 0+ (2n} (Majoron, SI=2)
DBDMODE_15 2nubb_bosonic_0 15 2nubb 0+ -> 0+ with bosonic neutrinos
DBDMODE_16 2nubb_bosonic_2 16 2nubb 0+ -> 2+ with bosonic neutrinos
DBDMODE_17 0nubb_rhc_eta_s 17 0nubb(rhc-eta) 0+ -> 0+ with RHC eta simplified expression
DBDMODE_18 0nubb_rhc_eta_nmes 18 0nubb(rhc-eta) 0+ -> 0+ with RHC eta and specific NMEs
DBDMODE_19 2nub_lv 19 2nubb(LV) 0+ -> 0+ with Lorentz violation
DBDMODE_20 0nu4b 20 0nu4b 0+ -> 0+ Quadruple beta decay

Comments on specific modes:

  • The original Decay0 code has changed the so-called modebb index of some Majoron modes from versions above 2017 with respect to former versions. We thus use an index map to ensure backward compatibility and force the BxDecay0 interface stable with respect to the original C++ port in the Bayeux/genbb module.
  • 5, 6, 13, 14 are Majoron modes with spectral index SI:
    • SI=1 - old Gelmini-Roncadelli Majoron
    • SI=3 - double M, vector M, charged M
    • SI=7
    • SI=2 - bulk M of Mohapatra
  • 20: quadruple beta decay, possible only for Zr96, Xe136, Nd150

From the background_isotopes.lis resource file:

  • Ac228
  • Am241
  • Ar39
  • Ar42
  • As79 (for As79+Se79m)
  • Bi207 (for Bi207+Pb207m)
  • Bi208
  • Bi210
  • Bi212 (for Bi212+Po212)
  • Bi214 (for Bi214+Po214)
  • Ca48 (for Ca48+Sc48)
  • C14
  • Cd113
  • Co60
  • Cs136
  • Cs137 (for Cs137+Ba137m)
  • Eu147
  • Eu152
  • Eu154
  • Gd146
  • Hf182
  • I126
  • I133
  • I134
  • I135
  • K40
  • K42
  • Kr81
  • Kr85
  • Mn54
  • Na22
  • P32
  • Pa231 (added from Bayeux : 2013-09-06)
  • Pa234m
  • Pb210
  • Pb211
  • Pb212
  • Pb214
  • Ra226 (added from Bayeux : 2013-07-11)
  • Ra228
  • Rb87
  • Rh106
  • Sb125
  • Sb126
  • Sb133
  • Sr90
  • Ta180m-B-
  • Ta180m-EC
  • Ta182
  • Te133
  • Te133m
  • Te134
  • Th234
  • Tl207
  • Tl208
  • Xe129m
  • Xe131m
  • Xe133
  • Xe135
  • Y88
  • Y90
  • Zn95
  • Zr96 (for Zr96+Nb96)

Unlike the original Decay0 code, BxDecay0 does not support the generation of so-called artifical events (Compton, Moller scattering, e+e- pair). It should not be difficult to implement such generators by yourself independently of BxDecay0.