9 #include <framework/logging/Logger.h>
10 #include <framework/gearbox/Unit.h>
11 #include <generators/aafh/AAFHInterface.h>
13 #include <TDatabasePDG.h>
22 double aafhdiag36_rndm_(
int*)
24 return gRandom->Rndm();
27 double aafhdiag36_rnf100_(
int*)
29 return gRandom->Rndm();
32 void aafhdiag36_esft_maxed_(
double* weight)
34 B2ERROR(
"AAFH: Maximum weight to small, increase maxFinalWeight to at least "
38 void aafhdiag36_eswe_maxed_(
double* weight)
40 B2ERROR(
"AAFH: Maximum weight to small, increase maxSubGeneratorWeight to at least "
45 extern void aafhdiag36_init_();
47 extern void aafhdiag36_event_(
int*);
49 extern void aafhdiag36_finish_();
147 input2_.itot = tries;
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");
161 std::copy_n(weights.begin(), 4, input_.wap);
162 if (weights.size() >= 8) {
163 std::copy_n(weights.begin() + 4, 4, input_.wbp);
169 input_.eswe = subgeneratorWeight;
170 input_.esft = finalWeight;
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];
187 B2ERROR(
"AAFH: Wrong number of suppression limits: supply either one or 4 values");
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);
216 return {factor_.face, factor_.facm, factor_.facl, factor_.proc};
231 input_.eb = beamEnergy;
235 input2_.iproc = mode;
236 input2_.irejec = rejection;
239 auto* part = TDatabasePDG::Instance()->GetParticle(
m_particle.c_str());
241 B2ERROR(
"AAFH: Particle '" <<
m_particle <<
"' does not exist");
245 masses_.xml = part->Mass() / beamEnergy;
247 charge_.qcharg = -part->Charge() / 3.0;
258 TLorentzVector vec(-q[0]*input_.eb, -q[1]*input_.eb, -q[2]*input_.eb, q[3]*input_.eb);
262 p.setMass(fabs(q[4])*input_.eb);
274 aafhdiag36_event_(&status);
276 B2ERROR(
"Failed to generate event after " << status <<
" tries, giving up");
291 std::array<int, 4> pdg {{0}};
292 switch (input2_.iproc) {
297 pdg = {{13, -13, 13, -13}};
300 pdg = {{11, -11, 13, -13}};
306 pdg = {{11, -11, 11, -11}};
320 aafhdiag36_finish_();
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);
345 std::array<double, 4> idealWAP;
346 std::array<double, 4> idealWBP;
349 for (
int i = 0; i < 4; ++i) {
350 if (input_.wap[i] != 0) {
356 const double w0 = westat_.mwe[reference] * input_.wbp[reference];
357 const double n0 = westat_.iwe[reference] / input_.wap[reference];
361 for (
int i = 0; i < 4; ++i) {
362 if (input_.wap[i] == 0) {
363 idealWBP[i] = input_.wbp[i];
366 idealWBP[i] = (westat_.mwe[i] * input_.wbp[i]) / w0;
367 idealWAP[i] = n0 * input_.wap[i] / westat_.iwe[i];
369 idealWAP[i] /= idealWAP[reference];
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] <<
"],");
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]
GraphParticle & addParticle()
Add new particle to the graph.
Abstract base class for different kinds of events.