Belle II Software development
AAFHInterface.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/aafh/AAFHInterface.h>
12#include <TRandom3.h>
13#include <TDatabasePDG.h>
14#include <Math/Vector4D.h>
15#include <algorithm>
16#include <array>
17#include <iomanip>
18
19using namespace Belle2;
20
21extern "C" {
23 double aafhdiag36_rndm_(int*)
24 {
25 return gRandom->Rndm();
26 }
28 double aafhdiag36_rnf100_(int*)
29 {
30 return gRandom->Rndm();
31 }
33 void aafhdiag36_esft_maxed_(const double* weight)
34 {
35 B2ERROR("AAFH: Maximum weight to small, increase maxFinalWeight to at least "
36 << *weight);
37 }
39 void aafhdiag36_eswe_maxed_(const double* weight)
40 {
41 B2ERROR("AAFH: Maximum weight to small, increase maxSubGeneratorWeight to at least "
42 << *weight);
43 }
44
46 extern void aafhdiag36_init_();
48 extern void aafhdiag36_event_(int*);
50 extern void aafhdiag36_finish_();
51
54 extern struct {
55 double eb;
56 double eswe;
57 double esft;
58 double wap[4];
59 double wbp[4];
60 double vap[4];
61 } input_;
62
65 extern struct {
66 double w2min;
67 double w2max;
68 } bound_;
69
71 extern struct {
72 int iproc;
73 int irejec;
74 int idump[4];
75 int info;
76 int itot;
77 } input2_;
78
80 extern struct {
81 double face;
82 double facl;
83 double facm;
84 double proc;
85 } factor_;
86
90 extern struct {
91 double q1[5];
92 double q2[5];
93 double q3[5];
94 double q4[5];
95 double q5[5];
96 double q6[5];
97 } momenz_;
98
101 extern struct {
102 double crosec;
103 double ercros;
104 } output_;
105
107 extern struct {
108 double xm;
109 double xmu;
110 double xml;
111 double xm2;
112 double xmu2;
113 double xml2;
114 } masses_;
115
117 extern struct {
118 double qcharg;
119 double qchrg2;
120 double qchrg3;
121 double qchrg4;
122 } charge_;
123
125 extern struct {
126 double swe[4];
127 double swek[4];
128 double mwe[4];
129 double sum;
130 double sumk;
131 double maxwe;
132 int iwe[4];
133 int igen;
134 } westat_;
135
137 extern struct {
138 double sumft;
139 double sumft2;
140 double ftmax;
141 int ieen;
142 } ftstat_;
143};
144
145
147{
148 input2_.itot = tries;
149}
150
151void AAFHInterface::setGeneratorWeights(const std::vector<double>& weights)
152{
153 if (weights.size() > 8) {
154 B2WARNING("AAFH: more than 8 generator weights supplied, "
155 "ignoring the extra weights");
156 } else if (weights.size() > 4 && weights.size() < 8) {
157 B2WARNING("AAFH: more than 4 but less than 8 generator weights supplied, "
158 "ignoring everything afer the first four");
159 } else if (weights.size() < 4) {
160 B2ERROR("AAFH: not enough generator weights, need exactly four or eight values");
161 }
162 std::copy_n(weights.begin(), 4, input_.wap);
163 if (weights.size() >= 8) {
164 std::copy_n(weights.begin() + 4, 4, input_.wbp);
165 }
166}
167
168void AAFHInterface::setMaxWeights(double subgeneratorWeight, double finalWeight)
169{
170 input_.eswe = subgeneratorWeight;
171 input_.esft = finalWeight;
172}
173
174void AAFHInterface::setSupressionLimits(std::vector<double> limits)
175{
176 if (limits.size() == 1) {
177 B2WARNING("AAFH: Only one suppression limit supplied, using same value for all limits");
178 factor_.face = limits[0];
179 factor_.facm = limits[0];
180 factor_.facl = limits[0];
181 factor_.proc = limits[0];
182 } else if (limits.size() == 4) {
183 factor_.face = limits[0];
184 factor_.facm = limits[1];
185 factor_.facl = limits[2];
186 factor_.proc = limits[3];
187 } else {
188 B2ERROR("AAFH: Wrong number of suppression limits: supply either one or 4 values");
189 }
190}
191
193{
194 return input2_.itot;
195}
196
197std::vector<double> AAFHInterface::getGeneratorWeights() const
198{
199 std::vector<double> weights(8);
200 std::copy_n(input_.wap, 4, weights.begin());
201 std::copy_n(input_.wbp, 4, weights.begin() + 4);
202 return weights;
203}
204
206{
207 return input_.eswe;
208}
209
211{
212 return input_.esft;
213}
214
215std::vector<double> AAFHInterface::getSuppressionLimits() const
216{
217 return {factor_.face, factor_.facm, factor_.facl, factor_.proc};
218}
219
221{
222 return output_.crosec * Unit::nb;
223}
224
226{
227 return output_.ercros * Unit::nb;
228}
229
230void AAFHInterface::initialize(double beamEnergy, EMode mode, ERejection rejection)
231{
232 input_.eb = beamEnergy;
233 //minimum invariant mass is in units of beam energy and needs to be squared
234 bound_.w2min = (m_minMass / beamEnergy) * (m_minMass / beamEnergy);
235 //set mode and rejection scheme
236 input2_.iproc = mode;
237 input2_.irejec = rejection;
238
239 //initialize the L
240 auto* part = TDatabasePDG::Instance()->GetParticle(m_particle.c_str());
241 if (!part) {
242 B2ERROR("AAFH: Particle '" << m_particle << "' does not exist");
243 return;
244 }
245 //set the mass in units of beam energy
246 masses_.xml = part->Mass() / beamEnergy;
247 //set the charge but with opposite sign
248 charge_.qcharg = -part->Charge() / 3.0;
249 //and remember pdg code
250 m_particlePDG = part->PdgCode();
251
252 aafhdiag36_init_();
253}
254
255void AAFHInterface::updateParticle(MCParticleGraph::GraphParticle& p, double* q, int pdg, bool isInitial)
256{
257 // all values in q are in units of beam energy
258 // z-direction inversed since AAFH uses different convention
259 ROOT::Math::PxPyPzEVector vec(-q[0]*input_.eb, -q[1]*input_.eb, -q[2]*input_.eb, q[3]*input_.eb);
260 p.set4Vector(vec);
261
262 // mass is multiplied by charge sign so take the absolute
263 p.setMass(fabs(q[4])*input_.eb);
264 p.setPDG(pdg);
266
267 if (isInitial) {
268 p.addStatus(MCParticle::c_Initial);
269 }
270}
271
273{
274 int status {0};
275 aafhdiag36_event_(&status);
276 if (status != 0) {
277 B2ERROR("Failed to generate event after " << status << " tries, giving up");
278 return;
279 }
280
281 //initial beam particles
282 auto& p1 = mpg.addParticle();
283 auto& p2 = mpg.addParticle();
284
285 //generated particles
286 auto& p3 = mpg.addParticle();
287 auto& p4 = mpg.addParticle();
288 auto& p5 = mpg.addParticle();
289 auto& p6 = mpg.addParticle();
290
291 //array of generated particle pdg codes depends on generated mode
292 std::array<int, 4> pdg {{0}};
293 switch (input2_.iproc) {
294 case c_MuonParticle:
295 pdg = {{m_particlePDG, -m_particlePDG, 13, -13}};
296 break;
297 case c_MuonMuon:
298 pdg = {{13, -13, 13, -13}};
299 break;
300 case c_ElectonMuon:
301 pdg = {{11, -11, 13, -13}};
302 break;
304 pdg = {{11, -11, m_particlePDG, -m_particlePDG}};
305 break;
307 pdg = {{11, -11, 11, -11}};
308 break;
309 }
310
311 updateParticle(p1, momenz_.q1, 11, true);
312 updateParticle(p2, momenz_.q2, -11, true);
313 updateParticle(p3, momenz_.q3, pdg[0]);
314 updateParticle(p4, momenz_.q4, pdg[1]);
315 updateParticle(p5, momenz_.q5, pdg[2]);
316 updateParticle(p6, momenz_.q6, pdg[3]);
317}
318
320{
321 aafhdiag36_finish_();
322 // The following code just prints the results. If B2RESULT has been compiled
323 // out this is unneccesary and cppcheck complains so we enclose this in an
324 // ifdef to only be executed if B2RESULT might actually be shown
325#ifndef LOG_NO_B2RESULT
326 B2RESULT("AAFH Generation finished.");
327 B2RESULT("Final cross section: "
328 << std::setprecision(3) << output_.crosec << "+-"
329 << std::setprecision(3) << output_.ercros << " nb");
330 B2RESULT("Minimum invariant mass: "
331 << sqrt(bound_.w2min)*input_.eb << " GeV");
332 B2RESULT("Maximum invariant mass: "
333 << sqrt(bound_.w2max)*input_.eb << " GeV");
334 B2RESULT("Maximum subgenerator weight: " << westat_.maxwe << " ("
335 << westat_.mwe[0] << ", "
336 << westat_.mwe[1] << ", "
337 << westat_.mwe[2] << ", "
338 << westat_.mwe[3] << ")");
339 B2RESULT("Maximum final weight: " << ftstat_.ftmax);
340 //The fortran source says that it's advisable if all subgenerators have the
341 //same probability to account for all peaks in the differental cross section
342 //and that it runs with the highest efficiency when the maximum weights in
343 //all sub generators are equal. So let's calculate appropriate subgenerator
344 //weights by looking how many events were generated per sub generator or the
345 //maximum weights
346 std::array<double, 4> idealWAP;
347 std::array<double, 4> idealWBP;
348 //find the first non-zero weight
349 int reference = 0;
350 for (int i = 0; i < 4; ++i) {
351 if (input_.wap[i] != 0) {
352 reference = i;
353 break;
354 }
355 }
356 //and use that es reference
357 const double w0 = westat_.mwe[reference] * input_.wbp[reference];
358 const double n0 = westat_.iwe[reference] / input_.wap[reference];
359 //now calculate the relative weights we need to get
360 //1) same maximum weight which is inversely proportional to to wbp
361 //2) same number of generated events which is proportional to wap
362 for (int i = 0; i < 4; ++i) {
363 if (input_.wap[i] == 0) {
364 idealWBP[i] = input_.wbp[i];
365 idealWAP[i] = 0;
366 } else {
367 idealWBP[i] = (westat_.mwe[i] * input_.wbp[i]) / w0;
368 idealWAP[i] = n0 * input_.wap[i] / westat_.iwe[i];
369 //and normalize to 1.0
370 idealWAP[i] /= idealWAP[reference];
371 }
372 }
373
374 B2RESULT("To get the same maximum event weight and chance for each sub "
375 "generator use these parameters:");
376 B2RESULT(" 'maxFinalWeight': " << std::setprecision(3)
377 << ftstat_.ftmax * 1.1 << ",");
378 B2RESULT(" 'maxSubgeneratorWeight': " << std::setprecision(3)
379 << westat_.maxwe * 1.1 << ",");
380 B2RESULT(" 'subgeneratorWeights': [" << std::scientific
381 << std::setprecision(3) << idealWAP[0] << ", "
382 << std::setprecision(3) << idealWAP[1] << ", "
383 << std::setprecision(3) << idealWAP[2] << ", "
384 << std::setprecision(3) << idealWAP[3] << ", ");
385 B2RESULT(" " << std::scientific
386 << std::setprecision(3) << idealWBP[0] << ", "
387 << std::setprecision(3) << idealWBP[1] << ", "
388 << std::setprecision(3) << idealWBP[2] << ", "
389 << std::setprecision(3) << idealWBP[3] << "],");
390#endif
391}
int m_particlePDG
pdg of the particle for modes c_MuonParticle and c_ElectronParticle, gets set in initialize()
double m_minMass
minimum invariant mass
double getMaxSubGeneratorWeight() const
Get the maximum expected final weight for the rejection method.
ERejection
Rejection mode.
Definition: AAFHInterface.h:46
void updateParticle(MCParticleGraph::GraphParticle &p, double *q, int pdg, bool isInitial=false)
update particle with generated values
double getTotalCrossSectionError() const
Return error on the total cross section.
double getTotalCrossSection() const
Return total cross section.
std::string m_particle
name of the particle for modes c_MuonParticle and c_ElectronParticle
void finish()
calculate total cross section
void generateEvent(MCParticleGraph &mpg)
generate one event and add it to the graph in CMS
void initialize(double beamEnergy, EMode mode, ERejection rejection)
initialize the generator
EMode
Generation mode.
Definition: AAFHInterface.h:32
@ c_ElectonElectron
E+E- -> E+E-E+E-.
Definition: AAFHInterface.h:42
@ c_ElectonMuon
E+E- -> E+E-Mu+Mu-.
Definition: AAFHInterface.h:38
@ c_MuonParticle
E+E- -> Mu+Mu-L+L- where L can be anything, defaults to tau.
Definition: AAFHInterface.h:34
@ c_MuonMuon
E+E- -> Mu+Mu-Mu+Mu-.
Definition: AAFHInterface.h:36
@ c_ElectronParticle
E+E- -> E+E-L+L- where L can be anything, defaults to tau.
Definition: AAFHInterface.h:40
void setMaxWeights(double subgeneratorWeight, double finalWeight)
Set the maximum expected weights for the rejection method.
void setMaxTries(int tries)
Set the maximum number of tries to generate an event.
std::vector< double > getSuppressionLimits() const
Get suppression limits.
void setSupressionLimits(std::vector< double > limits)
Set the suppression limits when calculation the matrix elements.
int getMaxTries() const
Get the maximum number of tries to generate an event.
void setGeneratorWeights(const std::vector< double > &weights)
Set the relative weights for the sub generators.
double getMaxFinalWeight() const
Get the maximum expected subgenerator weight for the rejection method.
std::vector< double > getGeneratorWeights() const
Get the relative weights for the sub generators.
Class to represent Particle data in graph.
Class to build, validate and sort a particle decay chain.
@ c_Initial
bit 5: Particle is initial such as e+ or e- and not going to Geant4
Definition: MCParticle.h:57
@ 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
static const double nb
[nanobarn]
Definition: Unit.h:84
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.