LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
GENIEReweight.cxx
Go to the documentation of this file.
1 
8 // C/C++ includes
9 #include <math.h>
10 #include <map>
11 #include <fstream>
12 
13 //ROOT includes
14 #include "TVector3.h"
15 #include "TLorentzVector.h"
16 #include "TSystem.h"
17 
18 //GENIE includes
19 #ifdef GENIE_PRE_R3
20  #include "Conventions/Units.h"
21  #include "EVGCore/EventRecord.h"
22  #include "GHEP/GHepUtils.h"
23  #include "Messenger/Messenger.h"
24 
25  #include "ReWeight/GReWeightI.h"
26  #include "ReWeight/GSystSet.h"
27  #include "ReWeight/GSyst.h"
28  #include "ReWeight/GReWeight.h"
29  #include "ReWeight/GReWeightNuXSecNCEL.h"
30  #include "ReWeight/GReWeightNuXSecCCQE.h"
31  #include "ReWeight/GReWeightNuXSecCCRES.h"
32  #include "ReWeight/GReWeightNuXSecCOH.h"
33  #include "ReWeight/GReWeightNonResonanceBkg.h"
34  #include "ReWeight/GReWeightFGM.h"
35  #include "ReWeight/GReWeightDISNuclMod.h"
36  #include "ReWeight/GReWeightResonanceDecay.h"
37  #include "ReWeight/GReWeightFZone.h"
38  #include "ReWeight/GReWeightINuke.h"
39  #include "ReWeight/GReWeightAGKY.h"
40  #include "ReWeight/GReWeightNuXSecCCQEvec.h"
41  #include "ReWeight/GReWeightNuXSecNCRES.h"
42  #include "ReWeight/GReWeightNuXSecDIS.h"
43  #include "ReWeight/GReWeightNuXSecNC.h"
44  #include "ReWeight/GSystUncertainty.h"
45  #include "ReWeight/GReWeightUtils.h"
46 
47 //#include "Geo/ROOTGeomAnalyzer.h"
48 //#include "Geo/GeomVolSelectorFiducial.h"
49 //#include "Geo/GeomVolSelectorRockBox.h"
50  #include "Utils/StringUtils.h"
51  #include "Utils/XmlParserUtils.h"
52  #include "Interaction/InitialState.h"
53  #include "Interaction/Interaction.h"
54  #include "Interaction/Kinematics.h"
55  #include "Interaction/KPhaseSpace.h"
56  #include "Interaction/ProcessInfo.h"
57  #include "Interaction/XclsTag.h"
58  #include "GHEP/GHepParticle.h"
59  #include "PDG/PDGCodeList.h"
60  #include "Conventions/Constants.h" //for calculating event kinematics
61 
62 #else
63  // GENIE includes R-3 and beyond
64  #include "GENIE/Framework/Messenger/Messenger.h"
65  #include "GENIE/Framework/Conventions/Units.h"
66  #include "GENIE/Framework/Conventions/Constants.h"
67  #include "GENIE/Framework/GHEP/GHepUtils.h"
68  #include "GENIE/Framework/EventGen/EventRecord.h"
69 
70 // #include "GENIE/Framework/Interaction/InitialState.h"
71 // #include "GENIE/Framework/Interaction/Interaction.h"
72 // #include "GENIE/Framework/Interaction/Kinematics.h"
73 // #include "GENIE/Framework/Interaction/KPhaseSpace.h"
74 // #include "GENIE/Framework/Interaction/ProcessInfo.h"
75 // #include "GENIE/Framework/Interaction/XclsTag.h"
76 
77 // #include "GENIE/Framework/ParticleData/PDGCodes.h"
78 // #include "GENIE/Framework/ParticleData/PDGCodeList.h"
79 // #include "GENIE/Framework/ParticleData/PDGLibrary.h"
80 // #include "GENIE/Framework/GHEP/GHepUtils.h"
81 // #include "GENIE/Framework/GHEP/GHepParticle.h"
82 
83  #include "RwFramework/GReWeightI.h"
84  #include "RwFramework/GSystSet.h"
85  #include "RwFramework/GSyst.h"
86  #include "RwFramework/GReWeight.h"
87  #include "RwFramework/GSystUncertainty.h"
88  #include "RwCalculators/GReWeightNuXSecNCEL.h"
89  #include "RwCalculators/GReWeightNuXSecCCQE.h"
90  #include "RwCalculators/GReWeightNuXSecCCRES.h"
91  #include "RwCalculators/GReWeightNuXSecCOH.h"
92  #include "RwCalculators/GReWeightNonResonanceBkg.h"
93  #include "RwCalculators/GReWeightFGM.h"
94  #include "RwCalculators/GReWeightDISNuclMod.h"
95  #include "RwCalculators/GReWeightResonanceDecay.h"
96  #include "RwCalculators/GReWeightFZone.h"
97  #include "RwCalculators/GReWeightINuke.h"
98  #include "RwCalculators/GReWeightAGKY.h"
99  #include "RwCalculators/GReWeightNuXSecCCQEvec.h"
100  #include "RwCalculators/GReWeightNuXSecNCRES.h"
101  #include "RwCalculators/GReWeightNuXSecDIS.h"
102  #include "RwCalculators/GReWeightNuXSecNC.h"
103  #include "RwCalculators/GReWeightUtils.h"
104 
105 //#include "Geo/ROOTGeomAnalyzer.h"
106 //#include "Geo/GeomVolSelectorFiducial.h"
107 //#include "Geo/GeomVolSelectorRockBox.h"
108 //#include "Utils/StringUtils.h"
109 //#include "Utils/XmlParserUtils.h"
110 
111 #endif
112 // Necessary because the GENIE LOG_* macros don't fully qualify Messenger
113 using genie::Messenger;
114 
115 
116 // NuGen includes
122 
123 // Framework includes
124 //#include "messagefacility/MessageLogger/MessageLogger.h"
125 
126 
127 
128 namespace rwgt {
131  fReweightNCEL(false),
132  fReweightQEMA(false),
133  fReweightQEVec(false),
134  fReweightCCRes(false),
135  fReweightNCRes(false),
136  fReweightResBkg(false),
137  fReweightResDecay(false),
138  fReweightNC(false),
139  fReweightDIS(false),
140  fReweightCoh(false),
141  fReweightAGKY(false),
142  fReweightDISNucMod(false),
143  fReweightFGM(false),
144  fReweightFZone(false),
145  fReweightINuke(false),
146  fReweightZexp(false),
147  fReweightMEC(false),
148  fMaQEshape(false),
149  fMaCCResShape(false),
150  fMaNCResShape(false),
151  fDISshape(false),
152  fUseSigmaDef(true) {
153 
154  LOG_INFO("GENIEReweight") << "Create GENIEReweight object";
155 
156  fWcalc = new genie::rew::GReWeight();
157  this->SetNominalValues();
158  }
159 
162  delete fWcalc;
163  }
164 
167  //CCQE Nominal Values
169 
171 
172  //CCQE Nominal Values
174 
176 
178 
181 
182  //Resonance Nominal Values
188 
194 
195  //Coherent pion nominal values
198 
199 
200  // Non-resonance background tweaking parameters:
217 
218 
219  // DIS tweaking parameters - applied for DIS events with [Q2>Q2o, W>Wo], typically Q2okReweight =1GeV^2, WokReweight =1.7-2.0GeV
224 
229 
231 
232 
233  fNominalParameters[(int)rwgt::fReweightRnubarnuCC] = 0.0;//v to vbar ratio reweighting is not working in GENIE at the moment
234  fNominalParameters[(int)rwgt::fReweightDISNuclMod] = 0.0;//The DIS nuclear model switch is not working in GENIE at the moment
235  //
236 
238 
239  //
240  // Hadronization [free nucleon target]
241  //
244 
245  //
246  // Medium-effects to hadronization
247  //
249 
250  //
251  // Intranuclear rescattering systematics.
255 #ifdef GENIE_PRE_R3
256  fNominalParameters[(int)rwgt::fReweightFrElas_pi] = 1.0;
257 #endif
262 #ifdef GENIE_PRE_R3
263  fNominalParameters[(int)rwgt::fReweightFrElas_N] = 1.0;
264 #endif
268 
269 
270  //
271  //RFG Nuclear model
272  //
275  //Continous "switch" at 0.0 full FG model. At 1.0 full spectral function model.
276  //From genie code it looks like weird stuff may happen for <0 and >1.
277  //This parameter does not have an "uncertainty" value associated with it.
278  //The tweaked dial value gets passed all the way through unchanged to the weight calculator
279 
280  //
281  // Resonance decays
282  //
286  //Continous "switch" at 0.0 full isotropic pion angular distribution. At 1.0 full R/S pion angular distribtuion.
287  //This parameter does not have an "uncertainty" value associated with it.
288  //The tweaked dial value gets passed all the way through unchanged to the weight calculator
289  }
290 
293  double nominal_value;
294  nominal_value = fNominalParameters[rLabel];
295  return nominal_value;
296  }
297 
300  int label = (int)rLabel;
301  bool in_loop = false;
302  bool found_par = false;
303  double parameter = -10000;
304  for(unsigned int i = 0; i < fReWgtParameterName.size(); i++) {
305  in_loop = true;
306  if(label == fReWgtParameterName[i]) {
307  parameter = fReWgtParameterValue[i];
308  found_par = true;
309  }
310  }
311  if(in_loop) {
312  if(found_par) return parameter;
313  else {
314  LOG_WARN("GENIEReweight") << "Parameter " << label << " not set yet or does not exist";
315  return parameter;
316  }
317  }
318  else {
319  LOG_WARN("GENIEReweight") << "Vector of parameters has not been initialized yet (Size is 0).";
320  return parameter;
321  }
322  }
323 
326  int label = (int)rLabel;
327  LOG_INFO("GENIEReweight") << "Adding parameter: " << genie::rew::GSyst::AsString(genie::rew::EGSyst(label)) << ". With value: " << value;
328  fReWgtParameterName.push_back(label);
329  fReWgtParameterValue.push_back(value);
330 
331  }
332 
335  int label = (int)rLabel;
336  bool found = false;
337  for(unsigned int i = 0; i < fReWgtParameterName.size(); i++) {
338  if(fReWgtParameterName[i] == label) {
340  found = true;
341  }
342  }
343  if(!found) {
344  this->AddReweightValue(rLabel, value);
345  }
346  }
347 
350  LOG_INFO("GENIEReweight") << "Configure weight calculator";
351 
352  for(unsigned int i = 0; i < fReWgtParameterName.size(); i++) {
353 
354  switch (fReWgtParameterName[i])
355  {
356  //NC Elastic parameters
359  fReweightNCEL = true;
360  break;
361 
362  //CC QE Axial parameters
367  fReweightQEMA = true;
368  break;
369 
370  //CC QE Vector parameters
372  fReweightQEVec = true;
373  break;
374 
375  //CC Resonance parameters
381  fReweightCCRes = true;
382  break;
383 
384  //NC Resonance parameters
390  fReweightNCRes = true;
391  break;
392 
393  //Coherent parameters
396  fReweightCoh = true;
397  break;
398 
399  //Non-resonance background KNO parameters
416  fReweightResBkg = true;
417  break;
418 
419  //DIS parameters
430  fReweightDIS = true;
431  break;
432 
433  //DIS nuclear model parameters
435  fReweightDISNucMod = true;
436  break;
437 
438  //NC cross section
439  case rwgt::fReweightNC:
440  fReweightNC = true;
441  break;
442 
443  //Hadronization parameters
446  fReweightAGKY = true;
447  break;
448 
449  //Elastic (and QE) nuclear model parameters
452  fReweightFGM = true;
453  break;
454 
455  //Formation Zone
457  fReweightFZone = true;
458  break;
459 
460  //Intranuke Parameters
464 #ifdef GENIE_PRE_R3
465  case rwgt::fReweightFrElas_pi:
466 #endif
471 #ifdef GENIE_PRE_R3
472  case rwgt::fReweightFrElas_N:
473 #endif
477  fReweightINuke = true;
478  break;
479 
480  //Resonance Decay parameters
484  fReweightResDecay = true;
485  break;
486 
487  //Z-expansion parameters
494  fReweightZexp = true;
495  break;
496 
497  } // switch(fReWgtParameterName[i])
498 
499  } //end for loop
500 
501  //configure the individual weight calculators
502  if(fReweightNCEL) this->ConfigureNCEL();
504  if(fReweightQEVec) this->ConfigureQEVec();
505  if(fReweightCCRes) this->ConfigureCCRes();
506  if(fReweightNCRes) this->ConfigureNCRes();
507  if(fReweightResBkg) this->ConfigureResBkg();
509  if(fReweightNC) this->ConfigureNC();
510  if(fReweightDIS) this->ConfigureDIS();
511  if(fReweightCoh) this->ConfigureCoh();
512  if(fReweightAGKY) this->ConfigureAGKY();
514  if(fReweightFGM) this->ConfigureFGM();
515  if(fReweightFZone) this->ConfigureFZone();
516  if(fReweightINuke) this->ConfigureINuke();
517  this->ConfigureParameters();
518 
519  }
520 
523  delete fWcalc;
524  fWcalc = new genie::rew::GReWeight();
525  this->Configure();
526  }
527 
531  void GENIEReweight::ReweightNCEL(double ma, double eta) {
532  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight for NC Elastic Reweighting";
533  if(ma!=0.0) {
535  }
536  if(eta!=0.0) {
538  }
539  this->Configure();
540  }
543  void GENIEReweight::ReweightQEMA(double ma) {
544  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight for QE Axial Mass Reweighting";
545  fMaQEshape=false;
547  this->Configure();
548  }
551  void GENIEReweight::ReweightQEVec(double mv) {
552  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight for QE Vector Mass Reweighting";
554  this->Configure();
555  }
556 
557  void GENIEReweight::ReweightQEZExp(double norm, double a1, double a2, double a3, double a4)
558  {
559  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight for Z-expansion QE Reweighting";
565  this->Configure();
566  }
569  void GENIEReweight::ReweightCCRes(double ma, double mv) {
570  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight for CC Resonance Reweighting";
571  fMaCCResShape=false;
573  if(mv!=0.0) {
575  }
576  this->Configure();
577  }
580  void GENIEReweight::ReweightNCRes(double ma, double mv) {
581  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight for NC Resonance Reweighting";
582  fMaNCResShape=false;
584  if(mv!=0.0) {
586  }
587  this->Configure();
588  }
591  void GENIEReweight::ReweightResGanged(double ma, double mv) {
592  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight for CC and NC Resonance Reweighting";
593  fMaCCResShape=false;
594  fMaNCResShape=false;
597  if(mv!=0.0) {
600  }
601  this->Configure();
602  }
605  void GENIEReweight::ReweightCoh(double ma, double r0) {
606  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight for Coherant Reweighting";
609  this->Configure();
610  }
611 
613  //Here it is being configured for v+p and vbar + n (1 pi) type interactions
614  void GENIEReweight::ReweightNonResRvp1pi(double sigma) {
615  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight for Non-Resonance Background Reweighting (Neutrino Single Pion)";
620  this->Configure();
621  }
622 
624  //Here it is being configured for v+n and vbar + p (1 pi) type interactions
625  void GENIEReweight::ReweightNonResRvbarp1pi(double sigma) {
626  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight for Non-Resonance Background Reweighting (Anti-Neutrino Single Pion)";
631  this->Configure();
632  }
635  void GENIEReweight::ReweightNonResRvp2pi(double sigma) {
636  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight for Non-Resonance Background Reweighting (Neutrino Two Pion)";
641  this->Configure();
642  }
643 
645  // Here it is being configured for v+n and vbar + p (2 pi) type interactions
646  void GENIEReweight::ReweightNonResRvbarp2pi(double sigma) {
647  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight for Non-Resonance Background Reweighting (Anti-Neutrino Two Pion)";
652  this->Configure();
653  }
656  void GENIEReweight::ReweightResDecay(double gamma, double eta, double theta) {
657  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight for Resoncance Decay Parameters";
658  if(gamma!=0.0) {
660  }
661  if(eta!=0.0) {
663  }
664  if(theta!=0.0) {
666  }
667  this->Configure();
668  }
671  void GENIEReweight::ReweightNC(double norm) {
672  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight for NC Cross Section Scale";
673  this->AddReweightValue(rwgt::fReweightNC, norm);
674  this->Configure();
675  }
678  void GENIEReweight::ReweightDIS(double aht, double bht, double cv1u, double cv2u) {
679  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight for DIS Form Factor Model Reweighting";
680  fDISshape = false;
681  if(aht != 0.0) {
683  }
684  if(bht != 0.0) {
686  }
687  if(cv1u != 0.0) {
689  }
690  if(cv2u != 0.0) {
692  }
693  this->Configure();
694  }
697  void GENIEReweight::ReweightDISnucl(bool mode) {
698  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight for DIS Nuclear Model";
700  this->Configure();
701  }
704  void GENIEReweight::ReweightAGKY(double xF, double pT) {
705  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight for DIS AGKY Hadronization Model Reweighting";
706  if(xF==0.0) {
708  }
709  if(pT==0.0) {
711  }
712  this->Configure();
713  }
716  void GENIEReweight::ReweightIntraNuke(ReweightLabel_t name, double sigma) {
717  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight for Intranuke Model Reweighting";
718  if ( (name==rwgt::fReweightMFP_pi) ||
719  (name==rwgt::fReweightMFP_N) ||
720  (name==rwgt::fReweightFrCEx_pi) ||
721 #ifdef GENIE_PRE_R3
722  (name==rwgt::fReweightFrElas_pi) ||
723 #endif
724  (name==rwgt::fReweightFrInel_pi) ||
725  (name==rwgt::fReweightFrAbs_pi) ||
726  (name==rwgt::fReweightFrPiProd_pi) ||
727  (name==rwgt::fReweightFrCEx_N) ||
728 #ifdef GENIE_PRE_R3
729  (name==rwgt::fReweightFrElas_N) ||
730 #endif
731  (name==rwgt::fReweightFrInel_N ) ||
732  (name==rwgt::fReweightFrAbs_N ) ||
733  (name==rwgt::fReweightFrPiProd_N) ) {
734  this->AddReweightValue(name, sigma);
735  }
736  else {
737  LOG_WARN("GENIEReweight") << "That is not a valid Intranuke parameter Intranuke not configured";
738  }
739  this->Configure();
740  }
743  void GENIEReweight::ReweightFormZone(double sigma) {
744  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight for Formation Zone Reweighting";
746  this->Configure();
747  }
750  void GENIEReweight::ReweightFGM(double kF, double sf) {
751  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight for Fermi Gas Model Reweighting";
754  this->Configure();
755  }
756 
758 
762  LOG_INFO("GENIEReweight") << "Adding NC elastic weight calculator";
763  fWcalc->AdoptWghtCalc( "xsec_ncel", new GReWeightNuXSecNCEL );
764  }
765 
768  LOG_INFO("GENIEReweight") << "Adding CCQE axial FF weight calculator ";
769  fWcalc->AdoptWghtCalc( "xsec_ccqe", new GReWeightNuXSecCCQE );
770  GReWeightNuXSecCCQE *rwccqe = dynamic_cast <GReWeightNuXSecCCQE*> (fWcalc->WghtCalc("xsec_ccqe"));
771  if (fReweightZexp)
772  {
773  LOG_INFO("GENIEReweight") << "in z-expansion mode";
774  rwccqe->SetMode(GReWeightNuXSecCCQE::kModeZExp);
775  }
776  else if(!fMaQEshape) {
777  LOG_INFO("GENIEReweight") << "in axial mass (QE) rate+shape mode";
778  rwccqe->SetMode(GReWeightNuXSecCCQE::kModeMa);
779  }
780  else {
781  LOG_INFO("GENIEReweight") << "in axial mass (QE) shape only mode";
782  }
783  }
784 
787  LOG_INFO("GENIEReweight") << "Adding CCQE vector FF weight calculator";
788  fWcalc->AdoptWghtCalc( "xsec_ccqe_vec", new GReWeightNuXSecCCQEvec );
789  }
790 
793  LOG_INFO("GENIEReweight") << "Adding CC resonance weight calculator";
794  fWcalc->AdoptWghtCalc( "xsec_ccres", new GReWeightNuXSecCCRES );
795  if(!fMaCCResShape) {
796  LOG_INFO("GENIEReweight") << "in axial mass (Res) rate+shape mode";
797  GReWeightNuXSecCCRES * rwccres = dynamic_cast<GReWeightNuXSecCCRES *> (fWcalc->WghtCalc("xsec_ccres"));
798  rwccres->SetMode(GReWeightNuXSecCCRES::kModeMaMv);
799  }
800  else {
801  LOG_INFO("GENIEReweight") << "in axial mass (Res) shape only mode";
802  }
803  }
804 
807  LOG_INFO("GENIEReweight") << "Adding NC resonance weight calculator";
808  fWcalc->AdoptWghtCalc( "xsec_ncres", new GReWeightNuXSecNCRES );
809  if(!fMaNCResShape) {
810  LOG_INFO("GENIEReweight") << "in axial mass (Res) rate+shape mode";
811  GReWeightNuXSecNCRES * rwncres = dynamic_cast<GReWeightNuXSecNCRES *> (fWcalc->WghtCalc("xsec_ncres"));
812  rwncres->SetMode(GReWeightNuXSecNCRES::kModeMaMv);
813  }
814  else {
815  LOG_INFO("GENIEReweight") << "in axial mass (Res) shape only mode";
816  }
817  }
818 
821  LOG_INFO("GENIEReweight") << "Adding low Q^2 DIS (KNO) weight calculator";
822  fWcalc->AdoptWghtCalc( "xsec_nonresbkg", new GReWeightNonResonanceBkg );
823  }
824 
827  LOG_INFO("GENIEReweight") << "Adding resonance decay weight calculator";
828  fWcalc->AdoptWghtCalc( "hadro_res_decay", new GReWeightResonanceDecay );
829  }
830 
833  LOG_INFO("GENIEReweight") << "Adding NC total cross section weight calculator";
834  fWcalc->AdoptWghtCalc( "xsec_nc", new GReWeightNuXSecNC );
835  }
836 
839  LOG_INFO("GENIEReweight") << "Adding DIS (Bodek-Yang) weight calculator";
840  fWcalc->AdoptWghtCalc( "xsec_dis", new GReWeightNuXSecDIS );
841  if(!fDISshape) {
842  LOG_INFO("GENIEReweight") << "in shape+rate mode";
843  GReWeightNuXSecDIS * rwdis = dynamic_cast<GReWeightNuXSecDIS *> (fWcalc->WghtCalc("xsec_dis"));
844  rwdis->SetMode(GReWeightNuXSecDIS::kModeABCV12u);
845  }
846  else {
847  LOG_INFO("GENIEReweight") << "in shape only mode";
848  }
849  }
850 
853  LOG_INFO("GENIEReweight") << "Adding coherant interaction model weight calculator";
854  fWcalc->AdoptWghtCalc( "xsec_coh", new GReWeightNuXSecCOH );
855  }
856 
859  LOG_INFO("GENIEReweight") << "Adding hadronization (AGKY) model weight calculator";
860  fWcalc->AdoptWghtCalc( "hadro_agky", new GReWeightAGKY );
861  }
862 
865  LOG_INFO("GENIEReweight") << "Adding DIS nuclear model weight calculator";
866  fWcalc->AdoptWghtCalc( "nuclear_dis", new GReWeightDISNuclMod );
867  }
868 
871  LOG_INFO("GENIEReweight") << "Adding Fermi Gas Model (FGM) weight calculator";
872  fWcalc->AdoptWghtCalc( "nuclear_qe", new GReWeightFGM );
873  }
874 
877  LOG_INFO("GENIEReweight") << "Adding Formation Zone weight calculator";
878  fWcalc->AdoptWghtCalc( "hadro_fzone", new GReWeightFZone );
879  }
880 
883  LOG_INFO("GENIEReweight") << "Adding the Intra-Nuke weight calculator";
884  fWcalc->AdoptWghtCalc( "hadro_intranuke", new GReWeightINuke );
885  }
886 
889  GSystSet & syst = fWcalc->Systematics();
890  for(unsigned int i = 0; i < fReWgtParameterName.size(); i++) {
891  LOG_INFO("GENIEReweight") << "Configuring GENIEReweight parameter: " << genie::rew::GSyst::AsString(genie::rew::EGSyst(fReWgtParameterName[i])) << " with value: " << fReWgtParameterValue[i];
892  if(fUseSigmaDef) {
893  syst.Set( (GSyst_t)fReWgtParameterName[i], fReWgtParameterValue[i]);
894  }
895  else {
897  syst.Set( (GSyst_t)fReWgtParameterName[i], parameter);
898  }
899  }
900  fWcalc->Reconfigure();
901  }
902 
905  double GENIEReweight::CalculateSigma(ReweightLabel_t label, double value) {
906  //double GENIEReweight::CalculateSigma(int label, double value) {
907  int iLabel = (int) label;
908  double nominal_def;
909  double frac_err;
910  double sigma;
911  int sign;
912  GSystUncertainty * gsysterr = GSystUncertainty::Instance();
915  label==rwgt::fReweightDISNuclMod) {
916  //These parameters don't use the sigma definition just pass them through the function unchanged
917  sigma = value;
918  }
919  else {
920  nominal_def = this->NominalParameterValue(label);
921  sign = genie::utils::rew::Sign(value-nominal_def);
922  frac_err = gsysterr->OneSigmaErr( (GSyst_t)iLabel, sign);
923  sigma = (value - nominal_def)/(frac_err*nominal_def);
924  }
925  return sigma;
926  }
927 
929  //double GENIEReweight::CalcWeight(simb::MCTruth truth, simb::GTruth gtruth) {
930  double GENIEReweight::CalculateWeight(const genie::EventRecord& evr) const {
931  //genie::EventRecord evr = this->RetrieveGHEP(truth, gtruth);
932  double wgt = fWcalc->CalcWeight(evr);
933  //mf::LogVerbatim("GENIEReweight") << "New Event Weight is: " << wgt;
934  return wgt;
935  }
936 
937  /*
939  genie::EventRecord GENIEReweight::RetrieveGHEP(simb::MCTruth truth, simb::GTruth gtruth) {
940 
941  genie::EventRecord newEvent;
942  newEvent.SetWeight(gtruth.fweight);
943  newEvent.SetProbability(gtruth.fprobability);
944  newEvent.SetXSec(gtruth.fXsec);
945  newEvent.SetDiffXSec(gtruth.fDiffXsec);
946  TLorentzVector vtx = gtruth.fVertex;
947  newEvent.SetVertex(vtx);
948 
949  for(int i = 0; i < truth.NParticles(); i++) {
950  simb::MCParticle mcpart = truth.GetParticle(i);
951 
952  int gmid = mcpart.PdgCode();
953  genie::GHepStatus_t gmst = (genie::GHepStatus_t)mcpart.StatusCode();
954  int gmmo = mcpart.Mother();
955  int ndaughters = mcpart.NumberDaughters();
956  //find the track ID of the first and last daughter particles
957  int fdtrkid = 0;
958  int ldtrkid = 0;
959  if(ndaughters !=0) {
960  fdtrkid = mcpart.Daughter(0);
961  if(ndaughters == 1) {
962  ldtrkid = 1;
963  }
964  else if(ndaughters >1) {
965  fdtrkid = mcpart.Daughter(ndaughters-1);
966  }
967  }
968  int gmfd = -1;
969  int gmld = -1;
970  //Genie uses the index in the particle array to reference the daughter particles.
971  //MCTruth keeps the particles in the same order so use the track ID to find the proper index.
972  for(int j = 0; j < truth.NParticles(); j++) {
973  simb::MCParticle temp = truth.GetParticle(i);
974  if(temp.TrackId() == fdtrkid) {
975  gmfd = j;
976  }
977  if(temp.TrackId() == ldtrkid) {
978  gmld = j;
979  }
980  }
981 
982  double gmpx = mcpart.Px(0);
983  double gmpy = mcpart.Py(0);
984  double gmpz = mcpart.Pz(0);
985  double gme = mcpart.E(0);
986  double gmvx = mcpart.Gvx();
987  double gmvy = mcpart.Gvy();
988  double gmvz = mcpart.Gvz();
989  double gmvt = mcpart.Gvt();
990  int gmri = mcpart.Rescatter();
991 
992  genie::GHepParticle gpart(gmid, gmst, gmmo, -1, gmfd, gmld, gmpx, gmpy, gmpz, gme, gmvx, gmvy, gmvz, gmvt);
993  gpart.SetRescatterCode(gmri);
994  TVector3 polz = mcpart.Polarization();
995  if(polz.x() !=0 || polz.y() !=0 || polz.z() !=0) {
996  gpart.SetPolarization(polz);
997  }
998  newEvent.AddParticle(gpart);
999 
1000  }
1001 
1002  genie::ProcessInfo proc_info;
1003  genie::ScatteringType_t gscty = (genie::ScatteringType_t)gtruth.fGscatter;
1004  genie::InteractionType_t ginty = (genie::InteractionType_t)gtruth.fGint;
1005 
1006  proc_info.Set(gscty,ginty);
1007 
1008  genie::XclsTag gxt;
1009 
1010  //Set Exclusive Final State particle numbers
1011  genie::Resonance_t gres = (genie::Resonance_t)gtruth.fResNum;
1012  gxt.SetResonance(gres);
1013  gxt.SetNPions(gtruth.fNumPiPlus, gtruth.fNumPi0, gtruth.fNumPiMinus);
1014  gxt.SetNNucleons(gtruth.fNumProton, gtruth.fNumNeutron);
1015 
1016  if(gtruth.fIsCharm) {
1017  gxt.SetCharm(0);
1018  }
1019  else {
1020  gxt.UnsetCharm();
1021  }
1022 
1023  //Set the GENIE kinematic info
1024  genie::Kinematics gkin;
1025  gkin.Setx(gtruth.fgX, true);
1026  gkin.Sety(gtruth.fgY, true);
1027  gkin.Sett(gtruth.fgT, true);
1028  gkin.SetW(gtruth.fgW, true);
1029  gkin.SetQ2(gtruth.fgQ2, true);
1030  gkin.Setq2(gtruth.fgq2, true);
1031  simb::MCNeutrino nu = truth.GetNeutrino();
1032  simb::MCParticle lep = nu.Lepton();
1033  gkin.SetFSLeptonP4(lep.Px(), lep.Py(), lep.Pz(), lep.E());
1034  gkin.SetHadSystP4(gtruth.fFShadSystP4.Px(), gtruth.fFShadSystP4.Py(), gtruth.fFShadSystP4.Pz(), gtruth.fFShadSystP4.E());
1035 
1036  //Set the GENIE final state interaction info
1037  genie::Interaction * p_gint = new genie::Interaction;
1038  genie::InitialState * p_ginstate = p_gint->InitStatePtr();
1039  //int Z = gtruth.ftgtZ;
1040  //int A = gtruth.ftgtA;
1041  int targetNucleon = nu.HitNuc();
1042  int struckQuark = nu.HitQuark();
1043  int incoming = gtruth.fProbePDG;
1044  p_ginstate->SetProbePdg(incoming);
1045 
1046  genie::Target* target123 = p_ginstate->TgtPtr();
1047 
1048  target123->SetId(gtruth.ftgtPDG);
1049  //target123->SetId(Z,A);
1050 
1051  //int pdg_code = pdg::IonPdgCode(A, Z);
1052  //TParticlePDG * p = PDGLibrary::Instance()->Find(pdg_code);
1053 
1054  //mf::LogWarning("GENIEReweight") << "Setting Target Z and A";
1055  //mf::LogWarning("GENIEReweight") << "Saved PDG: " << gtruth.ftgtPDG;
1056  //mf::LogWarning("GENIEReweight") << "Target PDG: " << target123->Pdg();
1057  target123->SetHitNucPdg(targetNucleon);
1058  target123->SetHitQrkPdg(struckQuark);
1059  target123->SetHitSeaQrk(gtruth.fIsSeaQuark);
1060 
1061  if(newEvent.HitNucleonPosition()> 0) {
1062  genie::GHepParticle * hitnucleon = newEvent.HitNucleon();
1063  std::auto_ptr<TLorentzVector> p4hitnucleon(hitnucleon->GetP4());
1064  target123->SetHitNucP4(*p4hitnucleon);
1065  }
1066  else {
1067  TLorentzVector dummy(0.,0.,0.,0.);
1068  target123->SetHitNucP4(dummy);
1069  }
1070 
1071  genie::GHepParticle * probe = newEvent.Probe();
1072  std::auto_ptr<TLorentzVector> p4probe(probe->GetP4());
1073  p_ginstate->SetProbeP4(*p4probe);
1074  if(newEvent.TargetNucleusPosition()> 0) {
1075  genie::GHepParticle * target = newEvent.TargetNucleus();
1076  std::auto_ptr<TLorentzVector> p4target(target->GetP4());
1077  p_ginstate->SetTgtP4(*p4target);
1078  } else {
1079  TLorentzVector dummy(0.,0.,0.,0.);
1080  p_ginstate->SetTgtP4(dummy);
1081  }
1082  p_gint->SetProcInfo(proc_info);
1083  p_gint->SetKine(gkin);
1084  p_gint->SetExclTag(gxt);
1085  newEvent.AttachSummary(p_gint);
1086 
1087 
1088  //For temporary debugging purposes
1089  //genie::Interaction *inter = newEvent.Summary();
1090  //const genie::InitialState &initState = inter->InitState();
1091  //const genie::Target &tgt = initState.Tgt();
1092  //std::cout << "TargetPDG as Recorded: " << gtruth.ftgtPDG << std::endl;
1093  //std::cout << "TargetZ as Recorded: " << gtruth.ftgtZ << std::endl;
1094  //std::cout << "TargetA as Recorded: " << gtruth.ftgtA << std::endl;
1095  //std::cout << "TargetPDG as Recreated: " << tgt.Pdg() << std::endl;
1096  //std::cout << "TargetZ as Recreated: " << tgt.Z() << std::endl;
1097  //std::cout << "TargetA as Recreated: " << tgt.A() << std::endl;
1098 
1099  return newEvent;
1100 
1101  }
1102  */
1103 
1104 
1105 }
void AddReweightValue(ReweightLabel_t rLabel, double value)
Change a reweight parameter. If it hasn&#39;t been added yet add it.
tweak inelastic probability for pions, for given total rescattering probability
void ReweightIntraNuke(ReweightLabel_t name, double sigma)
Simple Configuration of the Formation Zone reweight calculator.
tweak absorption probability for nucleons, for given total rescattering probability ...
void ReweightQEVec(double mv)
enum rwgt::EReweightLabel ReweightLabel_t
void ConfigureINuke()
configure the weight parameters being used
tweak the 2pi non-RES bkg in the RES region, for v+n CC
tweak the 1pi non-RES bkg in the RES region, for vbar+p NC
tweak Mv CCRES, affects d2sigma(CCRES)/dWdQ2 both in shape and normalization
tweak the Bodek-Yang model parameter B_{ht} - incl. both shape and normalization effect ...
tweak axial nucleon form factors (dipole -> z-expansion) - shape only effect of dsigma(CCQE)/dQ2 ...
void ConfigureCCRes()
Configure the NCRES calculator.
void ConfigureResBkg()
Configure the ResDecay weight calculator.
void ConfigureDIS()
Configure the Coherant model weight calculator.
void Configure()
Reconfigure the weight calculators.
tweak the Bodek-Yang model parameter CV2u - incl. both shape and normalization effect ...
tweak Ma CCRES, affects d2sigma(CCRES)/dWdQ2 in shape only (normalized to constant integral) ...
void ConfigureQEVec()
Configure the CCRES calculator.
void ReweightAGKY(double xF, double pT)
Simple Configuration of the Intranuke Nuclear model.
tweak the 2pi non-RES bkg in the RES region, for v+p CC
void ReweightResDecay(double gamma, double eta, double theta)
Simple Configuration of the Total NC cross section.
tweak inelastic probability for nucleons, for given total rescattering probability ...
void ReweightCCRes(double ma, double mv=0.0)
Simple Configurtion of the NC Resonance weight calculator.
void ReweightDIS(double aht, double bht, double cv1u, double cv2u)
Simple Configuration of the DIS nuclear model.
std::map< int, double > fNominalParameters
tweak Z-expansion CCQE normalization (energy independent)
void ConfigureNC()
Configure the DIS (Bodek-Yang) weight calculator.
void ReweightFGM(double kF, double sf)
End of Simple Reweight Configurations.
void ConfigureFGM()
Configure the Formation Zone weight calculator.
double CalculateWeight(const genie::EventRecord &evr) const
tweak Resonance -> X + gamma branching ratio, eg Delta+(1232) -> p gamma
void ReweightQEMA(double ma)
Simple Configuration of the CCQE vector weight calculator.
tweak the 1pi non-RES bkg in the RES region, for vbar+p CC
tweak pion production probability for nucleons, for given total rescattering probability ...
tweak xF distribution for low multiplicity (N + pi) DIS f/s produced by AGKY
std::vector< int > fReWgtParameterName
void ConfigureFZone()
Configure the intranuke weight calculator.
tweak the 2pi non-RES bkg in the RES region, for v+p NC
tweak CCQE normalization (maintains dependence on neutrino energy)
void ReweightQEZExp(double norm, double a1, double a2, double a3, double a4)
Simple Configuration of the CC Resonance weight calculator.
Particle class.
tweak the 2pi non-RES bkg in the RES region, for vbar+n CC
tweak CCQE normalization (energy independent)
tweak the 2pi non-RES bkg in the RES region, for vbar+n NC
tweak the 2pi non-RES bkg in the RES region, for vbar+p NC
void ReweightNonResRvbarp1pi(double sigma)
Simple Configuration of the Non-Resonance Background weight calculator. Here it is being configured f...
void ConfigureQEMA()
Configure the QE vector FF weight calculator.
void ReweightCoh(double ma, double r0)
Simple Configuration of the Non-Resonance Background weight calculator.
tweak formation zone
void ConfgureResDecay()
Configure the NC weight calculator.
#define a2
tweak pT distribution for low multiplicity (N + pi) DIS f/s produced by AGKY
#define a3
tweak the 2pi non-RES bkg in the RES region, for v+n NC
void SetNominalValues()
Return the nominal value for the given parameter.
tweak Ma NCRES, affects d2sigma(NCRES)/dWdQ2 both in shape and normalization
tweak Z-expansion coefficient 1, affects dsigma(CCQE)/dQ2 both in shape and normalization ...
GENIEReweight()
<constructor
tweak Mv NCRES, affects d2sigma(NCRES)/dWdQ2 both in shape and normalization
tweak the 1pi non-RES bkg in the RES region, for v+n CC
genie::rew::GReWeight * fWcalc
void ReweightResGanged(double ma, double mv=0.0)
Simple Configuration of the Coherant weight calculator.
void ConfigureNCRes()
Configure the ResBkg (kno) weight calculator.
tweak absorption probability for pions, for given total rescattering probability
tweak NCRES normalization
tweak Z-expansion coefficient 2, affects dsigma(CCQE)/dQ2 both in shape and normalization ...
void Reconfigure()
Simple Configuration functions for configuring a single weight calculator.
Wrapper for reweight neutrino interactions with GENIE base class.
void ChangeParameterValue(ReweightLabel_t rLabel, double value)
Configure the weight calculators.
void ReweightNonResRvbarp2pi(double sigma)
Simple Configuration of the Resonance decay model weight calculator.
tweak the Bodek-Yang model parameter CV1u - incl. both shape and normalization effect ...
tweak the 1pi non-RES bkg in the RES region, for vbar+n NC
void ReweightNCRes(double ma, double mv=0.0)
Simple Configuration of the NC and CC Resonance weight calculator with the axial mass parameter for N...
tweak the 1pi non-RES bkg in the RES region, for vbar+n CC
tweak mean free path for pions
tweak the 1pi non-RES bkg in the RES region, for v+p NC
tweak the Bodek-Yang model parameter A_{ht} - shape only effect to d2sigma(DIS)/dxdy ...
void ConfigureAGKY()
Configure the DIS nuclear model weight calculator.
double value
Definition: spectrum.C:18
void ConfigureNCEL()
Configure the MaQE weight calculator.
distort pi angular distribution in Delta -> N + pi
double NominalParameterValue(ReweightLabel_t rLabel)
Return the configured value of the given parameter.
tweak DIS nuclear modification (shadowing, anti-shadowing, EMC). Does not appear to be working in GEN...
~GENIEReweight()
Set the nominal values for the reweight parameters.
tweak the inclusive DIS CC normalization (not currently working in genie)
tweak the Bodek-Yang model parameter CV2u - shape only effect to d2sigma(DIS)/dxdy ...
tweak Z-expansion coefficient 3, affects dsigma(CCQE)/dQ2 both in shape and normalization ...
int sign(double val)
Definition: UtilFunc.cxx:104
tweak Mv CCRES, affects d2sigma(CCRES)/dWdQ2 in shape only (normalized to constant integral) ...
tweak Resonance -> X + eta branching ratio, eg N+(1440) -> p eta
std::vector< double > fReWgtParameterValue
double ReweightParameterValue(ReweightLabel_t rLabel)
Add reweight parameters to the list.
tweak the Bodek-Yang model parameter B_{ht} - shape only effect to d2sigma(DIS)/dxdy ...
Float_t norm
tweak Ma NCEL, affects dsigma(NCEL)/dQ2 both in shape and normalization
void ReweightNCEL(double ma, double eta)
Simple Configurtion of the CCQE axial weight calculator.
void ReweightNonResRvp1pi(double sigma)
Simple Configuration of the Non-Resonance Background weight calculator.
tweak charge exchange probability for nucleons, for given total rescattering probability ...
tweak the Bodek-Yang model parameter A_{ht} - incl. both shape and normalization effect ...
double CalculateSigma(ReweightLabel_t label, double value)
Calculate the weights.
tweak Z-expansion coefficient 4, affects dsigma(CCQE)/dQ2 both in shape and normalization ...
tweak Ma NCRES, affects d2sigma(NCRES)/dWdQ2 in shape only (normalized to constant integral) ...
tweak Ma CCRES, affects d2sigma(CCRES)/dWdQ2 both in shape and normalization
tweak NCEL strange axial form factor eta, affects dsigma(NCEL)/dQ2 both in shape and normalization ...
#define a4
tweak the 2pi non-RES bkg in the RES region, for vbar+p CC
tweak Ma CCQE, affects dsigma(CCQE)/dQ2 in shape only (normalized to constant integral) ...
tweak pion production probability for pions, for given total rescattering probability ...
tweak Ma CCQE, affects dsigma(CCQE)/dQ2 both in shape and normalization
tweak Mv NCRES, affects d2sigma(NCRES)/dWdQ2 in shape only (normalized to constant integral) ...
tweak charge exchange probability for pions, for given total rescattering probability ...
tweak the Bodek-Yang model parameter CV1u - shape only effect to d2sigma(DIS)/dxdy ...
tweak elastic nucleon form factors (BBA/default -> dipole) - shape only effect of dsigma(CCQE)/dQ2 ...
tweak mean free path for nucleons
void ReweightNC(double norm)
Simple Configuration of the DIS FF model weight calculator.
void ReweightNonResRvp2pi(double sigma)
Simple Configuration of the Non-Resonance Background weight calculator.
tweak the 1pi non-RES bkg in the RES region, for v+p CC
void ReweightFormZone(double sigma)
Simple Configuration of the Fermigas model reweight calculator.
tweak Ma for COH pion production
void ConfigureDISNucMod()
Configure the FG model weight calculator.
tweak CCRES normalization
tweak the 1pi non-RES bkg in the RES region, for v+n NC
tweak the ratio of ( CC) / ( CC) (not currently working in genie)
void ReweightDISnucl(bool mode)
Simple Configuration of the DIS AGKY hadronization model.
tweak R0 for COH pion production
#define a1
void ConfigureCoh()
Configure the hadronization (AGKY) weight calculator.