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

miniaod: cannot read RunIIFall17 MINIAOD with CMSSW_9_4_11_cand1 #25416

Closed
cbernet opened this issue Dec 4, 2018 · 15 comments
Closed

miniaod: cannot read RunIIFall17 MINIAOD with CMSSW_9_4_11_cand1 #25416

cbernet opened this issue Dec 4, 2018 · 15 comments

Comments

@cbernet
Copy link
Contributor

cbernet commented Dec 4, 2018

Hello, we encountered a very puzzling issue with the script below.

To test it, install a vanilla CMSSW_9_4_11_cand1, initialize your environment, and just run the script.

The problem is that when conversions are read, the gen particles get corrupted. In particular, all tau daughters cannot be accessed. After a segmentation fault, the script prints:

number of daughters 3
daughter 2 is not there

Are we doing anything wrong? This seems to point to a possible incompatibility between RunIIFall17 MINIAOD and CMSSW_9_4_11_cand1. If this incompatibility is known, what is the last release that can read RunIIFall17 MINIAOD?

Many thanks,

Colin, @lucastorterotot

ccing
@arizzi for mini AOD,
@guitargeek for conversions,
@Dr15Jones since the problem occurs in FWLite, maybe we are not initializing it properly

from ROOT import gSystem
from DataFormats.FWLite import Events, Handle

print('open input file...')

events = Events('root://cms-xrd-global.cern.ch//store/mc/RunIIFall17MiniAOD/DYJetsToLL_M-50_TuneCP5_13TeV-madgraphMLM-pythia8/MINIAODSIM/RECOSIMstep_94X_mc2017_realistic_v10-v1/00000/0293A280-B5F3-E711-8303-3417EBE33927.root')

conv_handle = Handle('std::vector<reco::Conversion>')
gen_ptc_handle  = Handle('std::vector<reco::GenParticle>')

n = 100000

print('start processing')

accepted = 0

for i,event in enumerate(events): 

    nEvent = event._event.id().event()
    print("processing event {0}: {1}...".format(i, nEvent))

    # comment the following 2 lines to run without crashing
    event.getByLabel(('reducedEgamma:reducedConversions'), conv_handle)
    convs     = conv_handle.product()
    
    event.getByLabel(('prunedGenParticles'), gen_ptc_handle)
    gen_ptcs = gen_ptc_handle.product()

    for p in gen_ptcs: 
        if abs(p.pdgId()) != 15 or not p.statusFlags().isPrompt():
            continue 
        tau = p
        for idau in xrange(tau.numberOfDaughters()):
            print 'number of daughters', tau.numberOfDaughters()    
            if not tau.daughter(idau) :
                print 'number of daughters', tau.numberOfDaughters()
                print 'daughter', idau, 'is not there'
                import pdb; pdb.set_trace() # to stop here if the error occurs

    accepted += 1
    if accepted==n:
        break
@cmsbuild
Copy link
Contributor

cmsbuild commented Dec 4, 2018

A new Issue was created by @cbernet Colin Bernet.

@davidlange6, @Dr15Jones, @smuzaffar, @fabiocos, @kpedro88 can you please review it and eventually sign/assign? Thanks.

cms-bot commands are listed here

@kpedro88
Copy link
Contributor

kpedro88 commented Dec 4, 2018

@cbernet is there a reason you are using the miniAODv1 rather than the miniAODv2?

@cbernet
Copy link
Contributor Author

cbernet commented Dec 4, 2018

@kpedro88 no particular reason (actually this file was provided by somebody else, I had a v2 file there before and I had not noticed the change). I just tried with v2, and the crash still occurs, unfortunately.

You can reproduce with this file:

root://cms-xrd-global.cern.ch//store/mc/RunIIFall17MiniAODv2/DYJetsToLL_M-50_TuneCP5_13TeV-madgraphMLM-pythia8/MINIAODSIM/PU2017RECOSIMstep_12Apr2018_94X_mc2017_realistic_v14-v1/90000/CAC74842-3348-E811-80B0-24BE05C68671.root 

@kpedro88
Copy link
Contributor

kpedro88 commented Dec 4, 2018

assign analysis

@cmsbuild
Copy link
Contributor

cmsbuild commented Dec 4, 2018

New categories assigned: analysis

@fgolf,@monttj,@peruzzim you have been requested to review this Pull request/Issue and eventually sign? Thanks

@Dr15Jones
Copy link
Contributor

The following code does work

