Belle II Software development
ParticleGun.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 <framework/logging/Logger.h>
10#include <framework/gearbox/Unit.h>
11#include <generators/particlegun/ParticleGun.h>
12
13#include <TRandom.h>
14
15#include <cmath>
16#include <limits>
17
18using namespace std;
19using namespace Belle2;
20
21namespace {
37 double randomPolyline(size_t n, const double* x, const double* y)
38 {
39 //Calculate the area of the step function which envelops the polyline and
40 //the area in each segment.
41 std::vector<double> weights(n - 1);
42 double sumw(0);
43 for (size_t i = 0; i < n - 1; ++i) {
44 weights[i] = (x[i + 1] - x[i]) * max(y[i], y[i + 1]);
45 sumw += weights[i];
46 }
47 while (true) {
48 //Generate x distributed according to the step function. This can be done
49 //with 100% efficiency
50 double weight = gRandom->Uniform(0, sumw);
51 size_t segment(0);
52 for (; segment < n - 1; ++segment) {
53 weight -= weights[segment];
54 if (weight <= 0) break;
55 }
56 const double x1 = x[segment];
57 const double x2 = x[segment + 1];
58 const double xx = gRandom->Uniform(x1, x2);
59 //Now calculate y(x) using line between (x1,y1) and (x2,y2)
60 const double y1 = y[segment];
61 const double y2 = y[segment + 1];
62 const double yy = (y2 == y1) ? y1 : y1 + (xx - x1) / (x2 - x1) * (y2 - y1);
63 //Generate a random y inside the step function
64 const double randY = gRandom->Uniform(0, max(y1, y2));
65 //And accept the x value if randY lies below the line
66 if (randY < yy) return xx;
67 //Otherwise repeat generation of x and y
68 }
69 }
70}
71
72
73double ParticleGun::generateValue(EDistribution dist, const vector<double>& params)
74{
75 double rand(0);
76 switch (dist) {
77 case c_fixedValue:
78 return params[0];
81 return gRandom->Uniform(params[0], params[1]);
83 return 1. / gRandom->Uniform(1. / params[1], 1. / params[0]);
85 return acos(gRandom->Uniform(cos(params[0]), cos(params[1])));
88 return exp(gRandom->Uniform(log(params[0]), log(params[1])));
91 return gRandom->Gaus(params[0], params[1]);
93 return acos(gRandom->Gaus(params[0], params[1]));
96 //Create weighted discrete distribution
97 //First we sum all the weights (second half of the array)
98 for (size_t i = params.size() / 2; i < params.size(); ++i) rand += params[i];
99 //Then generate a random variable between 0 and sum of weights
100 rand = gRandom->Uniform(0, rand);
101 //Go over the weights again and subtract from generated value
102 for (size_t i = params.size() / 2; i < params.size(); ++i) {
103 rand -= params[i];
104 //If the value is negative we found the correct bin, return that value
105 if (rand <= 0) return params[i - (params.size() / 2)];
106 }
107 B2FATAL("Something wrong with picking fixed spectra values");
110 return randomPolyline(params.size() / 2, params.data(), params.data() + params.size() / 2);
112 return acos(randomPolyline(params.size() / 2, params.data(), params.data() + params.size() / 2));
113 default:
114 B2FATAL("Unknown distribution");
115 }
116 return 0;
117}
118
120{
121 //generate the event vertex (possible the same for all particles in event)
125
126 // Time offset
127 double timeOffset = generateValue(m_params.timeDist, m_params.timeParams);
128
129 //Determine number of tracks
130 int nTracks = static_cast<int>(m_params.nTracks);
131 if (m_params.nTracks <= 0) {
132 nTracks = m_params.pdgCodes.size();
133 } else if (m_params.varyNumberOfTracks) {
134 nTracks = gRandom->Poisson(m_params.nTracks);
135 }
136
137 //Make list of particles
138 double firstMomentum;
139 for (int i = 0; i < nTracks; ++i) {
142 if (m_params.pdgCodes.size() == 1) {
143 //only one PDGcode available, always take this one
144 p.setPDG(m_params.pdgCodes[0]);
145 } else if (m_params.nTracks <= 0) {
146 //0 or negative nTracks, take the ids sequentially
147 p.setPDG(m_params.pdgCodes[i]);
148 } else {
149 //else choose randomly one of the available codes
150 int index = static_cast<int>(gRandom->Uniform(m_params.pdgCodes.size()));
151 p.setPDG(m_params.pdgCodes[index]);
152 }
153 p.setMassFromPDG();
154 p.setFirstDaughter(0);
155 p.setLastDaughter(0);
156
157 //lets generate the momentum vector:
159 // if the fixed momentum option enabled, either store the generated momentum for the first
160 // particle or reuse for all others
162 i == 0 ? firstMomentum = momentum : momentum = firstMomentum;
163 }
166
167 double pt = momentum * sin(theta);
171 //this means we are actually generating the Pt and not the P, so exchange values
172 pt = momentum;
173 momentum = (sin(theta) > 0) ? (pt / sin(theta)) : numeric_limits<double>::max();
174 }
175
176 double pz = momentum * cos(theta);
177 double px = pt * cos(phi);
178 double py = pt * sin(phi);
179 double m = p.getMass();
180 double e = sqrt(momentum * momentum + m * m);
181
182 p.setMomentum(px, py, pz);
183 p.setEnergy(e);
184 p.setProductionVertex(vx, vy, vz);
186 //Particle is stable in generator. We could use MCParticleGraph options to
187 //do this automatically but setting it here makes the particle correct
188 //independent of the options
189 p.setDecayTime(numeric_limits<double>::infinity());
190
191 // set time offset to check fit bias in e.g. the ECL waveform fits
192 p.setProductionTime(timeOffset);
193
195 //If we have independent vertices, generate new vertex for next particle
199 }
200 }// end loop over particles in event
201
202 return true;
203}
204
205bool ParticleGun::setParameters(const Parameters& p)
206{
207 //checking that all parameters are useful is the most complex part of the
208 //whole particle gun. We need to disallow distributions for certain variables
209 //and check the number of parameters and their validity
210
211 //As long as this is true, no error has been found. We could just return but
212 //it's so much nicer to display all possible errors at once so the user does
213 //not need to iterate over this all the time
214 bool ok(true);
215
216 //Make an enum -> name mapping for nice error messages
217 std::map<EDistribution, std::string> distributionNames = {
218 {c_fixedValue, "fixed"},
219 {c_uniformDistribution, "uniform"},
220 {c_uniformPtDistribution, "uniformPt"},
221 {c_uniformCosDistribution, "uniformCos"},
222 {c_uniformLogDistribution, "uniformLog"},
223 {c_uniformLogPtDistribution, "uniformLogPt"},
224 {c_normalDistribution, "normal"},
225 {c_normalPtDistribution, "normalPt"},
226 {c_normalCosDistribution, "normalCos"},
227 {c_inversePtDistribution, "inversePt"},
228 {c_polylineDistribution, "polyline"},
229 {c_polylinePtDistribution, "polylinePt"},
230 {c_polylineCosDistribution, "polylineCos"},
231 {c_discreteSpectrum, "discrete"},
232 {c_discretePtSpectrum, "discretePt"}
233 };
234 //Small helper lambda to get the distribution by name
235 auto getDist = [&p](const std::string & dist) -> EDistribution {
236 if (dist == "momentum") return p.momentumDist;
237 if (dist == "xVertex") return p.xVertexDist;
238 if (dist == "yVertex") return p.yVertexDist;
239 if (dist == "zVertex") return p.zVertexDist;
240 if (dist == "theta") return p.thetaDist;
241 if (dist == "phi") return p.phiDist;
242 if (dist == "timeOffset") return p.timeDist;
243 throw std::runtime_error("wrong parameter");
244 };
245 //Small helper lambda to get the parameters by name
246 auto getPars = [&p](const std::string & dist) -> const std::vector<double>& {
247 if (dist == "momentum") return p.momentumParams;
248 if (dist == "xVertex") return p.xVertexParams;
249 if (dist == "yVertex") return p.yVertexParams;
250 if (dist == "zVertex") return p.zVertexParams;
251 if (dist == "theta") return p.thetaParams;
252 if (dist == "phi") return p.phiParams;
253 if (dist == "timeOffset") return p.timeParams;
254 throw std::runtime_error("wrong parameter");
255 };
256 //Small helper lambda to produce a nice error message and set the error flag
257 //if the distribution is excluded.
258 auto excludeDistribution = [&](const std::string & dist, EDistribution excluded) {
259 if (getDist(dist) == excluded) {
260 B2ERROR(distributionNames[excluded] << " is not allowed for " << dist << " generation");
261 ok = false;
262 }
263 };
264
265 //Exclude some distributions
266 for (auto dist : {"momentum", "xVertex", "yVertex", "zVertex"}) {
267 excludeDistribution(dist, c_uniformCosDistribution);
268 excludeDistribution(dist, c_normalCosDistribution);
269 excludeDistribution(dist, c_polylineCosDistribution);
270 }
271 for (auto dist : {"xVertex", "yVertex", "zVertex", "phi", "theta", "timeOffset"}) {
272 excludeDistribution(dist, c_uniformPtDistribution);
273 excludeDistribution(dist, c_uniformLogPtDistribution);
274 excludeDistribution(dist, c_normalPtDistribution);
275 excludeDistribution(dist, c_inversePtDistribution);
276 excludeDistribution(dist, c_polylinePtDistribution);
277 }
278
279 //Check that we have some particle ids to generate
280 if (p.pdgCodes.empty()) {
281 B2ERROR("No pdgCodes specified for generation");
282 ok = false;
283 }
284
285 //Check minimum numbers of parameters
286 for (auto par : {"momentum", "xVertex", "yVertex", "zVertex", "theta", "phi", "timeOffset"}) {
287 const EDistribution dist = getDist(par);
288 size_t minParams = (dist == c_fixedValue) ? 1 : 2;
289 if (dist == c_polylineDistribution || dist == c_polylinePtDistribution || dist == c_polylineCosDistribution) minParams = 4;
290 const size_t hasParams = getPars(par).size();
291 if (hasParams < minParams) {
292 B2ERROR(par << " generation: " << distributionNames[dist]
293 << " distribution requires at least " << minParams
294 << " parameters, " << hasParams << " given");
295 ok = false;
296 }
297 //Check even number of parameters for discrete or polyline distributions
298 if (dist == c_discreteSpectrum || dist == c_discretePtSpectrum ||
300 if ((hasParams % 2) != 0) {
301 B2ERROR(par << " generation: " << distributionNames[dist] << " requires an even number of parameters");
302 ok = false;
303 } else if (dist == c_discreteSpectrum || dist == c_discretePtSpectrum || hasParams >= 4) {
304 //Check wellformedness of polyline pdf: ascending x coordinates and positive y coordinates with at least one nonzero value
305 //Discrete spectrum only requires positive weights, no sorting needed
306 const std::vector<double>& pvec = getPars(par);
307 const std::string parname = (dist == c_discreteSpectrum || dist == c_discretePtSpectrum)
308 ? "weight" : "y coordinate";
309 //Check for sorting for polylines
310 if ((dist != c_discreteSpectrum) && (dist != c_discretePtSpectrum)) {
311 for (size_t i = 0; i < (hasParams / 2) - 1; ++i) {
312 if (pvec[i] > pvec[i + 1]) {
313 B2ERROR(par << " generation: " << distributionNames[dist] << " requires x coordinates in ascending order");
314 ok = false;
315 }
316 }
317 }
318 //Check that values are non-negative and at least one value is positive
319 auto minmax = std::minmax_element(pvec.begin() + hasParams / 2, pvec.end());
320 if (*minmax.first < 0) {
321 B2ERROR(par << " generation: " << distributionNames[dist] << " requires "
322 << parname << "s to be non-negative");
323 ok = false;
324 }
325 if (*minmax.second <= 0) {
326 B2ERROR(par << " generation: " << distributionNames[dist] << " requires at least one "
327 << parname << " to be positive");
328 ok = false;
329 }
330 }
331 }
332 }
333
334 //Finally check that we do not have any problems with the inverse
335 if (p.momentumDist == c_inversePtDistribution && p.momentumParams.size() >= 2) {
336 if (p.momentumParams[0] == 0 || p.momentumParams[1] == 0) {
337 B2ERROR("inversePt distribution does not allow zero momentum");
338 ok = false;
339 }
340 }
341
342 //If everything is ok, set the new parameters, else return false
343 if (ok) { m_params = p; }
344 return ok;
345}
Class to represent Particle data in graph.
Class to build, validate and sort a particle decay chain.
@ c_PrimaryParticle
bit 0: Particle is primary particle.
Definition: MCParticle.h:47
@ c_StableInGenerator
bit 1: Particle is stable, i.e., not decaying in the generator.
Definition: MCParticle.h:49
EDistribution
enum containing all known distributions available for generation of values
Definition: ParticleGun.h:29
@ c_polylineDistribution
Distribution given as list of (x,y) points.
Definition: ParticleGun.h:57
@ c_uniformLogPtDistribution
Uniform distribution of the logarithm of the Pt, parameters are min and max value.
Definition: ParticleGun.h:41
@ c_normalDistribution
Normal distribution, parameters are mean and sigma.
Definition: ParticleGun.h:43
@ c_normalCosDistribution
Normal distribution of the cosinge, parameters are mean and sigma.
Definition: ParticleGun.h:47
@ c_uniformPtDistribution
Uniform distribution of Pt, parameters are min and max value.
Definition: ParticleGun.h:35
@ c_polylinePtDistribution
Same as polylineDistribution but for the transverse momentum.
Definition: ParticleGun.h:59
@ c_discretePtSpectrum
Discrete pt spectrum, parameters are first the values and then the weights (non-negative) for each va...
Definition: ParticleGun.h:51
@ c_uniformDistribution
Uniform distribution, parameters are min and max value.
Definition: ParticleGun.h:33
@ c_discreteSpectrum
Discrete spectrum, parameters are first the values and then the weights (non-negative) for each value...
Definition: ParticleGun.h:49
@ c_fixedValue
Fixed value, no random generation at all, 1 parameter.
Definition: ParticleGun.h:31
@ c_normalPtDistribution
Normal distribution of Pt, parameters are mean and sigma.
Definition: ParticleGun.h:45
@ c_uniformCosDistribution
Uniform distribution of the cosine of the values, parameters are min and max value.
Definition: ParticleGun.h:37
@ c_inversePtDistribution
Distribution uniform in the inverse pt distribution, that is uniform in track curvature.
Definition: ParticleGun.h:53
@ c_polylineCosDistribution
Same as polylineDistribution but for the cosine of the angle.
Definition: ParticleGun.h:61
@ c_uniformLogDistribution
Uniform distribution of the logarithm of the values, parameters are min and max value.
Definition: ParticleGun.h:39
Parameters m_params
All relevant parameters.
Definition: ParticleGun.h:133
bool generateEvent(MCParticleGraph &graph)
Generate the next event and store the result in the given MCParticle graph.
bool setParameters(const Parameters &parameters)
Set the parameters for generating the Particles.
double generateValue(EDistribution dist, const std::vector< double > &params)
Generate a value according to the given distribution with the given parameters.
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
GraphParticle & addParticle()
Add new particle to the graph.
Abstract base class for different kinds of events.
STL namespace.
std::vector< int > pdgCodes
List of PDG particle codes to pick from when generating particles.
Definition: ParticleGun.h:83
EDistribution xVertexDist
Distribution to use for x vertex generation.
Definition: ParticleGun.h:73
std::vector< double > thetaParams
Parameters for the polar angle generation, meaning depends on chosen distribution.
Definition: ParticleGun.h:89
std::vector< double > momentumParams
Parameters for the momentum generation, meaning depends on chosen distribution.
Definition: ParticleGun.h:85
EDistribution yVertexDist
Distribution to use for y vertex generation.
Definition: ParticleGun.h:75
std::vector< double > xVertexParams
Parameters for the x vertex generation, meaning depends on chosen distribution.
Definition: ParticleGun.h:91
EDistribution thetaDist
Distribution to use for polar angle generation.
Definition: ParticleGun.h:71
EDistribution zVertexDist
Distribution to use for z vertex generation.
Definition: ParticleGun.h:77
EDistribution timeDist
Distribution to use for time generation.
Definition: ParticleGun.h:79
EDistribution momentumDist
Distribution to use for momentum generation.
Definition: ParticleGun.h:67
EDistribution phiDist
Distribution to use for azimuth angle generation.
Definition: ParticleGun.h:69
double nTracks
Number of tracks to generate per event.
Definition: ParticleGun.h:81
bool independentVertices
If false, all particles of one event will have the same vertex, if true the vertex of each particle w...
Definition: ParticleGun.h:102
std::vector< double > zVertexParams
Parameters for the z vertex generation, meaning depends on chosen distribution.
Definition: ParticleGun.h:95
std::vector< double > phiParams
Parameters for the azimuth angle generation, meaning depends on chosen distribution.
Definition: ParticleGun.h:87
bool varyNumberOfTracks
If true, the number of tracks per event will fluctuate according to Poisson distribution.
Definition: ParticleGun.h:105
std::vector< double > timeParams
Parameters for the time generation, meaning depends on chosen distribution.
Definition: ParticleGun.h:97
bool fixedMomentumPerEvent
generates particle momentum according to the specified distribution and assigns this momentum to all ...
Definition: ParticleGun.h:109
std::vector< double > yVertexParams
Parameters for the y vertex generation, meaning depends on chosen distribution.
Definition: ParticleGun.h:93