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

Reweighting CCQE GENIE neutrino cross sections. More...

#include <GReWeightNuXSecCCQE.h>

Inheritance diagram for genie::rew::GReWeightNuXSecCCQE:
genie::rew::GReWeightI

Public Member Functions

 GReWeightNuXSecCCQE ()
 
 ~GReWeightNuXSecCCQE ()
 
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 SetMode (int mode)
 
void RewNue (bool tf)
 
void RewNuebar (bool tf)
 
void RewNumu (bool tf)
 
void RewNumubar (bool tf)
 
void SetMaPath (string p)
 
void SetZExpPath (string p)
 
- Public Member Functions inherited from genie::rew::GReWeightI
virtual ~GReWeightI ()
 

Static Public Attributes

static const int kModeMa = 0
 
static const int kModeNormAndMaShape = 1
 
static const int kModeZExp = 2
 
static const int fZExpMaxSyst = 4
 maximum number of systematics More...
 

Private Member Functions

void Init (void)
 
double CalcWeightNorm (const EventRecord &event)
 
double CalcWeightMaShape (const EventRecord &event)
 
double CalcWeightMa (const EventRecord &event)
 
double CalcWeightZExp (const EventRecord &event)
 

Private Attributes

XSecAlgorithmIfXSecModelDef
 default model More...
 
XSecAlgorithmIfXSecModel
 tweaked model More...
 
RegistryfXSecModelConfig
 config in tweaked model More...
 
string fFFModel
 
int fMode
 0: Ma, 1: Norm and MaShape, 2: Z-Expansion More...
 
bool fRewNue
 reweight nu_e CC? More...
 
bool fRewNuebar
 reweight nu_e_bar CC? More...
 
bool fRewNumu
 reweight nu_mu CC? More...
 
bool fRewNumubar
 reweight nu_mu_bar CC? More...
 
string fMaPath
 M_{A} path in config Registry. More...
 
double fNormTwkDial
 
double fNormDef
 
double fNormCurr
 
double fMaTwkDial
 
double fMaDef
 
double fMaCurr
 
int fZExpCurrIdx
 current coefficient index More...
 
int fZExpMaxCoef
 max number of coefficients to use More...
 
string fZExpPath
 algorithm path to get coefficients More...
 
double fZExpTwkDial [fZExpMaxSyst]
 
double fZExpDef [fZExpMaxSyst]
 
double fZExpCurr [fZExpMaxSyst]
 array of current parameter values More...
 

Additional Inherited Members

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

Detailed Description

Reweighting CCQE GENIE neutrino cross sections.

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

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

Aug 1, 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 45 of file GReWeightNuXSecCCQE.h.

Constructor & Destructor Documentation

GReWeightNuXSecCCQE::GReWeightNuXSecCCQE ( )

Definition at line 73 of file GReWeightNuXSecCCQE.cxx.

74 {
75  this->Init();
76 }
GReWeightNuXSecCCQE::~GReWeightNuXSecCCQE ( )

Definition at line 78 of file GReWeightNuXSecCCQE.cxx.

79 {
80 #ifdef _G_REWEIGHT_CCQE_DEBUG_
81  fTestFile->cd();
82  fTestNtp ->Write();
83  fTestFile->Close();
84  delete fTestFile;
85 #endif
86 }

Member Function Documentation

double GReWeightNuXSecCCQE::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 270 of file GReWeightNuXSecCCQE.cxx.