#include "FWCore/FWLite/interface/FWLiteEnabler.h"
#include "DataFormats/HepMCCandidate/interface/GenParticle.h"
#include "TFile.h"
#include "TTree.h"

#include <vector>
#include <memory>
#include <cassert>
#include <iostream>
#include <cmath>

int main() {

  FWLiteEnabler::enable();

  //std::unique_ptr<TFile> f{ TFile::Open("root://cms-xrd-global.cern.ch//store/mc/RunIIFall17MiniAOD/DYJetsToLL_M-50_TuneCP5_13TeV-madgraphMLM-pythia8/MINIAODSIM/RECOSIMstep_94X_mc2017_realistic_v10-v1/00000/0293A280-B5F3-E711-8303-3417EBE33927.root") };                    
  std::unique_ptr<TFile> f{ TFile::Open("0293A280-B5F3-E711-8303-3417EBE33927.root") };

  auto events = dynamic_cast<TTree*>(f->Get("Events"));
  assert(events != nullptr);
  auto branch = events->GetBranch("recoGenParticles_prunedGenParticles__PAT.obj");


  std::vector<reco::GenParticle> particles;
  auto* pParticles = &particles;

  int count = 0;
  for(long e = 0; e< events->GetEntries(); ++e) {
    if( (e % 100) == 0 ) {
      std::cout <<"event "<<e<<std::endl;
    }
    events->GetEntry(e);
    branch->SetAddress(pParticles);
    branch->GetEntry(e);

    for(auto const& p: particles) {
      if ( std::abs(p.pdgId()) != 15 or not p.statusFlags().isPrompt() ) {
        continue ;
      }
      //std::cout <<"found tau "<<std::endl;                                                                                                                                                                                                                                       
      for( int i=0; i< p.numberOfDaughters(); ++i) {
        auto const d = p.daughter(i);
        if( nullptr == d) {
          std::cout <<"***missing daughter "<<i<<std::endl;
        } else {
          auto v = d->pdgId();
          if( (count % 100) == 0 ) {
            std::cout<< "  daughter "<<i<<" "<<v<<std::endl;
          }
        }
      }
      ++count;
    }
  }

  return 0;
}

where I compiled it in a CMSSW_9_4_11_cand1 workarea using

g++ -g -I /cvmfs/cms.cern.ch/slc6_amd64_gcc630/lcg/root/6.10.08-elfike2/include -I $CMSSW_RELEASE_BASE/src test.cc -L $CMSSW_RELEASE_BASE/lib/slc6_amd64_gcc630 -l FWCoreFWLite -L/cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/external/slc6_amd64_gcc630/lib/ -l Tree -l RIO -lCore

@Dr15Jones
Copy link
Contributor

For the C++ code, If I just tried to call TBranch::GetEntry without first calling TTree::GetEntry I got the following exception message

terminate called after throwing an instance of 'cms::Exception'
  what():  An exception of category 'GetEntryNotCalled' occurred.
Exception Message:
please call GetEntry for the 'Events' TTree for each event in order to make edm::Ref's work.
 Also be sure to call 'SetAddress' for all Branches after calling the GetEntry.

That comes from the 'bare' FWLite handler, which is different from using fwlite::Event.

@Dr15Jones
Copy link
Contributor

The following also works

#include "DataFormats/FWLite/interface/Event.h"
#include "DataFormats/FWLite/interface/Handle.h"
#include "DataFormats/HepMCCandidate/interface/GenParticle.h"
#include "TFile.h"
#include "TTree.h"

#include <vector>
#include <memory>
#include <cassert>
#include <iostream>
#include <cmath>

int main() {

  //std::unique_ptr<TFile> f{ TFile::Open("root://cms-xrd-global.cern.ch//store/mc/RunIIFall17MiniAOD/DYJetsToLL_M-50_TuneCP5_13TeV-madgraphMLM-pythia8/MINIAODSIM/RECOSIMstep_94X_mc2017_realistic_v10-v1/00000/0293A280-B5F3-E711-8303-3417EBE33927.root") };                    
  std::unique_ptr<TFile> f{ TFile::Open("0293A280-B5F3-E711-8303-3417EBE33927.root") };

  fwlite::Event event(f.get());

  int count = 0;
  int e = 0;
  for(event.toBegin(); ! event.atEnd(); ++event) {
    if( (e % 100) == 0 ) {
      std::cout <<"event "<<e<<std::endl;
    }
    ++e;

    fwlite::Handle<std::vector<reco::GenParticle>> particles;
    particles.getByLabel(event, "prunedGenParticles");

    for(auto const& p: *particles) {
      if ( std::abs(p.pdgId()) != 15 or not p.statusFlags().isPrompt() ) {
        continue ;
      }
      //std::cout <<"found tau "<<std::endl;                                                                                                                                                                                                                                       
      for( int i=0; i< p.numberOfDaughters(); ++i) {
        auto const d = p.daughter(i);
        if( nullptr == d) {
          std::cout <<"***missing daughter "<<i<<std::endl;
        } else {
          auto v = d->pdgId();
          if( (count % 100) == 0 ) {
            std::cout<< "  daughter "<<i<<" "<<v<<std::endl;
          }
        }
      }
      ++count;
    }
  }

  return 0;
}

