Public Member Functions | Private Member Functions | Private Attributes | List of all members
genie::rew::GReWeightAGKY Class Reference

Reweighting the GENIE AGKY (free-nucleon) hadronization model. More...

#include <GReWeightAGKY.h>

Inheritance diagram for genie::rew::GReWeightAGKY:
genie::rew::GReWeightI

Public Member Functions

 GReWeightAGKY ()
 
 ~GReWeightAGKY ()
 
bool IsHandled (GSyst_t syst)
 does the current weight calculator handle the input nuisance param? More...
 
void SetSystematic (GSyst_t syst, double val)
 update the value for the specified nuisance param More...
 
void Reset (void)
 set all nuisance parameters to default values More...
 
void Reconfigure (void)
 propagate updated nuisance parameter values to actual MC, etc More...
 
double CalcWeight (const EventRecord &event)
 calculate a weight for the input event using the current nuisance param values More...
 
void RewNue (bool tf)
 
void RewNuebar (bool tf)
 
void RewNumu (bool tf)
 
void RewNumubar (bool tf)
 
void RewCC (bool tf)
 
void RewNC (bool tf)
 
- Public Member Functions inherited from genie::rew::GReWeightI
virtual ~GReWeightI ()
 

Private Member Functions

void Init (void)
 
double RewxFpT1pi (const EventRecord &event)
 

Private Attributes

bool fRewNue
 reweight nu_e? More...
 
bool fRewNuebar
 reweight nu_e_bar? More...
 
bool fRewNumu
 reweight nu_mu? More...
 
bool fRewNumubar
 reweight nu_mu_bar? More...
 
bool fRewCC
 reweight CC? More...
 
bool fRewNC
 reweight NC? More...
 
double fXFmin
 
double fXFmax
 
double fPT2min
 
double fPT2max
 
TF1 * fBaryonXFpdf
 
TF1 * fBaryonPT2pdf
 
TF1 * fBaryonXFpdfTwk
 
TF1 * fBaryonPT2pdfTwk
 
double fDefPeakBaryonXF
 
double fDefAvgPT2
 
double fPeakBaryonXFTwkDial
 
double fAvgPT2TwkDial
 
double fI0XFpdf
 
double fI0PT2pdf
 

Additional Inherited Members

- Protected Member Functions inherited from genie::rew::GReWeightI
 GReWeightI ()
 

Detailed Description

Reweighting the GENIE AGKY (free-nucleon) hadronization model.

Author
Jim Dobson <J.Dobson07 imperial.ac.uk> Imperial College London

Costas Andreopoulos <costas.andreopoulos stfc.ac.uk> University of Liverpool & STFC Rutherford Appleton Lab

Sep 10, 2009

Copyright (c) 2003-2017, GENIE Neutrino MC Generator Collaboration For the full text of the license visit http://copyright.genie-mc.org or see $GENIE/LICENSE

Definition at line 39 of file GReWeightAGKY.h.

Constructor & Destructor Documentation

GReWeightAGKY::GReWeightAGKY ( )

Definition at line 44 of file GReWeightAGKY.cxx.

44  :
45 GReWeightI()
46 {
47  this->Init();
48 }
GReWeightAGKY::~GReWeightAGKY ( )

Definition at line 50 of file GReWeightAGKY.cxx.

51 {
52  delete fBaryonXFpdf;
53  delete fBaryonPT2pdf;
54  delete fBaryonXFpdfTwk;
55  delete fBaryonPT2pdfTwk;
56 
57 #ifdef _G_REWEIGHT_AGKY_DEBUG_
58  fTestFile->cd();
59  fTestNtp ->Write();
60  fTestFile->Close();
61  delete fTestFile;
62 #endif
63 }

Member Function Documentation

double GReWeightAGKY::CalcWeight ( const EventRecord event)
virtual

calculate a weight for the input event using the current nuisance param values

