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 pointers double deletions from PoolDBOutputService in CondTools subsystem #35731

Merged
merged 2 commits into from
Oct 20, 2021
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
9 changes: 5 additions & 4 deletions CondTools/BeamSpot/plugins/BeamSpotOnlineRecordsWriter.cc
Original file line number Diff line number Diff line change
Expand Up @@ -207,17 +207,18 @@ void BeamSpotOnlineRecordsWriter::endJob() {
edm::LogPrint("BeamSpotOnlineRecordsWriter") << "new tag requested";
if (fuseNewSince) {
edm::LogPrint("BeamSpotOnlineRecordsWriter") << "Using a new Since: " << fnewSince;
poolDbService->createNewIOV<BeamSpotOnlineObjects>(abeam.get(), fnewSince, poolDbService->endOfTime(), fLabel);
poolDbService->createNewIOV<BeamSpotOnlineObjects>(
abeam.release(), fnewSince, poolDbService->endOfTime(), fLabel);
} else
poolDbService->createNewIOV<BeamSpotOnlineObjects>(
abeam.get(), poolDbService->beginOfTime(), poolDbService->endOfTime(), fLabel);
abeam.release(), poolDbService->beginOfTime(), poolDbService->endOfTime(), fLabel);
} else {
edm::LogPrint("BeamSpotOnlineRecordsWriter") << "no new tag requested";
if (fuseNewSince) {
edm::LogPrint("BeamSpotOnlineRecordsWriter") << "Using a new Since: " << fnewSince;
poolDbService->appendSinceTime<BeamSpotOnlineObjects>(abeam.get(), fnewSince, fLabel);
poolDbService->appendSinceTime<BeamSpotOnlineObjects>(abeam.release(), fnewSince, fLabel);
} else
poolDbService->appendSinceTime<BeamSpotOnlineObjects>(abeam.get(), poolDbService->currentTime(), fLabel);
poolDbService->appendSinceTime<BeamSpotOnlineObjects>(abeam.release(), poolDbService->currentTime(), fLabel);
}
}
edm::LogPrint("BeamSpotOnlineRecordsWriter") << "[BeamSpotOnlineRecordsWriter] endJob done \n";
Expand Down
Original file line number Diff line number Diff line change
@@ -0,0 +1,64 @@
// user include files
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/global/EDAnalyzer.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/ServiceRegistry/interface/Service.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "CondFormats/DataRecord/interface/SiPhase2OuterTrackerLorentzAngleRcd.h"
#include "CondFormats/SiPhase2TrackerObjects/interface/SiPhase2OuterTrackerLorentzAngle.h"

//
//
// class decleration
//
class SiPhase2OuterTrackerLorentzAngleReader : public edm::global::EDAnalyzer<> {
public:
explicit SiPhase2OuterTrackerLorentzAngleReader(const edm::ParameterSet&);
~SiPhase2OuterTrackerLorentzAngleReader() override = default;
void analyze(edm::StreamID, edm::Event const&, edm::EventSetup const&) const override;

static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);

private:
const uint32_t printdebug_;
const std::string label_;
const edm::ESGetToken<SiPhase2OuterTrackerLorentzAngle, SiPhase2OuterTrackerLorentzAngleRcd> laToken_;
};

SiPhase2OuterTrackerLorentzAngleReader::SiPhase2OuterTrackerLorentzAngleReader(const edm::ParameterSet& iConfig)
: printdebug_(iConfig.getUntrackedParameter<uint32_t>("printDebug", 5)),
label_(iConfig.getUntrackedParameter<std::string>("label", "")),
laToken_(esConsumes(edm::ESInputTag{"", label_})) {}