271 {
272  bool is_qe = event.Summary()->ProcInfo().IsQuasiElastic();
273  bool is_cc = event.Summary()->ProcInfo().IsWeakCC();
274  if(!is_qe || !is_cc) return 1.;
275 
276  bool charm = event.Summary()->ExclTag().IsCharmEvent(); // skip CCQE charm
277  if(charm) return 1.;
278 
279  int nupdg = event.Probe()->Pdg();
280  if(nupdg==kPdgNuMu && !fRewNumu ) return 1.;
281  if(nupdg==kPdgAntiNuMu && !fRewNumubar) return 1.;
282  if(nupdg==kPdgNuE && !fRewNue ) return 1.;
283  if(nupdg==kPdgAntiNuE && !fRewNuebar ) return 1.;
284 
285  if(fMode==kModeMa && strcmp(fFFModel.c_str(),kModelDipole) == 0) {
286  double wght = this->CalcWeightMa(event);
287  return wght;
288  }
289  else
290  if(fMode==kModeNormAndMaShape && strcmp(fFFModel.c_str(),kModelDipole) == 0) {
291  double wght =
292  this->CalcWeightNorm (event) *
293  this->CalcWeightMaShape (event);
294  return wght;
295  }
296  else
297  if(fMode==kModeZExp && strcmp(fFFModel.c_str(),kModelZExp) == 0) {
298  double wght =
299  this->CalcWeightNorm(event) *
300  this->CalcWeightZExp(event);
301  return wght;
302  }
303 
304  return 1.;
305 }
bool fRewNuebar
reweight nu_e_bar CC?
const int kPdgNuE
Definition: PDGCodes.h:25
const int kPdgAntiNuE
Definition: PDGCodes.h:26
const int kPdgNuMu
Definition: PDGCodes.h:27
double CalcWeightZExp(const EventRecord &event)
bool fRewNumubar
reweight nu_mu_bar CC?
double CalcWeightMaShape(const EventRecord &event)
static const char * kModelZExp
double CalcWeightMa(const EventRecord &event)
const int kPdgAntiNuMu
Definition: PDGCodes.h:28
int fMode
0: Ma, 1: Norm and MaShape, 2: Z-Expansion
bool fRewNumu
reweight nu_mu CC?
double CalcWeightNorm(const EventRecord &event)
static const char * kModelDipole
double GReWeightNuXSecCCQE::CalcWeightMa ( const EventRecord event)
private

Definition at line 382 of file GReWeightNuXSecCCQE.cxx.

383 {
384  bool tweaked = (TMath::Abs(fMaTwkDial) > controls::kASmallNum);
385  if(!tweaked) return 1.0;
386 
387  Interaction * interaction = event.Summary();
388 
389  interaction->KinePtr()->UseSelectedKinematics();
390  interaction->SetBit(kIAssumeFreeNucleon);
391 
392  double old_xsec = event.DiffXSec();
393  double old_weight = event.Weight();
394  double new_xsec = fXSecModel->XSec(interaction, kPSQ2fE);
395  double new_weight = old_weight * (new_xsec/old_xsec);
396 
397 //LOG("ReW", pDEBUG) << "differential cross section (old) = " << old_xsec;
398 //LOG("ReW", pDEBUG) << "differential cross section (new) = " << new_xsec;
399 //LOG("ReW", pDEBUG) << "event generation weight = " << old_weight;
400 //LOG("ReW", pDEBUG) << "new weight = " << new_weight;
401 
402  interaction->KinePtr()->ClearRunningValues();
403  interaction->ResetBit(kIAssumeFreeNucleon);
404 
405  return new_weight;
406 }
Kinematics * KinePtr(void) const
Definition: Interaction.h:73
virtual double XSec(const Interaction *i, KinePhaseSpace_t k=kPSfE) const =0
Compute the cross section for the input interaction.
Summary information for an interaction.
Definition: Interaction.h:53
void UseSelectedKinematics(void)
Definition: Kinematics.cxx:369
static const double kASmallNum
Definition: Controls.h:40
const UInt_t kIAssumeFreeNucleon
Definition: Interaction.h:46
void ClearRunningValues(void)
Definition: Kinematics.cxx:357
XSecAlgorithmI * fXSecModel
tweaked model
double GReWeightNuXSecCCQE::CalcWeightMaShape ( const EventRecord event)
private

Definition at line 408 of file GReWeightNuXSecCCQE.cxx.

