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