void SiPhase2OuterTrackerLorentzAngleReader::analyze(edm::StreamID,
edm::Event const& iEvent,
edm::EventSetup const& iSetup) const {
const auto& lorentzAngles = iSetup.getData(laToken_);
edm::LogInfo("SiPhase2OuterTrackerLorentzAngleReader")
<< "[SiPhase2OuterTrackerLorentzAngleReader::analyze] End Reading SiPhase2OuterTrackerLorentzAngle with label "
<< label_ << std::endl;

const auto& detid_la = lorentzAngles.getLorentzAngles();
std::unordered_map<unsigned int, float>::const_iterator it;
size_t count = 0;
for (it = detid_la.begin(); it != detid_la.end() && count < printdebug_; it++) {
edm::LogInfo("SiPhase2OuterTrackerLorentzAngleReader") << "detid " << it->first << " \t"
<< " Lorentz angle " << it->second;
count++;
}
}

// ------------ method fills 'descriptions' with the allowed parameters for the module ------------
void SiPhase2OuterTrackerLorentzAngleReader::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.setComment("Module to read SiPhase2OuterTrackerLorentzAngle Payloads");
desc.addUntracked<uint32_t>("printDebug", 5)->setComment("maximum amount of print-outs");
desc.addUntracked<std::string>("label", "")->setComment("label from which to read the payload");
descriptions.addWithDefaultLabel(desc);
}

//define this as a plug-in
DEFINE_FWK_MODULE(SiPhase2OuterTrackerLorentzAngleReader);
Original file line number Diff line number Diff line change
Expand Up @@ -51,9 +51,7 @@ class SiPhase2OuterTrackerLorentzAngleWriter : public edm::one::EDAnalyzer<> {
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);

private:
void beginJob() override;
void analyze(const edm::Event&, const edm::EventSetup&) override;
void endJob() override;

// ----------member data ---------------------------
const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> topoToken_;
Expand Down Expand Up @@ -131,27 +129,21 @@ void SiPhase2OuterTrackerLorentzAngleWriter::analyze(const edm::Event& iEvent, c
edm::LogInfo("SiPhase2OuterTrackerLorentzAngleWriter")
<< " There are " << detsLAtoDB.size() << " OT Lorentz Angle values assigned" << std::endl;

// SiStripLorentzAngle object

// SiPhase2OuterTrackerLorentzAngle object
auto lorentzAngle = std::make_unique<SiPhase2OuterTrackerLorentzAngle>();
lorentzAngle->putLorentzAngles(detsLAtoDB);
edm::LogInfo("SiPhase2OuterTrackerLorentzAngleWriter") << "currentTime " << mydbservice->currentTime() << std::endl;
mydbservice->writeOne(lorentzAngle.get(), mydbservice->currentTime(), m_record);
mydbservice->writeOne(lorentzAngle.release(), mydbservice->currentTime(), m_record);
}

// ------------ method called once each job just before starting event loop ------------
void SiPhase2OuterTrackerLorentzAngleWriter::beginJob() {}

// ------------ method called once each job just after ending the event loop ------------
void SiPhase2OuterTrackerLorentzAngleWriter::endJob() {}

// ------------ method fills 'descriptions' with the allowed parameters for the module ------------
void SiPhase2OuterTrackerLorentzAngleWriter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.add<std::string>("record", "SiPhase2OuterTrackerLorentzAngleRcd");
desc.add<std::string>("tag", "SiPhase2OuterTrackerLorentzAngle");
desc.add<double>("value", 0.07);
descriptions.addDefault(desc);
desc.setComment("Module to write SiPhase2OuterTrackerLorentzAngle Payloads");
desc.add<std::string>("record", "SiPhase2OuterTrackerLorentzAngleRcd")->setComment("record to write");
desc.add<std::string>("tag", "SiPhase2OuterTrackerLorentzAngle")->setComment("tag to write");
desc.add<double>("value", 0.07)->setComment("value to be put in the payload");
descriptions.addWithDefaultLabel(desc);
}

//define this as a plug-in
Expand Down
Original file line number Diff line number Diff line change
@@ -0,0 +1,78 @@
#! /usr/bin/env cmsRun
# Author: Marco Musich (October 2021)
import FWCore.ParameterSet.Config as cms
process = cms.Process("TEST")

