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

Add SiStripApproximateClusterCollection as a simple format for RAW (re-vamped) #42495

Merged
merged 7 commits into from
Aug 15, 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
Original file line number Diff line number Diff line change
Expand Up @@ -132,7 +132,6 @@
workflows[140.5611] = ['RunHI2018reMINIAOD',['RunHI2018AOD','REMINIAODHID18','HARVESTHI18MINIAOD']]
workflows[140.57] = ['RunHI2018Reduced',['RunHI2018Reduced','RECOHID18','HARVESTDHI18']]
workflows[140.58] = ['',['RunHI2018','RAWPRIMEHI18','RECOHID18APPROXCLUSTERS','HARVESTDHI18']]
workflows[140.60] = ['',['RunHI2022','RECOHID22APPROXCLUSTERS','HARVESTDHI22']]
workflows[140.61] = ['',['RunHI2022FullFormat','RECOHID22','HARVESTDHI22']]

### run2 2015B 50ns ###
Expand Down
21 changes: 0 additions & 21 deletions Configuration/PyReleaseValidation/python/relval_steps.py
Original file line number Diff line number Diff line change
Expand Up @@ -2192,27 +2192,6 @@ def lhegensim2018ml(fragment,howMuch):
'-n':'100'
}

steps['RAWPRIMEHI22']={ '--scenario':'pp',
'--conditions':'auto:run3_data_prompt',
'-s':'REPACK:DigiToApproxClusterRaw',
'--datatier':'GEN-SIM-DIGI-RAW-HLTDEBUG',
'--eventcontent':'REPACKRAW',
'--era':'Run3_pp_on_PbPb_approxSiStripClusters',
'-n':'10',
'--customise_commands':'\"process.hltSiStripRawToDigi.ProductLabel=\'rawDataCollector\';process.hltScalersRawToDigi.scalersInputTag=\'rawDataCollector\'\"',
'--process':'REHLT'
}

steps['RECOHID22APPROXCLUSTERS']=merge([{ '--scenario':'pp',
'--conditions':'auto:run3_data_prompt',
'-s':'RAW2DIGI,L1Reco,RECO,DQM:@commonFakeHLT+@standardDQMFakeHLT',
'--datatier':'AOD,DQMIO',
'--eventcontent':'AOD,DQM',
'--era':'Run3_pp_on_PbPb_approxSiStripClusters',
'--repacked':'',
'-n':'100'
},steps['RECOHID15']])