409 {
410  bool tweaked = (TMath::Abs(fMaTwkDial) > controls::kASmallNum);
411  if(!tweaked) return 1.0;
412 
413  Interaction * interaction = event.Summary();
414 
415  interaction->KinePtr()->UseSelectedKinematics();
416  interaction->SetBit(kIAssumeFreeNucleon);
417 
418  double old_xsec = event.DiffXSec();
419  double old_weight = event.Weight();
420  double new_xsec = fXSecModel->XSec(interaction, kPSQ2fE);
421  double new_weight = old_weight * (new_xsec/old_xsec);
422 
423 //LOG("ReW", pDEBUG) << "differential cross section (old) = " << old_xsec;
424 //LOG("ReW", pDEBUG) << "differential cross section (new) = " << new_xsec;
425 //LOG("ReW", pDEBUG) << "event generation weight = " << old_weight;
426 //LOG("ReW", pDEBUG) << "new weight = " << new_weight;
427 
428 //double old_integrated_xsec = event.XSec();
429  double old_integrated_xsec = fXSecModelDef -> Integral(interaction);
430  double new_integrated_xsec = fXSecModel -> Integral(interaction);
431  assert(new_integrated_xsec > 0);
432  new_weight *= (old_integrated_xsec/new_integrated_xsec);
433 
434 //LOG("ReW", pDEBUG) << "integrated cross section (old) = " << old_integrated_xsec;
435 //LOG("ReW", pDEBUG) << "integrated cross section (new) = " << new_integrated_xsec;
436 //LOG("ReW", pDEBUG) << "new weight (normalized to const integral) = " << new_weight;
437 
438  interaction->KinePtr()->ClearRunningValues();
439  interaction->ResetBit(kIAssumeFreeNucleon);
440 
441 #ifdef _G_REWEIGHT_CCQE_DEBUG_
442  double E = interaction->InitState().ProbeE(kRfHitNucRest);
443  double Q2 = interaction->Kine().Q2(true);
444  fTestNtp->Fill(E,Q2,new_weight);
445 #endif
446 
447  return new_weight;
448 }
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:642
Kinematics * KinePtr(void) const
Definition: Interaction.h:73
virtual double XSec(const Interaction *i, KinePhaseSpace_t k=kPSfE) const =0
Compute the cross section for the input interaction.
XSecAlgorithmI * fXSecModelDef
default model
Summary information for an interaction.
Definition: Interaction.h:53
void UseSelectedKinematics(void)
Definition: Kinematics.cxx:369
const Kinematics & Kine(void) const
Definition: Interaction.h:68
static const double kASmallNum
Definition: Controls.h:40
const UInt_t kIAssumeFreeNucleon
Definition: Interaction.h:46
const InitialState & InitState(void) const
Definition: Interaction.h:66
double Q2(bool selected=false) const
Definition: Kinematics.cxx:135
void ClearRunningValues(void)
Definition: Kinematics.cxx:357
XSecAlgorithmI * fXSecModel
tweaked model
double ProbeE(RefFrame_t rf) const
double GReWeightNuXSecCCQE::CalcWeightNorm ( const EventRecord event)
private

Definition at line 373 of file GReWeightNuXSecCCQE.cxx.

374 {
375  bool tweaked = (TMath::Abs(fNormTwkDial) > controls::kASmallNum);
376  if(!tweaked) return 1.0;
377 
378  double wght = fNormCurr;
379  return wght;
380 }
static const double kASmallNum
Definition: Controls.h:40
double GReWeightNuXSecCCQE::CalcWeightZExp ( const EventRecord event)
private

Definition at line 450 of file GReWeightNuXSecCCQE.cxx.

451 {
452  // very similar to CalcWeightMa
453  bool tweaked = false;
454  for (int i=0;i<fZExpMaxCoef;i++)
455  {
456  tweaked = tweaked || (TMath::Abs(fZExpTwkDial[i]) > controls::kASmallNum);
457  }
458  if(!tweaked) { return 1.0; }
459 
460  Interaction * interaction = event.Summary();
461 
462  interaction->KinePtr()->UseSelectedKinematics();
463  interaction->SetBit(kIAssumeFreeNucleon);
464 
465  double old_xsec = event.DiffXSec();
466  double old_weight = event.Weight();
467  double new_xsec = fXSecModel->XSec(interaction, kPSQ2fE);
468  double new_weight = old_weight * (new_xsec/old_xsec);
469 
470  interaction->KinePtr()->ClearRunningValues();
471  interaction->ResetBit(kIAssumeFreeNucleon);
472 
473  return new_weight;
474 }
Kinematics * KinePtr(void) const
Definition: Interaction.h:73
virtual double XSec(const Interaction *i, KinePhaseSpace_t k=kPSfE) const =0
Compute the cross section for the input interaction.
Summary information for an interaction.
Definition: Interaction.h:53
void UseSelectedKinematics(void)
Definition: Kinematics.cxx:369
static const double kASmallNum
Definition: Controls.h:40
const UInt_t kIAssumeFreeNucleon
Definition: Interaction.h:46
void ClearRunningValues(void)
Definition: Kinematics.cxx:357
XSecAlgorithmI * fXSecModel
tweaked model
int fZExpMaxCoef
max number of coefficients to use
void GReWeightNuXSecCCQE::Init ( void  )
private

Definition at line 307 of file GReWeightNuXSecCCQE.cxx.