Implements genie::rew::GReWeightI.

Definition at line 106 of file GReWeightAGKY.cxx.

107 {
108  Interaction * interaction = event.Summary();
109 
110  bool is_cc = interaction->ProcInfo().IsWeakCC();
111  bool is_nc = interaction->ProcInfo().IsWeakNC();
112  if(is_cc && !fRewCC) return 1.;
113  if(is_nc && !fRewNC) return 1.;
114 
115  int nupdg = interaction->InitState().ProbePdg();
116  if(nupdg==kPdgNuMu && !fRewNumu ) return 1.;
117  if(nupdg==kPdgAntiNuMu && !fRewNumubar) return 1.;
118  if(nupdg==kPdgNuE && !fRewNue ) return 1.;
119  if(nupdg==kPdgAntiNuE && !fRewNuebar ) return 1.;
120 
121  // Skip events not handled by the AGKY hadronization model
122  if(! utils::rew::HadronizedByAGKY(event)) return 1.0;
123 
124  // Skip high-W events hadronized by AGKY/PYTHIA
125  if(utils::rew::HadronizedByAGKYPythia(event)) return 1.0;
126 
127  //
128  // Reweight events hadronized by AGKY/KNO
129  //
130 
131  double wght =
132  this->RewxFpT1pi(event); /* reweight pT2 and xF for `nucleon+pion' states */
133 
134  return wght;
135 }
bool fRewNue
reweight nu_e?
Definition: GReWeightAGKY.h:65
bool IsWeakCC(void) const
const int kPdgNuE
Definition: PDGCodes.h:25
bool fRewNumu
reweight nu_mu?
Definition: GReWeightAGKY.h:67
const int kPdgAntiNuE
Definition: PDGCodes.h:26
const int kPdgNuMu
Definition: PDGCodes.h:27
bool fRewNuebar
reweight nu_e_bar?
Definition: GReWeightAGKY.h:66
bool fRewNumubar
reweight nu_mu_bar?
Definition: GReWeightAGKY.h:68
Summary information for an interaction.
Definition: Interaction.h:53
bool IsWeakNC(void) const
bool HadronizedByAGKY(const EventRecord &event)
bool fRewNC
reweight NC?
Definition: GReWeightAGKY.h:70
bool fRewCC
reweight CC?
Definition: GReWeightAGKY.h:69
int ProbePdg(void) const
Definition: InitialState.h:54
double RewxFpT1pi(const EventRecord &event)
const int kPdgAntiNuMu
Definition: PDGCodes.h:28
bool HadronizedByAGKYPythia(const EventRecord &event)
const InitialState & InitState(void) const
Definition: Interaction.h:66
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:67
void GReWeightAGKY::Init ( void  )
private

Definition at line 324 of file GReWeightAGKY.cxx.

