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

Reweighter of thin target K production. More...

#include <ThinTargetpCKaonReweighter.h>

Inheritance diagram for NeutrinoFluxReweight::ThinTargetpCKaonReweighter:
NeutrinoFluxReweight::IInteractionReweighting

Public Member Functions

 ThinTargetpCKaonReweighter (int iuniv, const ParameterTable &cv_pars, const ParameterTable &univ_pars)
 
virtual ~ThinTargetpCKaonReweighter ()
 
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_kap
 
std::vector< float > vbin_data_kam
 
std::vector< float > mipp_vbin_data_kap_pip
 
std::vector< float > mipp_vbin_data_kam_pim
 

Private Attributes

int iUniv
 
const ParameterTablecvPars
 
const ParameterTableunivPars
 
ThinTargetpCPionReweightertt_pCPionRew
 

Detailed Description

Reweighter of thin target K production.

Definition at line 16 of file ThinTargetpCKaonReweighter.h.

Constructor & Destructor Documentation

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

Definition at line 10 of file ThinTargetpCKaonReweighter.cpp.

10  :iUniv(iuniv),cvPars(cv_pars),univPars(univ_pars){
11 
12  ThinTargetBins* Thinbins = ThinTargetBins::getInstance();
13  vbin_data_kap.reserve(Thinbins->GetNbins_pC_KX_NA49());
14  vbin_data_kam.reserve(Thinbins->GetNbins_pC_KX_NA49());
15  mipp_vbin_data_kap_pip.reserve(Thinbins->GetNbins_pC_KX_MIPP());
16  mipp_vbin_data_kam_pim.reserve(Thinbins->GetNbins_pC_KX_MIPP());
17 
18  // const boost::interprocess::flat_map<std::string, double>& cv_table = cvPars.getMap();
19  // const boost::interprocess::flat_map<std::string, double>& univ_table = univPars.getMap();
20 
21  data_prod_xs = univPars.getParameterValue("prod_prtC_xsec");
22  char namepar[100];
23  for(int ii=0;ii<Thinbins->GetNbins_pC_KX_NA49();ii++){
24 
25  sprintf(namepar,"ThinTargetLowxF_pC_%s_sys_%d","kap",ii);
26  double data_cv = cvPars.getParameterValue(std::string(namepar));
27  double data_sys = univPars.getParameterValue(std::string(namepar));
28  sprintf(namepar,"ThinTargetLowxF_pC_%s_stats_%d","kap",ii);
29  double data_sta = univPars.getParameterValue(std::string(namepar));
30  vbin_data_kap.push_back(data_sta + data_sys - data_cv);
31 
32  sprintf(namepar,"ThinTargetLowxF_pC_%s_sys_%d","kam",ii);
33  data_cv = cvPars.getParameterValue(std::string(namepar));
34  data_sys = univPars.getParameterValue(std::string(namepar));
35  sprintf(namepar,"ThinTargetLowxF_pC_%s_stats_%d","kam",ii);
36  data_sta = univPars.getParameterValue(std::string(namepar));
37  vbin_data_kam.push_back(data_sta + data_sys - data_cv);
38 
39  }
40 
41  for(int ii=0;ii<Thinbins->GetNbins_pC_KX_MIPP();ii++){
42 
43  sprintf(namepar,"ThinTarget_kap_pip_sys_%d",ii);
44  double data_cv = cvPars.getParameterValue(std::string(namepar));
45  double data_sys = univPars.getParameterValue(std::string(namepar));
46  sprintf(namepar,"ThinTarget_kap_pip_stats_%d",ii);
47  double data_sta = univPars.getParameterValue(std::string(namepar));
48  mipp_vbin_data_kap_pip.push_back(data_sta + data_sys - data_cv);
49 
50  sprintf(namepar,"ThinTarget_kam_pim_sys_%d",ii);
51  data_cv = cvPars.getParameterValue(std::string(namepar));
52  data_sys = univPars.getParameterValue(std::string(namepar));
53  sprintf(namepar,"ThinTarget_kam_pim_stats_%d",ii);
54  data_sta = univPars.getParameterValue(std::string(namepar));
55  mipp_vbin_data_kam_pim.push_back(data_sta + data_sys - data_cv);
56 
57  }
58 
59  }
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::ThinTargetpCKaonReweighter::~ThinTargetpCKaonReweighter ( )
virtual

Definition at line 61 of file ThinTargetpCKaonReweighter.cpp.