308 {
309  AlgId id("genie::LwlynSmithQELCCPXSec","Default");
310 
311  AlgFactory * algf = AlgFactory::Instance();
312 
313  Algorithm * algdef = algf->AdoptAlgorithm(id);
314  fXSecModelDef = dynamic_cast<XSecAlgorithmI*>(algdef);
316 
317  Algorithm * alg = algf->AdoptAlgorithm(id);
318  fXSecModel = dynamic_cast<XSecAlgorithmI*>(alg);
320 
322  fFFModel = fXSecModelConfig->GetAlg("FormFactorsAlg/AxialFormFactorModel").name;
323 
324  this->RewNue (true);
325  this->RewNuebar (true);
326  this->RewNumu (true);
327  this->RewNumubar(true);
328 
329  //this->SetMaPath("FormFactorsAlg/Ma");
330  this->SetMaPath("FormFactorsAlg/AxialFormFactorModel/QEL-Ma");
331  this->SetZExpPath("FormFactorsAlg/AxialFormFactorModel/");
332 
333  if (strcmp(fFFModel.c_str(),kModelDipole) == 0)
334  {
337  fZExpMaxCoef = 0;
338  } else
339  if (strcmp(fFFModel.c_str(),kModelZExp) == 0)
340  {
341  this->SetMode(kModeZExp);
342  fMaDef = 0.;
343  fZExpMaxCoef =
344  TMath::Min(fXSecModelConfig->GetInt(fZExpPath + "QEL-Kmax"),
345  this->fZExpMaxSyst);
346  }
347 
348  fNormTwkDial = 0.;
349  fNormDef = 1.;
351  fMaTwkDial = 0.;
352  fMaCurr = fMaDef;
353 
354  ostringstream alg_key;
355  for (int i=0;i<fZExpMaxSyst;i++)
356  {
357  alg_key.str("");
358  alg_key << fZExpPath << "QEL-Z_A" << i+1;
359  if (strcmp(fFFModel.c_str(),kModelZExp) == 0 && i < fZExpMaxCoef)
360  { fZExpDef[i] = fXSecModelConfig->GetDouble(alg_key.str()); }
361  else
362  { fZExpDef[i] = 0.; }
363  fZExpTwkDial[i] = 0.;
364  fZExpCurr [i] = fZExpDef[i];
365  }
366 
367 #ifdef _G_REWEIGHT_CCQE_DEBUG_
368  fTestFile = new TFile("./ccqe_reweight_test.root","recreate");
369  fTestNtp = new TNtupleD("testntp","","E:Q2:wght");
370 #endif
371 }
virtual const Registry & GetConfig(void) const
Get configuration registry.
Definition: Algorithm.h:63
Cross Section Calculation Interface.
string fMaPath
M_{A} path in config Registry.
Algorithm abstract base class.
Definition: Algorithm.h:48
RgDbl GetDouble(RgKey key) const
Definition: Registry.cxx:488
double fZExpCurr[fZExpMaxSyst]
array of current parameter values
RgInt GetInt(RgKey key) const
Definition: Registry.cxx:481
XSecAlgorithmI * fXSecModelDef
default model
Registry * fXSecModelConfig
config in tweaked model
void AdoptSubstructure(void)
Definition: Algorithm.cxx:287
static const char * kModelZExp
Algorithm * AdoptAlgorithm(const AlgId &algid) const
Definition: AlgFactory.cxx:127
Algorithm ID (algorithm name + configuration set name)
Definition: AlgId.h:35
static AlgFactory * Instance()
Definition: AlgFactory.cxx:75
A registry. Provides the container for algorithm configuration parameters.
Definition: Registry.h:66
string fZExpPath
algorithm path to get coefficients
The GENIE Algorithm Factory.
Definition: AlgFactory.h:40
XSecAlgorithmI * fXSecModel
tweaked model
int fZExpMaxCoef
max number of coefficients to use
RgAlg GetAlg(RgKey key) const
Definition: Registry.cxx:502
static const int fZExpMaxSyst
maximum number of systematics
static const char * kModelDipole
bool GReWeightNuXSecCCQE::IsHandled ( GSyst_t  syst)
virtual

does the current weight calculator handle the input nuisance param?

Implements genie::rew::GReWeightI.

Definition at line 88 of file GReWeightNuXSecCCQE.cxx.

