00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include <TCanvas.h>
00011 #include <TPostScript.h>
00012 #include <TLegend.h>
00013 #include <TH2D.h>
00014 #include <TGraphErrors.h>
00015 #include <TGraphAsymmErrors.h>
00016 #include <TStyle.h>
00017 #include <TLatex.h>
00018 #include <TMarker.h>
00019 #include <TArrow.h>
00020 #include <TGaxis.h>
00021 #include <TF1.h>
00022 #include <TLine.h>
00023
00024 #include <iostream>
00025 #include <fstream>
00026
00027 #include "BAT/BCModel.h"
00028 #include "BAT/BCSummaryPriorModel.h"
00029 #include "BAT/BCH1D.h"
00030 #include "BAT/BCH2D.h"
00031 #include "BAT/BCLog.h"
00032 #include "BAT/BCMath.h"
00033
00034 #include "BAT/BCSummaryTool.h"
00035
00036 unsigned int BCSummaryTool::fHCounter=0;
00037
00038
00039 BCSummaryTool::BCSummaryTool()
00040 : fModel(0)
00041 , fPriorModel(0)
00042 , fFlagInfoMarg(false)
00043 , fFlagInfoOpt(false)
00044 {
00045
00046 fSumProb.push_back(0.05);
00047 fSumProb.push_back(0.10);
00048 fSumProb.push_back(0.1587);
00049 fSumProb.push_back(0.50);
00050 fSumProb.push_back(0.8413);
00051 fSumProb.push_back(0.90);
00052 fSumProb.push_back(0.95);
00053
00054
00055 gStyle->SetPaintTextFormat(".2g");
00056 }
00057
00058
00059 BCSummaryTool::BCSummaryTool(BCModel * model)
00060 : fModel(model)
00061 , fPriorModel(0)
00062 , fFlagInfoMarg(false)
00063 , fFlagInfoOpt(false)
00064 {
00065
00066 fSumProb.push_back(0.05);
00067 fSumProb.push_back(0.10);
00068 fSumProb.push_back(0.1587);
00069 fSumProb.push_back(0.50);
00070 fSumProb.push_back(0.8413);
00071 fSumProb.push_back(0.90);
00072 fSumProb.push_back(0.95);
00073
00074
00075 gStyle->SetPaintTextFormat(".2g");
00076 }
00077
00078
00079 BCSummaryTool::~BCSummaryTool()
00080 {
00081 delete fPriorModel;
00082 }
00083
00084
00085 int BCSummaryTool::CopySummaryData()
00086 {
00087
00088 if (!fModel)
00089 return 0;
00090
00091
00092 fParName.clear();
00093 fParMin.clear();
00094 fParMax.clear();
00095 fMean.clear();
00096 fMargMode.clear();
00097 fGlobalMode.clear();
00098 fQuantiles.clear();
00099 fSmallInt.clear();
00100 fRMS.clear();
00101 fCorrCoeff.clear();
00102
00103
00104 int npar = fModel->GetNParameters();
00105 int nquantiles = int( fSumProb.size() );
00106
00107
00108 for (int i = 0; i < npar; ++i) {
00109
00110
00111 fParName.push_back( (fModel->GetParameter(i)->GetName()) );
00112 fParMin.push_back( fModel->GetParameter(i)->GetLowerLimit() );
00113 fParMax.push_back( fModel->GetParameter(i)->GetUpperLimit() );
00114
00115
00116 if (fModel->MCMCGetFlagRun()) {
00117 fFlagInfoMarg = true;
00118 BCH1D * bch1d_temp = fModel->GetMarginalized( fModel->GetParameter(i) );
00119 if (bch1d_temp) {
00120 fMean.push_back( bch1d_temp->GetMean() );
00121 fRMS.push_back( bch1d_temp->GetRMS() );
00122 fMargMode.push_back( bch1d_temp->GetMode() );
00123 for (int j = 0; j < nquantiles; ++j)
00124 fQuantiles.push_back( bch1d_temp->GetQuantile( fSumProb.at(j) ) );
00125 std::vector <double> intervals = bch1d_temp->GetSmallestIntervals();
00126 int nintervals = int(intervals.size() / 5);
00127 fSmallInt.push_back(nintervals);
00128 fSmallInt.insert( fSmallInt.end(), intervals.begin(), intervals.end() );
00129 }
00130 else {
00131 double tmpval = fModel->GetParameter(i)->GetUpperLimit() - fModel->GetParameter(i)->GetLowerLimit();
00132 tmpval = fModel->GetParameter(i)->GetLowerLimit() - 2. * tmpval;
00133 fMean.push_back( tmpval );
00134 fRMS.push_back( tmpval );
00135 fMargMode.push_back( tmpval );
00136 for (int j = 0; j < nquantiles; ++j)
00137 fQuantiles.push_back( tmpval );
00138 fSmallInt.push_back( 0 );
00139
00140 }
00141
00142
00143 for (int j = 0; j < npar; ++j) {
00144 if (i!=j) {
00145 BCH2D * bch2d_temp = fModel->GetMarginalized(fModel->GetParameter(i),fModel->GetParameter(j));
00146 if ( bch2d_temp )
00147 fCorrCoeff.push_back( bch2d_temp->GetHistogram()->GetCorrelationFactor() );
00148 else
00149 fCorrCoeff.push_back( 0. );
00150 }
00151 else
00152 fCorrCoeff.push_back(1.);
00153 }
00154 }
00155 else {
00156
00157 }
00158
00159
00160 if ((fModel->GetBestFitParameters()).size() > 0) {
00161 fFlagInfoOpt = true;
00162 fGlobalMode.push_back ( (fModel->GetBestFitParameters()).at(i) );
00163 }
00164 else {
00165
00166 }
00167 }
00168
00169
00170 return 1;
00171 };
00172
00173
00174 int BCSummaryTool::PrintParameterPlot(const char * filename)
00175 {
00176
00177 if (!CopySummaryData())
00178 return 0;
00179
00180
00181 int npar = fModel->GetNParameters();
00182 int nquantiles = int( fSumProb.size() );
00183
00184
00185 TH1D * hist_axes = new TH1D(
00186 TString::Format("hist_axes_par_%d",getNextIndex()),
00187 ";;Scaled parameter range [a.u.]",npar, -0.5, npar-0.5);
00188 hist_axes->SetStats(kFALSE);
00189 for (int i = 0; i < npar; ++i)
00190 hist_axes->GetXaxis()->SetBinLabel( i+1, fParName.at(i).c_str() );
00191
00192 hist_axes->GetXaxis()->SetLabelSize(0.06);
00193 hist_axes->GetXaxis()->SetTickLength(0.0);
00194 hist_axes->GetYaxis()->SetRangeUser(-0.1, 1.1);
00195 hist_axes->GetYaxis()->SetTickLength(0.0);
00196
00197
00198 TGraphErrors * graph_quantiles = new TGraphErrors(npar*nquantiles);
00199 graph_quantiles->SetMarkerSize(0);
00200 graph_quantiles->SetLineColor(38);
00201 graph_quantiles->SetLineStyle(2);
00202
00203 TGraphErrors * graph_mean = new TGraphErrors(npar);
00204 graph_mean->SetMarkerColor(kBlack);
00205 graph_mean->SetMarkerStyle(20);
00206
00207 TGraphErrors * graph_mode = new TGraphErrors(npar);
00208 graph_mode->SetMarkerColor(kRed);
00209 graph_mode->SetMarkerStyle(20);
00210
00211 TGraphAsymmErrors * graph_intervals = new TGraphAsymmErrors(0);
00212 graph_intervals->SetFillColor(kYellow);
00213 graph_intervals->SetLineStyle(2);
00214 graph_intervals->SetLineColor(kRed);
00215 graph_intervals->SetMarkerSize(0);
00216
00217
00218 int indexintervals = 0;
00219
00220
00221 if (fFlagInfoMarg) {
00222 for (int i = 0; i < npar; ++i) {
00223 for (int j = 0; j < nquantiles; ++j) {
00224 graph_quantiles->SetPoint(i*nquantiles+j,double(i),
00225 (fQuantiles.at(i*nquantiles+j) - fParMin.at(i))/(fParMax.at(i)-fParMin.at(i)));
00226 graph_quantiles->SetPointError(i*nquantiles+j, 0.5, 0.0);
00227 }
00228 }
00229 }
00230
00231
00232 if (fFlagInfoMarg) {
00233 for (int i = 0; i < npar; ++i) {
00234
00235 graph_mean->SetPoint(i, double(i), (fMean.at(i) - fParMin.at(i))/(fParMax.at(i)-fParMin.at(i)));
00236 graph_mean->SetPointError(i, 0.0, fRMS.at(i)/(fParMax.at(i)-fParMin.at(i)));
00237 }
00238 }
00239
00240
00241 if (fFlagInfoOpt) {
00242 for (int i = 0; i < npar; ++i)
00243 graph_mode->SetPoint(i, double(i), (fGlobalMode.at(i) - fParMin.at(i))/(fParMax.at(i)-fParMin.at(i)));
00244 }
00245
00246
00247 if (fFlagInfoMarg) {
00248 for (int i = 0; i < npar; ++i) {
00249 int nintervals = int(fSmallInt.at(indexintervals++));
00250 for (int j = 0; j < nintervals; ++j) {
00251 double xmin = fSmallInt.at(indexintervals++);
00252 double xmax = fSmallInt.at(indexintervals++);
00253 indexintervals++;
00254 double xlocalmaxpos = fSmallInt.at(indexintervals++);
00255 indexintervals++;
00256 int npoints = graph_intervals->GetN();
00257 graph_intervals->SetPoint(npoints,double(i),
00258 (xlocalmaxpos - fParMin.at(i))/(fParMax.at(i)-fParMin.at(i)));
00259 graph_intervals->SetPointError(npoints,0.5, 0.5,
00260 (xlocalmaxpos - xmin)/(fParMax.at(i)-fParMin.at(i)),
00261 (xmax - xlocalmaxpos)/(fParMax.at(i)-fParMin.at(i)));
00262 }
00263 }
00264 }
00265
00266
00267 TLegend * legend = new TLegend(0.15, 0.88, 0.85, 0.99);
00268 legend->SetBorderSize(0);
00269 legend->SetFillColor(0);
00270
00271
00272 TLatex * latex = new TLatex();
00273 latex->SetTextSize(0.02);
00274
00275
00276 TLine * line_top = new TLine(-0.5, 1.0, npar-0.5, 1.0);
00277 line_top->SetLineColor(kBlack);
00278 line_top->SetLineStyle(1);
00279 line_top->SetLineWidth(2);
00280
00281 TLine * line_bot = new TLine(-0.5, 0.0, npar-0.5, 0.0);
00282 line_bot->SetLineColor(kBlack);
00283 line_bot->SetLineStyle(1);
00284 line_bot->SetLineWidth(2);
00285
00286
00287 TCanvas * c_par = new TCanvas(TString::Format("c_par_%d",getNextIndex()));
00288 c_par->cd();
00289 hist_axes->Draw();
00290 line_top->Draw();
00291 line_bot->Draw();
00292 if (fFlagInfoMarg) {
00293 graph_intervals->DrawClone("SAME2");
00294 for (int i = 0; i < graph_intervals->GetN(); ++i)
00295 graph_intervals->SetPointError(i, 0.5, 0.5, 0.0, 0.0);
00296 graph_intervals->Draw("SAMEPZ");
00297 graph_quantiles->Draw("SAMEPZ");
00298 graph_mean->Draw("SAMEP");
00299 legend->AddEntry(graph_quantiles, "Quantiles (5%, 10%, 16%, 50%, 84%, 90%, 95%)", "L");
00300 legend->AddEntry(graph_mean, "Mean and RMS", "LEP");
00301 legend->AddEntry(graph_intervals, "Smallest 68% intervals and local modes", "FL");
00302 }
00303 if (fFlagInfoOpt) {
00304 graph_mode->Draw("SAMEP");
00305 legend->AddEntry(graph_mode, "Global mode", "P");
00306 }
00307 for (int i = 0; i < npar;++i) {
00308
00309
00310 latex->DrawLatex(double(i)-0.1, 0.010-0.07, Form("%+3.3f", fParMin.at(i)));
00311 latex->DrawLatex(double(i)-0.1, 0.965+0.07, Form("%+3.3f", fParMax.at(i)));
00312 }
00313 latex->SetNDC();
00314 latex->DrawLatex(0.9, 0.175, "Par. min.");
00315 latex->DrawLatex(0.9, 0.83, "Par. max.");
00316 legend->Draw("SAME");
00317 gPad->RedrawAxis();
00318 c_par->Print(filename);
00319
00320
00321 return 1;
00322 }
00323
00324
00325 int BCSummaryTool::PrintCorrelationMatrix(const char * filename)
00326 {
00327
00328 if (!CopySummaryData())
00329 return 0;
00330
00331
00332 if (!fFlagInfoMarg)
00333 return 0;
00334
00335
00336 int npar = fModel->GetNParameters();
00337
00338
00339 TH2D * hist_corr = new TH2D(
00340 TString::Format("hist_corr_%d",getNextIndex()),
00341 ";;",npar, -0.5, npar-0.5,npar, -0.5, npar-0.5);
00342 hist_corr->SetStats(kFALSE);
00343 hist_corr->GetXaxis()->SetTickLength(0.0);
00344
00345 hist_corr->GetYaxis()->SetTickLength(0.0);
00346
00347 hist_corr->GetZaxis()->SetRangeUser(-1.0, 1.0);
00348
00349 for (int i = 0; i < npar; ++i) {
00350 hist_corr->GetXaxis()->SetLabelSize(0.06);
00351 hist_corr->GetYaxis()->SetLabelSize(0.06);
00352 if (npar < 5) {
00353 hist_corr->GetXaxis()->SetBinLabel( i+1, fParName.at(i).c_str() );
00354 hist_corr->GetYaxis()->SetBinLabel( npar-i, fParName.at(i).c_str() );
00355 }
00356 else {
00357 hist_corr->GetXaxis()->SetBinLabel( i+1, TString::Format("%d",i) );
00358 hist_corr->GetYaxis()->SetBinLabel( npar-i, TString::Format("%d",i) );
00359 }
00360 }
00361
00362
00363 for (int i = 0; i < npar; ++i)
00364 for (int j = 0; j < npar; ++j) {
00365 int index = i * npar + j;
00366 double corr = fCorrCoeff.at(index);
00367 hist_corr->SetBinContent(i+1, npar-j, corr);
00368 }
00369
00370
00371 TCanvas * c_corr = new TCanvas(TString::Format("c_corr_matrix_%d",getNextIndex()));
00372 c_corr->cd();
00373 hist_corr->Draw("colz text");
00374
00375 TF1 * f = new TF1("fUp","x",-0.5,npar-0.5);
00376 TGaxis * A1 = new TGaxis(-0.5,npar-0.5,npar-0.5,npar-0.5,"fUp",100,"-");
00377 A1->ImportAxisAttributes(hist_corr->GetXaxis());
00378 A1->Draw();
00379
00380
00381
00382 hist_corr->GetXaxis()->SetLabelSize(0.);
00383 hist_corr->Draw("colz text same");
00384
00385 for (int i = 0; i < npar; ++i)
00386 for (int j = 0; j < npar; ++j) {
00387 BCH2D * bch2d_temp = fModel->GetMarginalized(fModel->GetParameter(i),fModel->GetParameter(j));
00388 if ( bch2d_temp || i==j )
00389 continue;
00390
00391 TBox * bempty = new TBox(
00392 hist_corr->GetXaxis()->GetBinLowEdge(i+1),
00393 hist_corr->GetYaxis()->GetBinLowEdge(npar-j),
00394 hist_corr->GetXaxis()->GetBinLowEdge(i+2),
00395 hist_corr->GetYaxis()->GetBinLowEdge(npar-j+1)
00396 );
00397 bempty->SetLineStyle(0);
00398 bempty->SetLineWidth(0);
00399 bempty->SetFillColor(kWhite);
00400 bempty->Draw();
00401 }
00402
00403
00404 TLine * lA1 = new TLine(-0.5,npar-0.5,npar-0.5,npar-0.5);
00405 lA1->Draw("same");
00406 TLine * lA2 = new TLine(npar-0.5,npar-0.5,npar-0.5,-0.5);
00407 lA2->Draw("same");
00408
00409 gPad->RedrawAxis();
00410 c_corr->Print(filename);
00411
00412 delete f;
00413 delete A1;
00414 delete lA1;
00415 delete lA2;
00416 delete hist_corr;
00417 delete c_corr;
00418
00419
00420 return 1;
00421 }
00422
00423
00424 int BCSummaryTool::PrintCorrelationPlot(const char * filename)
00425 {
00426
00427 if (!fModel->MCMCGetFlagRun()) {
00428 BCLog::OutError("BCSummaryTool::PrintCorrelationPlot : MCMC was not run. Marginalized distributions not available.");
00429 return 0;
00430 }
00431
00432
00433 int npar = fModel->GetNParameters();
00434
00435 TCanvas * c = new TCanvas(TString::Format("c_corr_%d",getNextIndex()));
00436 c->cd();
00437
00438 double margin = .1;
00439 double padsize = (1. - 2.*margin) / (double)npar;
00440
00441
00442 TPad ** pad = new TPad*[npar*npar];
00443
00444
00445 double xlow, xup, ylow, yup;
00446 double marginleft, marginright, margintop, marginbottom;
00447 marginleft=marginright=margintop=marginbottom=0.01;
00448
00449
00450 for (int i=npar-1;i>=0;--i) {
00451 xlow = (double)i*padsize + margin;
00452 xup = xlow+padsize;
00453
00454 for (int j=i;j<=npar-1;j++) {
00455 yup = 1. - (double)j*padsize - margin;
00456 ylow = yup-padsize;
00457
00458
00459 int ipad=i*npar+j;
00460 pad[ipad]=new TPad(TString::Format("pad_%d_%d",ipad,getNextIndex()), "", xlow, ylow, xup, yup);
00461 pad[ipad]->SetTopMargin(margintop);
00462 pad[ipad]->SetBottomMargin(marginbottom);
00463 pad[ipad]->SetLeftMargin(marginleft);
00464 pad[ipad]->SetRightMargin(marginright);
00465 pad[ipad]->SetFillColor(kWhite);
00466
00467 pad[ipad]->Draw();
00468 pad[ipad]->cd();
00469
00470
00471 TH1 * hh = 0;
00472
00473 if(i==j) {
00474 BCH1D * bh1 = fModel->GetMarginalized(fModel->GetParameter(i));
00475 if (bh1)
00476 hh = (TH1D*)bh1->GetHistogram()->Clone();
00477 }
00478 else {
00479 BCH2D * bh2 = fModel->GetMarginalized(fModel->GetParameter(i),fModel->GetParameter(j));
00480 if (bh2)
00481 hh = (TH2D*)bh2->GetHistogram()->Clone();
00482 }
00483
00484
00485 if (!hh) {
00486 pad[ipad]->SetFillColor(kGray);
00487 TBox * box = new TBox(marginleft,marginbottom,1.-marginright,1.-margintop);
00488 box->SetLineWidth(1);
00489 box->SetLineColor(kGray+1);
00490 box->SetFillColor(kWhite);
00491 box->Draw();
00492
00493 TText * text_na = new TText(.5,.5,"N/A");
00494 text_na->SetTextFont(42);
00495 text_na->SetTextAlign(22);
00496 text_na->SetTextSize(.8/(double)npar);
00497 text_na->SetTextColor(kGray+1);
00498 text_na->Draw();
00499 }
00500
00501 else {
00502 if(i==j) {
00503 hh->SetFillStyle(1001);
00504 hh->SetFillColor(kYellow);
00505 hh->Draw("hist");
00506 }
00507 else {
00508 hh->SetContour(20);
00509 hh->Draw("col");
00510 }
00511
00512 hh->GetXaxis()->SetLabelOffset(999);
00513 hh->GetYaxis()->SetLabelOffset(999);
00514
00515 hh->GetXaxis()->SetTitleSize(0.);
00516 hh->GetYaxis()->SetTitleSize(0.);
00517 }
00518
00519 c->cd();
00520
00521
00522 double labelsize = .8/(double)npar/5.;
00523 double xtext, ytext;
00524
00525
00526 if(i==0) {
00527 TText * label = new TText;
00528 label->SetTextFont(62);
00529 label->SetTextSize(labelsize);
00530 label->SetTextAlign(22);
00531 label->SetNDC();
00532
00533 label->SetTextAngle(90);
00534
00535 xtext = margin * (1. - 8. * labelsize);
00536 ytext = yup - padsize / 2.;
00537
00538 label->DrawText(xtext,ytext,fModel->GetParameter(j)->GetName().c_str());
00539 }
00540
00541
00542 if(j==npar-1) {
00543 TText * label = new TText;
00544 label->SetTextFont(62);
00545 label->SetTextSize(labelsize);
00546 label->SetTextAlign(22);
00547 label->SetNDC();
00548
00549 xtext = xlow + padsize / 2.;
00550 ytext = margin * (1. - 6. * labelsize);
00551
00552 label->DrawText(xtext,ytext,fModel->GetParameter(i)->GetName().c_str());
00553 }
00554 }
00555 }
00556
00557 gPad->RedrawAxis();
00558 c->Print(filename);
00559
00560 return 1;
00561 }
00562
00563
00564 int BCSummaryTool::PrintKnowledgeUpdatePlots(const char * filename)
00565 {
00566
00567 CalculatePriorModel();
00568
00569
00570 TPostScript * ps = new TPostScript(filename);
00571
00572
00573 TCanvas * c_update = new TCanvas(TString::Format("c_update_%d",getNextIndex()));
00574
00575 c_update->Update();
00576 ps->NewPage();
00577 c_update->cd();
00578
00579
00580 TLegend * legend1d = new TLegend(0.15, 0.88, 0.85, 0.94);
00581 legend1d->SetBorderSize(0);
00582 legend1d->SetFillColor(0);
00583
00584
00585 int npar = fModel->GetNParameters();
00586 for (int i = 0; i < npar; ++i) {
00587
00588 c_update->Update();
00589 ps->NewPage();
00590 c_update->cd();
00591
00592
00593 BCParameter * par = fModel->GetParameter(i);
00594 TH1D * hist_prior = fPriorModel->GetMarginalized(par)->GetHistogram();
00595 hist_prior->SetLineColor(kRed);
00596 TH1D * hist_posterior = fModel->GetMarginalized(par)->GetHistogram();
00597
00598
00599 if (!i) {
00600 legend1d->AddEntry(hist_prior, "Prior probability", "L");
00601 legend1d->AddEntry(hist_posterior, "Posterior probability", "L");
00602 }
00603
00604
00605 hist_prior->Scale(hist_posterior->Integral()/hist_prior->Integral());
00606
00607
00608 double max_prior = hist_prior->GetMaximum();
00609 double max_posterior = hist_posterior->GetMaximum();
00610 double max = 1.1 * TMath::Max(max_prior, max_posterior);
00611
00612
00613 c_update->cd();
00614 hist_prior->GetXaxis()->SetNdivisions(508);
00615 hist_posterior->GetXaxis()->SetNdivisions(508);
00616 hist_prior->Draw();
00617 hist_posterior->Draw("SAME");
00618 legend1d->Draw("SAME");
00619
00620
00621 hist_prior->GetYaxis()->SetRangeUser(0.0, max);
00622 hist_posterior->GetYaxis()->SetRangeUser(0.0, max);
00623 }
00624
00625
00626 TLegend * legend2d = new TLegend(0.15, 0.88, 0.85, 0.99);
00627 legend2d->SetBorderSize(0);
00628 legend2d->SetFillColor(0);
00629
00630
00631 TMarker * marker_prior = new TMarker();
00632 marker_prior->SetMarkerStyle(24);
00633 marker_prior->SetMarkerColor(kRed);
00634
00635 TMarker * marker_posterior = new TMarker();
00636 marker_posterior->SetMarkerStyle(24);
00637 marker_posterior->SetMarkerColor(kBlack);
00638
00639 TArrow * arrow = new TArrow();
00640 arrow->SetArrowSize(0.02);
00641 arrow->SetLineColor(kBlue);
00642 arrow->SetLineStyle(2);
00643
00644
00645 for (int i = 0; i < npar; ++i) {
00646 for (int j = 0; j < i; ++j) {
00647
00648 c_update->Update();
00649 ps->NewPage();
00650 c_update->cd();
00651
00652
00653 BCParameter * par1 = fModel->GetParameter(i);
00654 BCParameter * par2 = fModel->GetParameter(j);
00655
00656
00657 BCH2D * bch2d_2dprior = fPriorModel->GetMarginalized(par1, par2);
00658 BCH2D * bch2d_2dposterior = fModel->GetMarginalized(par1, par2);
00659
00660
00661 TH2D * hist_2dprior = bch2d_2dprior->GetHistogram();
00662 hist_2dprior->SetLineColor(kRed);
00663 TH2D * hist_2dposterior = bch2d_2dposterior->GetHistogram();
00664
00665
00666 hist_2dprior->Scale(1.0/hist_2dprior->Integral("width"));
00667 hist_2dposterior->Scale(1.0/hist_2dposterior->Integral("width"));
00668
00669
00670 bch2d_2dprior->CalculateIntegratedHistogram();
00671 bch2d_2dposterior->CalculateIntegratedHistogram();
00672
00673 double level[1] = {bch2d_2dprior->GetLevel(0.32)};
00674 hist_2dprior->SetContour(1, level);
00675 hist_2dprior->Draw("CONT3");
00676 level[0] = bch2d_2dposterior->GetLevel(0.32);
00677 hist_2dposterior->SetContour(1, level);
00678 hist_2dposterior->Draw("CONT3 SAME");
00679
00680 std::vector <double> mode_prior = fPriorModel->GetBestFitParameters();
00681 std::vector <double> mode_posterior = fModel->GetBestFitParameters();
00682
00683 marker_prior->DrawMarker(mode_prior.at(j), mode_prior.at(i));
00684 marker_posterior->DrawMarker(mode_posterior.at(j), mode_posterior.at(i));
00685 arrow->DrawArrow(mode_prior.at(j), mode_prior.at(i), mode_posterior.at(j), mode_posterior.at(i));
00686
00687 if (i==1 && j == 0) {
00688 legend2d->AddEntry(hist_2dprior, "68% prior contour", "L");
00689 legend2d->AddEntry(hist_2dposterior, "68% posterior contour", "L");
00690 legend2d->AddEntry(marker_prior, "Prior mode", "P");
00691 legend2d->AddEntry(marker_posterior, "Posterior mode", "P");
00692 legend2d->AddEntry(arrow, "Change in mode", "L");
00693 }
00694 legend2d->Draw();
00695 }
00696 }
00697
00698
00699 c_update->Update();
00700 ps->Close();
00701
00702
00703 delete legend1d;
00704 delete legend2d;
00705 delete marker_prior;
00706 delete marker_posterior;
00707 delete arrow;
00708 delete ps;
00709 delete c_update;
00710
00711
00712 return 1;
00713 }
00714
00715
00716
00717
00718
00719
00720
00721
00722
00723
00724
00725
00726
00727
00728
00729
00730
00731
00732
00733
00734
00735
00736
00737
00738
00739
00740
00741
00742
00743
00744
00745
00746
00747
00748
00749
00750
00751
00752
00753
00754
00755
00756 int BCSummaryTool::PrintParameterLatex(const char * filename)
00757 {
00758
00759 std::ofstream ofi(filename);
00760 ofi.precision(3);
00761
00762
00763 if(!ofi.is_open()) {
00764 std::cerr << "Couldn't open file " << filename <<std::endl;
00765 return 0;
00766 }
00767
00768
00769 int npar = fModel->GetNParameters();
00770
00771
00772 ofi
00773 << "\\documentclass[11pt, a4paper]{article}" << std::endl
00774 << std::endl
00775 << "\\begin{document}" << std::endl
00776 << std::endl
00777 << "\\begin{table}[ht!]" << std::endl
00778 << "\\begin{center}" << std::endl
00779 <<"\\begin{tabular}{llllllll}" << std::endl
00780 << "\\hline" << std::endl
00781 << "Parameter & Mean & RMS & Gl. mode & Mode & Median & 16\\% quant. & 84\\% quant. \\\\" << std::endl
00782 << "\\hline" << std::endl;
00783
00784 for (int i = 0; i < npar; ++i) {
00785 BCParameter * par = fModel->GetParameter(i);
00786 BCH1D * bch1d = fModel->GetMarginalized(par);
00787 ofi
00788 << par->GetName() << " & "
00789 << bch1d->GetMean() << " & "
00790 << bch1d->GetRMS() << " & "
00791 << fModel->GetBestFitParameters().at(i) << " & "
00792 << bch1d->GetMode() << " & "
00793 << bch1d->GetMedian() << " & "
00794 << bch1d->GetQuantile(0.16) << " & "
00795 << bch1d->GetQuantile(0.84) << " \\\\" << std::endl;
00796 }
00797 ofi
00798 << "\\hline" << std::endl
00799 << "\\end{tabular}" << std::endl
00800 << "\\caption{Summary of the parameter estimates.}" << std::endl
00801 << "\\end{center}" << std::endl
00802 << "\\end{table}" << std::endl
00803 << std::endl
00804 << "\\end{document}" << std::endl;
00805
00806
00807 ofi.close();
00808
00809
00810 return 1;
00811 }
00812
00813
00814 int BCSummaryTool::CalculatePriorModel()
00815 {
00816
00817 delete fPriorModel;
00818
00819 fPriorModel = new BCSummaryPriorModel();
00820
00821
00822 fPriorModel->SetModel(fModel);
00823
00824
00825 fPriorModel->MarginalizeAll();
00826
00827
00828 fPriorModel->FindMode( fPriorModel->GetBestFitParameters() );
00829
00830
00831 return 1;
00832 }
00833
00834