Public Member Functions | Public Attributes | Private Attributes | List of all members
NeutrinoFluxReweight::ThinTargetpCNucleonReweighter Class Reference

Reweighter of thin target p,n production. More...

#include <ThinTargetpCNucleonReweighter.h>

Inheritance diagram for NeutrinoFluxReweight::ThinTargetpCNucleonReweighter:
NeutrinoFluxReweight::IInteractionReweighting

Public Member Functions

 ThinTargetpCNucleonReweighter (int iuniv, const ParameterTable &cv_pars, const ParameterTable &univ_pars)
 
virtual ~ThinTargetpCNucleonReweighter ()
 
virtual bool canReweight (const InteractionData &aa)
 can the particular instance of this class reweight this interaction? More...
 
virtual double calculateWeight (const InteractionData &aa)
 calculate a weight for this interaction given the central value parameters and the parameters for this universe. The weight is something like: f(cv)/f(MC) * f(univ)/f(cv) where cv in this case corresponds to the best value of the parameter, given the data. If univ_pars=cv_pars then we are calculating a central value weight More...
 
double calculateDataScale (int inc_pdg, double inc_mom, int prod_pdg, double xf, double pt)
 
- Public Member Functions inherited from NeutrinoFluxReweight::IInteractionReweighting
virtual ~IInteractionReweighting ()
 

Public Attributes

double data_prod_xs
 
std::vector< float > vbin_data_p
 
std::vector< float > vbin_data_n
 

Private Attributes

int iUniv
 
const ParameterTablecvPars
 
const ParameterTableunivPars
 

Detailed Description

Reweighter of thin target p,n production.

Definition at line 15 of file ThinTargetpCNucleonReweighter.h.

Constructor & Destructor Documentation

NeutrinoFluxReweight::ThinTargetpCNucleonReweighter::ThinTargetpCNucleonReweighter ( int  iuniv,
const ParameterTable cv_pars,
const ParameterTable univ_pars 
)

Definition at line 10 of file ThinTargetpCNucleonReweighter.cpp.

10  :iUniv(iuniv),cvPars(cv_pars),univPars(univ_pars){
11 
12  ThinTargetBins* Thinbins = ThinTargetBins::getInstance();
13 
14  vbin_data_p.reserve(Thinbins->GetNbins_pC_pX_NA49());
15  vbin_data_n.reserve(Thinbins->GetNbins_pC_nX_NA49());
16  // const boost::interprocess::flat_map<std::string, double>& cv_table = cvPars.getMap();
17  // const boost::interprocess::flat_map<std::string, double>& univ_table = univPars.getMap();
18 
19  data_prod_xs = univPars.getParameterValue("prod_prtC_xsec");
20 
21  //the number of bins needs to be written from the xmls files
22  char namepar[100];
23  for(int ii=0;ii<Thinbins->GetNbins_pC_pX_NA49();ii++){
24 
25  sprintf(namepar,"ThinTarget_pC_%s_sys_%d","p",ii);
26  double data_cv = cvPars.getParameterValue(std::string(namepar));
27  double data_sys = univPars.getParameterValue(std::string(namepar));
28  sprintf(namepar,"ThinTarget_pC_%s_stats_%d","p",ii);
29  double data_sta = univPars.getParameterValue(std::string(namepar));
30  vbin_data_p.push_back(data_sta + data_sys - data_cv);
31 
32  }
33 
34  for(int ii=0;ii<Thinbins->GetNbins_pC_nX_NA49();ii++){
35 
36  sprintf(namepar,"ThinTarget_pC_%s_sys_%d","n",ii);
37  double data_cv = cvPars.getParameterValue(std::string(namepar));
38  double data_sys = univPars.getParameterValue(std::string(namepar));
39  sprintf(namepar,"ThinTarget_pC_%s_stats_%d","n",ii);
40  double data_sta = univPars.getParameterValue(std::string(namepar));
41  vbin_data_n.push_back(data_sta + data_sys - data_cv);
42 
43  }
44 
45 
46  }
std::string string
Definition: nybbler.cc:12
double getParameterValue(const std::string &name) const
get the value of a parameter. throw an exception of a well defined type if we don&#39;t have it ...
static ThinTargetBins * getInstance()
NeutrinoFluxReweight::ThinTargetpCNucleonReweighter::~ThinTargetpCNucleonReweighter ( )
virtual

Definition at line 48 of file ThinTargetpCNucleonReweighter.cpp.

48  {
49 
50  }

Member Function Documentation

double NeutrinoFluxReweight::ThinTargetpCNucleonReweighter::calculateDataScale ( int  inc_pdg,
double  inc_mom,
int  prod_pdg,
double  xf,
double  pt 
)

Definition at line 119 of file ThinTargetpCNucleonReweighter.cpp.