89 {
90 // read form factor model and compare to mode
91  bool handle;
92 
93  switch(syst) {
94 
95  case ( kXSecTwkDial_NormCCQE ) :
96  if(fMode==kModeNormAndMaShape && strcmp(fFFModel.c_str(),kModelDipole) == 0)
97  {
98  handle = true;
99  } else {
100  handle = false;
101  }
102  break;
103 
104  case ( kXSecTwkDial_MaCCQEshape ) :
105  if(fMode==kModeNormAndMaShape && strcmp(fFFModel.c_str(),kModelDipole) == 0)
106  {
107  handle = true;
108  } else {
109  handle = false;
110  }
111  break;
112 
113  case ( kXSecTwkDial_MaCCQE ) :
114  if(fMode==kModeMa && strcmp(fFFModel.c_str(),kModelDipole) == 0)
115  {
116  handle = true;
117  } else {
118  handle = false;
119  }
120  break;
121 
122  case ( kXSecTwkDial_ZNormCCQE ) :
123  if(fMode==kModeZExp && strcmp(fFFModel.c_str(),kModelZExp) == 0)
124  {
125  handle = true;
126  } else {
127  handle = false;
128  }
129  break;
130 
131  case ( kXSecTwkDial_ZExpA1CCQE ):
132  case ( kXSecTwkDial_ZExpA2CCQE ):
133  case ( kXSecTwkDial_ZExpA3CCQE ):
134  case ( kXSecTwkDial_ZExpA4CCQE ):
135  if(fMode==kModeZExp && strcmp(fFFModel.c_str(),kModelZExp) == 0)
136  {
137  handle = true;
138  } else {
139  handle = false;
140  }
141  break;
142 
143  default:
144  handle = false;
145  break;
146  }
147 
148  return handle;
149 }
tweak Z-expansion coefficient 1, affects dsigma(CCQE)/dQ2 both in shape and normalization ...
Definition: GSyst.h:157
tweak CCQE normalization (energy independent)
Definition: GSyst.h:49
tweak Z-expansion coefficient 4, affects dsigma(CCQE)/dQ2 both in shape and normalization ...
Definition: GSyst.h:160
tweak Z-expansion coefficient 3, affects dsigma(CCQE)/dQ2 both in shape and normalization ...
Definition: GSyst.h:159
tweak Z-expansion coefficient 2, affects dsigma(CCQE)/dQ2 both in shape and normalization ...
Definition: GSyst.h:158
tweak Z-expansion CCQE normalization (energy independent)
Definition: GSyst.h:156
static const char * kModelZExp
tweak Ma CCQE, affects dsigma(CCQE)/dQ2 both in shape and normalization
Definition: GSyst.h:52
int fMode
0: Ma, 1: Norm and MaShape, 2: Z-Expansion
tweak Ma CCQE, affects dsigma(CCQE)/dQ2 in shape only (normalized to constant integral) ...
Definition: GSyst.h:51
static const char * kModelDipole
void GReWeightNuXSecCCQE::Reconfigure ( void  )
virtual

propagate updated nuisance parameter values to actual MC, etc

Implements genie::rew::GReWeightI.

Definition at line 202 of file GReWeightNuXSecCCQE.cxx.

