Belle II Software  release-05-02-19
TrepsInputModule.cc
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 #include <generators/modules/trepsinput/TrepsInputModule.h>
12 #include <string.h>
13 #include <string>
14 #include <boost/filesystem.hpp>
15 #include <TFile.h>
16 #include <TRandom.h>
17 
18 #include <map>
19 
20 using namespace Belle2;
21 
22 //-----------------------------------------------------------------
23 // Register the Module
24 //-----------------------------------------------------------------
25 REG_MODULE(TrepsInput)
26 
27 //-----------------------------------------------------------------
28 // Implementation
29 //-----------------------------------------------------------------
30 
31 TrepsInputModule::TrepsInputModule() : Module(), m_initial(BeamParameters::c_smearVertex)
32 {
33  // Set module properties
34  setDescription("Input from TREPS generator (No-tag), Input from TREPS generator for ee->ee hadrons");
35 
36  // Parameter definitions
37  addParam("ParameterFile", m_parameterFile,
38  "parameter file for TREPS input", std::string("treps_par.dat"));
39  addParam("DifferentialCrossSectionFile", m_differentialCrossSectionFile,
40  "file name for differential cross section table input. If UseDiscreteAndSortedW is true, the file is used",
41  std::string("pipidcs.dat"));
42  addParam("WListTableFile", m_wListTableFile,
43  "file name for W-List table input. If UseDiscreteAndSortedW is false (default), the file is used", std::string("wlist_table.dat"));
44  addParam("UseDiscreteAndSortedW", m_useDiscreteAndSortedW,
45  "if true, use WListTable for discrete and sorted W. if false (default), use DifferentialCrossSection", false);
46 
47  addParam("MaximalQ2", m_maximalQ2,
48  "Maximal :math:`Q^2 = -q^2`, where q is the difference between the initial "
49  "and final electron or positron momenta. Negative means no cut.",
50  -1.0);
51  addParam("MaximalAbsCosTheta", m_maximalAbsCosTheta,
52  "Maximal :math:`|\\cos(\\theta)|`, where :math:`\\theta` is the final-state particle "
53  "polar angle.", 1.01);
54  addParam("ApplyCosThetaCutCharged", m_applyCosThetaCutCharged,
55  "Whether to apply cut on :math:`|cos(theta)|` for charged particles only.",
56  true);
57  addParam("MinimalTransverseMomentum", m_minimalTransverseMomentum,
58  "Minimal transverse momentum of the final-state particles.", 0.0);
59  addParam("ApplyTransverseMomentumCutCharged",
60  m_applyTransverseMomentumCutCharged,
61  "Whether to apply cut on the minimal transverse momentum for "
62  "charged particles only.", true);
63 
64 }
65 
67 {
69  m_mcparticles.registerInDataStore();
70 }
71 
73 {
74  /*
75  * Check if the BeamParameters have changed (if they do, abort the job;
76  * otherwise cross section calculation will be a nightmare).
77  */
78  if (m_beamParams.hasChanged()) {
79  if (!m_initialized) {
81  B2INFO("Initialized the TREPS generator!");
82  } else {
83  B2FATAL("TrepsInputModule::event(): BeamParameters have changed within "
84  "a job, this is not supported for TREPS!");
85  }
86  }
87 
88  /* Generation of the initial particle from beam parameters. */
90  TVector3 vertex = initial.getVertex();
91 
93  if (m_generator.inmode != 0) return;
94 
95  m_generator.wtcount++;
96  m_generator.wf = (float)(m_generator.wtable(1));
97 
98  if (abs((double)m_generator.wf - m_generator.w) >= 0.001 && m_generator.wf > 0.01) {
99  B2INFO(" W value changed. " << m_generator.w << " to " << m_generator.wf);
100  m_generator.w = (double)m_generator.wf;
101  m_generator.updateW();
102  }
103  } else {
104  m_generator.w = simulateW();
105  m_generator.updateW();
106  }
107 
108  int idummy = 0;
109  int iret = m_generator.event_gen(idummy);
110  m_mpg.clear();
111 
112  // fill data of the initial electron and positron
115  electron.setPDG(11);
116  electron.setMomentum(m_generator.getElectronMomentum());
117  electron.setMass(m_generator.me);
118  electron.setEnergy(sqrt(m_generator.me * m_generator.me + m_generator.getElectronMomentum().Mag2()));
119 
122  positron.setPDG(-11);
123  positron.setMomentum(m_generator.getPositronMomentum());
124  positron.setMass(m_generator.me);
125  positron.setEnergy(sqrt(m_generator.me * m_generator.me + m_generator.getPositronMomentum().Mag2()));
126 
127  if (iret >= 1) {
128 
129  const Part_gen* part = m_generator.partgen;
130 
131  // fill data of the final-state particles
132  for (int i = 0; i < m_generator.npart ; i++) {
133  auto& p = m_mpg.addParticle();
134  p.setPDG(part[i].part_prop.icode);
135  p.set4Vector(part[i].p);
136  p.setMass(part[i].part_prop.pmass);
138  p.setProductionVertex(vertex);
139  p.setValidVertex(true);
140  }
141  // fill data of the recoil electron and positron
142  auto& p1 = m_mpg.addParticle();
143  p1.setPDG(11);
144  p1.set4Vector(m_generator.pe);
145  p1.setMass(m_generator.me);
147  p1.setProductionVertex(vertex);
148  p1.setValidVertex(true);
149 
150  auto& p2 = m_mpg.addParticle();
151  p2.setPDG(-11);
152  p2.set4Vector(m_generator.pp);
153  p2.setMass(m_generator.me);
155  p2.setProductionVertex(vertex);
156  p2.setValidVertex(true);
157 
158  }
159  //Fill MCParticle List
161 
162 }
163 
165 {
166 }
167 
169 {
170  if (m_generator.diffCrossSectionOfW.size() == 0) {
171  B2FATAL("Cross Section Table is empty !!!");
172  return 0.;
173  }
174 
175  // lower_bound returns first iterator which meets >=W condition. --> upper side
176  auto it_upper = m_generator.diffCrossSectionOfW.lower_bound(W);
177  auto it_lower = it_upper;
178  it_lower--;
179 
180  return (it_upper->second - it_lower->second) / (it_upper->first - it_lower->first) * (W - it_lower->first) + it_lower->second;
181 }
182 
184 {
185 
186  while (1) {
187  double crossSectionForMC = gRandom->Uniform(0.0, m_generator.totalCrossSectionForMC);
188 
189  auto it_upper = m_generator.WOfCrossSectionForMC.lower_bound(crossSectionForMC);
190  auto it_lower = it_upper;
191  it_lower--;
192 
193  double diffCrossSectionAtUpper = m_generator.diffCrossSectionOfW.at(it_upper->second);
194  double diffCrossSectionAtLower = m_generator.diffCrossSectionOfW.at(it_lower->second);
195  double limit = (diffCrossSectionAtUpper > diffCrossSectionAtLower) ? diffCrossSectionAtUpper * 1.01 : diffCrossSectionAtLower *
196  1.01;
197  // Higher diffCrossSection * 1.01 is set as limit. Same definition of limit is used in TrepsB::wtable().
198 
199  double W = (crossSectionForMC - it_lower->first) / limit + it_lower->second ;
200  if (W < m_generator.diffCrossSectionOfW.begin()->first or W > m_generator.diffCrossSectionOfW.rbegin()->first)
201  B2FATAL("W has to be in [" << m_generator.diffCrossSectionOfW.begin()->first << ", "
202  << m_generator.diffCrossSectionOfW.rbegin()->first
203  << "] !!! W = " << W << ", crossSectionForMC = " << crossSectionForMC);
204 
205  double trial = gRandom->Uniform(0.0, limit);
206  double crossSection = getCrossSection(W);
207 
208  if (trial < crossSection)
209  return W;
210 
211  }
212 
213  return 0;
214 }
215 
216 
218 {
219  // Set parameter files
220  m_generator.setParameterFile(m_parameterFile);
221  m_generator.setDiffcrosssectionFile(m_differentialCrossSectionFile);
222  m_generator.setWlistFile(m_wListTableFile);
223 
224  // Set cut parameters
225  m_generator.setMaximalQ2(m_maximalQ2);
226  m_generator.setMaximalAbsCosTheta(m_maximalAbsCosTheta);
227  m_generator.applyCosThetaCutCharged(m_applyCosThetaCutCharged);
228  m_generator.setMinimalTransverseMomentum(m_minimalTransverseMomentum);
229  m_generator.applyTransverseMomentumCutCharged(m_applyTransverseMomentumCutCharged);
230 
231  // Initialize the initial particle information
232  TVector3 p3;
233  const BeamParameters& nominalBeam = m_initial.getBeamParameters();
234  m_generator.setBeamEnergy(nominalBeam.getMass() / 2.);
235  p3 = nominalBeam.getHER().Vect();
236  m_generator.setElectronMomentum(p3);
237  p3 = nominalBeam.getLER().Vect();
238  m_generator.setPositronMomentum(p3);
239 
240  // Initialize generator;
241  m_generator.initp();
242 
243  m_generator.create_hist();
244  m_generator.initg();
245 
247  // Initialize wtable with WListTable
248  B2INFO("Discrete W-list is used !!!");
249 
250  m_generator.wtcount = 0;
251  m_generator.wtable(0);
252 
253  m_generator.w = (double)m_generator.wf;
254 
255  m_generator.updateW();
256  } else {
257  // Initialize wtable with DifferentialCrossSection
258  m_generator.wtable();
259  }
260 
261  m_initialized = true;
262 }
Belle2::MCParticleGraph::generateList
void generateList(const std::string &name="", int options=c_setNothing)
Generates the MCParticle list and stores it in the StoreArray with the given name.
Definition: MCParticleGraph.cc:211
Belle2::MCInitialParticles::getLER
const TLorentzVector & getLER() const
Get 4vector of the low energy beam.
Definition: MCInitialParticles.h:140
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
Input from TREPS generator for ee->eeff.
Definition: TrepsInputModule.h:43
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
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
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::InitialParticleGeneration::getBeamParameters
const BeamParameters & getBeamParameters() const
Return reference to nominal beam parameters.
Definition: InitialParticleGeneration.h:66
Belle2::MCInitialParticles::getHER
const TLorentzVector & getHER() const
Get 4vector of the high energy beam.
Definition: MCInitialParticles.h:137
Belle2::TrepsInputModule::m_parameterFile
std::string m_parameterFile
Parameter file which configures the setting of beam, production particle etc.
Definition: TrepsInputModule.h:75
Belle2::Part_gen
Definition: Particle_array.h:50
Belle2::MCParticleGraph::c_setDecayInfo
@ c_setDecayInfo
Set decay time and vertex.
Definition: MCParticleGraph.h:74
Belle2::Module
Base class for Modules.
Definition: Module.h:74
Belle2::TrepsInputModule::m_initialized
bool m_initialized
True if generator has been initialized.
Definition: TrepsInputModule.h:117
Belle2::MCInitialParticles
This class contains the initial state for the given event.
Definition: MCInitialParticles.h:35
Belle2::TrepsInputModule::m_applyCosThetaCutCharged
bool m_applyCosThetaCutCharged
Whether to apply cut on |cos(theta)| for charged particles only.
Definition: TrepsInputModule.h:106
Belle2::MCInitialParticles::getMass
double getMass() const
Get the invariant mass of the collision (= energy in CMS)
Definition: MCInitialParticles.h:152
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::MCParticle::setPDG
void setPDG(int pdg)
Set PDG code of the particle.
Definition: MCParticle.h:343
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::MCParticle::c_Initial
@ c_Initial
bit 5: Particle is initial such as e+ or e- and not going to Geant4
Definition: MCParticle.h:68
Belle2::MCParticleGraph::addParticle
GraphParticle & addParticle()
Add new particle to the graph.
Definition: MCParticleGraph.h:305
Belle2::MCParticle::addStatus
void addStatus(unsigned short int bitmask)
Add bitmask to current status.
Definition: MCParticle.h:361
Belle2::MCInitialParticles::getVertex
const TVector3 & getVertex() const
Get the position of the collision.
Definition: MCInitialParticles.h:143
Belle2::MCParticle::c_StableInGenerator
@ c_StableInGenerator
bit 1: Particle is stable, i.e., not decaying in the generator.
Definition: MCParticle.h:60
Belle2::MCParticle::setEnergy
void setEnergy(float energy)
Set energy.
Definition: MCParticle.h:380
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::InitialParticleGeneration::initialize
void initialize()
function to be executed on initialize()
Definition: InitialParticleGeneration.cc:29
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::MCParticle::setMass
void setMass(float mass)
Set particle mass.
Definition: MCParticle.h:374
Belle2::TrepsInputModule::m_mcparticles
StoreArray< MCParticle > m_mcparticles
MCParticle collection.
Definition: TrepsInputModule.h:89
Belle2::UtrepsB::initg
void initg()
initialization of Pparametrization of pi+pi- partial waves
Definition: UtrepsB.cc:21
Belle2::TrepsInputModule::m_beamParams
DBObjPtr< BeamParameters > m_beamParams
BeamParameter.
Definition: TrepsInputModule.h:85
Belle2::InitialParticleGeneration::generate
MCInitialParticles & generate()
Generate a new event.
Definition: InitialParticleGeneration.cc:84
Belle2::MCParticleGraph::clear
void clear()
Reset particles and decay information to make the class reusable.
Definition: MCParticleGraph.h:298
Belle2::MCParticle::setMomentum
void setMomentum(const TVector3 &momentum)
Set particle momentum.
Definition: MCParticle.h:425
Belle2::BeamParameters
This class contains the nominal beam parameters and the parameters used for smearing of the primary v...
Definition: BeamParameters.h:33
Belle2::TrepsInputModule::getCrossSection
double getCrossSection(double W)
Returns the differential cross section for given W [GeV].
Definition: TrepsInputModule.cc:168
Belle2::MCParticleGraph::GraphParticle
Class to represent Particle data in graph.
Definition: MCParticleGraph.h:86
Belle2::MCParticle::c_PrimaryParticle
@ c_PrimaryParticle
bit 0: Particle is primary particle.
Definition: MCParticle.h:58
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