LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
PPFXCVWeightCalc.cxx
Go to the documentation of this file.
1 
4 
5 #include "CLHEP/Random/RandGaussQ.h"
6 
7 #include "MakeReweight.h"
8 #include "cetlib/search_path.h"
9 #include "dk2nu/tree/dk2nu.h"
10 #include "dk2nu/tree/dkmeta.h"
12 
13 #include "TSystem.h"
14 
15 namespace evwgh {
16  class PPFXCVWeightCalc : public WeightCalc {
17  public:
19  void Configure(fhicl::ParameterSet const& p, CLHEP::HepRandomEngine& engine) override;
20  std::vector<std::vector<double>> GetWeight(art::Event& e) override;
21 
22  private:
23  std::string fGenieModuleLabel;
24 
25  std::vector<std::string> fInputLabels;
26  std::string fPPFXMode;
27  std::string fMode;
28  std::string fHorn;
29  std::string fTarget;
30  int fSeed;
31  int fVerbose;
32  NeutrinoFluxReweight::MakeReweight* fPPFXrw;
33 
35  };
36 
38 
39  void PPFXCVWeightCalc::Configure(fhicl::ParameterSet const& p, CLHEP::HepRandomEngine& engine)
40  {
41  //get configuration for this function
43  fGenieModuleLabel = p.get<std::string>("genie_module_label");
44 
45  //ppfx setup
46  fInputLabels = pset.get<std::vector<std::string>>("input_labels");
47  fPPFXMode = pset.get<std::string>("ppfx_mode");
48  fVerbose = pset.get<int>("verbose");
49  fMode = pset.get<std::string>("mode");
50  fHorn = pset.get<std::string>("horn_curr");
51  fTarget = pset.get<std::string>("target_config");
52  fSeed = pset.get<int>("random_seed", -1);
53 
54  gSystem->Setenv("MODE", fPPFXMode.c_str());
55 
56  fPPFXrw = NeutrinoFluxReweight::MakeReweight::getInstance();
57 
58  std::string inputOptions; // Full path.
59  std::string mode_file = "inputs_" + fPPFXMode + ".xml"; // Just file name.
60  cet::search_path sp("FW_SEARCH_PATH");
61  sp.find_file(mode_file, inputOptions);
62 
63  if (fSeed != -1) fPPFXrw->setBaseSeed(fSeed); // Set the random seed
64  std::cout << "is PPFX setup : " << fPPFXrw->AlreadyInitialized() << std::endl;
65  std::cout << "Setting PPFX, inputs: " << inputOptions << std::endl;
66  std::cout << "Setting Horn Current Configuration to: " << fHorn << std::endl;
67  std::cout << "Setting Target Configuration to: " << fTarget << std::endl;
68  if (!(fPPFXrw->AlreadyInitialized())) { fPPFXrw->SetOptions(inputOptions); }
69  std::cout << "PPFX just set with mode: " << fPPFXMode << std::endl;
70  }
71 
72  std::vector<std::vector<double>> PPFXCVWeightCalc::GetWeight(art::Event& e)
73  {
74  std::vector<std::vector<double>> weight;
76  //calculate weight(s) here
77 
78  int nmctruth = 0, nmatched = 0;
79  bool flag = true;
80  int ievt = -1;
81  std::vector<art::Handle<std::vector<bsim::Dk2Nu>>> dk2nus2 =
82  e.getMany<std::vector<bsim::Dk2Nu>>();
83  for (size_t dk2 = 0; dk2 < dk2nus2.size(); ++dk2) {
84  art::Handle<std::vector<bsim::Dk2Nu>> dk2nus = dk2nus2[dk2];
85  }
86  while ((flag = mcitr.Next())) {
87  std::string label = mcitr.GetLabel();
88  const simb::MCTruth* pmctruth = mcitr.GetMCTruth();
89  // const simb::GTruth* pgtruth = mcitr.GetGTruth();
90  //not-used//const simb::MCFlux* pmcflux = mcitr.GetMCFlux();
91  const bsim::Dk2Nu* pdk2nu = mcitr.GetDk2Nu();
92  //not-used//const bsim::NuChoice* pnuchoice = mcitr.GetNuChoice();
93  // art::Ptr<simb::MCTruth> mctruthptr = mcitr.GetMCTruthPtr();
94 
95  ++ievt;
96  ++nmctruth;
97  if (fVerbose > 0) {
98  std::cout << "FluxWeightCalculator [" << std::setw(4) << ievt << "] "
99  << " label \"" << label << "\" MCTruth@ " << pmctruth << " Dk2Nu@ " << pdk2nu
100  << std::endl;
101  }
102  if (!pdk2nu) continue;
103  ++nmatched;
104 
105  //RWH//bsim::Dk2Nu* tmp_dk2nu = fTmpDK2NUConversorAlg.construct_toy_dk2nu( &mctruth,&mcflux);
106  //RWH//bsim::DkMeta* tmp_dkmeta = fTmpDK2NUConversorAlg.construct_toy_dkmeta(&mctruth,&mcflux);
107  //RWH// those appear to have been memory leaks
108  //RWH// herein begins the replacment for the "construct_toy_dkmeta"
109 
110  static bsim::DkMeta
111  dkmeta_obj; //RWH// create this on stack (destroyed when out-of-scope) ... or static
112  dkmeta_obj.tgtcfg = fTarget;
113  dkmeta_obj.horncfg = fHorn;
114  (dkmeta_obj.vintnames).push_back("Index_Tar_In_Ancestry");
115  (dkmeta_obj.vintnames).push_back("Playlist");
116  bsim::DkMeta* tmp_dkmeta = &dkmeta_obj;
117  //RWH// herein ends this block
118 
119  // sigh ....
120  //RWH// this is the signature in NeutrinoFluxReweight::MakeReweight :
121  // //! create an interaction chain from the new dk2nu(dkmeta) format
122  // void calculateWeights(bsim::Dk2Nu* nu, bsim::DkMeta* meta);
123  //RWH// these _should_ be "const <class>*" because we don't need to change them
124  //RWH// and the pointers we get out of the ART record are going to be const.
125  bsim::Dk2Nu* tmp_dk2nu = const_cast<bsim::Dk2Nu*>(pdk2nu); // remove const-ness
126  try {
127  fPPFXrw->calculateWeights(tmp_dk2nu, tmp_dkmeta);
128  }
129  catch (...) {
130  std::cerr << "Failed to calcualte weight" << std::endl;
131  continue;
132  }
133  //Get weights:
134  if (fMode == "reweight") {
135  double ppfx_cv_wgt = fPPFXrw->GetCVWeight();
136  std::vector<double> wvec = {ppfx_cv_wgt};
137  weight.push_back(wvec);
138  }
139  else {
140  std::vector<double> vmipppion = fPPFXrw->GetWeights("MIPPNumiPionYields");
141  std::vector<double> vmippkaon = fPPFXrw->GetWeights("MIPPNumiKaonYields");
142  std::vector<double> vtgtatt = fPPFXrw->GetWeights("TargetAttenuation");
143  std::vector<double> vabsorp = fPPFXrw->GetWeights("TotalAbsorption");
144  std::vector<double> vttpcpion = fPPFXrw->GetWeights("ThinTargetpCPion");
145  std::vector<double> vttpckaon = fPPFXrw->GetWeights("ThinTargetpCKaon");
146  std::vector<double> vttpcnucleon = fPPFXrw->GetWeights("ThinTargetpCNucleon");
147  std::vector<double> vttncpion = fPPFXrw->GetWeights("ThinTargetnCPion");
148  std::vector<double> vttnucleona = fPPFXrw->GetWeights("ThinTargetnucleonA");
149  std::vector<double> vttmesoninc = fPPFXrw->GetWeights("ThinTargetMesonIncident");
150  std::vector<double> vothers = fPPFXrw->GetWeights("Other");
151 
152  std::vector<double> tmp_vhptot;
153  for (unsigned int iuniv = 0; iuniv < vtgtatt.size(); iuniv++) {
154  tmp_vhptot.push_back(float(vmipppion[iuniv] * vmippkaon[iuniv] * vtgtatt[iuniv] *
155  vabsorp[iuniv] * vttpcpion[iuniv] * vttpckaon[iuniv] *
156  vttpcnucleon[iuniv] * vttncpion[iuniv] * vttnucleona[iuniv] *
157  vttmesoninc[iuniv] * vothers[iuniv]));
158  }
159  weight.push_back(tmp_vhptot);
160  }
161  }
162 
163  return weight;
164  }
166 }
NeutrinoFluxReweight::MakeReweight * fPPFXrw
void Configure(fhicl::ParameterSet const &p, CLHEP::HepRandomEngine &engine) override
const simb::MCTruth * GetMCTruth() const
std::string GetName()
Definition: WeightCalc.h:26
#define DECLARE_WEIGHTCALC(wghcalc)
T get(std::string const &key) const
Definition: ParameterSet.h:314
std::vector< std::string > fInputLabels
double weight
Definition: plottest35.C:25
Event generator information.
Definition: MCTruth.h:32
Float_t e
Definition: plot.C:35
std::vector< std::vector< double > > GetWeight(art::Event &e) override
std::string GetLabel() const
const bsim::Dk2Nu * GetDk2Nu() const
std::vector< Handle< PROD > > getMany(SelectorBase const &selector=MatchAllSelector{}) const
#define REGISTER_WEIGHTCALC(wghcalc)