GiBUUData.h
Go to the documentation of this file.
1 //____________________________________________________________________________
2 /*!
3 
4 \class genie::GiBUUData
5 
6 \brief Singleton to load and serve data tables provided by the GiBUU group
7 
8 \ref http://gibuu.physik.uni-giessen.de/GiBUU
9  Specific references for each piece of data included in given below.
10 
11 \author Costas Andreopoulos <costas.andreopoulos \at stfc.ac.uk>
12  University of Liverpool & STFC Rutherford Appleton Lab
13 
14 \created May 30, 2009
15 
16 */
17 //____________________________________________________________________________
18 
19 #ifndef _GIBUU_DATA_H_
20 #define _GIBUU_DATA_H_
21 
24 
25 namespace genie {
26 
27 class Spline;
28 
29 class GiBUUData
30 {
31 public:
32 
33  class FormFactors;
34 
35  // access GiBUUData singleton instance
36  static GiBUUData * Instance (void);
37 
38  // access form factor data
39  const FormFactors & FF(void) const;
40 
41  //
42  // Resonance form factors.
43  // Details given in Phys. Rev. C 79, 034601 (2009).
44  //
45  class FormFactors {
46 
47  public:
48 
49  FormFactors();
50  ~FormFactors();
51 
52  // the following is non-zero for I=1/2 (N) resonances
53  double C3V (double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const;
54  double C4V (double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const;
55  double C5V (double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const;
56  double C6V (double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const;
57  double C3A (double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const;
58  double C4A (double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const;
59  double C5A (double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const;
60  double C6A (double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const;
61 
62  // the following is non-zero for I=3/2 (Delta) resonances
63  double F1V (double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const;
64  double F2V (double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const;
65  double FA (double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const;
66  double FP (double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const;
67 
68  double Q2min (void) const { return fMinQ2; }
69  double Q2max (void) const { return fMaxQ2; }
70 
71  private:
72 
73  static double fMinQ2; ///< min Q2 for which resonance f/f data are given
74  static double fMaxQ2; ///< max Q2 for which resonance f/f data are given
75 
76  // The first array index is the resonance id.
77  // Tina provided GiBUU form factor data for 13 resonances given below along with
78  // the corresponding GENIE resonance ids
79  // GENIE Resonance_t as integer
80  // P33(1232) -> kP33_1232 0
81  // S11(1535) -> kS11_1535 1
82  // D13(1520) -> kD13_1520 2
83  // S11(1650) -> kS11_1650 3
84  // D15(1675) -> kD15_1675 5
85  // S31(1620) -> kS31_1620 6
86  // D33(1700) -> kD33_1700 7
87  // P11(1440) -> kP11_1440 8
88  // P13(1720) -> kP13_1720 10
89  // F15(1680) -> kF15_1680 11
90  // P31(1910) -> kP31_1910 12
91  // F35(1905) -> kF35_1905 14
92  // F37(1950) -> kF37_1950 15
93  // The remaining 3 array indices are:
94  // 0 1 2 0 1 0 1 2 3 4 5 6 7 8 9 10 11
95  // [CC,NC,EM][n,p][F1V,F2V,FA,FP,C3V,C4V,C5V,C6V,C3A,C4A,C5A,C6A]
96  // |-------------| for I=1/2 resonances
97  // |-------------------------------| for I=3/2 resonances
98 
99  static const int kNRes = 18;
100  static const int kNCurr = 3;
101  static const int kNHitNuc = 2;
102  static const int kNFFRes = 12;
103 
104  //! actual form factor data = f(Q2)
106 
107  //! func to retrieve interpolated form factor values
108  double FFRes (double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it, int ffid) const;
109 
110  //! load all form factor data tables
111  void LoadTables(void);
112 
113  friend class GiBUUData;
114 
115  }; // FormFactors nested class
116 
117 private:
118  GiBUUData();
119  GiBUUData(const GiBUUData & gibuu_data);
120  ~GiBUUData();
121 
122  // load all data tables
123  void LoadTables(void);
124 
125  // singleton 'self'
127 
128  // form factor data
130 
131  // singleton cleaner
132  struct Cleaner {
135  if (GiBUUData::fInstance !=0) {
136  delete GiBUUData::fInstance;
138  }
139  }
140  };
141  friend struct Cleaner;
142 
143 }; // GiBUUData class
144 
145 } // genie namespace
146 
147 #endif // _GIBUU_DATA_H_
148 
Singleton to load and serve data tables provided by the GiBUU group.
Definition: GiBUUData.h:29
double FA(double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const
Definition: GiBUUData.cxx:310
#include "Numerical/GSFunc.h"
Definition: AlgCmp.h:26
static const int kNCurr
Definition: GiBUUData.h:100
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:642
double C6V(double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const
Definition: GiBUUData.cxx:261
A numeric analysis tool class for interpolating 1-D functions.
Definition: Spline.h:47
void DummyMethodAndSilentCompiler()
Definition: GiBUUData.h:133
static GiBUUData * fInstance
Definition: GiBUUData.h:126
void LoadTables(void)
load all form factor data tables
Definition: GiBUUData.cxx:115
double C4A(double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const
Definition: GiBUUData.cxx:275
double Q2max(void) const
Definition: GiBUUData.h:69
enum genie::EResonance Resonance_t
static const int kNFFRes
Definition: GiBUUData.h:102
double C4V(double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const
Definition: GiBUUData.cxx:247
FormFactors * fFormFactors
Definition: GiBUUData.h:129
static const int kNRes
Definition: GiBUUData.h:99
static double fMaxQ2
max Q2 for which resonance f/f data are given
Definition: GiBUUData.h:74
double C5V(double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const
Definition: GiBUUData.cxx:254
double C5A(double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const
Definition: GiBUUData.cxx:282
static const int kNHitNuc
Definition: GiBUUData.h:101
double FFRes(double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it, int ffid) const
func to retrieve interpolated form factor values
Definition: GiBUUData.cxx:324
double Q2min(void) const
Definition: GiBUUData.h:68
double C3A(double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const
Definition: GiBUUData.cxx:268
Spline * fFFRes[kNRes][kNCurr][kNHitNuc][kNFFRes]
actual form factor data = f(Q2)
Definition: GiBUUData.h:105
double FP(double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const
Definition: GiBUUData.cxx:317
static double fMinQ2
min Q2 for which resonance f/f data are given
Definition: GiBUUData.h:73
double C6A(double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const
Definition: GiBUUData.cxx:289
const FormFactors & FF(void) const
Definition: GiBUUData.cxx:77
double F1V(double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const
Definition: GiBUUData.cxx:296
double C3V(double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const
Definition: GiBUUData.cxx:240
double F2V(double Q2, Resonance_t res, int nucleon_pdg, InteractionType_t it) const
Definition: GiBUUData.cxx:303
enum genie::EInteractionType InteractionType_t
static GiBUUData * Instance(void)
Definition: GiBUUData.cxx:60