compiled using

g++ -g -I /cvmfs/cms.cern.ch/slc6_amd64_gcc630/lcg/root/6.10.08-elfike2/include -I $CMSSW_RELEASE_BASE/src test_fwlite.cc -L $CMSSW_RELEASE_BASE/lib/slc6_amd64_gcc630 -l FWCoreFWLite -L/cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/external/slc6_amd64_gcc630/lib/ -l Tree -l RIO -lCore -l DataFormatsCommon -l DataFormatsFWLite

Also it runs substantially faster than the 'bare' case I first made :).

@Dr15Jones
Copy link
Contributor

For the fwlite::Event case, adding the lines

    fwlite::Handle<std::vector<reco::Conversion>> conversions;
    conversions.getByLabel(event, "reducedEgamma","reducedConversions");
    assert(conversions.isValid());

before getting the reco::GenParticles did not cause a problem.

@Dr15Jones
Copy link
Contributor

If I add reco::Conversion access to the 'bare' FWLite, it crashes with

===========================================================
#5  0x00007f48fccab339 in edm::PtrVectorBase::PtrVectorBase() () from /cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/lib/slc6_amd64_gcc630/libDataFormatsCommon.so
#6  0x00007f48f24c29ba in reco::Conversion::Conversion() () from /cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/lib/slc6_amd64_gcc630/libDataFormatsEgammaCandidates.so
#7  0x00007f48f24b0cb4 in reco::Conversion* std::__uninitialized_default_n_1<false>::__uninit_default_n<reco::Conversion*, unsigned long>(reco::Conversion*, unsigned long) () from /cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/lib/slc6_amd64_gcc630/libDataFormatsEgammaCandidates.so
#8  0x00007f48f24b19fc in std::vector<reco::Conversion, std::allocator<reco::Conversion> >::_M_default_append(unsigned long) () from /cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/lib/slc6_amd64_gcc630/libDataFormatsEgammaCandidates.so
#9  0x00007f48fd208d04 in TGenCollectionProxy::Allocate(unsigned int, bool) () from /cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/external/slc6_amd64_gcc630/lib/libRIO.so
#10 0x00007f48fd2a35f4 in int TStreamerInfoActions::ReadSTL<&TStreamerInfoActions::ReadSTLMemberWiseSameClass, &TStreamerInfoActions::ReadSTLObjectWiseFastArray>(TBuffer&, void*, TStreamerInfoActions::TConfiguration const*) () from /cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/external/slc6_amd64_gcc630/lib/libRIO.so
#11 0x00007f48fd1ada75 in TBufferFile::ApplySequence(TStreamerInfoActions::TActionSequence const&, void*) () from /cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/external/slc6_amd64_gcc630/lib/libRIO.so
#12 0x00007f48fd46fe1d in TBranchElement::ReadLeavesMember(TBuffer&) () from /cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/external/slc6_amd64_gcc630/lib/libTree.so
#13 0x00007f48fd46307a in TBranch::GetEntry(long long, int) () from /cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/external/slc6_amd64_gcc630/lib/libTree.so
#14 0x00007f48fd47b011 in TBranchElement::GetEntry(long long, int) () from /cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/external/slc6_amd64_gcc630/lib/libTree.so
#15 0x00007f48fd47afca in TBranchElement::GetEntry(long long, int) () from /cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/external/slc6_amd64_gcc630/lib/libTree.so
#16 0x00007f48fd4bdcc1 in TTree::GetEntry(long long, int) () from /cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/external/slc6_amd64_gcc630/lib/libTree.so
#17 0x00000000004012fe in main () at test.cc:36
===========================================================

The updated test program is

#include "FWCore/FWLite/interface/FWLiteEnabler.h"
#include "DataFormats/HepMCCandidate/interface/GenParticle.h"
#include "DataFormats/EgammaCandidates/interface/Conversion.h"
#include "TFile.h"
#include "TTree.h"