61  {
62 
63  }

Member Function Documentation

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

Definition at line 204 of file ThinTargetpCKaonReweighter.cpp.

204  {
205  double scaling_violation = 1.0;
206  ThinTargetMC* dtH = ThinTargetMC::getInstance();
207  //temporary:
208  const int Nscl = 11;
209  const int moms[Nscl] = {12,20,31,40,50,60,70,80,100,120,158};
210 
211  int idx_part = -1;
212  if(prod_pdg == 321)idx_part = 2;
213  if(prod_pdg ==-321)idx_part = 3;
214  if(idx_part<0){
215  std::cout<<"Error in the prod particle"<<std::endl;
216  return 1.0;
217  }
218 
219  int binid = dtH->hTTScl[idx_part][Nscl-1]->FindBin(xf,pt);
220  double scl_ref158 = dtH->hTTScl[idx_part][Nscl-1]->GetBinContent(binid);
221 
222  int idx_lowp = -1;
223  int idx_hip = -1;
224  for(int i=0;i<Nscl-1;i++){
225  if(inc_mom>=double(moms[i]) && inc_mom<double(moms[i+1])){
226  idx_lowp=i;
227  idx_hip =i+1;
228  }
229  }
230  if(idx_lowp<0 || idx_hip<0){
231  std::cout<<"Error calculating the scaling"<<std::endl;
232  return 1.0;
233  }
234  double scl_low = dtH->hTTScl[idx_part][idx_lowp]->GetBinContent(binid);
235  double scl_hi = dtH->hTTScl[idx_part][idx_hip]->GetBinContent(binid);
236  double scl_m = scl_low + (inc_mom-double(moms[idx_lowp]))*(scl_hi-scl_low)/(double(moms[idx_hip])-double(moms[idx_lowp]));
237  if(scl_ref158<1.e-10){
238  // std::cout<<"ref158 zero!!! "<<scl_ref158<<std::endl;
239  return 1.0;
240  }
241  scaling_violation = scl_m/scl_ref158;
242  return scaling_violation;
243  }
const double e
QTextStream & endl(QTextStream &s)
static ThinTargetMC * getInstance()
double NeutrinoFluxReweight::ThinTargetpCKaonReweighter::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 109 of file ThinTargetpCKaonReweighter.cpp.

