65int main(
int argc,
char ** argv )
70 std::string isample =
"";
78 if ( isample.empty() )
81 std::cout <<
" Missing GENIE event file on input " << std::endl;
90 genie::rew::GSyst_t param_to_tweak = genie::rew::kXSecTwkDial_MaCCQE ;
94 genie::rew::GReWeight rw;
99 rw.AdoptWghtCalc(
"xsec_ccqe",
new genie::rew::GReWeightNuXSecCCQE() );
104 genie::rew::GSystSet& syst = rw.Systematics();
106 syst.Init( param_to_tweak );
112 genie::rew::GReWeightNuXSecCCQE* rwccqe =
dynamic_cast<genie::rew::GReWeightNuXSecCCQE*
>( rw.WghtCalc(
"xsec_ccqe") );
113 rwccqe->SetMode( genie::rew::GReWeightNuXSecCCQE::kModeMa );
117 TFile* file =
new TFile( isample.c_str(),
"UPDATE" );
123 std::cout <<
" Can NOT open input GENIE file " << isample << std::endl;
131 rwtree =
dynamic_cast<TTree*
>( file->Get(
"reweighting") );
134 rwtree =
new TTree(
"reweighting",
"GENIE weights tree" );
135 TTree::SetBranchStyle(1);
136 rwtree->SetAutoSave( 200000000 );
146 genie::rew::GReWeightIORecord* rwrec = 0;
147 std::string param_name = genie::rew::GSyst::AsString( param_to_tweak );
148 TBranch* rwbr = rwtree->Branch( param_name.c_str(),
149 "genie::rew::GReWeightIORecord",
152 rwbr->SetAutoDelete(kFALSE);
165 genie::rew::GSystUncertainty* syser = genie::rew::GSystUncertainty::Instance();
166 double sigpls = syser->OneSigmaErr( param_to_tweak, 1 );
167 double sigmin = syser->OneSigmaErr( param_to_tweak, -1 );
169 rwtree->GetUserInfo()->Add(
new genie::rew::GReWeightIOBranchDesc( param_name, 0.99, sigpls, sigmin ) );
173 TTree* evtree =
dynamic_cast<TTree*
>( file->Get(
"gtree") );
178 evtree->SetBranchAddress(
"gmcrec", &mcrec );
182 evtree->AddFriend( rwtree );
186 int nevt_total = evtree->GetEntries();
191 for (
int iev=0; iev<nevt_total; ++iev )
194 evtree->GetEntry(iev);
210 rwrec =
new genie::rew::GReWeightIORecord();
211 rwrec->SetOriginalEvtNumber(iev);
218 rwrec =
new genie::rew::GReWeightIORecord();
220 rwrec->SetOriginalEvtNumber( iev );
223 syst.Set( param_to_tweak, twk );
225 wt = rw.CalcWeight(evt);
226 rwrec->Insert( twk, wt );
229 syst.Set( param_to_tweak, twk );
231 wt = rw.CalcWeight(evt);
232 rwrec->Insert( twk, wt );
239 if ( rwrec )
delete rwrec;
245 rwtree->Write(
"",TObject::kOverwrite);
260 TFile* tfile =
new TFile( isample.c_str(),
"READ" );
265 TTree* rwtree_test =
dynamic_cast<TTree*
>( tfile->Get(
"reweighting") );
267 TList* hdr = rwtree_test->GetUserInfo();
270 int nentries = rwtree_test->GetEntries();
272 std::cout <<
" num of entries of re-test: " << nentries << std::endl;
273 genie::rew::GReWeightIORecord* rwrec_test = 0;
274 rwtree_test->SetBranchAddress( param_name.c_str(), &rwrec_test );
275 for (
int i=0; i<nentries; ++i )
277 rwtree_test->GetEntry(i);
279 int nres = rwrec_test->GetNumOfRWResults();
280 if ( nres <= 0 )
continue;
281 for (
int ir=0; ir<nres; ++ir )
283 double twk_test = rwrec_test->GetTweak( ir );
284 double wt_test = rwrec_test->GetWeight( ir );
285 std::cout <<
" twk = " << twk_test <<
" wt = " << wt_test << std::endl;
290 std::cout <<
" num of re-weighted results of re-test: " << nrw << std::endl;
294 int nhdr = hdr->GetEntries();
295 for (
int i=0; i<nhdr; ++i )
297 genie::rew::GReWeightIOBranchDesc* brdesc =
dynamic_cast<genie::rew::GReWeightIOBranchDesc*
>( hdr->At(i) );
298 std::cout <<
" branch name: " << brdesc->GetParameterName() << std::endl;
299 std::cout <<
" parameter: " << brdesc->GetParameterMean() <<
" "
300 << brdesc->GetParameterSigmaPlus() <<
" " << brdesc->GetParameterSigmaMinus() << std::endl;