325 {
326  this->RewNue (true);
327  this->RewNuebar (true);
328  this->RewNumu (true);
329  this->RewNumubar(true);
330  this->RewCC (true);
331  this->RewNC (true);
332 
333  // Baryon pT^2 and xF parameterizations used as PDFs in AGKY.
334  // Code from KNOHadronization.cxx:
335  //
336 
337  fXFmin = -1.0;
338  fXFmax = 0.5;
339  fPT2min = 0;
340  fPT2max = 0.6;
341 
342  fBaryonXFpdf = new TF1("fBaryonXFpdf",
343  "0.083*exp(-0.5*pow(x+0.385,2.)/0.131)", fXFmin, fXFmax);
344  fBaryonPT2pdf = new TF1("fBaryonPT2pdf",
345  "exp(-0.214-6.625*x)", fPT2min, fPT2max);
346 
347  fI0XFpdf = fBaryonXFpdf ->Integral(fXFmin, fXFmax);
348  fI0PT2pdf = fBaryonPT2pdf->Integral(fPT2min,fPT2max);
349 
350  //
351  // Now, define same function as above, but insert tweaking dials.
352  // - For the xF PDF add a parameter to shift the peak of the distribution.
353  // - For the pT2 PDF add a parameter to shift the <pT2>
354  // Include parameter to allow re-normalizing tweaked PDF to 1.
355 
356  fDefPeakBaryonXF = -0.385;
357  fBaryonXFpdfTwk = new TF1("fBaryonXFpdf",
358  "[0]*0.083*exp(-0.5*pow(x-[1],2.)/0.131)", fXFmin, fXFmax);
359  fBaryonXFpdfTwk->SetParameter(0, 1.); // norm
360  fBaryonXFpdfTwk->SetParameter(1, fDefPeakBaryonXF);
361 
362  fDefAvgPT2 = 1./6.625;
363  fBaryonPT2pdfTwk = new TF1("fBaryonPT2pdf",
364  "[0]*exp(-0.214-x/[1])", fPT2min, fPT2max);
365  fBaryonPT2pdfTwk->SetParameter(0, 1.); // norm
366  fBaryonPT2pdfTwk->SetParameter(1,fDefAvgPT2);
367 
368  // init tweaking dials
370  fAvgPT2TwkDial = 0.;
371 
372 #ifdef _G_REWEIGHT_AGKY_DEBUG_
373  fTestFile = new TFile("./agky_reweight_test.root","recreate");
374  fTestNtp = new TNtupleD("testntp","","W:xF:pT2:xFtwkdial:pT2twkdial:wght");
375 #endif
376 }
bool GReWeightAGKY::IsHandled ( GSyst_t  syst)
virtual

does the current weight calculator handle the input nuisance param?

Implements genie::rew::GReWeightI.

Definition at line 65 of file GReWeightAGKY.cxx.

66 {
67  switch(syst) {
68  case ( kHadrAGKYTwkDial_xF1pi ) :
69  case ( kHadrAGKYTwkDial_pT1pi ) :
70  return true;
71  break;
72  default:
73  return false;
74  break;
75  }
76 
77  return false;
78 }
tweak xF distribution for low multiplicity (N + pi) DIS f/s produced by AGKY
Definition: GSyst.h:105
tweak pT distribution for low multiplicity (N + pi) DIS f/s produced by AGKY
Definition: GSyst.h:106
void GReWeightAGKY::Reconfigure ( void  )
virtual

propagate updated nuisance parameter values to actual MC, etc

Implements genie::rew::GReWeightI.

Definition at line 101 of file GReWeightAGKY.cxx.

102 {
103 
104 }
void GReWeightAGKY::Reset ( void  )
virtual

set all nuisance parameters to default values

Implements genie::rew::GReWeightI.

Definition at line 95 of file GReWeightAGKY.cxx.

96 {
98  fAvgPT2TwkDial = 0.;
99 }
void genie::rew::GReWeightAGKY::RewCC ( bool  tf)
inline

Definition at line 57 of file GReWeightAGKY.h.

57 { fRewCC = tf; }
Definition: tf_graph.h:23
bool fRewCC
reweight CC?
Definition: GReWeightAGKY.h:69
void genie::rew::GReWeightAGKY::RewNC ( bool  tf)
inline

Definition at line 58 of file GReWeightAGKY.h.

58 { fRewNC = tf; }
Definition: tf_graph.h:23
bool fRewNC
reweight NC?
Definition: GReWeightAGKY.h:70
void genie::rew::GReWeightAGKY::RewNue ( bool  tf)
inline

Definition at line 53 of file GReWeightAGKY.h.

53 { fRewNue = tf; }
bool fRewNue
reweight nu_e?
Definition: GReWeightAGKY.h:65
Definition: tf_graph.h:23
void genie::rew::GReWeightAGKY::RewNuebar ( bool  tf)
inline

Definition at line 54 of file GReWeightAGKY.h.