119  {
120 
121  double scaling_violation = 1.0;
122  ThinTargetMC* dtH = ThinTargetMC::getInstance();
123  //temporary:
124  const int Nscl = 11;
125  const int moms[Nscl] = {12,20,31,40,50,60,70,80,100,120,158};
126 
127  int idx_part = -1;
128  if(prod_pdg == 2212)idx_part = 4;
129  else if(prod_pdg == 2112)idx_part = 5;
130  else{
131  std::cout<<"Error in the prod particle"<<std::endl;
132  return 1.0;
133  }
134 
135  int idx_lowp = -1;
136  int idx_hip = -1;
137  for(int i=0;i<Nscl-1;i++){
138  if(inc_mom>=double(moms[i]) && inc_mom<double(moms[i+1])){
139  idx_lowp=i;
140  idx_hip =i+1;
141  }
142  }
143  if(idx_lowp<0 || idx_hip<0){
144  std::cout<<"Error calculating the scaling"<<std::endl;
145  return 1.0;
146  }
147 
148  double scl_ref158 = -1.0;
149  double scl_m = 0.0;
150 
151  if(idx_part==4){
152  int binid = dtH->hTTScl[idx_part][Nscl-1]->FindBin(xf,pt);
153  scl_ref158 = dtH->hTTScl[idx_part][Nscl-1]->GetBinContent(binid);
154  //just provisional... the scaling just reach up to xF=0.8975... consering a close xf value
155  if(xf>0.8975){
156  binid = dtH->hTTScl[idx_part][Nscl-1]->FindBin(0.89,pt);
157  scl_ref158 = dtH->hTTScl[idx_part][Nscl-1]->GetBinContent(binid);
158  }
159  double scl_low = dtH->hTTScl[idx_part][idx_lowp]->GetBinContent(binid);
160  double scl_hi = dtH->hTTScl[idx_part][idx_hip]->GetBinContent(binid);
161  scl_m = scl_low + (inc_mom-double(moms[idx_lowp]))*(scl_hi-scl_low)/(double(moms[idx_hip])-double(moms[idx_lowp]));
162  }
163  else if(idx_part==5){
164  int binid = dtH->hTTScl_n[Nscl-1]->FindBin(xf);
165  scl_ref158 = (double)dtH->hTTScl_n[Nscl-1]->GetBinContent(binid);
166  double scl_low = dtH->hTTScl_n[idx_lowp]->GetBinContent(binid);
167  double scl_hi = dtH->hTTScl_n[idx_hip]->GetBinContent(binid);
168  scl_m = scl_low + (inc_mom-double(moms[idx_lowp]))*(scl_hi-scl_low)/(double(moms[idx_hip])-double(moms[idx_lowp]));
169  }
170  else{
171  std::cout<<"still error, not expected here!!"<<std::endl;
172  return 1.0;
173  }
174 
175  if(scl_ref158 < 1.e-10 || scl_m<1.e-10){
176  std::cout<<"scale problems: "<<scl_ref158<<" "<<scl_m<<" "<<inc_pdg<<" "<<inc_mom<<" "<<prod_pdg<<" "<<xf<<" "<<pt<<std::endl;
177  return 1.0;
178  }
179  scaling_violation = scl_m/scl_ref158;
180  return scaling_violation;
181  }
const double e
QTextStream & endl(QTextStream &s)
static ThinTargetMC * getInstance()
double NeutrinoFluxReweight::ThinTargetpCNucleonReweighter::calculateWeight ( const InteractionData inter_data)
virtual

calculate a weight for this interaction given the central value parameters and the parameters for this universe. The weight is something like: f(cv)/f(MC) * f(univ)/f(cv) where cv in this case corresponds to the best value of the parameter, given the data. If univ_pars=cv_pars then we are calculating a central value weight

Implements NeutrinoFluxReweight::IInteractionReweighting.

Definition at line 74 of file ThinTargetpCNucleonReweighter.cpp.

