GENIEGenerator
Loading...
Searching...
No Matches
AGCharmPythiaBaseHadro2023.h
Go to the documentation of this file.
1//____________________________________________________________________________
2/*!
3
4\class genie::AGCharmPythiaBaseHadro2023
5
6\brief Andreopoulos - Gallagher (AG) GENIE Charm Hadronization model.
7
8 The model relies on empirical charm fragmentation and pT functions,
9 as well as on experimentally-determined charm fractions, to produce
10 the ID and 4-momentum of charmed hadron in charm production events.
11
12 The remnant (non-charm) system is hadronised by a call to PYTHIA.
13
14 Is a concrete implementation of the EventRecordVisitorI interface.
15
16\author Costas Andreopoulos <c.andreopoulos \at cern.ch>
17 University of Liverpool
18
19 Hugh Gallagher <gallag@minos.phy.tufts.edu>
20 Tufts University
21
22\created August 17, 2004
23
24\cpright Copyright (c) 2003-2025, The GENIE Collaboration
25 For the full text of the license visit http://copyright.genie-mc.org
26*/
27//____________________________________________________________________________
28
29#ifndef _CHARM_HADRONIZATION_H_
30#define _CHARM_HADRONIZATION_H_
31
32#include <TGenPhaseSpace.h>
33
34#include "Framework/Conventions/GBuild.h"
37
38class TF1;
39
40namespace genie {
41
42class Spline;
44
46
47public:
48
49 // Implement the EventRecordVisitorI interface
50 void ProcessEventRecord(GHepRecord * event) const;
51
52protected:
53
55 AGCharmPythiaBaseHadro2023(string config);
56 AGCharmPythiaBaseHadro2023(string name, string config);
58
59protected:
60 void LoadConfig (void);
61 void Initialize (void) const ;
62
63 virtual bool HadronizeRemnant(int qrkSyst1, int qrkSyst2, double WR, TLorentzVector p4R,
64 unsigned int& rpos, TClonesArray * particle_list) const = 0;
65 // needs to append to TClonesArray of GHepParticle starting at rpos
66
67 // Overload the Algorithm::Configure() methods to load private data
68 // members from configuration options
69 //
70 void Configure(const Registry & config);
71 void Configure(string config);
72
73private:
74
75 TClonesArray * Hadronize (const Interaction* ) const ;
76 // returns TClonesArray of GHepParticle
77
78 int GenerateCharmHadron (int nupdg, double EvLab) const ;
79
80 double Weight (void) const ;
81
82 mutable TGenPhaseSpace fPhaseSpaceGenerator; ///< a phase space generator
83
84 // Configuration parameters
85 //
86 bool fCharmOnly; ///< don't hadronize non-charm blob
87 TF1 * fCharmPT2pdf; ///< charm hadron pT^2 pdf
88 const FragmentationFunctionI * fFragmFunc; ///< charm hadron fragmentation func
89
90 double fFracMaxEnergy ; ///< Maximum energy available for the Meson fractions
91
92 Spline * fD0FracSpl; ///< nu charm fraction vs Ev: D0
93 Spline * fDpFracSpl; ///< nu charm fraction vs Ev: D+
94 Spline * fDsFracSpl; ///< nu charm fraction vs Ev: Ds+
95 double fD0BarFrac; ///< nubar \bar{D0} charm fraction
96 double fDmFrac; ///< nubar D- charm fraction
97};
98
99} // genie namespace
100
101#endif // _CHARM_HADRONIZATION__H_
double fD0BarFrac
nubar \bar{D0} charm fraction
virtual bool HadronizeRemnant(int qrkSyst1, int qrkSyst2, double WR, TLorentzVector p4R, unsigned int &rpos, TClonesArray *particle_list) const =0
bool fCharmOnly
don't hadronize non-charm blob
Spline * fD0FracSpl
nu charm fraction vs Ev: D0
TGenPhaseSpace fPhaseSpaceGenerator
a phase space generator
TClonesArray * Hadronize(const Interaction *) const
Spline * fDsFracSpl
nu charm fraction vs Ev: Ds+
Spline * fDpFracSpl
nu charm fraction vs Ev: D+
void ProcessEventRecord(GHepRecord *event) const
int GenerateCharmHadron(int nupdg, double EvLab) const
double fFracMaxEnergy
Maximum energy available for the Meson fractions.
const FragmentationFunctionI * fFragmFunc
charm hadron fragmentation func
Pure abstract base class. Defines the FragmentationFunctionI interface to be implemented by any algor...
GENIE's GHEP MC event record.
Definition GHepRecord.h:45
Summary information for an interaction.
Definition Interaction.h:56
A registry. Provides the container for algorithm configuration parameters.
Definition Registry.h:65
A numeric analysis tool class for interpolating 1-D functions.
Definition Spline.h:58
THE MAIN GENIE PROJECT NAMESPACE
Definition AlgCmp.h:25