54 { fRewNuebar = tf; }
bool fRewNuebar
reweight nu_e_bar?
Definition: GReWeightAGKY.h:66
Definition: tf_graph.h:23
void genie::rew::GReWeightAGKY::RewNumu ( bool  tf)
inline

Definition at line 55 of file GReWeightAGKY.h.

55 { fRewNumu = tf; }
bool fRewNumu
reweight nu_mu?
Definition: GReWeightAGKY.h:67
Definition: tf_graph.h:23
void genie::rew::GReWeightAGKY::RewNumubar ( bool  tf)
inline

Definition at line 56 of file GReWeightAGKY.h.

56 { fRewNumubar = tf; }
Definition: tf_graph.h:23
bool fRewNumubar
reweight nu_mu_bar?
Definition: GReWeightAGKY.h:68
double GReWeightAGKY::RewxFpT1pi ( const EventRecord event)
private

Definition at line 137 of file GReWeightAGKY.cxx.

138 {
139  bool PT2tweaked = (TMath::Abs(fAvgPT2TwkDial) > controls::kASmallNum);
140  bool XFtweaked = (TMath::Abs(fPeakBaryonXFTwkDial) > controls::kASmallNum);
141 
142  bool tweaked = (PT2tweaked || XFtweaked);
143  if(!tweaked) return 1.;
144 
145  //
146  // Did the hadronization code produced a `nucleon+pion' system?
147  // If yes, keep the nucleon xF and pT (in HCM) for event reweighting.
148  //
149 
150  GHepParticle * p = 0;
151  TIter event_iter(&event);
152  int i=-1, ihadsyst = -1;
153  while ( (p = dynamic_cast<GHepParticle *>(event_iter.Next())) ) {
154  i++;
155  if(p->Pdg() != kPdgHadronicSyst) continue;
156  ihadsyst = i;
157  break;
158  }
159  if(ihadsyst<0) return 1.;
160 
161  LOG("ReW", pDEBUG)
162  << "Found HadronicSystem pseudo-particle at position = " << ihadsyst;
163 
164  int fd = event.Particle(ihadsyst)->FirstDaughter();
165  int ld = event.Particle(ihadsyst)->LastDaughter();
166  int nd = ld-fd+1;
167 
168  LOG("ReW", pDEBUG)
169  << "HadronicSystem pseudo-particle's num. of daughters = " << nd
170  << " at positions (" << fd << ", " << ld << ")";
171 
172  bool is_Npi = false;
173  int iN = -1;
174  if(nd==2) {
175  int fdpdg = event.Particle(fd)->Pdg();
176  int ldpdg = event.Particle(ld)->Pdg();
177  if( pdg::IsNucleon(fdpdg) && pdg::IsPion(ldpdg) ) {
178  is_Npi = true;
179  iN = fd;
180  }
181  else
182  if( pdg::IsNucleon(ldpdg) && pdg::IsPion(fdpdg) ) {
183  is_Npi = true;
184  iN = ld;
185  }
186  }
187  if(!is_Npi) return 1.;
188 
189  LOG("ReW", pDEBUG)
190  << "A DIS event with a 'nucleon+pion' primary hadronic state";
191 
192  // Nucleon 4-momentum at LAB
193  GHepParticle * N = event.Particle(fd);
194  TLorentzVector p4N(N->Px(), N->Py(), N->Pz(), N->E());
195  TVector3 p3N = p4N.Vect();
196  LOG("ReW", pDEBUG)
197  << "4-p @ LAB: px = " << p4N.Px() << ", py = " << p4N.Py()
198  << ", pz = " << p4N.Pz() << ", E = " << p4N.Energy();
199 
200  // Hadronic system 4-momentum
201  TLorentzVector p4had = utils::rew::Hadronic4pLAB(event);
202  TVector3 p3had = p4had.Vect();
203 
204  // Nucleon 4-momentum at LAB' (LAB but with z' rotated
205  // along the hadron shower direction)
206  double pTlabp = p3N.Pt(p3had);
207  double pLlabp = TMath::Sqrt(p3N.Mag2()-pTlabp*pTlabp);
208  double Elabp = p4N.Energy();
209  TLorentzVector p4Nr(0,pTlabp,pLlabp,Elabp);
210 
211  // Boost velocity (LAB' -> HCM)
212  TVector3 bv(0,0,-1.*p4had.P()/p4had.Energy());
213 
214  // Nucleon 4-momentum at HCM
215  p4Nr.Boost(bv);
216 
217  // Get baryon xF and pT2 in HCM
218  bool selected = true;
219  double W = event.Summary()->Kine().W(selected);
220  double PZmax = W/2.;
221  double XF = p4Nr.Pz() / PZmax;
222  double PT2 = p4Nr.Perp2();
223 
224  LOG("ReW", pDEBUG)
225  << "@ HCM: pT2 = " << PT2 << ", xF = " << XF;
226 
227  bool XFinrange = (XF > fXFmin && XF < fXFmax);
228  bool PT2inrange = (PT2 > fPT2min && PT2 < fPT2max);
229  bool inrange = XFinrange && PT2inrange;
230  if(!inrange) return 1.;
231 
232  // Get tweaked weighting functions
233 
235  if(XFtweaked) {
236  double frcerr = uncertainty->OneSigmaErr(kHadrAGKYTwkDial_xF1pi);
237  double XFpeak = fDefPeakBaryonXF * (1. + fPeakBaryonXFTwkDial * frcerr);
238  fBaryonXFpdfTwk->SetParameter(1, XFpeak);
239  double I = fBaryonXFpdfTwk->Integral(fXFmin,fXFmax);
240  if(I>0.) {
241  double norm = fI0XFpdf/I;
242  fBaryonXFpdfTwk->SetParameter(0, norm);
243  }
244  }
245  if(PT2tweaked) {
246  double frcerr = uncertainty->OneSigmaErr(kHadrAGKYTwkDial_pT1pi);
247  double PT2avg = fDefAvgPT2* (1. + fAvgPT2TwkDial * frcerr);
248  fBaryonPT2pdfTwk->SetParameter(1, PT2avg);
249  double I = fBaryonPT2pdfTwk->Integral(fPT2min,fPT2max);
250  if(I>0.) {
251  double norm = fI0PT2pdf/I;
252  fBaryonPT2pdfTwk->SetParameter(0, norm);
253  }
254  }
255 
256  // Default and tweaked nucleon xF:pT2 distribution at given W and for
257  // given tweaking dials.
258  Double_t masses[2] = { kNucleonMass, kPionMass } ;
259  int ndec=20000;
260  int nbin=20;
261  RandomGen * rnd = RandomGen::Instance();
262  TLorentzVector ph(0,0,0,W);
263  TGenPhaseSpace phspgen;
264  phspgen.SetDecay(ph, 2, masses);
265  TH2F hdef("hdef","def xF:pT2 at given W", nbin, fXFmin, fXFmax, nbin, fPT2min, fPT2max);
266  TH2F htwk("htwk","twk xF:pT2 at given W", nbin, fXFmin, fXFmax, nbin, fPT2min, fPT2max);
267  for (Int_t n=0;n<ndec;n++) {
268  double dec_weight = phspgen.Generate();
269  TLorentzVector * nucp4 = phspgen.GetDecay(0);
270  double dec_px=nucp4->Px();
271  double dec_py=nucp4->Py();
272  double dec_pz=nucp4->Pz();
273  double dec_pT2 = dec_px*dec_px+dec_py*dec_py;
274  double dec_xF = dec_pz/(W/2.);
275 
276  double fpT2max = 1.1 * fBaryonXFpdf ->GetMaximum(fXFmin, fXFmax );
277  double fxFmax = 1.1 * fBaryonPT2pdf->GetMaximum(fPT2min,fPT2max);
278  double fpT2rnd = fpT2max * rnd->RndHadro().Rndm();
279  double fxFrnd = fxFmax * rnd->RndHadro().Rndm();
280  double fpT2pdf = fBaryonPT2pdf->Eval(dec_pT2);
281  double fxFpdf = fBaryonXFpdf ->Eval(dec_xF );
282  if(fxFrnd < fxFpdf && fpT2rnd < fpT2pdf) {
283  hdef.Fill(dec_xF, dec_pT2, dec_weight);
284  }
285 
286  fpT2max = 1.1 * fBaryonXFpdfTwk ->GetMaximum(fXFmin, fXFmax );
287  fxFmax = 1.1 * fBaryonPT2pdfTwk->GetMaximum(fPT2min,fPT2max);
288  fpT2rnd = fpT2max * rnd->RndHadro().Rndm();
289  fxFrnd = fxFmax * rnd->RndHadro().Rndm();
290  fpT2pdf = fBaryonPT2pdfTwk->Eval(dec_pT2);
291  fxFpdf = fBaryonXFpdfTwk ->Eval(dec_xF );
292  if(fxFrnd < fxFpdf && fpT2rnd < fpT2pdf) {
293  htwk.Fill(dec_xF, dec_pT2, dec_weight);
294  }
295  }//ndec
296 
297  double Idef = hdef.Integral("width");
298  double Itwk = htwk.Integral("width");
299  if(Idef > 0 && Itwk > 0) {
300  hdef.Scale(1./hdef.Integral("width"));
301  htwk.Scale(1./htwk.Integral("width"));
302  } else {
303  return 1.;
304  }
305 
306  int xfbin = hdef.GetXaxis()->FindBin(XF);
307  int pt2bin = hdef.GetYaxis()->FindBin(PT2);
308 
309  double prob_def = hdef.GetBinContent(xfbin,pt2bin);
310  double prob_twk = htwk.GetBinContent(xfbin,pt2bin);
311  if(prob_def <= 0 || prob_twk < 0) {
312  return 1.;
313  }
314 
315  double wght = prob_twk/prob_def;
316 
317 #ifdef _G_REWEIGHT_AGKY_DEBUG_
318  fTestNtp->Fill(W,XF,PT2,fPeakBaryonXFTwkDial,fAvgPT2TwkDial,wght);
319 #endif
320 
321  return wght;
322 }
bool IsPion(int pdgc)
Definition: PDGUtils.cxx:296
double E(void) const
Get energy.
Definition: GHepParticle.h:90
static const double kNucleonMass
Definition: Constants.h:78
static RandomGen * Instance()
Access instance.
Definition: RandomGen.cxx:79
bool IsNucleon(int pdgc)
Definition: PDGUtils.cxx:316
tweak xF distribution for low multiplicity (N + pi) DIS f/s produced by AGKY
Definition: GSyst.h:105
TLorentzVector Hadronic4pLAB(const EventRecord &event)
A singleton holding random number generator classes. All random number generation in GENIE should tak...
Definition: RandomGen.h:30
double Pz(void) const
Get Pz.
Definition: GHepParticle.h:89
tweak pT distribution for low multiplicity (N + pi) DIS f/s produced by AGKY
Definition: GSyst.h:106
double Px(void) const
Get Px.
Definition: GHepParticle.h:87
double OneSigmaErr(GSyst_t syst, int sign=0) const
int Pdg(void) const
Definition: GHepParticle.h:64
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:87
static const double kASmallNum
Definition: Controls.h:40
auto norm(Vector const &v)
Return norm of the specified vector.
p
Definition: test.py:228
TRandom3 & RndHadro(void) const
rnd number generator used by hadronization models
Definition: RandomGen.h:54
static const double kPionMass
Definition: Constants.h:74
static GSystUncertainty * Instance(void)
const int kPdgHadronicSyst
Definition: PDGCodes.h:181
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:40
#define pDEBUG
Definition: Messenger.h:54
double Py(void) const
Get Py.
Definition: GHepParticle.h:88
void GReWeightAGKY::SetSystematic ( GSyst_t  syst,
double  val 
)
virtual

