Skip to content

[RF] Uninitialized value issue reported by Valgrind #18323

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

Open
1 task done
TomasDado opened this issue Apr 9, 2025 · 5 comments
Open
1 task done

[RF] Uninitialized value issue reported by Valgrind #18323

TomasDado opened this issue Apr 9, 2025 · 5 comments
Assignees
Labels

Comments

@TomasDado
Copy link
Contributor

Check duplicate issues.

  • Checked for duplicates

Description

Using the following valgrind setup:

valgrind --leak-check=full          --show-leak-kinds=all          --track-origins=yes          --verbose          --log-file=valgrind-out.txt --suppressions=$ROOTSYS/etc/valgrind-root.supp

Gives

==42205== 112 errors in context 5 of 27:
==42205== Conditional jump or move depends on uninitialised value(s)
==42205==    at 0x5B355D0: RooAbsReal::getValV(RooArgSet const*) const (RooAbsReal.cxx:301)
==42205==    by 0x5CD978A: getVal (RooAbsReal.h:117)
==42205==    by 0x5CD978A: RooProduct::evaluate() const (RooProduct.cxx:360)
==42205==    by 0x5B35415: RooAbsReal::traceEval(RooArgSet const*) const (RooAbsReal.cxx:320)
==42205==    by 0x5B355F6: RooAbsReal::getValV(RooArgSet const*) const (RooAbsReal.cxx:307)
==42205==    by 0x5CF6D72: getVal (RooAbsReal.h:117)
==42205==    by 0x5CF6D72: RooRealSumPdf::evaluate(RooAbsReal const&, RooArgList const&, RooArgList const&, bool, bool&) (RooRealSumPdf.cxx:229)
==42205==    by 0x5B14CEA: RooAbsPdf::getValV(RooArgSet const*) const (RooAbsPdf.cxx:324)
==42205==    by 0x5B721D5: integrand (RooAbsIntegrator.h:35)
==42205==    by 0x5B721D5: RooBinIntegrator::recursive_integration(unsigned int, double, ROOT::Math::KahanSum<double, 1u>&) (RooBinIntegrator.cxx:190)
==42205==    by 0x5B72277: integral (RooBinIntegrator.cxx:170)
==42205==    by 0x5B72277: RooBinIntegrator::integral(double const*) (RooBinIntegrator.cxx:163)
==42205==    by 0x5B0A48A: RooAbsIntegrator::calculate(double const*) (RooAbsIntegrator.cxx:54)
==42205==    by 0x5CEC819: integrate (RooRealIntegral.cxx:972)
==42205==    by 0x5CEC819: RooRealIntegral::sum() const (RooRealIntegral.cxx:958)
==42205==    by 0x5CEF2DC: RooRealIntegral::evaluate() const (RooRealIntegral.cxx:841)
==42205==    by 0x5B35415: RooAbsReal::traceEval(RooArgSet const*) const (RooAbsReal.cxx:320)
==42205==  Uninitialised value was created by a stack allocation
==42205==    at 0x677A031: xRooNode::BuildHistogram(RooAbsLValue*, bool, bool, int, int, xRooNode const&, bool, bool, int, TH1*, bool, bool) const (xRooNode.cxx:8055)

Reproducer