#include <vector>
#include <memory>
#include <cassert>
#include <iostream>
#include <cmath>

int main() {

  FWLiteEnabler::enable();

  //std::unique_ptr<TFile> f{ TFile::Open("root://cms-xrd-global.cern.ch//store/mc/RunIIFall17MiniAOD/DYJetsToLL_M-50_TuneCP5_13TeV-madgraphMLM-pythia8/MINIAODSIM/RECOSIMstep_94X_mc2017_realistic_v10-v1/00000/0293A280-B5F3-E711-8303-3417EBE33927.root") };                    
  std::unique_ptr<TFile> f{ TFile::Open("0293A280-B5F3-E711-8303-3417EBE33927.root") };

  auto events = dynamic_cast<TTree*>(f->Get("Events"));
  assert(events != nullptr);
  auto branch = events->GetBranch("recoGenParticles_prunedGenParticles__PAT.obj");

  std::vector<reco::GenParticle> particles;
  auto* pParticles = &particles;

  auto conversionBranch = events->GetBranch("recoConversions_reducedEgamma_reducedConversions_PAT.obj");
  std::vector<reco::Conversion> conversions;
  auto* pConversions = &conversions;

  int count = 0;
  for(long e = 0; e< events->GetEntries(); ++e) {
    if( (e % 100) == 0 ) {
      std::cout <<"event "<<e<<std::endl;
    }
    events->GetEntry(e);
    conversionBranch->SetAddress(pConversions);
    conversionBranch->GetEntry(e);

    branch->SetAddress(pParticles);
    branch->GetEntry(e);

    for(auto const& p: particles) {
      if ( std::abs(p.pdgId()) != 15 or not p.statusFlags().isPrompt() ) {
        continue ;
      }
      //std::cout <<"found tau "<<std::endl;                                                                                                                                                                                                                                       
      for( int i=0; i< p.numberOfDaughters(); ++i) {
        auto const d = p.daughter(i);
        if( nullptr == d) {
          std::cout <<"***missing daughter "<<i<<std::endl;
        } else {
          auto v = d->pdgId();
          if( (count % 100) == 0 ) {
            std::cout<< "  daughter "<<i<<" "<<v<<std::endl;
          }
        }
      }
      ++count;
    }
  }

  return 0;
}

compiled with

g++ -g -I /cvmfs/cms.cern.ch/slc6_amd64_gcc630/lcg/root/6.10.08-elfike2/include -I $CMSSW_RELEASE_BASE/src test.cc -L $CMSSW_RELEASE_BASE/lib/slc6_amd64_gcc630 -l FWCoreFWLite -L/cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/external/slc6_amd64_gcc630/lib/ -l Tree -l RIO -lCore -I /cvmfs/cms.cern.ch/slc6_amd64_gcc630/external/clhep/2.3.4.2-fmblme/include -l DataFormatsCommon

@Dr15Jones
Copy link
Contributor

Running valgrind on the failing 'bare' FWLite give

