Belle II Software  release-08-01-10
BBBremInputModule.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/bbbreminput/BBBremInputModule.h>
10 
11 #include <framework/logging/Logger.h>
12 #include <framework/datastore/DataStore.h>
13 #include <framework/datastore/StoreArray.h>
14 #include <framework/datastore/StoreObjPtr.h>
15 #include <framework/dataobjects/EventMetaData.h>
16 
17 #include <Math/Vector3D.h>
18 
19 using namespace std;
20 using namespace Belle2;
21 
22 //-----------------------------------------------------------------
23 // Register the Module
24 //-----------------------------------------------------------------
25 REG_MODULE(BBBremInput);
26 
27 
28 //-----------------------------------------------------------------
29 // Implementation
30 //-----------------------------------------------------------------
31 
32 BBBremInputModule::BBBremInputModule() : Module(), m_initial(BeamParameters::c_smearVertex)
33 {
34  //Set module properties
35  setDescription("Generates low scattering angle radiative Bhabha events (Beam-Beam Bremsstrahlung).");
36 
37  //Parameter definition
38  addParam("MinPhotonEnergyFraction", m_photonEFrac, "Fraction of the minimum photon energy.", 0.000001);
39  addParam("Unweighted", m_unweighted, "Produce unweighted or weighted events.", true);
40  addParam("MaxWeight", m_maxWeight, "The max weight (only for Unweighted=True).", 2000.0);
41  addParam("DensityCorrectionMode", m_densityCorrectionMode, "Mode for bunch density correction (none=0, hard=1 (default), soft=2)",
42  1);
43  addParam("DensityCorrectionParameter", m_DensityCorrectionParameter, "Density correction parameter tc (=(hbarc/sigma_y)^2)",
44  1.68e-17);
45 }
46 
47 
49 {
50 
51 }
52 
53 
55 {
56  StoreArray<MCParticle> mcparticle;
57  mcparticle.registerInDataStore();
58 
59  //Beam Parameters, initial particle - BBBREM cannot handle beam energy spread
61 
62 }
63 
64 
66 {
67  // Check if the BeamParameters have changed (if they do, abort the job! otherwise cross section calculation will be a nightmare.)
68  if (m_beamParams.hasChanged()) {
69  if (!m_initialized) {
71  } else {
72  B2FATAL("BBBremInputModule::event(): BeamParameters have changed within a job, this is not supported for BBBREM!");
73  }
74  }
75 
76  StoreObjPtr<EventMetaData> evtMetaData("EventMetaData", DataStore::c_Event);
77 
78  // initial particle from beam parameters
79  const MCInitialParticles& initial = m_initial.generate();
80 
81  // true boost
82  ROOT::Math::LorentzRotation boost = initial.getCMSToLab();
83 
84  // vertex
85  ROOT::Math::XYZVector vertex = initial.getVertex();
86 
87  m_mcGraph.clear();
88  double weight = m_generator.generateEvent(m_mcGraph, vertex, boost);
89 
91 
92  //store the weight (1.0 for unweighted events)
93  evtMetaData->setGeneratedWeight(weight);
94 }
95 
96 
98 {
99  m_generator.term();
100 
101  B2RESULT("Cross-section (weighted): " << m_generator.getCrossSection() << " +/- " <<
102  m_generator.getCrossSectionError() << " [mb]");
103  B2RESULT("Maximum weight delivered: " << m_generator.getMaxWeightDelivered());
104  B2RESULT("Overweight bias cross-section (unweighted): " << m_generator.getCrossSectionOver() << " +/- " <<
106 }
107 
109 {
110 
111  const BeamParameters& nominal = m_initial.getBeamParameters();
112  double centerOfMassEnergy = nominal.getMass();
113 
116 
117  m_initialized = true;
118 
119 }
DBObjPtr< BeamParameters > m_beamParams
BeamParameter.
bool m_initialized
True if generator has been initialized.
double m_DensityCorrectionParameter
Density correction parameter tc.
virtual void initialize() override
Initializes the module.
BBBrem m_generator
Variables.
bool m_unweighted
True if BBBrem should produce unweighted events.
virtual void event() override
Method is called for each event.
virtual void terminate() override
Method is called at the end of the event processing.
double m_maxWeight
The maximum weight.
MCParticleGraph m_mcGraph
The MCParticle graph object.
virtual ~BBBremInputModule()
Destructor.
InitialParticleGeneration m_initial
initial particle used by BeamParameter class
double m_photonEFrac
Module parameters.
int m_densityCorrectionMode
Mode for bunch density correction.
void initializeGenerator()
Method is called to initialize the generator.
double getCrossSectionErrorOver()
Returns the error on the total overweight bias cross section of the generated process in millibarn.
Definition: BBBrem.h:116
void term()
Terminates the generator.
Definition: BBBrem.cc:135
double generateEvent(MCParticleGraph &mcGraph, ROOT::Math::XYZVector vertex, ROOT::Math::LorentzRotation boost)
Generates one single event.
Definition: BBBrem.cc:80
void init(double cmsEnergy, double minPhotonEFrac, bool unweighted=true, double maxWeight=2000.0, int densitymode=1, double densityparameter=1.68e-17)
Initializes the generator.
Definition: BBBrem.cc:21
double getCrossSectionOver()
Returns the total overweight bias cross section of the generated process in millibarn.
Definition: BBBrem.h:111
double getCrossSection()
Returns the total cross section of the generated process in millibarn.
Definition: BBBrem.h:101
double getCrossSectionError()
Returns the error on the total cross section of the generated process in millibarn.
Definition: BBBrem.h:106
double getMaxWeightDelivered()
Returns the maximum weight given by the BBBrem generation.
Definition: BBBrem.h:121
This class contains the nominal beam parameters and the parameters used for smearing of the primary v...
@ c_Event
Different object in each event, all objects/arrays are invalidated after event() function has been ca...
Definition: DataStore.h:59
const BeamParameters & getBeamParameters() const
Return reference to nominal beam parameters.
This class contains the initial state for the given event.
const ROOT::Math::XYZVector & getVertex() const
Get the position of the collision.
const ROOT::Math::LorentzRotation & getCMSToLab() const
Return the LorentzRotation to convert from CMS to lab frame.
@ c_checkCyclic
Check for cyclic dependencies.
@ c_setDecayInfo
Set decay time and vertex.
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
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
bool registerInDataStore(DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut)
Register the object/array in the DataStore.
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#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.
void clear()
Reset particles and decay information to make the class reusable.
Abstract base class for different kinds of events.