###################################################################
# Messages
###################################################################
process.load('FWCore.MessageService.MessageLogger_cfi')
process.MessageLogger.cerr.enable = False
process.MessageLogger.SiPhase2OuterTrackerLorentzAngleReader=dict()
process.MessageLogger.SiPhase2OuterTrackerLorentzAngle=dict()
process.MessageLogger.cout = cms.untracked.PSet(
enable = cms.untracked.bool(True),
enableStatistics = cms.untracked.bool(True),
threshold = cms.untracked.string("INFO"),
default = cms.untracked.PSet(limit = cms.untracked.int32(0)),
FwkReport = cms.untracked.PSet(limit = cms.untracked.int32(-1),
reportEvery = cms.untracked.int32(1000)
),
SiPhase2OuterTrackerLorentzAngleReader = cms.untracked.PSet( limit = cms.untracked.int32(-1)),
SiPhase2OuterTrackerLorentzAngle = cms.untracked.PSet( limit = cms.untracked.int32(-1)),
)

###################################################################
# A data source must always be defined.
# We don't need it, so here's a dummy one.
###################################################################
process.source = cms.Source("EmptyIOVSource",
timetype = cms.string('runnumber'),
firstValue = cms.uint64(1),
lastValue = cms.uint64(1),
interval = cms.uint64(1)
)

###################################################################
# Input data
###################################################################
tag = 'SiPhase2OuterTrackerLorentzAngle_T15'
suffix = 'v0'
inFile = tag+'_'+suffix+'.db'
inDB = 'sqlite_file:'+inFile

process.load("CondCore.CondDB.CondDB_cfi")
# input database (in this case the local sqlite file)
process.CondDB.connect = inDB

process.PoolDBESSource = cms.ESSource("PoolDBESSource",
process.CondDB,
DumpStat=cms.untracked.bool(True),
toGet = cms.VPSet(cms.PSet(record = cms.string("SiPhase2OuterTrackerLorentzAngleRcd"),
tag = cms.string(tag))
)
)

###################################################################
# check the ES data getter
###################################################################
process.get = cms.EDAnalyzer("EventSetupRecordDataGetter",
toGet = cms.VPSet(cms.PSet(
record = cms.string(' SiPhase2OuterTrackerLorentzAngleRcd'),
data = cms.vstring('SiPhase2OuterTrackerLorentzAngle')
)),
verbose = cms.untracked.bool(True)
)

###################################################################
# Payload reader
###################################################################
import CondTools.SiPhase2Tracker.siPhase2OuterTrackerLorentzAngleReader_cfi as _mod
process.LAPayloadReader = _mod.siPhase2OuterTrackerLorentzAngleReader.clone(printDebug = 10,
label = "")

###################################################################
# Path
###################################################################
process.p = cms.Path(process.LAPayloadReader)

Original file line number Diff line number Diff line change
Expand Up @@ -2,13 +2,18 @@
TEST_DIR=$CMSSW_BASE/src/CondTools/SiPhase2Tracker/test
echo "test dir: $TEST_DIR"

printf "testing writing Phase2 Lorentz Angle \n\n"
## need to be in order (don't read before writing)
cmsRun ${TEST_DIR}/SiPhase2OuterTrackerLorentzAngleWriter_cfg.py
cmsRun ${TEST_DIR}/SiPhase2OuterTrackerLorentzAngleReader_cfg.py

printf "testing writing Phase2 Tracker Cabling Map (test) \n\n"
## need to be in order (don't read before writing)
cmsRun ${TEST_DIR}/DTCCablingMapTestProducer_write.py
cmsRun ${TEST_DIR}/DTCCablingMapTestProducer_retrieve.py
cmsRun ${TEST_DIR}/DTCCablingMapTestProducer_dump.py

printf "testing writing Phase2 Tracker Cabling Map \n\n"
## need to be in order (don't read before writing)
cmsRun ${TEST_DIR}/DTCCablingMapProducer_write.py
cmsRun ${TEST_DIR}/DTCCablingMapProducer_retrieve.py
Expand Down