steps['RECOHID22']=merge([{ '--scenario':'pp',
'--conditions':'auto:run3_data_prompt',
'-s':'RAW2DIGI,L1Reco,RECO,DQM:@commonFakeHLT+@standardDQMFakeHLT',
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -21,6 +21,7 @@
#include "DQMServices/Core/interface/MonitorElement.h"
#include "DataFormats/Common/interface/DetSet.h"
#include "DataFormats/Common/interface/DetSetVectorNew.h"
#include "DataFormats/SiStripCluster/interface/SiStripApproximateClusterCollection.h"
#include "DataFormats/SiStripCluster/interface/SiStripApproximateCluster.h"
#include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
#include "FWCore/Framework/interface/Event.h"
Expand Down Expand Up @@ -102,7 +103,7 @@ class SiStripMonitorApproximateCluster : public DQMEDAnalyzer {
MonitorElement* h_deltaEndStrip_{nullptr};

// Event Data
edm::EDGetTokenT<edmNew::DetSetVector<SiStripApproximateCluster>> approxClustersToken_;
edm::EDGetTokenT<SiStripApproximateClusterCollection> approxClustersToken_;
edm::EDGetTokenT<edmNew::DetSetVector<SiStripCluster>> stripClustersToken_;
const edmNew::DetSetVector<SiStripCluster>* stripClusterCollection_;

Expand All @@ -117,7 +118,7 @@ SiStripMonitorApproximateCluster::SiStripMonitorApproximateCluster(const edm::Pa
: folder_(iConfig.getParameter<std::string>("folder")),
compareClusters_(iConfig.getParameter<bool>("compareClusters")),
// Poducer name of input StripClusterCollection
approxClustersToken_(consumes<edmNew::DetSetVector<SiStripApproximateCluster>>(
approxClustersToken_(consumes<SiStripApproximateClusterCollection>(
iConfig.getParameter<edm::InputTag>("ApproxClustersProducer"))) {
tkGeomToken_ = esConsumes();
if (compareClusters_) {
Expand All @@ -139,7 +140,7 @@ void SiStripMonitorApproximateCluster::analyze(const edm::Event& iEvent, const e
const auto tkDets = tkGeom->dets();

// get collection of DetSetVector of clusters from Event
edm::Handle<edmNew::DetSetVector<SiStripApproximateCluster>> approx_cluster_detsetvector;
edm::Handle<SiStripApproximateClusterCollection> approx_cluster_detsetvector;
iEvent.getByToken(approxClustersToken_, approx_cluster_detsetvector);
if (!approx_cluster_detsetvector.isValid()) {
edm::LogError("SiStripMonitorApproximateCluster")
Expand All @@ -164,11 +165,11 @@ void SiStripMonitorApproximateCluster::analyze(const edm::Event& iEvent, const e
}

int nApproxClusters{0};
const edmNew::DetSetVector<SiStripApproximateCluster>* clusterCollection = approx_cluster_detsetvector.product();
const SiStripApproximateClusterCollection* clusterCollection = approx_cluster_detsetvector.product();

for (const auto& detClusters : *clusterCollection) {
edmNew::DetSet<SiStripCluster> strip_clusters_detset;
const auto& detid = detClusters.detId(); // get the detid of the current detset
const auto& detid = detClusters.id(); // get the detid of the current detset

// starts here comaparison with regular clusters
if (compareClusters_) {
Expand Down Expand Up @@ -233,6 +234,7 @@ void SiStripMonitorApproximateCluster::analyze(const edm::Event& iEvent, const e

} // loop on clusters in a detset
} // loop on the detset vector

h_nclusters_->Fill(nApproxClusters);
}

Expand Down
Original file line number Diff line number Diff line change
@@ -0,0 +1,107 @@
#ifndef DataFormats_SiStripCluster_SiStripApproximateClusterCollection_h
#define DataFormats_SiStripCluster_SiStripApproximateClusterCollection_h

#include <vector>

#include "DataFormats/SiStripCluster/interface/SiStripApproximateCluster.h"

/**
* This class provides a minimal interface that resembles
* edmNew::DetSetVector, but is crafted such that we are comfortable
* to provide an infinite backwards compatibility guarantee for it
* (like all RAW data). Any modifications need to be made with care.
* Please consult core software group if in doubt.
**/
class SiStripApproximateClusterCollection {
public:
// Helper classes to make creation and iteration easier
class Filler {
public:
void push_back(SiStripApproximateCluster const& cluster) { clusters_.push_back(cluster); }

private:
friend SiStripApproximateClusterCollection;
Filler(std::vector<SiStripApproximateCluster>& clusters) : clusters_(clusters) {}

std::vector<SiStripApproximateCluster>& clusters_;
};

class const_iterator;
class DetSet {
public:
using const_iterator = std::vector<SiStripApproximateCluster>::const_iterator;

unsigned int id() const { return coll_->detIds_[detIndex_]; }

const_iterator begin() const { return coll_->clusters_.begin() + clusBegin_; }
const_iterator cbegin() const { return begin(); }
const_iterator end() const { return coll_->clusters_.begin() + clusEnd_; }
const_iterator cend() const { return end(); }

private:
friend SiStripApproximateClusterCollection::const_iterator;
DetSet(SiStripApproximateClusterCollection const* coll, unsigned int detIndex)
: coll_(coll),
detIndex_(detIndex),
clusBegin_(coll_->beginIndices_[detIndex]),
clusEnd_(detIndex == coll_->beginIndices_.size() - 1 ? coll->clusters_.size()
: coll_->beginIndices_[detIndex + 1]) {}

SiStripApproximateClusterCollection const* const coll_;
unsigned int const detIndex_;
unsigned int const clusBegin_;
unsigned int const clusEnd_;
};

class const_iterator {
public:
DetSet operator*() const { return DetSet(coll_, index_); }

const_iterator& operator++() {
++index_;
if (index_ == coll_->detIds_.size()) {
*this = const_iterator();
}
return *this;
}

const_iterator operator++(int) {
const_iterator clone = *this;
++(*this);
return clone;
}

bool operator==(const_iterator const& other) const { return coll_ == other.coll_ and index_ == other.index_; }
bool operator!=(const_iterator const& other) const { return not operator==(other); }

private:
friend SiStripApproximateClusterCollection;
// default-constructed object acts as the sentinel
const_iterator() = default;
const_iterator(SiStripApproximateClusterCollection const* coll) : coll_(coll) {}

SiStripApproximateClusterCollection const* coll_ = nullptr;
unsigned int index_ = 0;
};

// Actual public interface
SiStripApproximateClusterCollection() = default;

void reserve(std::size_t dets, std::size_t clusters);
Filler beginDet(unsigned int detId);

const_iterator begin() const { return const_iterator(this); }
const_iterator cbegin() const { return begin(); }
const_iterator end() const { return const_iterator(); }
const_iterator cend() const { return end(); }

private:
// The detIds_ and beginIndices_ have one element for each Det. An
// element of beginIndices_ points to the first cluster of the Det
// in clusters_.
std::vector<unsigned int> detIds_; // DetId for the Det
std::vector<unsigned int> beginIndices_;
std::vector<SiStripApproximateCluster> clusters_;
};

#endif
Original file line number Diff line number Diff line change
@@ -0,0 +1,13 @@
#include "DataFormats/SiStripCluster/interface/SiStripApproximateClusterCollection.h"

void SiStripApproximateClusterCollection::reserve(std::size_t dets, std::size_t clusters) {
detIds_.reserve(dets);
beginIndices_.reserve(dets);
clusters_.reserve(clusters);
}

SiStripApproximateClusterCollection::Filler SiStripApproximateClusterCollection::beginDet(unsigned int detId) {
detIds_.push_back(detId);
beginIndices_.push_back(clusters_.size());
return Filler(clusters_);
}
1 change: 1 addition & 0 deletions DataFormats/SiStripCluster/src/classes.h
Original file line number Diff line number Diff line change
Expand Up @@ -7,6 +7,7 @@
#include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
#include "DataFormats/SiStripCluster/interface/SiStripClustersSOA.h"
#include "DataFormats/SiStripCluster/interface/SiStripApproximateCluster.h"
#include "DataFormats/SiStripCluster/interface/SiStripApproximateClusterCollection.h"
#include "DataFormats/Common/interface/ContainerMask.h"

#endif // SISTRIPCLUSTER_CLASSES_H
4 changes: 4 additions & 0 deletions DataFormats/SiStripCluster/src/classes_def.xml
Original file line number Diff line number Diff line change
Expand Up @@ -31,6 +31,10 @@
<version ClassVersion="4" checksum="2854791577"/>
<version ClassVersion="3" checksum="2041370183"/>
</class>
<class name="SiStripApproximateClusterCollection" ClassVersion="3">
<version ClassVersion="3" checksum="3101417750"/>
</class>
<class name="edm::Wrapper<SiStripApproximateClusterCollection>"/>

<class name="edmNew::DetSetVector<SiStripApproximateCluster>"/>
<class name="edm::Wrapper<edmNew::DetSetVector<SiStripApproximateCluster>>"/>
Expand Down
Original file line number Diff line number Diff line change
@@ -1,6 +1,6 @@
#include "DataFormats/Common/interface/DetSetVector.h"
#include "DataFormats/Common/interface/DetSetVectorNew.h"
#include "DataFormats/SiStripCluster/interface/SiStripApproximateCluster.h"
#include "DataFormats/SiStripCluster/interface/SiStripApproximateClusterCollection.h"
#include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
Expand All @@ -25,13 +25,12 @@ class SiStripApprox2Clusters : public edm::global::EDProducer<> {
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);

private:
edm::EDGetTokenT<edmNew::DetSetVector<SiStripApproximateCluster>> clusterToken_;
edm::EDGetTokenT<SiStripApproximateClusterCollection> clusterToken_;
edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> tkGeomToken_;
};

SiStripApprox2Clusters::SiStripApprox2Clusters(const edm::ParameterSet& conf) {
clusterToken_ = consumes<edmNew::DetSetVector<SiStripApproximateCluster>>(
conf.getParameter<edm::InputTag>("inputApproxClusters"));
clusterToken_ = consumes(conf.getParameter<edm::InputTag>("inputApproxClusters"));
tkGeomToken_ = esConsumes();
produces<edmNew::DetSetVector<SiStripCluster>>();
}
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -9,10 +9,11 @@
#include "DataFormats/GeometryVector/interface/GlobalPoint.h"
#include "DataFormats/GeometryVector/interface/LocalPoint.h"
#include "DataFormats/SiStripCluster/interface/SiStripApproximateCluster.h"
#include "DataFormats/SiStripCluster/interface/SiStripApproximateClusterCollection.h"
#include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
#include "DataFormats/SiStripCommon/interface/ConstantsForHardwareSystems.h"
#include "DataFormats/TrackReco/interface/Track.h"
#include "DataFormats/TrackReco/interface/TrackBase.h"
#include "DataFormats/SiStripCommon/interface/ConstantsForHardwareSystems.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/Framework/interface/stream/EDProducer.h"
Expand Down Expand Up @@ -83,13 +84,13 @@ SiStripClusters2ApproxClusters::SiStripClusters2ApproxClusters(const edm::Parame
csfToken_ = esConsumes(edm::ESInputTag("", csfLabel_));

stripNoiseToken_ = esConsumes();

produces<edmNew::DetSetVector<SiStripApproximateCluster> >();
produces<SiStripApproximateClusterCollection>();
}

void SiStripClusters2ApproxClusters::produce(edm::Event& event, edm::EventSetup const& iSetup) {
auto result = std::make_unique<edmNew::DetSetVector<SiStripApproximateCluster> >();
const auto& clusterCollection = event.get(clusterToken);
auto result = std::make_unique<SiStripApproximateClusterCollection>();
result->reserve(clusterCollection.size(), clusterCollection.dataSize());

auto const beamSpotHandle = event.getHandle(beamSpotToken_);
auto const& bs = beamSpotHandle.isValid() ? *beamSpotHandle : reco::BeamSpot();
Expand All @@ -103,7 +104,7 @@ void SiStripClusters2ApproxClusters::produce(edm::Event& event, edm::EventSetup
const auto& theNoise_ = &iSetup.getData(stripNoiseToken_);

for (const auto& detClusters : clusterCollection) {
edmNew::DetSetVector<SiStripApproximateCluster>::FastFiller ff{*result, detClusters.id()};
auto ff = result->beginDet(detClusters.id());

unsigned int detId = detClusters.id();
const GeomDet* det = tkGeom->idToDet(detId);
Expand Down