Skip to content

Commit 1e835b2

Browse files
committed
[RF] Avoid warning about duplicate command arguments
Avoid warning about duplicate command arguments by removing old arguments when overriding them. (cherry picked from commit d4b150c)
1 parent 63be8d3 commit 1e835b2

File tree

3 files changed

+25
-7
lines changed

3 files changed

+25
-7
lines changed

roofit/roofitcore/src/RooAbsPdf.cxx

Lines changed: 10 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -210,6 +210,15 @@ inline double getLog(double prob, RooAbsReal const *caller)
210210
return std::log(prob);
211211
}
212212

213+
void replaceOrAdd(RooLinkedList &lst, TObject &obj)
214+
{
215+
TObject *old = lst.FindObject(obj.GetName());
216+
if (old)
217+
lst.Replace(old, &obj);
218+
else
219+
lst.Add(&obj);
220+
}
221+
213222
} // namespace
214223

215224
using std::endl, std::string, std::ostream, std::vector, std::pair, std::make_pair;
@@ -2173,7 +2182,7 @@ RooPlot* RooAbsPdf::plotOn(RooPlot* frame, RooLinkedList& cmdList) const
21732182
// Append overriding scale factor command at end of original command list
21742183
RooCmdArg tmp = RooFit::Normalization(scaleFactor,Raw) ;
21752184
tmp.setInt(1,1) ; // Flag this normalization command as created for internal use (so that VisualizeError can strip it)
2176-
cmdList.Add(&tmp) ;
2185+
replaceOrAdd(cmdList, tmp);
21772186

21782187
// Was a component selected requested
21792188
if (haveCompSel) {

roofit/roofitcore/src/RooSimultaneous.cxx

Lines changed: 14 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -87,6 +87,15 @@ std::map<std::string, RooAbsPdf *> createPdfMap(const RooArgList &inPdfList, Roo
8787
return pdfMap;
8888
}
8989

90+
void replaceOrAdd(RooLinkedList &lst, TObject &obj)
91+
{
92+
TObject *old = lst.FindObject(obj.GetName());
93+
if (old)
94+
lst.Replace(old, &obj);
95+
else
96+
lst.Add(&obj);
97+
}
98+
9099
} // namespace
91100

92101
RooSimultaneous::InitializationOutput::~InitializationOutput() = default;
@@ -774,9 +783,9 @@ RooPlot* RooSimultaneous::plotOn(RooPlot *frame, RooLinkedList& cmdList) const
774783
// WVE -- do not adjust normalization for asymmetry plots
775784
RooLinkedList cmdList2(cmdList) ;
776785
if (!cmdList.find("Asymmetry")) {
777-
cmdList2.Add(&tmp1) ;
786+
replaceOrAdd(cmdList2, tmp1);
778787
}
779-
cmdList2.Add(&tmp2) ;
788+
replaceOrAdd(cmdList2, tmp2);
780789

781790
// Plot single component
782791
RooPlot* retFrame = getPdf(idxCatClone->getCurrentLabel())->plotOn(frame,cmdList2);
@@ -875,15 +884,15 @@ RooPlot* RooSimultaneous::plotOn(RooPlot *frame, RooLinkedList& cmdList) const
875884
RooCmdArg tmp2 = RooFit::ProjWData(*projDataSet,*projDataTmp) ;
876885
// WVE -- do not adjust normalization for asymmetry plots
877886
if (!cmdList.find("Asymmetry")) {
878-
cmdList2.Add(&tmp1) ;
887+
replaceOrAdd(cmdList2, tmp1);
879888
}
880-
cmdList2.Add(&tmp2) ;
889+
replaceOrAdd(cmdList2, tmp2);
881890

882891
RooPlot* frame2 ;
883892
if (!projSetTmp.empty()) {
884893
// Plot temporary function
885894
RooCmdArg tmp3 = RooFit::Project(projSetTmp) ;
886-
cmdList2.Add(&tmp3) ;
895+
replaceOrAdd(cmdList2, tmp3);
887896
frame2 = plotVar.plotOn(frame,cmdList2) ;
888897
} else {
889898
// Plot temporary function

tutorials/roofit/roofit/rf501_simultaneouspdf.C

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -92,7 +92,7 @@ void rf501_simultaneouspdf()
9292
// ---------------------------------------------------
9393

9494
// Perform simultaneous fit of model to data and model_ctl to data_ctl
95-
std::unique_ptr<RooFitResult> fitResult{simPdf.fitTo(combData, PrintLevel(-1), Save(), PrintLevel(-1))};
95+
std::unique_ptr<RooFitResult> fitResult{simPdf.fitTo(combData, Save(), PrintLevel(-1))};
9696
fitResult->Print();
9797

9898
// P l o t m o d e l s l i c e s o n d a t a s l i c e s

0 commit comments

Comments
 (0)