9#include <framework/logging/Logger.h>
10#include <framework/gearbox/Unit.h>
11#include <generators/aafh/AAFHInterface.h>
13#include <TDatabasePDG.h>
14#include <Math/Vector4D.h>
23 double aafhdiag36_rndm_(
int*)
25 return gRandom->Rndm();
28 double aafhdiag36_rnf100_(
int*)
30 return gRandom->Rndm();
33 void aafhdiag36_esft_maxed_(
const double* weight)
35 B2ERROR(
"AAFH: Maximum weight to small, increase maxFinalWeight to at least "
39 void aafhdiag36_eswe_maxed_(
const double* weight)
41 B2ERROR(
"AAFH: Maximum weight to small, increase maxSubGeneratorWeight to at least "
46 extern void aafhdiag36_init_();
48 extern void aafhdiag36_event_(
int*);
50 extern void aafhdiag36_finish_();
148 input2_.itot = tries;
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");
162 std::copy_n(weights.begin(), 4, input_.wap);
163 if (weights.size() >= 8) {
164 std::copy_n(weights.begin() + 4, 4, input_.wbp);
170 input_.eswe = subgeneratorWeight;
171 input_.esft = finalWeight;
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];
188 B2ERROR(
"AAFH: Wrong number of suppression limits: supply either one or 4 values");
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);
217 return {factor_.face, factor_.facm, factor_.facl, factor_.proc};
232 input_.eb = beamEnergy;
236 input2_.iproc = mode;
237 input2_.irejec = rejection;
240 auto* part = TDatabasePDG::Instance()->GetParticle(
m_particle.c_str());
242 B2ERROR(
"AAFH: Particle '" <<
m_particle <<
"' does not exist");
246 masses_.xml = part->Mass() / beamEnergy;
248 charge_.qcharg = -part->Charge() / 3.0;
259 ROOT::Math::PxPyPzEVector vec(-q[0]*input_.eb, -q[1]*input_.eb, -q[2]*input_.eb, q[3]*input_.eb);
263 p.setMass(fabs(q[4])*input_.eb);
275 aafhdiag36_event_(&status);
277 B2ERROR(
"Failed to generate event after " << status <<
" tries, giving up");
292 std::array<int, 4> pdg {{0}};
293 switch (input2_.iproc) {
298 pdg = {{13, -13, 13, -13}};
301 pdg = {{11, -11, 13, -13}};
307 pdg = {{11, -11, 11, -11}};
321 aafhdiag36_finish_();
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);
346 std::array<double, 4> idealWAP;
347 std::array<double, 4> idealWBP;
350 for (
int i = 0; i < 4; ++i) {
351 if (input_.wap[i] != 0) {
357 const double w0 = westat_.mwe[reference] * input_.wbp[reference];
358 const double n0 = westat_.iwe[reference] / input_.wap[reference];
362 for (
int i = 0; i < 4; ++i) {
363 if (input_.wap[i] == 0) {
364 idealWBP[i] = input_.wbp[i];
367 idealWBP[i] = (westat_.mwe[i] * input_.wbp[i]) / w0;
368 idealWAP[i] = n0 * input_.wap[i] / westat_.iwe[i];
370 idealWAP[i] /= idealWAP[reference];
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] <<
"],");
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.
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
@ c_ElectonElectron
E+E- -> E+E-E+E-.
@ c_ElectonMuon
E+E- -> E+E-Mu+Mu-.
@ c_MuonParticle
E+E- -> Mu+Mu-L+L- where L can be anything, defaults to tau.
@ c_MuonMuon
E+E- -> Mu+Mu-Mu+Mu-.
@ c_ElectronParticle
E+E- -> E+E-L+L- where L can be anything, defaults to tau.
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
@ c_PrimaryParticle
bit 0: Particle is primary particle.
@ c_StableInGenerator
bit 1: Particle is stable, i.e., not decaying in the generator.
static const double nb
[nanobarn]
double sqrt(double a)
sqrt for double
GraphParticle & addParticle()
Add new particle to the graph.
Abstract base class for different kinds of events.