Belle II Software  release-06-01-15
PhokharaInputModule.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/phokharainput/PhokharaInputModule.h>
10 
11 #include <framework/datastore/StoreArray.h>
12 #include <framework/utilities/FileSystem.h>
13 #include <framework/utilities/IOIntercept.h>
14 
15 using namespace std;
16 using namespace Belle2;
17 
18 //-----------------------------------------------------------------
19 // Register the Module
20 //-----------------------------------------------------------------
21 REG_MODULE(PhokharaInput)
22 
23 //-----------------------------------------------------------------
24 // Implementation
25 //-----------------------------------------------------------------
27 {
28  //Set module properties
29  setDescription("Generates radiative return events with PHOKHARA 10.");
30 
31  //Parameter definition
32  addParam("FinalState", m_finalState,
33  "Final state: mu+mu-(0), pi+pi-(1), 2pi0pi+pi-(2), 2pi+2pi-(3), ppbar(4), nnbar(5), K+K-(6), K0K0bar(7), pi+pi-pi0(8), lamb(->pi-p)lambbar(->pi+pbar)(9), eta pi+ pi- (10), chi_c1 production (11), chi_c2 production (12), pi0 g (13), eta g (14), etaP g (15)",
34  1);
35  addParam("ReplaceMuonsByVirtualPhoton", m_replaceMuonsByVirtualPhoton,
36  "Replace muons by a virtual photon (for FinalState == 0 only).",
37  false);
38  addParam("SearchMax", m_nSearchMax, "Number of events used to search for maximum of differential cross section", 100000);
39  addParam("Epsilon", m_epsilon, "Soft/hard photon separator", 0.0001);
40  addParam("nMaxTrials", m_nMaxTrials, "Maximum trials per event", 10000);
41  addParam("Weighted", m_weighted, "generate weighted events", 0);
42 
43  addParam("LO", m_LO, "ph0 Born: 1ph(0), Born: 0ph(1), only Born: 0ph(-1)", 0);
44  addParam("NLO", m_NLO, "1 photon : Born(0), NLO(1)", 0);
45  addParam("FullNLO", m_fullNLO, "full NLO : No(0), Yes(1). Allowed only for ph0=1, nlo=1, fsr=2, fsrnlo=1", 0);
46  addParam("QED", m_QED, "ISR only(0, default), ISR+FSR(1), ISR+INT+FSR(2)", 0);
47  addParam("IFSNLO", m_IFSNLO, "IFSNLO: no(0), yes(1)", 0);
48  addParam("Alpha", m_alpha, "Vacuum polarization switch: off (0), on (1,[by Fred Jegerlehner], default), on (2,[by Thomas Teubner])",
49  1);
50  addParam("PionFF", m_pionff, "Pion FF: KS PionFormFactor(0, default), GS old 1), GS new(2)", 0);
51  addParam("PionStructure", m_pionstructure,
52  "For pi+pi- only: f0+f0(600): K+K- model(0, default), no structure model(1), no f0+f0(600)(2), f0 KLOE(3)", 0);
53  addParam("KaonFF", m_kaonff,
54  "Kaon FF: KaonFormFactor constrained(0, default), KaonFormFactor unconstrained(1) KaonFormFactor old(2)", 0);
55  addParam("NarrowRes", m_narres, "Only for m_finalState = 0,1,6,7: No narrow resonances (0, default), J/Psi (1) and Psi(2S) (2)", 0);
56  addParam("ProtonFF", m_protonff, "ProtonFormFactor old(0), ProtonFormFactor new(1)", 1);
57  addParam("ChiSW", m_chi_sw,
58  "chi_sw: Radiative return(0), Chi production(1), Radiative return + Chi production (2). Works only for pion=11 and pion=12.", 0);
59  addParam("SwitchBeamResolution", m_be_r,
60  "be_r: without beam resolution(0), with beam resolution(1). Works only for pion=11 and pion=12; (0) assumes exact CMS-Energy, (1) each beam ennergy (=CMS-Energy/2) is smeared with Gaussian distribution of the given variance=BeamResolution**2.",
61  0);
62  addParam("BeamResolution", m_beamres, "beamres - beam resolution for pion==11 and pion==12 only", 0.);
63 
64  addParam("ScatteringAngleRangePhoton", m_ScatteringAngleRangePhoton,
65  "Min [0] and Max [1] value for the scattering angle of photons [deg], default (0, 180)", make_vector(0.0, 180.0));
66  addParam("ScatteringAngleRangeFinalStates", m_ScatteringAngleRangeFinalStates,
67  "Min [0] and Max [1] value for the scattering angle of pions(muons,nucleons,kaons) [deg], default (0, 180)", make_vector(0.0,
68  180.0));
69 
70  addParam("MinInvMassHadronsGamma", m_MinInvMassHadronsGamma, "Minimal hadrons/muons-gamma invariant mass squared [GeV^2]", 0.0);
71  addParam("MinInvMassHadrons", m_MinInvMassHadrons, "Minimal hadrons/muons invariant mass squared [GeV^2]", 0.2);
72  addParam("ForceMinInvMassHadronsCut", m_ForceMinInvMassHadronsCut,
73  "Force application of the MinInvMassHadrons cut "
74  "It is ignored by PHOKHARA with LO = 1, NLO = 1.",
75  false);
76  addParam("MaxInvMassHadrons", m_MaxInvMassHadrons, "Maximal hadrons/muons invariant mass squared [GeV^2]", 112.0);
77  addParam("MinEnergyGamma", m_MinEnergyGamma, "Minimal photon energy/missing energy, must be greater than 0.0098 * CMS energy [GeV]",
78  0.15);
79 
80  const std::string defaultParameterFile =
81  FileSystem::findFile("/data/generators/phokhara/const_and_model_paramall10.0.dat");
82  addParam("ParameterFile", m_ParameterFile,
83  "File that contains particle properties.",
84  defaultParameterFile);
85  addParam("BeamEnergySpread", m_BeamEnergySpread,
86  "Simulate beam-energy spread (initializes PHOKHARA for every "
87  "event - very slow).", false);
88 // addParam("InputFile", m_InputFile, "File that contain input configuration", FileSystem::findFile("/data/generators/phokhara/input_9.1.dat"));
89 
90 }
91 
92 //-----------------------------------------------------------------
93 // Destructor
94 //-----------------------------------------------------------------
95 PhokharaInputModule::~PhokharaInputModule()
96 {
97 
98 }
99 
100 //-----------------------------------------------------------------
101 // Initialize
102 //-----------------------------------------------------------------
103 void PhokharaInputModule::initialize()
104 {
105  StoreArray<MCParticle> mcparticle;
106  mcparticle.registerInDataStore();
107 
108  if (m_replaceMuonsByVirtualPhoton) {
109  if (m_finalState != 0) {
110  B2FATAL("You requested to replace muons by a virtual photon, but the "
111  "final state is not mu+ mu-.");
112  }
113  if (m_QED != 0) {
114  B2FATAL("You requested to replace muons by a virtual photon. In this "
115  "mode, PHOKHARA works as an ISR generator. The parameter QED "
116  "should be set to 0 (ISR only). If FSR is taken into account "
117  "(QED = 1 or 2), the results will be incorrect.");
118  }
119  if (m_IFSNLO != 0) {
120  B2FATAL("You requested to replace muons by a virtual photon. In this "
121  "mode, PHOKHARA works as an ISR generator. The parameter IFSNLO "
122  "should be set to 0 (off). If simultaneous emission of initial "
123  "and final-state photons is taken into account (IFSNLO = 1), "
124  "the results will be incorrect.");
125  }
126  }
127  if (m_finalState == 11 || m_finalState == 12) {
128  if (m_narres != 1 || m_LO != 0 || m_NLO != 0) {
129  B2FATAL("Generation of the chi_c1 or chi_c2 (final states 11 and 12, "
130  "respectively) requires to turn on the J/psi (NarrowRes = 1) "
131  "and use LO radiative-return mode (LO = 0, NLO = 0)");
132  }
133  }
134  //Beam Parameters, initial particle - PHOKHARA cannot handle beam energy spread
135  if (m_BeamEnergySpread) {
136  m_initial.setAllowedFlags(BeamParameters::c_smearVertex |
137  BeamParameters::c_smearBeam);
138  }
139  m_initial.initialize();
140 
141 }
142 
143 //-----------------------------------------------------------------
144 // Event
145 //-----------------------------------------------------------------
146 void PhokharaInputModule::event()
147 {
148 
149  // Check if the BeamParameters have changed (if they do, abort the job! otherwise cross section calculation will be a nightmare.)
150  if (m_beamParams.hasChanged()) {
151  if (!m_initialized) {
152  initializeGenerator();
153  } else {
154  B2FATAL("PhokharaInputModule::event(): BeamParameters have changed within a job, this is not supported for PHOKHARA!");
155  }
156  }
157 
158  StoreObjPtr<EventMetaData> evtMetaData("EventMetaData", DataStore::c_Event);
159 
160  // initial particle from beam parameters
161  const MCInitialParticles& initial = m_initial.generate();
162 
163  // true boost
164  TLorentzRotation boost = initial.getCMSToLab();
165 
166  // vertex
167  TVector3 vertex = initial.getVertex();
168 
169  m_mcGraph.clear();
170  if (m_BeamEnergySpread) {
171  // PHOKHARA does not support beam-energy spread. To simulate it,
172  // the generator is initialized with the new energy for every event.
173  IOIntercept::OutputToLogMessages initLogCapture(
174  "PHOKHARA", LogConfig::c_Debug, LogConfig::c_Debug, 100, 100);
175  initLogCapture.start();
176  m_generator.setCMSEnergy(initial.getMass());
177  m_generator.init(m_ParameterFile);
178  initLogCapture.finish();
179  }
180  double weight = m_generator.generateEvent(m_mcGraph, vertex, boost);
181  m_mcGraph.generateList("", MCParticleGraph::c_setDecayInfo | MCParticleGraph::c_checkCyclic);
182 
183  //store the weight (1.0 for unweighted events)
184  evtMetaData->setGeneratedWeight(weight);
185 }
186 
187 
188 //-----------------------------------------------------------------
189 // Terminate
190 //-----------------------------------------------------------------
191 void PhokharaInputModule::terminate()
192 {
193 
194  m_generator.term();
195 }
196 
197 void PhokharaInputModule::initializeGenerator()
198 {
199 
200  const BeamParameters& nominal = m_initial.getBeamParameters();
201  m_cmsEnergy = nominal.getMass();
202 
203  m_generator.setCMSEnergy(m_cmsEnergy);
204 
205  m_generator.setNSearchMax(m_nSearchMax);
206  m_generator.setWeighted(m_weighted);
207  m_generator.setFinalState(m_finalState);
208  m_generator.setReplaceMuonsByVirtualPhoton(m_replaceMuonsByVirtualPhoton);
209  m_generator.setNMaxTrials(m_nMaxTrials);
210  m_generator.setEpsilon(m_epsilon);
211  m_generator.setLO(m_LO);
212  m_generator.setNLO(m_NLO);
213  m_generator.setFullNLO(m_fullNLO);
214  m_generator.setQED(m_QED);
215  m_generator.setIFSNLO(m_IFSNLO);
216  m_generator.setAlpha(m_alpha);
217  m_generator.setPionFF(m_pionff);
218  m_generator.setKaonFF(m_kaonff);
219  m_generator.setPionStructure(m_pionstructure);
220  m_generator.setNarrowRes(m_narres);
221  m_generator.setProtonFF(m_protonff);
222  m_generator.setChiSW(m_chi_sw);
223  m_generator.setSwitchBeamResolution(m_be_r);
224  m_generator.setBeamResolution(m_beamres);
225 
226 
227  m_generator.setMinInvMassHadronsGamma(m_MinInvMassHadronsGamma);
228  m_generator.setm_MinInvMassHadrons(m_MinInvMassHadrons);
229  m_generator.setForceMinInvMassHadronsCut(m_ForceMinInvMassHadronsCut);
230  m_generator.setm_MaxInvMassHadrons(m_MaxInvMassHadrons);
231  m_generator.setMinEnergyGamma(m_MinEnergyGamma);
232  m_generator.setScatteringAngleRangePhoton(vectorToPair(m_ScatteringAngleRangePhoton, "ScatteringAngleRangePhoton"));
233  m_generator.setScatteringAngleRangeFinalStates(vectorToPair(m_ScatteringAngleRangeFinalStates, "ScatteringAngleRangeFinalStates"));
234 
235  m_generator.init(m_ParameterFile);
236 
237  m_initialized = true;
238 
239 }
This class contains the nominal beam parameters and the parameters used for smearing of the primary v...
bool start()
Start intercepting the output.
Definition: IOIntercept.h:155
Capture stdout and stderr and convert into log messages.
Definition: IOIntercept.h:226
bool finish()
Finish the capture and emit the message if output has appeard on stdout or stderr.
Definition: IOIntercept.cc:235
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.
double getMass() const
Get the invariant mass of the collision (= energy in CMS)
Base class for Modules.
Definition: Module.h:72
The Phokhara Generator module.vectorToPair.
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:95
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.