74  {
75 
76  //quick check:
77  double wgt = 1.0;
78  ThinTargetBins* Thinbins = ThinTargetBins::getInstance();
79  int bin_p = Thinbins->BinID_pC_p(aa.xF,aa.Pt,aa.Prod_pdg);
80  int bin_n = Thinbins->BinID_pC_n(aa.xF,aa.Prod_pdg);
81  if(bin_p < 0 && bin_n < 0){
82  //std::cout<<"Not bin found "<<std::endl;
83  return wgt;
84  }
85 
86  //Calculating the scale:
87  double data_scale = calculateDataScale(aa.Inc_pdg,aa.Inc_P,aa.Prod_pdg,aa.xF,aa.Pt);
88  double dataval = -1;
89  if(aa.Prod_pdg==2212)dataval = vbin_data_p[bin_p];
90  if(aa.Prod_pdg==2112)dataval = vbin_data_n[bin_n];
91 
92  if(dataval<1.e-12){
93  //std::cout<<"Not data found "<<std::endl;
94  return wgt;
95  }
96 
97  //checking if this is the first interaction:
98  if(aa.gen==0 && aa.Prod_pdg==2212)dataval /= data_prod_xs;
99  if(aa.gen>0 && aa.Prod_pdg==2212)dataval /= 1.0;
100  if(aa.gen==0 && aa.Prod_pdg==2112)dataval /= 1.0;
101  if(aa.gen>0 && aa.Prod_pdg==2112)dataval *= data_prod_xs;
102 
103  dataval *= data_scale;
104 
105  ThinTargetMC* mc = ThinTargetMC::getInstance();
106  double mc_cv = mc->getMCval_pC_X(aa.Inc_P,aa.xF,aa.Pt,aa.Prod_pdg);
107  double mc_prod = mc->getMCxs_pC_nucleon(aa.gen,aa.Prod_pdg,aa.Inc_P);
108  mc_cv /= mc_prod;
109  if(mc_cv<1.e-12)return wgt;
110  wgt = dataval/mc_cv;
111  if(wgt<0){
112  std::cout<<"TTPCNU check wgt(<0) "<<iUniv<<" "<<wgt<<" "<<aa.Inc_P<<" "<<aa.xF<<" "<<aa.Pt<<" "<<aa.Prod_pdg<<" "<<(mc->getMCxs_pC_nucleon(aa.gen,aa.Prod_pdg,aa.Inc_P))<<std::endl;
113  return 1.0;
114  }
115 
116  return wgt;
117  }
double calculateDataScale(int inc_pdg, double inc_mom, int prod_pdg, double xf, double pt)
const double e
static ThinTargetBins * getInstance()
QTextStream & endl(QTextStream &s)
static ThinTargetMC * getInstance()
bool NeutrinoFluxReweight::ThinTargetpCNucleonReweighter::canReweight ( const InteractionData aa)
virtual

can the particular instance of this class reweight this interaction?

Implements NeutrinoFluxReweight::IInteractionReweighting.

Definition at line 51 of file ThinTargetpCNucleonReweighter.cpp.

51  {
52 
53  //checking:
54  std::string mode(getenv("MODE"));
55  if(aa.Inc_pdg != 2212)return false;
56  if(aa.Inc_P < 12.0)return false;
57  //volume check:
58  bool is_wrong_volume = aa.Vol != "TGT1" && aa.Vol != "BudalMonitor" && aa.Vol != "Budal_HFVS" && aa.Vol != "Budal_VFHS";
59  if( (mode=="REF") || (mode=="OPT") ){
60  is_wrong_volume = aa.Vol != "TargetFinHorizontal" && aa.Vol != "TargetNoSplitSegment" && aa.Vol!="tCoreLog";
61  }
62  if(is_wrong_volume)return false;
63  //
64  if(aa.Prod_pdg != 2212 && aa.Prod_pdg != 2112)return false;
65 
66  ThinTargetBins* Thinbins = ThinTargetBins::getInstance();
67  int bin_p = Thinbins->BinID_pC_p(aa.xF,aa.Pt,aa.Prod_pdg);
68  int bin_n = Thinbins->BinID_pC_n(aa.xF,aa.Prod_pdg);
69  if(bin_p < 0 && bin_n < 0)return false;
70  else return true;
71 
72  }
std::string string
Definition: nybbler.cc:12
std::string getenv(std::string const &name)
Definition: getenv.cc:15
static ThinTargetBins * getInstance()

Member Data Documentation

const ParameterTable& NeutrinoFluxReweight::ThinTargetpCNucleonReweighter::cvPars
private

Definition at line 29 of file ThinTargetpCNucleonReweighter.h.

double NeutrinoFluxReweight::ThinTargetpCNucleonReweighter::data_prod_xs

Definition at line 24 of file ThinTargetpCNucleonReweighter.h.

int NeutrinoFluxReweight::ThinTargetpCNucleonReweighter::iUniv
private

Definition at line 28 of file ThinTargetpCNucleonReweighter.h.

const ParameterTable& NeutrinoFluxReweight::ThinTargetpCNucleonReweighter::univPars
private

Definition at line 30 of file ThinTargetpCNucleonReweighter.h.

std::vector<float> NeutrinoFluxReweight::ThinTargetpCNucleonReweighter::vbin_data_n

Definition at line 25 of file ThinTargetpCNucleonReweighter.h.

std::vector<float> NeutrinoFluxReweight::ThinTargetpCNucleonReweighter::vbin_data_p

Definition at line 25 of file ThinTargetpCNucleonReweighter.h.


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