Belle II Software  release-05-02-19
TrepsInputModule.h
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2019 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Kiyoshi Hayasaka, Yo Sato *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #pragma once
12 
13 #include <framework/core/Module.h>
14 #include <framework/datastore/StoreArray.h>
15 
16 #include <generators/treps/Treps3B.h>
17 #include <generators/treps/UtrepsB.h>
18 
19 #include <mdst/dataobjects/MCParticle.h>
20 #include <mdst/dataobjects/MCParticleGraph.h>
21 #include <generators/utilities/InitialParticleGeneration.h>
22 #include <string>
23 
24 namespace Belle2 {
35  class TrepsInputModule : public Module {
36 
37  public:
42 
44  virtual void initialize() override;
45 
47  virtual void terminate() override;
48 
50  virtual void event() override;
51 
53  double simulateW();
54 
56  double getCrossSection(double W);
57 
59  void initializeGenerator();
60 
61  private:
62 
64  UtrepsB m_generator;
65 
67  std::string m_parameterFile;
68 
70  std::string m_wListTableFile;
71 
74 
89  double m_maximalQ2;
90 
95  double m_maximalAbsCosTheta;
96 
99 
102 
108 
109  bool m_initialized{false};
111  };
113 }
114 
Belle2::TrepsInputModule::terminate
virtual void terminate() override
initialization for trepsinput
Definition: TrepsInputModule.cc:164
Belle2::TrepsInputModule::m_generator
UtrepsB m_generator
Treps generator.
Definition: TrepsInputModule.h:72
Belle2::TrepsInputModule::m_initial
InitialParticleGeneration m_initial
initial particle used by BeamParameter class
Definition: TrepsInputModule.h:87
Belle2::TrepsInputModule::initialize
virtual void initialize() override
initialization for trepsinput
Definition: TrepsInputModule.cc:66
Belle2::MCParticleGraph
Class to build, validate and sort a particle decay chain.
Definition: MCParticleGraph.h:48
Belle2::TrepsInputModule::event
virtual void event() override
input event from TREPS
Definition: TrepsInputModule.cc:72
Belle2::TrepsInputModule::m_maximalQ2
double m_maximalQ2
Maximal Q^2 = -q^2, where q is the difference between the initial and final electron or positron mome...
Definition: TrepsInputModule.h:97
Belle2::TrepsInputModule::simulateW
double simulateW()
Simulate W distribution according to given input file of cross section table.
Definition: TrepsInputModule.cc:183
Belle2::TrepsInputModule::m_differentialCrossSectionFile
std::string m_differentialCrossSectionFile
Differential cross section table file.
Definition: TrepsInputModule.h:81
Belle2::TrepsInputModule::initializeGenerator
void initializeGenerator()
Initialize the TREPS generator
Definition: TrepsInputModule.cc:217
Belle2::TrepsInputModule::m_parameterFile
std::string m_parameterFile
Parameter file which configures the setting of beam, production particle etc.
Definition: TrepsInputModule.h:75
Belle2::InitialParticleGeneration
Generate Collision.
Definition: InitialParticleGeneration.h:35
Belle2::DBObjPtr
Class for accessing objects in the database.
Definition: DBObjPtr.h:31
Belle2::TrepsInputModule::m_initialized
bool m_initialized
True if generator has been initialized.
Definition: TrepsInputModule.h:117
Belle2::TrepsInputModule::m_applyCosThetaCutCharged
bool m_applyCosThetaCutCharged
Whether to apply cut on |cos(theta)| for charged particles only.
Definition: TrepsInputModule.h:106
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::TrepsInputModule::m_mpg
MCParticleGraph m_mpg
An instance of the MCParticle graph.
Definition: TrepsInputModule.h:83
Belle2::TrepsInputModule::m_wListTableFile
std::string m_wListTableFile
W-List table file.
Definition: TrepsInputModule.h:78
Belle2::TrepsInputModule::TrepsInputModule
TrepsInputModule()
Constructor: Sets the description, the properties and the parameters of the module.
Definition: TrepsInputModule.cc:31
Belle2::TrepsInputModule::m_minimalTransverseMomentum
double m_minimalTransverseMomentum
Minimal transverse momentum of the final-state particles.
Definition: TrepsInputModule.h:109
Belle2::TrepsInputModule::m_maximalAbsCosTheta
double m_maximalAbsCosTheta
Maximal |cos(theta)|, where theta is the final-state particle polar angle.
Definition: TrepsInputModule.h:103
Belle2::TrepsInputModule::m_applyTransverseMomentumCutCharged
bool m_applyTransverseMomentumCutCharged
Whether to apply cut on the minimal transverse momentum for charged particles only.
Definition: TrepsInputModule.h:115
Belle2::TrepsInputModule::m_mcparticles
StoreArray< MCParticle > m_mcparticles
MCParticle collection.
Definition: TrepsInputModule.h:89
Belle2::StoreArray< MCParticle >
Belle2::TrepsInputModule::m_beamParams
DBObjPtr< BeamParameters > m_beamParams
BeamParameter.
Definition: TrepsInputModule.h:85
Belle2::TrepsInputModule::getCrossSection
double getCrossSection(double W)
Returns the differential cross section for given W [GeV].
Definition: TrepsInputModule.cc:168
Belle2::TrepsInputModule::m_useDiscreteAndSortedW
bool m_useDiscreteAndSortedW
if it is true, W-list table is used for discrete and sorted W
Definition: TrepsInputModule.h:91