Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Fix usage of deprecated methods returning TIterator (replace with range loops or begin/end) #43373

Merged
merged 2 commits into from
Nov 27, 2023
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
26 changes: 11 additions & 15 deletions PhysicsTools/TagAndProbe/src/TagProbeFitter.cc
Original file line number Diff line number Diff line change
Expand Up @@ -329,9 +329,7 @@ string TagProbeFitter::calculateEfficiency(string dirName,
} else {
// disactive not needed branches
inputTree->SetBranchStatus("*", false);
TIterator* iter = dataVars.createIterator();
TObject* obj(nullptr);
while ((obj = iter->Next()))
for (TObject* obj : dataVars)
Copy link
Contributor

@guitargeek guitargeek Dec 5, 2023

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Thanks for these changes! This even fixed a memory leak of the TIterator !

inputTree->SetBranchStatus(obj->GetName(), true);
}

Expand Down Expand Up @@ -480,8 +478,8 @@ string TagProbeFitter::calculateEfficiency(string dirName,
if (data_bin->numEntries() > 0) {
//set the values of binnedVariables to the mean value in this data bin
RooArgSet meanOfVariables;
RooLinkedListIter vit = binnedVariables.iterator();
for (RooRealVar* v = (RooRealVar*)vit.Next(); v != nullptr; v = (RooRealVar*)vit.Next()) {
for (const RooAbsArg* vv : binnedVariables) {
const RooRealVar* v = dynamic_cast<const RooRealVar*>(vv);
meanOfVariables.addClone(*v);
double mean = w->data("data")->mean(*v);
RooBinning binning((RooBinning&)v->getBinning());
Expand Down Expand Up @@ -788,8 +786,7 @@ void TagProbeFitter::saveFitPlot(RooWorkspace* w) {
RooAbsPdf& pdf = *w->pdf("simPdf");
std::unique_ptr<RooArgSet> obs(pdf.getObservables(*dataAll));
RooRealVar* mass = nullptr;
RooLinkedListIter it = obs->iterator();
for (RooAbsArg* v = (RooAbsArg*)it.Next(); v != nullptr; v = (RooAbsArg*)it.Next()) {
for (RooAbsArg* v : *obs) {
if (!v->InheritsFrom("RooRealVar"))
continue;
mass = (RooRealVar*)v;
Expand Down Expand Up @@ -867,8 +864,7 @@ void TagProbeFitter::saveDistributionsPlot(RooWorkspace* w) {

const RooArgSet* vars = dataAll->get();
vector<RooRealVar*> reals;
RooLinkedListIter it = vars->iterator();
for (RooAbsArg* v = (RooAbsArg*)it.Next(); v != nullptr; v = (RooAbsArg*)it.Next()) {
for (RooAbsArg* v : *vars) {
if (!v->InheritsFrom("RooRealVar"))
continue;
reals.push_back((RooRealVar*)v);
Expand Down Expand Up @@ -909,15 +905,15 @@ void TagProbeFitter::saveEfficiencyPlots(RooDataSet& eff,
const TString& effName,
RooArgSet& binnedVariables,
RooArgSet& mappedCategories) {
RooLinkedListIter v1it = binnedVariables.iterator();
bool isOnePoint =
(eff.numEntries() == 1); // for datasets with > 1 entry, we don't make plots for variables with just one bin
for (RooRealVar* v1 = (RooRealVar*)v1it.Next(); v1 != nullptr; v1 = (RooRealVar*)v1it.Next()) {
for (auto it1 = binnedVariables.begin(); it1 != binnedVariables.end(); it1++) {
RooRealVar* v1 = dynamic_cast<RooRealVar*>(*it1);
RooArgSet binCategories1D;
if (v1->numBins() == 1 && !isOnePoint)
continue;
RooLinkedListIter v2it = binnedVariables.iterator();
for (RooRealVar* v2 = (RooRealVar*)v2it.Next(); v2 != nullptr; v2 = (RooRealVar*)v2it.Next()) {
for (auto it2 = binnedVariables.begin(); it2 != binnedVariables.end(); it2++) {
RooRealVar* v2 = dynamic_cast<RooRealVar*>(*it2);
if (v2 == v1)
continue;
if (v2->numBins() == 1 && !isOnePoint)
Expand All @@ -926,8 +922,8 @@ void TagProbeFitter::saveEfficiencyPlots(RooDataSet& eff,
RooBinningCategory(TString(v2->GetName()) + "_bins", TString(v2->GetName()) + "_bins", *v2));

RooArgSet binCategories2D;
RooLinkedListIter v3it = binnedVariables.iterator();
for (RooRealVar* v3 = (RooRealVar*)v3it.Next(); v3 != nullptr; v3 = (RooRealVar*)v3it.Next()) {
for (auto it3 = binnedVariables.begin(); it3 != binnedVariables.end(); it3++) {
RooRealVar* v3 = dynamic_cast<RooRealVar*>(*it3);
if (v3 == v1 || v3 == v2)
continue;
binCategories2D.addClone(
Expand Down
14 changes: 3 additions & 11 deletions PhysicsTools/Utilities/src/SideBandSubtraction.cc
Original file line number Diff line number Diff line change
Expand Up @@ -98,15 +98,14 @@ int SideBandSubtract::doSubtraction(RooRealVar* variable,
//out how to do this in one shot to avoid a loop
//O(N_vars*N_events)...

TIterator* iter = (TIterator*)Data->get()->createIterator();
RooAbsArg* var = nullptr;
RooRealVar* sep_var = nullptr;
while ((var = (RooAbsArg*)iter->Next())) {
for (const auto& var : *Data->get()) {
Copy link
Contributor

@guitargeek guitargeek Dec 5, 2023

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

And this also fixed a memory leak 🙂

if ((string)var->GetName() == (string)SeparationVariable->GetName()) {
sep_var = (RooRealVar*)var;
break;
}
}

for (int i = 0; i < Data->numEntries(); i++) {
Data->get(i);
Double_t value = variable->getVal();
Expand Down Expand Up @@ -431,21 +430,14 @@ int SideBandSubtract::doGlobalFit() {

//need to grab sbs objects after each global fit, because they get reset
resetSBSProducts();
TIterator* iter = (TIterator*)Data->get()->createIterator();
RooAbsArg* variable;
while ((variable = (RooAbsArg*)iter->Next())) {
for (const auto& variable : *Data->get()) {
for (unsigned int i = 0; i < BaseHistos.size(); i++) {
if ((string)variable->GetName() != (string)SeparationVariable->GetName() &&
(string)variable->GetName() == (string)BaseHistos[i]->GetName())
doSubtraction((RooRealVar*)variable, SignalSidebandRatio, i);
}
}

// clean up our memory...
if (variable)
delete variable;
if (iter)
delete iter;
return 0;
}
void SideBandSubtract::doFastSubtraction(TH1F& Total, TH1F& Result, SbsRegion& leftRegion, SbsRegion& rightRegion) {
Expand Down