203 {
205 
206  if(fMode==kModeMa && strcmp(fFFModel.c_str(),kModelDipole) == 0) {
207  int sign_matwk = utils::rew::Sign(fMaTwkDial);
208  double fracerr_ma = fracerr->OneSigmaErr(kXSecTwkDial_MaCCQE, sign_matwk);
209  fMaCurr = fMaDef * (1. + fMaTwkDial * fracerr_ma);
210  fMaCurr = TMath::Max(0., fMaCurr );
211  }
212  else
213  if(fMode==kModeNormAndMaShape && strcmp(fFFModel.c_str(),kModelDipole) == 0) {
214  int sign_normtwk = utils::rew::Sign(fNormTwkDial);
215  int sign_mashtwk = utils::rew::Sign(fMaTwkDial );
216  double fracerr_norm = fracerr->OneSigmaErr(kXSecTwkDial_NormCCQE, sign_normtwk);
217  double fracerr_mash = fracerr->OneSigmaErr(kXSecTwkDial_MaCCQEshape, sign_mashtwk);
218  fNormCurr = fNormDef * (1. + fNormTwkDial * fracerr_norm);
219  fMaCurr = fMaDef * (1. + fMaTwkDial * fracerr_mash);
220  fNormCurr = TMath::Max(0., fNormCurr);
221  fMaCurr = TMath::Max(0., fMaCurr );
222  }
223  else
224  if(fMode==kModeZExp && strcmp(fFFModel.c_str(),kModelZExp) == 0) {
225  int sign_twk = 0;
226  int sign_normtwk = utils::rew::Sign(fNormTwkDial);
227  double fracerr_zexp = 0.;
228  double fracerr_norm = fracerr->OneSigmaErr(kXSecTwkDial_ZNormCCQE, sign_normtwk);
229  GSyst_t syst;
230  // loop over all indices and update each
231  for (int i=0;i<fZExpMaxCoef;i++)
232  {
233  switch(i){
234  case 0: syst = kXSecTwkDial_ZExpA1CCQE; break;
235  case 1: syst = kXSecTwkDial_ZExpA2CCQE; break;
236  case 2: syst = kXSecTwkDial_ZExpA3CCQE; break;
237  case 3: syst = kXSecTwkDial_ZExpA4CCQE; break;
238  default: return; break;
239  }
240  sign_twk = utils::rew::Sign(fZExpTwkDial[i]);
241  fracerr_zexp = fracerr->OneSigmaErr(syst, sign_twk);
242  fZExpCurr[i] = fZExpDef[i] * (1. + fZExpTwkDial[i] * fracerr_zexp);
243  }
244  fNormCurr = fNormDef * (1. + fNormTwkDial * fracerr_norm);
245  fNormCurr = TMath::Max(0., fNormCurr);
246  }
247  else {
248  return;
249  }
250 
253  {
254  r.Set(fMaPath, fMaCurr);
255  }
256  else
257  if (fMode==kModeZExp)
258  {
259  ostringstream alg_key;
260  for (int i=0;i<fZExpMaxCoef;i++)
261  {
262  alg_key.str(""); // algorithm key for each coefficient
263  alg_key << fZExpPath << "QEL-Z_A" << i+1;
264  r.Set(alg_key.str(), fZExpCurr[i]);
265  }
266  }
267  fXSecModel->Configure(r);
268 }
tweak Z-expansion coefficient 1, affects dsigma(CCQE)/dQ2 both in shape and normalization ...
Definition: GSyst.h:157
tweak CCQE normalization (energy independent)
Definition: GSyst.h:49
string fMaPath
M_{A} path in config Registry.
double fZExpCurr[fZExpMaxSyst]
array of current parameter values
double OneSigmaErr(GSyst_t syst, int sign=0) const
tweak Z-expansion coefficient 4, affects dsigma(CCQE)/dQ2 both in shape and normalization ...
Definition: GSyst.h:160
Registry * fXSecModelConfig
config in tweaked model
tweak Z-expansion coefficient 3, affects dsigma(CCQE)/dQ2 both in shape and normalization ...
Definition: GSyst.h:159
tweak Z-expansion coefficient 2, affects dsigma(CCQE)/dQ2 both in shape and normalization ...
Definition: GSyst.h:158
virtual void Configure(const Registry &config)
Configure the algorithm.
Definition: Algorithm.cxx:70
tweak Z-expansion CCQE normalization (energy independent)
Definition: GSyst.h:156
An enumeration of systematic parameters.
static const char * kModelZExp
A registry. Provides the container for algorithm configuration parameters.
Definition: Registry.h:66
tweak Ma CCQE, affects dsigma(CCQE)/dQ2 both in shape and normalization
Definition: GSyst.h:52
string fZExpPath
algorithm path to get coefficients
int fMode
0: Ma, 1: Norm and MaShape, 2: Z-Expansion
static GSystUncertainty * Instance(void)
XSecAlgorithmI * fXSecModel
tweaked model
int fZExpMaxCoef
max number of coefficients to use
void Set(RgIMapPair entry)
Definition: Registry.cxx:281
tweak Ma CCQE, affects dsigma(CCQE)/dQ2 in shape only (normalized to constant integral) ...
Definition: GSyst.h:51
int Sign(double twkdial)
static const char * kModelDipole
void GReWeightNuXSecCCQE::Reset ( void  )
virtual

set all nuisance parameters to default values

Implements genie::rew::GReWeightI.

Definition at line 186 of file GReWeightNuXSecCCQE.cxx.

187 {
188  fNormTwkDial = 0.;
190  fMaTwkDial = 0.;
191  fMaCurr = fMaDef;
192 
193  for (int i=0;i<fZExpMaxSyst;i++)
194  {
195  fZExpTwkDial[i] = 0.;
196  fZExpCurr [i] = fZExpDef[i];
197  }
198 
199  this->Reconfigure();
200 }
double fZExpCurr[fZExpMaxSyst]
array of current parameter values
void Reconfigure(void)
propagate updated nuisance parameter values to actual MC, etc
static const int fZExpMaxSyst
maximum number of systematics
void genie::rew::GReWeightNuXSecCCQE::RewNue ( bool  tf)
inline

Definition at line 66 of file GReWeightNuXSecCCQE.h.

66 { fRewNue = tf; }
Definition: tf_graph.h:23
void genie::rew::GReWeightNuXSecCCQE::RewNuebar ( bool  tf)
inline

Definition at line 67 of file GReWeightNuXSecCCQE.h.