109  {
110 
111  double wgt = 1.0;
112 
113  ThinTargetBins* Thinbins = ThinTargetBins::getInstance();
114  int bin = Thinbins->BinID_pC_k(aa.xF,aa.Pt,aa.Prod_pdg);
115  int mipp_bin = Thinbins->mipp_BinID_pC_k(aa.Pz,aa.Pt,aa.Prod_pdg);
116  if(bin<0 && mipp_bin<0)return wgt; //double check
117 
118  double dataval = -1;
119 
120  if(bin>=0){
121  if(aa.Prod_pdg == 321)dataval = vbin_data_kap[bin];
122  else if(aa.Prod_pdg == -321) dataval = vbin_data_kam[bin];
123  else if(aa.Prod_pdg == 130 || aa.Prod_pdg == 310){
124  dataval = 0.25*(vbin_data_kap[bin]+3.*vbin_data_kam[bin]);
125  }
126  else{
127  std::cout<<"Something is wrong, pdg_prod: "<< aa.Prod_pdg <<std::endl;
128  return wgt;
129  }
130  }
131  else if(mipp_bin>=0){
132 
133  MakeReweight* makerew = MakeReweight::getInstance();
134  if(iUniv==-1)tt_pCPionRew = (makerew->cv_rw)->THINTARGET_PC_PION_Universe;
135  else tt_pCPionRew = (makerew->vec_rws[iUniv])->THINTARGET_PC_PION_Universe;
136 
137  double inc_mom[3] = {aa.Inc_P4[0], aa.Inc_P4[1], aa.Inc_P4[2]};
138  double prod_mom[3] = {aa.Prod_P4[0],aa.Prod_P4[1],aa.Prod_P4[2]};
139  double vtx_int[3] = {aa.Vtx[0],aa.Vtx[1],aa.Vtx[2]};
140  double pip_data = -1.0;
141  double pim_data = -1.0;
142  int bin_pi = -1;
143  InteractionData aux_aa(aa.gen, inc_mom,2212,prod_mom,211,aa.Vol,aa.Proc,vtx_int);
144  if(aux_aa.xF<=0.5){
145  bin_pi = Thinbins->BinID_pC_pi(aux_aa.xF,aux_aa.Pt,aux_aa.Prod_pdg);
146  if(bin_pi<0)return 1.0;
147  pip_data = tt_pCPionRew->vbin_data_pip[bin_pi];
148  pim_data = tt_pCPionRew->vbin_data_pim[bin_pi];
149  }
150  else if(aux_aa.xF>0.5){
151  bin_pi = Thinbins->barton_BinID_pC_pi(aux_aa.xF,aux_aa.Pt,aux_aa.Prod_pdg);
152  if(bin_pi<0)return 1.0;
153  pip_data = tt_pCPionRew->bart_vbin_data_pip[bin_pi];
154  pim_data = tt_pCPionRew->bart_vbin_data_pim[bin_pi];
155  }
156 
157  if(aa.Prod_pdg == 321){
158  dataval = mipp_vbin_data_kap_pip[mipp_bin]*pip_data;
159  }
160  else if(aa.Prod_pdg ==-321){
161  dataval = mipp_vbin_data_kam_pim[mipp_bin]*pim_data;
162  }
163  else if(aa.Prod_pdg == 130 || aa.Prod_pdg == 310){
164  dataval = 0.25*(mipp_vbin_data_kap_pip[mipp_bin]*pip_data +3.0*(mipp_vbin_data_kam_pim[mipp_bin]*pim_data));
165  }
166  }
167 
168  //checking if this is the first interaction:
169  if(aa.gen==0)dataval /= data_prod_xs;
170  else if(aa.gen>0)dataval /= 1.0;
171  else{
172  std::cout<<"Something is wrong with gen "<<std::endl;
173  return wgt;
174  }
175 
176  double data_scale = -1.0;
177  if(aa.Prod_pdg == 321 || aa.Prod_pdg == -321){
178  data_scale = calculateDataScale(aa.Inc_pdg,aa.Inc_P,aa.Prod_pdg,aa.xF,aa.Pt);
179  }
180  else if(aa.Prod_pdg == 130 || aa.Prod_pdg == 310){
181  data_scale = 0.25*(calculateDataScale(aa.Inc_pdg,aa.Inc_P,321,aa.xF,aa.Pt)+3.*calculateDataScale(aa.Inc_pdg,aa.Inc_P,-321,aa.xF,aa.Pt));
182 
183  }
184  dataval *= data_scale;
185 
186  ThinTargetMC* mc = ThinTargetMC::getInstance();
187  double mc_cv = mc->getMCval_pC_X(aa.Inc_P,aa.xF,aa.Pt,aa.Prod_pdg);
188  double mc_prod = mc->getMCxs_pC_piK(aa.gen,aa.Inc_P);
189 
190  mc_cv /= mc_prod;
191 
192  if(mc_cv<1.e-12)return wgt;
193  wgt = dataval/mc_cv;
194 
195  if(wgt<0){
196  //std::cout<<"TTPCK check wgt(<0) "<<iUniv<<" "<<aa.Inc_P<<" "<<aa.xF<<" "<<aa.Pt<<" "<<aa.Prod_pdg<<std::endl;
197  return 1.0;
198  }
199 
200  return wgt;
201 
202  }
static MakeReweight * getInstance()
const double e
static ThinTargetBins * getInstance()
QTextStream & bin(QTextStream &s)
double calculateDataScale(int inc_pdg, double inc_mom, int prod_pdg, double xf, double pt)
QTextStream & endl(QTextStream &s)
static ThinTargetMC * getInstance()
bool NeutrinoFluxReweight::ThinTargetpCKaonReweighter::canReweight ( const InteractionData aa)
virtual

can the particular instance of this class reweight this interaction?

Implements NeutrinoFluxReweight::IInteractionReweighting.

Definition at line 64 of file ThinTargetpCKaonReweighter.cpp.