update the value for the specified nuisance param

Implements genie::rew::GReWeightI.

Definition at line 80 of file GReWeightAGKY.cxx.

81 {
82  switch(syst) {
83  case ( kHadrAGKYTwkDial_xF1pi ) :
84  fPeakBaryonXFTwkDial = twk_dial;
85  break;
86  case ( kHadrAGKYTwkDial_pT1pi ) :
87  fAvgPT2TwkDial = twk_dial;
88  break;
89  default:
90  return;
91  break;
92  }
93 }
tweak xF distribution for low multiplicity (N + pi) DIS f/s produced by AGKY
Definition: GSyst.h:105
tweak pT distribution for low multiplicity (N + pi) DIS f/s produced by AGKY
Definition: GSyst.h:106

Member Data Documentation

double genie::rew::GReWeightAGKY::fAvgPT2TwkDial
private

Definition at line 82 of file GReWeightAGKY.h.

TF1* genie::rew::GReWeightAGKY::fBaryonPT2pdf
private

Definition at line 76 of file GReWeightAGKY.h.

TF1* genie::rew::GReWeightAGKY::fBaryonPT2pdfTwk
private

Definition at line 78 of file GReWeightAGKY.h.

TF1* genie::rew::GReWeightAGKY::fBaryonXFpdf
private