67 { fRewNuebar = tf; }
bool fRewNuebar
reweight nu_e_bar CC?
Definition: tf_graph.h:23
void genie::rew::GReWeightNuXSecCCQE::RewNumu ( bool  tf)
inline

Definition at line 68 of file GReWeightNuXSecCCQE.h.

68 { fRewNumu = tf; }
Definition: tf_graph.h:23
bool fRewNumu
reweight nu_mu CC?
void genie::rew::GReWeightNuXSecCCQE::RewNumubar ( bool  tf)
inline

Definition at line 69 of file GReWeightNuXSecCCQE.h.

69 { fRewNumubar = tf; }
Definition: tf_graph.h:23
bool fRewNumubar
reweight nu_mu_bar CC?
void genie::rew::GReWeightNuXSecCCQE::SetMaPath ( string  p)
inline

Definition at line 70 of file GReWeightNuXSecCCQE.h.

70 { fMaPath = p; }
string fMaPath
M_{A} path in config Registry.
p
Definition: test.py:228
void genie::rew::GReWeightNuXSecCCQE::SetMode ( int  mode)
inline

Definition at line 65 of file GReWeightNuXSecCCQE.h.

65 { fMode = mode; }
int fMode
0: Ma, 1: Norm and MaShape, 2: Z-Expansion
void GReWeightNuXSecCCQE::SetSystematic ( GSyst_t  syst,
double  val 
)
virtual

update the value for the specified nuisance param

Implements genie::rew::GReWeightI.

Definition at line 151 of file GReWeightNuXSecCCQE.cxx.

152 {
153  if(!this->IsHandled(syst))
154  {
155  LOG("ReW",pWARN) << "Systematic " << GSyst::AsString(syst) << " is not handled for algorithm "
156  << fFFModel << " and mode " << fMode;
157  return;
158  }
159 
160  switch(syst) {
161  case ( kXSecTwkDial_NormCCQE ) :
162  case ( kXSecTwkDial_ZNormCCQE ) :
163  fNormTwkDial = twk_dial;
164  break;
165  case ( kXSecTwkDial_MaCCQEshape ) :
166  case ( kXSecTwkDial_MaCCQE ) :
167  fMaTwkDial = twk_dial;
168  break;
169  case ( kXSecTwkDial_ZExpA1CCQE ) :
170  if(fZExpMaxCoef>0){ fZExpTwkDial[0] = twk_dial; }
171  break;
172  case ( kXSecTwkDial_ZExpA2CCQE ) :
173  if(fZExpMaxCoef>1){ fZExpTwkDial[1] = twk_dial; }
174  break;
175  case ( kXSecTwkDial_ZExpA3CCQE ) :
176  if(fZExpMaxCoef>2){ fZExpTwkDial[2] = twk_dial; }
177  break;
178  case ( kXSecTwkDial_ZExpA4CCQE ) :
179  if(fZExpMaxCoef>3){ fZExpTwkDial[3] = twk_dial; }
180  break;
181  default:
182  break;
183  }
184 }
tweak Z-expansion coefficient 1, affects dsigma(CCQE)/dQ2 both in shape and normalization ...
Definition: GSyst.h:157
tweak CCQE normalization (energy independent)
Definition: GSyst.h:49
tweak Z-expansion coefficient 4, affects dsigma(CCQE)/dQ2 both in shape and normalization ...
Definition: GSyst.h:160
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:87
tweak Z-expansion coefficient 3, affects dsigma(CCQE)/dQ2 both in shape and normalization ...
Definition: GSyst.h:159
tweak Z-expansion coefficient 2, affects dsigma(CCQE)/dQ2 both in shape and normalization ...
Definition: GSyst.h:158
tweak Z-expansion CCQE normalization (energy independent)
Definition: GSyst.h:156
#define pWARN
Definition: Messenger.h:51
tweak Ma CCQE, affects dsigma(CCQE)/dQ2 both in shape and normalization
Definition: GSyst.h:52
static string AsString(GSyst_t syst)
Definition: GSyst.h:175
int fMode
0: Ma, 1: Norm and MaShape, 2: Z-Expansion
int fZExpMaxCoef
max number of coefficients to use
tweak Ma CCQE, affects dsigma(CCQE)/dQ2 in shape only (normalized to constant integral) ...
Definition: GSyst.h:51
bool IsHandled(GSyst_t syst)
does the current weight calculator handle the input nuisance param?
void genie::rew::GReWeightNuXSecCCQE::SetZExpPath ( string  p)
inline

Definition at line 72 of file GReWeightNuXSecCCQE.h.