64  {
65  //checking:
66  std::string mode(getenv("MODE"));
67  if(aa.Inc_pdg != 2212)return false;
68  if(aa.Inc_P < 12.0)return false;
69  //volume check:
70  bool is_wrong_volume = aa.Vol != "TGT1" && aa.Vol != "BudalMonitor" && aa.Vol != "Budal_HFVS" && aa.Vol != "Budal_VFHS";
71  if( (mode=="REF") || (mode=="OPT") ){
72  is_wrong_volume = aa.Vol != "TargetFinHorizontal" && aa.Vol != "TargetNoSplitSegment" && aa.Vol!="tCoreLog";
73  }
74  if(is_wrong_volume)return false;
75  //
76  if(aa.Prod_pdg != 321 && aa.Prod_pdg != -321 && aa.Prod_pdg != 310 && aa.Prod_pdg != 130)return false;
77 
78  //Looking for low pz kaon:
79  ThinTargetBins* Thinbins = ThinTargetBins::getInstance();
80  int bin = Thinbins->BinID_pC_k(aa.xF,aa.Pt,aa.Prod_pdg);
81  if(bin>=0)return true; //found NA49 kaon
82 
83  //Looking for high pz kaon:
84  int mipp_bin = Thinbins->mipp_BinID_pC_k(aa.Pz,aa.Pt,aa.Prod_pdg);
85  if(mipp_bin<0)return false; //no mipp thin target kaon
86 
87  double inc_mom[3] = {aa.Inc_P4[0], aa.Inc_P4[1], aa.Inc_P4[2]};
88  double prod_mom[3] = {aa.Prod_P4[0],aa.Prod_P4[1],aa.Prod_P4[2]};
89  double vtx_int[3] = {aa.Vtx[0],aa.Vtx[1],aa.Vtx[2]};
90 
91  MakeReweight* makerew = MakeReweight::getInstance();
92  if(iUniv==-1)tt_pCPionRew = (makerew->cv_rw)->THINTARGET_PC_PION_Universe;
93  else tt_pCPionRew = (makerew->vec_rws[iUniv])->THINTARGET_PC_PION_Universe;
94 
95  InteractionData aux_intData(aa.gen, inc_mom,2212,prod_mom,211,aa.Vol,aa.Proc,vtx_int);
96  bool there_is_pip = tt_pCPionRew->canReweight(aux_intData);
97 
98  InteractionData aux_intData2(aa.gen, inc_mom,2212,prod_mom,-211,aa.Vol,aa.Proc,vtx_int);
99  bool there_is_pim = tt_pCPionRew->canReweight(aux_intData2);
100 
101  if(aa.Prod_pdg == 321)return there_is_pip;
102  else if(aa.Prod_pdg ==-321)return there_is_pim;
103  else if(aa.Prod_pdg == 130 || aa.Prod_pdg == 310)return there_is_pim && there_is_pip;
104  else{
105  return false;
106  }
107  }
std::string string
Definition: nybbler.cc:12
static MakeReweight * getInstance()
std::string getenv(std::string const &name)
Definition: getenv.cc:15
virtual bool canReweight(const InteractionData &aa)
can the particular instance of this class reweight this interaction?
static ThinTargetBins * getInstance()
QTextStream & bin(QTextStream &s)

Member Data Documentation

const ParameterTable& NeutrinoFluxReweight::ThinTargetpCKaonReweighter::cvPars
private

Definition at line 30 of file ThinTargetpCKaonReweighter.h.

double NeutrinoFluxReweight::ThinTargetpCKaonReweighter::data_prod_xs

Definition at line 23 of file ThinTargetpCKaonReweighter.h.

int NeutrinoFluxReweight::ThinTargetpCKaonReweighter::iUniv
private

Definition at line 29 of file ThinTargetpCKaonReweighter.h.

std::vector<float> NeutrinoFluxReweight::ThinTargetpCKaonReweighter::mipp_vbin_data_kam_pim

Definition at line 25 of file ThinTargetpCKaonReweighter.h.

std::vector<float> NeutrinoFluxReweight::ThinTargetpCKaonReweighter::mipp_vbin_data_kap_pip

Definition at line 25 of file ThinTargetpCKaonReweighter.h.

ThinTargetpCPionReweighter* NeutrinoFluxReweight::ThinTargetpCKaonReweighter::tt_pCPionRew
private

Definition at line 32 of file ThinTargetpCKaonReweighter.h.

const ParameterTable& NeutrinoFluxReweight::ThinTargetpCKaonReweighter::univPars
private

Definition at line 31 of file ThinTargetpCKaonReweighter.h.

std::vector<float> NeutrinoFluxReweight::ThinTargetpCKaonReweighter::vbin_data_kam

Definition at line 24 of file ThinTargetpCKaonReweighter.h.

std::vector<float> NeutrinoFluxReweight::ThinTargetpCKaonReweighter::vbin_data_kap

Definition at line 24 of file ThinTargetpCKaonReweighter.h.


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