==4122165== Use of uninitialised value of size 8
==4122165==    at 0x50B2339: edm::PtrVectorBase::PtrVectorBase() (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/lib/slc6_amd64_gcc630/libDataFormatsCommon.so)
==4122165==    by 0x140249B9: reco::Conversion::Conversion() (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/lib/slc6_amd64_gcc630/libDataFormatsEgammaCandidates.so)
==4122165==    by 0x14012CB3: reco::Conversion* std::__uninitialized_default_n_1<false>::__uninit_default_n<reco::Conversion*, unsigned long>(reco::Conversion*, unsigned long) (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/lib/slc6_amd64_gcc630/libDataFormatsEgammaCandidates.so)
==4122165==    by 0x140139FB: std::vector<reco::Conversion, std::allocator<reco::Conversion> >::_M_default_append(unsigned long) (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/lib/slc6_amd64_gcc630/libDataFormatsEgammaCandidates.so)
==4122165==    by 0x4A93D03: TGenCollectionProxy::Allocate(unsigned int, bool) (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/lcg/root/6.10.08-elfike2/lib/libRIO.so)
==4122165==    by 0x4B2E5F3: int TStreamerInfoActions::ReadSTL<&TStreamerInfoActions::ReadSTLMemberWiseSameClass, &TStreamerInfoActions::ReadSTLObjectWiseFastArray>(TBuffer&, void*, TStreamerInfoActions::TConfiguration const*) (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/lcg/root/6.10.08-elfike2/lib/libRIO.so)
==4122165==    by 0x4A38A74: TBufferFile::ApplySequence(TStreamerInfoActions::TActionSequence const&, void*) (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/lcg/root/6.10.08-elfike2/lib/libRIO.so)
==4122165==    by 0x48CAE1C: TBranchElement::ReadLeavesMember(TBuffer&) (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/lcg/root/6.10.08-elfike2/lib/libTree.so)
==4122165==    by 0x48BE079: TBranch::GetEntry(long long, int) (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/lcg/root/6.10.08-elfike2/lib/libTree.so)
==4122165==    by 0x48D6010: TBranchElement::GetEntry(long long, int) (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/lcg/root/6.10.08-elfike2/lib/libTree.so)
==4122165==    by 0x48D5FC9: TBranchElement::GetEntry(long long, int) (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/lcg/root/6.10.08-elfike2/lib/libTree.so)
==4122165==    by 0x4918CC0: TTree::GetEntry(long long, int) (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/lcg/root/6.10.08-elfike2/lib/libTree.so)
==4122165== 
==4122165== Invalid write of size 8
==4122165==    at 0x50B2339: edm::PtrVectorBase::PtrVectorBase() (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/lib/slc6_amd64_gcc630/libDataFormatsCommon.so)
==4122165==    by 0x140249B9: reco::Conversion::Conversion() (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/lib/slc6_amd64_gcc630/libDataFormatsEgammaCandidates.so)
==4122165==    by 0x14012CB3: reco::Conversion* std::__uninitialized_default_n_1<false>::__uninit_default_n<reco::Conversion*, unsigned long>(reco::Conversion*, unsigned long) (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/lib/slc6_amd64_gcc630/libDataFormatsEgammaCandidates.so)
==4122165==    by 0x140139FB: std::vector<reco::Conversion, std::allocator<reco::Conversion> >::_M_default_append(unsigned long) (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/cms/cmssw/CMSSW_9_4_11_cand1/lib/slc6_amd64_gcc630/libDataFormatsEgammaCandidates.so)
==4122165==    by 0x4A93D03: TGenCollectionProxy::Allocate(unsigned int, bool) (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/lcg/root/6.10.08-elfike2/lib/libRIO.so)
==4122165==    by 0x4B2E5F3: int TStreamerInfoActions::ReadSTL<&TStreamerInfoActions::ReadSTLMemberWiseSameClass, &TStreamerInfoActions::ReadSTLObjectWiseFastArray>(TBuffer&, void*, TStreamerInfoActions::TConfiguration const*) (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/lcg/root/6.10.08-elfike2/lib/libRIO.so)
==4122165==    by 0x4A38A74: TBufferFile::ApplySequence(TStreamerInfoActions::TActionSequence const&, void*) (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/lcg/root/6.10.08-elfike2/lib/libRIO.so)
==4122165==    by 0x48CAE1C: TBranchElement::ReadLeavesMember(TBuffer&) (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/lcg/root/6.10.08-elfike2/lib/libTree.so)
==4122165==    by 0x48BE079: TBranch::GetEntry(long long, int) (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/lcg/root/6.10.08-elfike2/lib/libTree.so)
==4122165==    by 0x48D6010: TBranchElement::GetEntry(long long, int) (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/lcg/root/6.10.08-elfike2/lib/libTree.so)
==4122165==    by 0x48D5FC9: TBranchElement::GetEntry(long long, int) (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/lcg/root/6.10.08-elfike2/lib/libTree.so)
==4122165==    by 0x4918CC0: TTree::GetEntry(long long, int) (in /cvmfs/cms.cern.ch/slc6_amd64_gcc630/lcg/root/6.10.08-elfike2/lib/libTree.so)
==4122165==  Address 0x10 is not stack'd, malloc'd or (recently) free'd

which isn't as useful as one could have hoped.

@Dr15Jones
Copy link
Contributor

#25420 should fix the problem. [Thanks to @pcanal for helping find the issue.]

@cbernet
Copy link
Contributor Author

cbernet commented Dec 5, 2018

Thanks a lot for the very fast fix, @Dr15Jones and @pcanal , indeed it works. now trying it in the context of our analysis.

@fabiocos
Copy link
Contributor

So is this issue solved? It seems so

@cbernet
Copy link
Contributor Author

cbernet commented Dec 18, 2018

Hi Fabio, yes, it works now.

@cbernet cbernet closed this as completed Dec 18, 2018
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Projects
None yet
Development

No branches or pull requests

5 participants