Tagging @will-cern who kindly said that if a reproducer is needed he can provide it```

ROOT version

ROOT Version: 6.35.01
Built for linuxx8664gcc on Apr 08 2025, 06:57:24
From heads/master@7f390ce5ef

Installation method

build from source

Operating system

Ubuntu 24.10

Additional context

No response

@TomasDado TomasDado added the bug label Apr 9, 2025
@ferdymercury ferdymercury changed the title [RF] Memory issue reported by Valgrind [RF] Uninitialized value issue reported by Valgrind Apr 9, 2025
@ferdymercury
Copy link
Collaborator

@will-cern could you provide a reproducer?

@TomasDado
Copy link
Contributor Author

Sorry it was a mistake on my side, I will try to provide the reproducer and not Will. He just wanted to be tagged.

@TomasDado
Copy link
Contributor Author

TomasDado commented Apr 10, 2025

Here is the reproducer:

#include "RooFit/xRooFit/xRooFit.h"
using namespace ROOT::Experimental::XRooFit;

int main(int argv, char** argc) {
  xRooNode node("FitExampleNtuple_combined_FitExampleNtuple_model.root");
  auto pdf = node["simPdf"];
  auto channel = pdf->find("ljets_5j3b_HT");
  auto samples = channel->find("samples");
  auto sample = samples->find("ttlight_ljets_5j3b_HT_shapes");

  const auto bins = sample->GetBinContents();

  const double central = channel->GetBinContent(1);
}
cmake_minimum_required(VERSION 3.27)
project(mwe)

set( CMAKE_CXX_FLAGS "-g" )

set( _xRooFitComp "RooFitXRooFit" )
 message(STATUS "Using xRooFit From ROOT")

find_package( ROOT COMPONENTS RooFitCore RooFit ${_xRooFitComp} )

add_executable(mwe mwe.cxx)
target_include_directories(mwe PUBLIC ${ROOT_INCLUDE_DIRS} ${XROOFIT_INCLUDE_DIRS} )
target_link_libraries(mwe PUBLIC ${ROOT_LIBRARIES} ${XROOFIT_LIBRARIES} )

Note that just calling const auto bins = sample->GetBinContents(); seems fine, the call to channel->GetBinContent(1) seems to lead to the problem:

==19638== 6 errors in context 1 of 23:
==19638== Conditional jump or move depends on uninitialised value(s)
==19638==    at 0x6249750: RooAbsReal::getValV(RooArgSet const*) const (RooAbsReal.cxx:301)
==19638==    by 0x63EDA4A: getVal (RooAbsReal.h:117)
==19638==    by 0x63EDA4A: RooProduct::evaluate() const (RooProduct.cxx:360)
==19638==    by 0x6249595: RooAbsReal::traceEval(RooArgSet const*) const (RooAbsReal.cxx:320)
==19638==    by 0x6249776: RooAbsReal::getValV(RooArgSet const*) const (RooAbsReal.cxx:307)
==19638==    by 0x640B032: getVal (RooAbsReal.h:117)
==19638==    by 0x640B032: RooRealSumPdf::evaluate(RooAbsReal const&, RooArgList const&, RooArgList const&, bool, bool&) (RooRealSumPdf.cxx:229)
==19638==    by 0x6228E6A: RooAbsPdf::getValV(RooArgSet const*) const (RooAbsPdf.cxx:324)
==19638==    by 0x6286355: integrand (RooAbsIntegrator.h:35)
==19638==    by 0x6286355: RooBinIntegrator::recursive_integration(unsigned int, double, ROOT::Math::KahanSum<double, 1u>&) (RooBinIntegrator.cxx:190)
==19638==    by 0x62863F7: integral (RooBinIntegrator.cxx:170)
==19638==    by 0x62863F7: RooBinIntegrator::integral(double const*) (RooBinIntegrator.cxx:163)
==19638==    by 0x621E60A: RooAbsIntegrator::calculate(double const*) (RooAbsIntegrator.cxx:54)
==19638==    by 0x6400AD9: integrate (RooRealIntegral.cxx:972)
==19638==    by 0x6400AD9: RooRealIntegral::sum() const (RooRealIntegral.cxx:958)
==19638==    by 0x640359C: RooRealIntegral::evaluate() const (RooRealIntegral.cxx:841)
==19638==    by 0x6249595: RooAbsReal::traceEval(RooArgSet const*) const (RooAbsReal.cxx:320)
==19638==  Uninitialised value was created by a stack allocation
==19638==    at 0x4E6078D: ROOT::Experimental::XRooFit::xRooNode::BuildHistogram(RooAbsLValue*, bool, bool, int, int, ROOT::Experimental::XRooFit::xRooNode const&, bool, bool, int, TH1*, bool, bool) const (xRooNode.cxx:8055)

The test file is available at /afs/cern.ch/user/t/tdado/public/Temp/FitExampleNtuple_combined_FitExampleNtuple_model.root

@ferdymercury
Copy link
Collaborator

Thanks Tomas!
Unrelated: Did you check the new force-pushed commit of #18182 wrt runtime speed ?

@TomasDado
Copy link
Contributor Author

Thanks Tomas! Unrelated: Did you check the new force-pushed commit of #18182 wrt runtime speed ?

I did, but it seems the problem must be something else since we actually disable the histogram registration globally in the code with TH1::AddDirectory(kFALSE);

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
Projects
None yet
Development

No branches or pull requests

3 participants