72 { fZExpPath = p; }
p
Definition: test.py:228
string fZExpPath
algorithm path to get coefficients

Member Data Documentation

string genie::rew::GReWeightNuXSecCCQE::fFFModel
private

Definition at line 85 of file GReWeightNuXSecCCQE.h.

double genie::rew::GReWeightNuXSecCCQE::fMaCurr
private

Definition at line 98 of file GReWeightNuXSecCCQE.h.

double genie::rew::GReWeightNuXSecCCQE::fMaDef
private

Definition at line 97 of file GReWeightNuXSecCCQE.h.

string genie::rew::GReWeightNuXSecCCQE::fMaPath
private

M_{A} path in config Registry.

Definition at line 92 of file GReWeightNuXSecCCQE.h.

double genie::rew::GReWeightNuXSecCCQE::fMaTwkDial
private

Definition at line 96 of file GReWeightNuXSecCCQE.h.

int genie::rew::GReWeightNuXSecCCQE::fMode
private

0: Ma, 1: Norm and MaShape, 2: Z-Expansion

Definition at line 87 of file GReWeightNuXSecCCQE.h.

double genie::rew::GReWeightNuXSecCCQE::fNormCurr
private

Definition at line 95 of file GReWeightNuXSecCCQE.h.

double genie::rew::GReWeightNuXSecCCQE::fNormDef
private

Definition at line 94 of file GReWeightNuXSecCCQE.h.

double genie::rew::GReWeightNuXSecCCQE::fNormTwkDial
private

Definition at line 93 of file GReWeightNuXSecCCQE.h.

bool genie::rew::GReWeightNuXSecCCQE::fRewNue
private

reweight nu_e CC?

Definition at line 88 of file GReWeightNuXSecCCQE.h.

bool genie::rew::GReWeightNuXSecCCQE::fRewNuebar
private

reweight nu_e_bar CC?

Definition at line 89 of file GReWeightNuXSecCCQE.h.

bool genie::rew::GReWeightNuXSecCCQE::fRewNumu
private

reweight nu_mu CC?

Definition at line 90 of file GReWeightNuXSecCCQE.h.

bool genie::rew::GReWeightNuXSecCCQE::fRewNumubar
private

reweight nu_mu_bar CC?

Definition at line 91 of file GReWeightNuXSecCCQE.h.

XSecAlgorithmI* genie::rew::GReWeightNuXSecCCQE::fXSecModel
private

tweaked model

Definition at line 83 of file GReWeightNuXSecCCQE.h.

Registry* genie::rew::GReWeightNuXSecCCQE::fXSecModelConfig
private

config in tweaked model

Definition at line 84 of file GReWeightNuXSecCCQE.h.

XSecAlgorithmI* genie::rew::GReWeightNuXSecCCQE::fXSecModelDef
private

default model

Definition at line 82 of file GReWeightNuXSecCCQE.h.

double genie::rew::GReWeightNuXSecCCQE::fZExpCurr[fZExpMaxSyst]
private

array of current parameter values

Definition at line 105 of file GReWeightNuXSecCCQE.h.

int genie::rew::GReWeightNuXSecCCQE::fZExpCurrIdx
private

current coefficient index

Definition at line 100 of file GReWeightNuXSecCCQE.h.

double genie::rew::GReWeightNuXSecCCQE::fZExpDef[fZExpMaxSyst]
private

Definition at line 104 of file GReWeightNuXSecCCQE.h.

int genie::rew::GReWeightNuXSecCCQE::fZExpMaxCoef
private

max number of coefficients to use

Definition at line 101 of file GReWeightNuXSecCCQE.h.

const int GReWeightNuXSecCCQE::fZExpMaxSyst = 4
static

maximum number of systematics

Definition at line 52 of file GReWeightNuXSecCCQE.h.

string genie::rew::GReWeightNuXSecCCQE::fZExpPath
private

algorithm path to get coefficients

Definition at line 102 of file GReWeightNuXSecCCQE.h.

double genie::rew::GReWeightNuXSecCCQE::fZExpTwkDial[fZExpMaxSyst]
private

Definition at line 103 of file GReWeightNuXSecCCQE.h.

const int GReWeightNuXSecCCQE::kModeMa = 0
static

Definition at line 48 of file GReWeightNuXSecCCQE.h.

const int GReWeightNuXSecCCQE::kModeNormAndMaShape = 1
static

Definition at line 49 of file GReWeightNuXSecCCQE.h.

const int GReWeightNuXSecCCQE::kModeZExp = 2
static

Definition at line 50 of file GReWeightNuXSecCCQE.h.


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