Definition at line 75 of file GReWeightAGKY.h.

TF1* genie::rew::GReWeightAGKY::fBaryonXFpdfTwk
private

Definition at line 77 of file GReWeightAGKY.h.

double genie::rew::GReWeightAGKY::fDefAvgPT2
private

Definition at line 80 of file GReWeightAGKY.h.

double genie::rew::GReWeightAGKY::fDefPeakBaryonXF
private

Definition at line 79 of file GReWeightAGKY.h.

double genie::rew::GReWeightAGKY::fI0PT2pdf
private

Definition at line 84 of file GReWeightAGKY.h.

double genie::rew::GReWeightAGKY::fI0XFpdf
private

Definition at line 83 of file GReWeightAGKY.h.

double genie::rew::GReWeightAGKY::fPeakBaryonXFTwkDial
private

Definition at line 81 of file GReWeightAGKY.h.

double genie::rew::GReWeightAGKY::fPT2max
private

Definition at line 74 of file GReWeightAGKY.h.

double genie::rew::GReWeightAGKY::fPT2min
private

Definition at line 73 of file GReWeightAGKY.h.

bool genie::rew::GReWeightAGKY::fRewCC
private

reweight CC?

Definition at line 69 of file GReWeightAGKY.h.

bool genie::rew::GReWeightAGKY::fRewNC
private

reweight NC?

Definition at line 70 of file GReWeightAGKY.h.

bool genie::rew::GReWeightAGKY::fRewNue
private

reweight nu_e?

Definition at line 65 of file GReWeightAGKY.h.

bool genie::rew::GReWeightAGKY::fRewNuebar
private

reweight nu_e_bar?

Definition at line 66 of file GReWeightAGKY.h.

bool genie::rew::GReWeightAGKY::fRewNumu
private

reweight nu_mu?

Definition at line 67 of file GReWeightAGKY.h.

bool genie::rew::GReWeightAGKY::fRewNumubar
private

reweight nu_mu_bar?

Definition at line 68 of file GReWeightAGKY.h.

double genie::rew::GReWeightAGKY::fXFmax
private

Definition at line 72 of file GReWeightAGKY.h.

double genie::rew::GReWeightAGKY::fXFmin
private

Definition at line 71 of file GReWeightAGKY.h.


The documentation for this class was generated from the following files: