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
Comments
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 |
@cbernet is there a reason you are using the miniAODv1 rather than the miniAODv2? |
@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:
|
assign analysis |
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
|
For the C++ code, If I just tried to call
That comes from the 'bare' FWLite handler, which is different from using |
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
Also it runs substantially faster than the 'bare' case I first made :). |
For the fwlite::Handle<std::vector<reco::Conversion>> conversions;
conversions.getByLabel(event, "reducedEgamma","reducedConversions");
assert(conversions.isValid()); before getting the |
If I add
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
|
Running valgrind on the failing 'bare' FWLite give
which isn't as useful as one could have hoped. |
Thanks a lot for the very fast fix, @Dr15Jones and @pcanal , indeed it works. now trying it in the context of our analysis. |
So is this issue solved? It seems so |
Hi Fabio, yes, it works now. |
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:
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
The text was updated successfully, but these errors were encountered: