@@ -356,6 +356,7 @@ void HFInvMassFitter::doFit()
356356 mTotalPdf ->plotOn (mInvMassFrame , Name (" Tot_c" ), LineColor (kBlue ));
357357 mSgnPdf ->plotOn (mInvMassFrame , Normalization (1.0 , RooAbsReal::RelativeExpected), DrawOption (" F" ), FillColor (TColor::GetColorTransparent (kBlue , 0.2 )), VLines ());
358358 mChiSquareOverNdfTotal = mInvMassFrame ->chiSquare (" Tot_c" , " data_c" ); // calculate reduced chi2 / DNF
359+
359360 // plot residual distribution
360361 mResidualFrame = mass->frame (Title (" Residual Distribution" ));
361362 RooHist* residualHistogram = mInvMassFrame ->residHist (" data_c" , " Bkg_c" );
@@ -372,6 +373,7 @@ void HFInvMassFitter::doFit()
372373 calculateSignal (mRawYield , mRawYieldErr );
373374 countSignal (mRawYieldCounted , mRawYieldCountedErr );
374375 calculateSignificance (mSignificance , mSignificanceErr );
376+ calculateFitToDataRatio ();
375377 }
376378}
377379
@@ -637,6 +639,37 @@ void HFInvMassFitter::drawResidual(TVirtualPad* pad)
637639 highlightPeakRegion (mResidualFrame );
638640}
639641
642+ // draw ratio on canvas
643+ void HFInvMassFitter::drawRatio (TVirtualPad* pad)
644+ {
645+ pad->cd ();
646+ mRatioFrame ->GetYaxis ()->SetTitle (" " );
647+ TPaveText* textInfo = new TPaveText (0.12 , 0.65 , 0.47 , .89 , " NDC" );
648+ textInfo->SetBorderSize (0 );
649+ textInfo->SetFillStyle (0 );
650+ textInfo->SetTextColor (kBlue );
651+ textInfo->AddText (Form (" S = %.0f #pm %.0f " , mRawYield , mRawYieldErr ));
652+ textInfo->AddText (Form (" S_{count} = %.0f #pm %.0f " , mRawYieldCounted , mRawYieldCountedErr ));
653+ textInfo->AddText (Form (" mean = %.3f #pm %.3f" , mRooMeanSgn ->getVal (), mRooMeanSgn ->getError ()));
654+ if (mTypeOfSgnPdf == DoubleGaus) {
655+ auto const & baseSigmaSgn = mWorkspace ->var (" sigma" );
656+ textInfo->AddText (Form (" sigma = %.3f #pm %.3f" , baseSigmaSgn->getVal (), baseSigmaSgn->getError ()));
657+ textInfo->AddText (Form (" sigma 2 = %.3f #pm %.3f" , mRooSigmaSgn ->getVal (), mRooSigmaSgn ->getError ()));
658+ } else {
659+ textInfo->AddText (Form (" sigma = %.3f #pm %.3f" , mRooSigmaSgn ->getVal (), mRooSigmaSgn ->getError ()));
660+ }
661+ mRatioFrame ->addObject (textInfo);
662+ double xMin = mRatioFrame ->GetXaxis ()->GetXmin ();
663+ double xMax = mRatioFrame ->GetXaxis ()->GetXmax ();
664+ TLine* line = new TLine (xMin, 1.0 , xMax, 1.0 );
665+ line->SetLineColor (kGray );
666+ line->SetLineStyle (2 );
667+ line->SetLineWidth (2 );
668+ mRatioFrame ->addObject (line);
669+ mRatioFrame ->Draw ();
670+ highlightPeakRegion (mRatioFrame );
671+ }
672+
640673// draw peak region with vertical lines
641674void HFInvMassFitter::highlightPeakRegion (const RooPlot* plot, Color_t color, Width_t width, Style_t style) const
642675{
@@ -818,6 +851,42 @@ void HFInvMassFitter::plotRefl(RooAbsPdf* pdf)
818851 pdf->plotOn (mInvMassFrame , Components (namesOfReflPdf.at (mTypeOfReflPdf ).c_str ()), Name (" Refl_c" ), LineColor (kGreen ));
819852}
820853
854+ // Calculate fit to data ratio
855+ void HFInvMassFitter::calculateFitToDataRatio ()
856+ {
857+ if (!mInvMassFrame )
858+ return ;
859+ // Create a new RooPlot for the ratio
860+ mRatioFrame = mWorkspace ->var (" mass" )->frame (Title (" Fit/Data Ratio" ));
861+
862+ // Get the data and fit curves from the frame
863+ auto * dataHist = dynamic_cast <RooHist*>(mInvMassFrame ->findObject (" data_c" ));
864+ auto * fitCurve = dynamic_cast <RooCurve*>(mInvMassFrame ->findObject (" Tot_c" )); // or the relevant fit curve
865+
866+ if (!dataHist || !fitCurve)
867+ return ;
868+
869+ RooHist* ratioHist = new RooHist ();
870+
871+ for (int i = 0 ; i < dataHist->GetN (); ++i) {
872+ double x, dataY, dataErr;
873+ dataHist->GetPoint (i, x, dataY);
874+ dataErr = dataHist->GetErrorY (i);
875+
876+ double fitY = fitCurve->Eval (x);
877+
878+ double ratio = dataY != 0 ? fitY / dataY : 0 ;
879+ double err = dataY != 0 ? ratio * dataErr / dataY : 0 ;
880+
881+ ratioHist->SetPoint (i, x, ratio);
882+ ratioHist->SetPointError (i, 0 , 0 , err, err);
883+ }
884+
885+ mRatioFrame ->SetMinimum (0.0 );
886+ mRatioFrame ->SetMaximum (2.0 );
887+ mRatioFrame ->addPlotable (ratioHist, " P" );
888+ }
889+
821890// Fix reflection pdf
822891void HFInvMassFitter::setReflFuncFixed ()
823892{
0 commit comments