Belle II Software  release-06-00-14
AafhInputModule.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #include <generators/modules/aafhinput/AafhInputModule.h>
10 #include <framework/logging/Logger.h>
11 
12 using namespace Belle2;
13 
14 //-----------------------------------------------------------------
15 // Register the Module
16 //-----------------------------------------------------------------
17 REG_MODULE(AafhInput)
18 
19 //-----------------------------------------------------------------
20 // Implementation
21 //-----------------------------------------------------------------
22 
23 AafhInputModule::AafhInputModule() : Module(), m_initial(BeamParameters::c_smearVertex)
24 {
25  // Set module properties
26  setDescription("AAFH Generator to generate non-radiative two-photon events like e+e- -> e+e-e+e-");
27 
28  // Parameter definitions
29  addParam("mode", m_mode,
30  "decay mode to generate."
31  " 1: e+e- -> mu+mu-L+L-,"
32  " 2: e+e- -> mu+mu-mu+mu-,"
33  " 3: e+e- -> e+e-mu+mu-,"
34  " 4: e+e- -> e+e-L+L-,"
35  " 5: e+e- -> e+e-e+e-. "
36  "L is a user defined particle (default: tau)",
37  5);
38  addParam("rejection", m_rejection,
39  "Rejection method."
40  " 1: apply once for the final event weight,"
41  " 2: apply for the subgenerator weight and for the final weight",
42  2);
43  addParam("maxTries", m_maxTries,
44  "Maximum tries for event generation (ITOT)",
45  m_generator.getMaxTries());
46  addParam("maxFinalWeight", m_maxFinalWeight,
47  "maximum expected final weight for rejection scheme (ESFT)",
48  m_generator.getMaxFinalWeight());
49  addParam("maxSubgeneratorWeight", m_maxSubgeneratorWeight,
50  "maximum expected subgenerator weight for rejection scheme (ESWE)",
51  m_generator.getMaxSubGeneratorWeight());
52  addParam("subgeneratorWeights", m_subgeneratorWeights,
53  "relative weights of the subgenerators: this must be a list of four "
54  "or eight values (first four are interpreted as WAP, rest as WBP) "
55  "which specify the relativ weights for each of the "
56  "four sub generators. The orginial code states that it the program "
57  "run most efficient when the maximum weight is equal in all sub "
58  "generators and that if one wants to be sure that all peaks in the "
59  "differential cross section are accounted the chance to enter each "
60  "sub generator should be equal. Values which try to fullfill both "
61  "conditions are printed at after generation when the output level "
62  "is set to INFO",
63  m_generator.getGeneratorWeights());
64  addParam("suppressionLimits", m_suppressionLimits,
65  "suppression limits for the matrix element calculations. This must "
66  "be a list of four values, [FACE, FACM, FACL, PROC]. For FACE, FACM "
67  "and FACL these specify the size of the propagator for which the "
68  "corresponding spin configurations will be omitted. A value of 1e3 "
69  "will calculate the dominat part correctly and a very large value "
70  "compared to the mass (i.e. 1e50) will calculate it exactly. PROC "
71  "specifies that feynman diagrams which contribute less than 1/PROC "
72  "of the biggest contribution are omitted. For exact calculation it "
73  "should be very big but 1e9 should be considerably faster without "
74  "affecting the result very much",
75  m_generator.getSuppressionLimits());
76  addParam("minMass", m_minMass,
77  "minimum invariant mass of generated fermion pair. Will "
78  "automatically be set to at least 2 times the generated particle "
79  "mass",
80  0.);
81  addParam("particle", m_particle,
82  "name of the generated particle for mode 1 and 4",
83  m_generator.getParticle());
84 }
85 
87 {
89 
90  //Beam Parameters, initial particle - AAFH cannot handle beam energy spread
92 }
93 
95 {
96 
97  // Check if the BeamParameters have changed (if they do, abort the job! otherwise cross section calculation will be a nightmare.)
98  if (m_beamParams.hasChanged()) {
99  if (!m_initialized) {
101  } else {
102  B2FATAL("AafhInputModule::event(): BeamParameters have changed within a job, this is not supported for AAFH!");
103  }
104  }
105 
106  // Initial particle from beam parameters (for random vertex)
107  const MCInitialParticles& initial = m_initial.generate();
108 
109  // True boost.
110  TLorentzRotation boost = initial.getCMSToLab();
111 
112  // vertex.
113  TVector3 vertex = initial.getVertex();
114 
115  MCParticleGraph mpg;
116 
117  //Generate event.
119 
120  //Boost to lab and set vertex.
121  for (size_t i = 0; i < mpg.size(); ++i) {
122  mpg[i].set4Vector(boost * mpg[i].get4Vector());
123 
124  TVector3 v3 = mpg[i].getProductionVertex();
125  v3 = v3 + vertex;
126  mpg[i].setProductionVertex(v3);
127  mpg[i].setValidVertex(true);
128  }
129 
130  //Fill MCParticle List
132 }
133 
135 {
137 }
138 
140 {
141  const BeamParameters& nominal = m_initial.getBeamParameters();
142  const double beamEnergy = nominal.getMass() / 2.;
143 
144  //Set Generator options
145  if (m_mode < 1 || m_mode > 5) {
146  B2ERROR("AafhInputModule::initializeGenerator: 'mode' must be a value between 1 and 5");
147  }
148  if (m_rejection < 1 || m_rejection > 2) {
149  B2ERROR("AafhInputModule::initializeGenerator: 'rejection' must be a value between 1 and 2");
150  }
159 
160  m_initialized = true;
161 
162 }
163 
164 
ERejection
Rejection mode.
Definition: AAFHInterface.h:46
void setMinimumMass(double minMass)
Set the minimum invariant mass for the generated event.
Definition: AAFHInterface.h:67
void finish()
calculate total cross section
void generateEvent(MCParticleGraph &mpg)
generate one event and add it to the graph in CMS
void initialize(double beamEnergy, EMode mode, ERejection rejection)
initialize the generator
EMode
Generation mode.
Definition: AAFHInterface.h:32
void setMaxWeights(double subgeneratorWeight, double finalWeight)
Set the maximum expected weights for the rejection method.
void setMaxTries(int tries)
Set the maximum number of tries to generate an event.
void setSupressionLimits(std::vector< double > limits)
Set the suppression limits when calculation the matrix elements.
void setParticle(const std::string &particle)
Set the particle type for modes c_MuonParticle and c_ElectronParticle.
Definition: AAFHInterface.h:69
void setGeneratorWeights(const std::vector< double > &weights)
Set the relative weights for the sub generators.
AAFH Generator to generate 2-fermion events like e+e- -> e+e-e+e-.
double m_minMass
minimum invariant mass of generated fermion pair
DBObjPtr< BeamParameters > m_beamParams
BeamParameter.
bool m_initialized
True if generator has been initialized.
int m_rejection
Rejection mode.
virtual void initialize() override
initialize generator
std::vector< double > m_suppressionLimits
suppression limits for the matrix element calculations
virtual void event() override
generate event
StoreArray< MCParticle > m_mcparticles
MCParticle collection.
int m_mode
Generator mode.
std::string m_particle
name of the generated particle for modes c_MuonParticle and c_ElectronParticle
virtual void terminate() override
calculate cross section
int m_maxTries
maximum number of tries for generating an event
std::vector< double > m_subgeneratorWeights
relative subgenerator weights
double m_maxSubgeneratorWeight
maximum expected subgenerator weight for rejection scheme
InitialParticleGeneration m_initial
initial particle used by BeamParameter class
double m_maxFinalWeight
maximum expected final weight for rejection scheme
AAFHInterface m_generator
interface to the generator
void initializeGenerator()
Method is called to initialize the generator.
This class contains the nominal beam parameters and the parameters used for smearing of the primary v...
const BeamParameters & getBeamParameters() const
Return reference to nominal beam parameters.
This class contains the initial state for the given event.
const TLorentzRotation & getCMSToLab() const
Return the LorentzRotation to convert from CMS to lab frame.
const TVector3 & getVertex() const
Get the position of the collision.
Class to build, validate and sort a particle decay chain.
@ c_setDecayInfo
Set decay time and vertex.
size_t size() const
Return the number of particles in the graph.
void generateList(const std::string &name="", int options=c_setNothing)
Generates the MCParticle list and stores it in the StoreArray with the given name.
Base class for Modules.
Definition: Module.h:72
const std::string & getName() const
Return name under which the object is saved in the DataStore.
bool registerInDataStore(DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut)
Register the object/array in the DataStore.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
void initialize()
function to be executed on initialize()
MCInitialParticles & generate()
Generate a new event.
Abstract base class for different kinds of events.