Project
Loading...
Searching...
No Matches
PrimaryVertexingSpec.cxx
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
13
14#include <vector>
15#include <TStopwatch.h>
37
38using namespace o2::framework;
41
42namespace o2
43{
44namespace vertexing
45{
46
47namespace o2d = o2::dataformats;
48
50{
51 public:
52 PrimaryVertexingSpec(std::shared_ptr<DataRequest> dr, std::shared_ptr<o2::base::GRPGeomRequest> gr, GTrackID::mask_t src, bool skip, bool validateWithIR, bool useMC)
53 : mDataRequest(dr), mGGCCDBRequest(gr), mTrackSrc(src), mSkip(skip), mUseMC(useMC), mValidateWithIR(validateWithIR) {}
54 ~PrimaryVertexingSpec() override = default;
55 void init(InitContext& ic) final;
56 void run(ProcessingContext& pc) final;
57 void endOfStream(EndOfStreamContext& ec) final;
58 void finaliseCCDB(ConcreteDataMatcher& matcher, void* obj) final;
59
60 private:
61 void updateTimeDependentParams(ProcessingContext& pc);
62 std::shared_ptr<DataRequest> mDataRequest;
63 std::shared_ptr<o2::base::GRPGeomRequest> mGGCCDBRequest;
65 GTrackID::mask_t mTrackSrc{};
66 bool mSkip{false};
67 bool mUseMC{false};
68 bool mValidateWithIR{false};
69 float mITSROFrameLengthMUS = 0.;
70 float mITSROFBiasMUS = 0.;
71 TStopwatch mTimer;
72};
73
75{
76 mTimer.Stop();
77 mTimer.Reset();
79 mVertexer.setValidateWithIR(mValidateWithIR);
80 auto dumpDir = ic.options().get<std::string>("pool-dumps-directory");
81 if (!(dumpDir.empty() || dumpDir == "/dev/null") && !o2::utils::Str::pathIsDirectory(dumpDir)) {
82 throw std::runtime_error(fmt::format("directory {} for raw data dumps does not exist", dumpDir));
83 }
84 mVertexer.setPoolDumpDirectory(dumpDir);
85 mVertexer.setTrackSources(mTrackSrc);
86}
87
89{
90 double timeCPU0 = mTimer.CpuTime(), timeReal0 = mTimer.RealTime();
91 mTimer.Start(false);
92 std::vector<PVertex> vertices;
93 std::vector<GIndex> vertexTrackIDs;
94 std::vector<V2TRef> v2tRefs;
95 std::vector<o2::MCEventLabel> lblVtx;
96
97 if (!mSkip) {
99 recoData.collectData(pc, *mDataRequest.get()); // select tracks of needed type, with minimal cuts, the real selected will be done in the vertexer
100 updateTimeDependentParams(pc); // Make sure this is called after recoData.collectData, which may load some conditions
101
102 std::vector<TrackWithTimeStamp> tracks;
103 std::vector<o2::MCCompLabel> tracksMCInfo;
104 std::vector<o2d::GlobalTrackID> gids;
105 auto maxTrackTimeError = PVertexerParams::Instance().maxTimeErrorMUS;
106 auto trackMaxX = PVertexerParams::Instance().trackMaxX;
107 auto minIBHits = PVertexerParams::Instance().minIBHits;
108 auto halfROFITS = 0.5 * mITSROFrameLengthMUS + mITSROFBiasMUS;
109 auto hw2ErrITS = 2.f / std::sqrt(12.f) * mITSROFrameLengthMUS; // conversion from half-width to error for ITS
110
111 auto creator = [maxTrackTimeError, hw2ErrITS, halfROFITS, trackMaxX, minIBHits, &tracks, &gids, &recoData](auto& _tr, GTrackID _origID, float t0, float terr) {
112 if constexpr (isBarrelTrack<decltype(_tr)>()) {
113 if (!_origID.includesDet(DetID::ITS) || _tr.getX() > trackMaxX) {
114 return true; // just in case this selection was not done on RecoContainer filling level
115 }
116 auto itsID = recoData.getITSContributorGID(_origID);
117 if ((itsID.getSource() == GTrackID::ITS && o2::math_utils::numberOfBitsSet(recoData.getITSTrack(itsID).getPattern() & 7) < minIBHits) ||
118 (itsID.getSource() == GTrackID::ITSAB && o2::math_utils::numberOfBitsSet(recoData.getITSABRef(itsID).pattern & 7) < minIBHits)) { // do not accept ITSAB tracklets
119 return true;
120 }
121 if constexpr (isITSTrack<decltype(_tr)>()) {
122 t0 += halfROFITS; // ITS time is supplied in \mus as beginning of ROF
123 terr *= hw2ErrITS; // error is supplied as a half-ROF duration, convert to \mus
124 }
125 // for all other barrel tracks the time is in \mus with gaussian error
126 if (terr < maxTrackTimeError) {
127 tracks.emplace_back(TrackWithTimeStamp{_tr, {t0, terr}});
128 gids.emplace_back(_origID);
129 }
130 }
131 return true;
132 };
133
134 recoData.createTracksVariadic(creator, mTrackSrc); // create track sample considered for vertexing
135
136 if (mUseMC) {
137 recoData.fillTrackMCLabels(gids, tracksMCInfo);
138 }
139 mVertexer.setStartIR(recoData.startIR);
140 static std::vector<InteractionCandidate> ft0Data;
141 if (mValidateWithIR) { // select BCs for validation
142 ft0Data.clear();
144 auto ft0all = recoData.getFT0RecPoints();
145 for (const auto& ftRP : ft0all) {
146 if (ft0Params.isSelected(ftRP)) {
147 ft0Data.emplace_back(InteractionCandidate{ftRP.getInteractionRecord(),
148 float(ftRP.getInteractionRecord().differenceInBC(recoData.startIR) * o2::constants::lhc::LHCBunchSpacingMUS),
149 float(ftRP.getTrigger().getAmplA() + ftRP.getTrigger().getAmplC()),
151 }
152 }
153 }
154 mVertexer.process(tracks, gids, ft0Data, vertices, vertexTrackIDs, v2tRefs, tracksMCInfo, lblVtx);
155
156 // flag vertices using UPC ITS mode
157 auto itsrofs = recoData.getITSTracksROFRecords();
158 std::vector<bool> itsTrUPC(recoData.getITSTracks().size());
159 for (auto& rof : itsrofs) {
160 if (rof.getFlag(o2::itsmft::ROFRecord::VtxUPCMode)) {
161 for (int i = rof.getFirstEntry(); i < rof.getFirstEntry() + rof.getNEntries(); i++) {
162 itsTrUPC[i] = true;
163 }
164 }
165 }
166 int nv = vertices.size();
167 for (int iv = 0; iv < nv; iv++) {
168 int idMin = v2tRefs[iv].getFirstEntry(), idMax = idMin + v2tRefs[iv].getEntries();
169 int nits = 0, nitsUPC = 0;
170 for (int id = idMin; id < idMax; id++) {
171 auto gid = recoData.getITSContributorGID(vertexTrackIDs[id]);
172 if (gid.getSource() == GIndex::ITS) {
173 nits++;
174 if (itsTrUPC[gid.getIndex()]) {
175 nitsUPC++;
176 }
177 }
178 }
179 if (nitsUPC > nits / 2) {
180 vertices[iv].setFlags(PVertex::UPCMode);
181 }
182 }
183 }
184
185 pc.outputs().snapshot(Output{"GLO", "PVTX", 0}, vertices);
186 pc.outputs().snapshot(Output{"GLO", "PVTX_CONTIDREFS"}, v2tRefs);
187 pc.outputs().snapshot(Output{"GLO", "PVTX_CONTID", 0}, vertexTrackIDs);
188
189 if (mUseMC) {
190 pc.outputs().snapshot(Output{"GLO", "PVTX_MCTR", 0}, lblVtx);
191 }
192
193 mTimer.Stop();
194 LOGP(info, "Found {} PVs, Time CPU/Real:{:.3f}/{:.3f} (DBScan: {:.4f}, Finder:{:.4f}, MADSel:{:.4f}, Rej.Debris:{:.4f}, Reattach:{:.4f}) | {} trials for {} TZ-clusters, max.trials: {}, Slowest TZ-cluster: {} ms of mult {} | NInitial:{}, Rejections: NoFilledBC:{}, NoIntCand:{}, Debris:{}, Quality:{}, ITSOnly:{}",
195 vertices.size(), mTimer.CpuTime() - timeCPU0, mTimer.RealTime() - timeReal0,
196 mVertexer.getTimeDBScan().CpuTime(), mVertexer.getTimeVertexing().CpuTime(), mVertexer.getTimeMADSel().CpuTime(), mVertexer.getTimeDebris().CpuTime(),
197 mVertexer.getTimeReAttach().CpuTime(), mVertexer.getTotTrials(), mVertexer.getNTZClusters(), mVertexer.getMaxTrialsPerCluster(),
198 mVertexer.getLongestClusterTimeMS(), mVertexer.getLongestClusterMult(), mVertexer.getNIniFound(),
199 mVertexer.getNKilledBCValid(), mVertexer.getNKilledIntCand(), mVertexer.getNKilledDebris(), mVertexer.getNKilledQuality(), mVertexer.getNKilledITSOnly());
200}
201
203{
204 mVertexer.end();
205 LOGF(info, "Primary vertexing total timing: Cpu: %.3e Real: %.3e s in %d slots",
206 mTimer.CpuTime(), mTimer.RealTime(), mTimer.Counter() - 1);
207}
208
210{
212 return;
213 }
214 // Note: strictly speaking, for Configurable params we don't need finaliseCCDB check, the singletons are updated at the CCDB fetcher level
215 if (matcher == ConcreteDataMatcher("ITS", "ALPIDEPARAM", 0)) {
216 LOG(info) << "ITS Alpide param updated";
218 par.printKeyValues();
219 return;
220 }
221 if (matcher == ConcreteDataMatcher("GLO", "MEANVERTEX", 0)) {
222 LOG(info) << "Imposing new MeanVertex: " << ((const o2::dataformats::MeanVertexObject*)obj)->asString();
224 return;
225 }
226}
227
228void PrimaryVertexingSpec::updateTimeDependentParams(ProcessingContext& pc)
229{
231 static bool initOnceDone = false;
232 if (!initOnceDone) { // this params need to be queried only once
233 initOnceDone = true;
234 // Note: reading of the ITS AlpideParam needed for ITS timing is done by the RecoContainer
237 if (!grp->isDetContinuousReadOut(DetID::ITS)) {
238 mITSROFrameLengthMUS = alpParams.roFrameLengthTrig / 1.e3; // ITS ROFrame duration in \mus
239 } else {
240 mITSROFrameLengthMUS = alpParams.roFrameLengthInBC * o2::constants::lhc::LHCBunchSpacingNS * 1e-3; // ITS ROFrame duration in \mus
241 }
242 mITSROFBiasMUS = alpParams.roFrameBiasInBC * o2::constants::lhc::LHCBunchSpacingNS * 1e-3;
243 if (o2::base::GRPGeomHelper::instance().getGRPECS()->getRunType() != o2::parameters::GRPECSObject::RunType::COSMICS) {
244 mVertexer.setBunchFilling(o2::base::GRPGeomHelper::instance().getGRPLHCIF()->getBunchFilling());
245 }
246 mVertexer.setITSROFrameLength(mITSROFrameLengthMUS);
247 mVertexer.init();
249 PVertexerParams::Instance().printKeyValues();
250 }
251 }
252 // we may have other params which need to be queried regularly
254}
255
256DataProcessorSpec getPrimaryVertexingSpec(GTrackID::mask_t src, bool skip, bool validateWithFT0, bool useMC, bool useGeom)
257{
258 std::vector<OutputSpec> outputs;
259 auto dataRequest = std::make_shared<DataRequest>();
260
261 dataRequest->requestTracks(src, useMC);
262 if (validateWithFT0 && src[GTrackID::FT0]) {
263 dataRequest->requestFT0RecPoints(false);
264 }
265
266 outputs.emplace_back("GLO", "PVTX", 0, Lifetime::Timeframe);
267 outputs.emplace_back("GLO", "PVTX_CONTID", 0, Lifetime::Timeframe);
268 outputs.emplace_back("GLO", "PVTX_CONTIDREFS", 0, Lifetime::Timeframe);
269
270 if (useMC) {
271 outputs.emplace_back("GLO", "PVTX_MCTR", 0, Lifetime::Timeframe);
272 }
273
274 auto ggRequest = std::make_shared<o2::base::GRPGeomRequest>(false, // orbitResetTime
275 true, // GRPECS=true
276 true, // GRPLHCIF
277 true, // GRPMagField
278 true, // askMatLUT
280 dataRequest->inputs,
281 true);
282 dataRequest->inputs.emplace_back("meanvtx", "GLO", "MEANVERTEX", 0, Lifetime::Condition, ccdbParamSpec("GLO/Calib/MeanVertex", {}, 1));
283
284 return DataProcessorSpec{
285 "primary-vertexing",
286 dataRequest->inputs,
287 outputs,
288 AlgorithmSpec{adaptFromTask<PrimaryVertexingSpec>(dataRequest, ggRequest, src, skip, validateWithFT0, useMC)},
289 Options{{"pool-dumps-directory", VariantType::String, "", {"Destination directory for the tracks pool dumps"}}}};
290}
291
292} // namespace vertexing
293} // namespace o2
Wrapper container for different reconstructed object types.
Definition of the GeometryManager class.
Definition of the FIT RecPoints class.
int32_t i
Helper for geometry and GRP related CCDB requests.
Global index for barrel track: provides provenance (detectors combination), index in respective array...
Definition of the Names Generator class.
Primary vertex finder.
Wrapper container for different reconstructed object types.
Result of refitting TPC-ITS matched track.
Reference on ITS/MFT clusters set.
void checkUpdates(o2::framework::ProcessingContext &pc)
static GRPGeomHelper & instance()
void setRequest(std::shared_ptr< GRPGeomRequest > req)
Static class with identifiers, bitmasks and names for ALICE detectors.
Definition DetID.h:58
static constexpr ID ITS
Definition DetID.h:63
void snapshot(const Output &spec, T const &object)
ConfigParamRegistry const & options()
Definition InitContext.h:33
decltype(auto) get(R binding, int part=0) const
DataAllocator & outputs()
The data allocator is used to allocate memory for the output data.
InputRecord & inputs()
The inputs associated with this processing context.
ServiceRegistryRef services()
The services registry associated with this processing context.
void setTrackSources(GTrackID::mask_t s)
auto getNIniFound() const
Definition PVertexer.h:127
int process(const TR &tracks, const gsl::span< o2d::GlobalTrackID > gids, const gsl::span< InteractionCandidate > intCand, std::vector< PVertex > &vertices, std::vector< o2d::VtxTrackIndex > &vertexTrackIDs, std::vector< V2TRef > &v2tRefs, const gsl::span< const o2::MCCompLabel > lblTracks, std::vector< o2::MCEventLabel > &lblVtx)
Definition PVertexer.h:333
auto getTotTrials() const
Definition PVertexer.h:118
void setStartIR(const o2::InteractionRecord &ir)
set InteractionRecods for the beginning of the TF
Definition PVertexer.h:75
auto getLongestClusterMult() const
Definition PVertexer.h:120
auto getLongestClusterTimeMS() const
Definition PVertexer.h:121
TStopwatch & getTimeDebris()
Definition PVertexer.h:131
auto getMaxTrialsPerCluster() const
Definition PVertexer.h:119
auto getNKilledITSOnly() const
Definition PVertexer.h:126
TStopwatch & getTimeVertexing()
Definition PVertexer.h:130
TStopwatch & getTimeMADSel()
Definition PVertexer.h:132
auto getNKilledQuality() const
Definition PVertexer.h:125
void setPoolDumpDirectory(const std::string &d)
Definition PVertexer.h:135
auto getNTZClusters() const
Definition PVertexer.h:117
auto getNKilledDebris() const
Definition PVertexer.h:124
void setBunchFilling(const o2::BunchFilling &bf)
void setValidateWithIR(bool v)
Definition PVertexer.h:88
TStopwatch & getTimeReAttach()
Definition PVertexer.h:133
TStopwatch & getTimeDBScan()
Definition PVertexer.h:129
auto getNKilledBCValid() const
Definition PVertexer.h:122
void setMeanVertex(const o2d::MeanVertexObject *v)
Definition PVertexer.h:96
void setITSROFrameLength(float v)
Definition PVertexer.h:106
auto getNKilledIntCand() const
Definition PVertexer.h:123
void run(ProcessingContext &pc) final
PrimaryVertexingSpec(std::shared_ptr< DataRequest > dr, std::shared_ptr< o2::base::GRPGeomRequest > gr, GTrackID::mask_t src, bool skip, bool validateWithIR, bool useMC)
~PrimaryVertexingSpec() override=default
void finaliseCCDB(ConcreteDataMatcher &matcher, void *obj) final
void endOfStream(EndOfStreamContext &ec) final
This is invoked whenever we have an EndOfStream event.
GLenum src
Definition glcorearb.h:1767
GLuint GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat t0
Definition glcorearb.h:5034
constexpr double LHCBunchSpacingMUS
constexpr double LHCBunchSpacingNS
Definition of a container to keep/associate and arbitrary number of labels associated to an index wit...
Defining PrimaryVertex explicitly as messageable.
Definition TFIDInfo.h:20
std::vector< ConfigParamSpec > ccdbParamSpec(std::string const &path, int runDependent, std::vector< CCDBMetadata > metadata={}, int qrate=0)
std::vector< ConfigParamSpec > Options
o2::framework::DataProcessorSpec getPrimaryVertexingSpec(o2::dataformats::GlobalTrackID::mask_t src, bool skip, bool validateWithFT0, bool useMC, bool useGeom)
create a processor spec
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
size_t inputTimesliceId
The time pipelining id of this particular device.
Definition DeviceSpec.h:68
bool isSelected(const RecPoints &rp) const
GTrackID getITSContributorGID(GTrackID source) const
void createTracksVariadic(T creator, GTrackID::mask_t srcSel=GTrackID::getSourcesMask("all")) const
const o2::itsmft::TrkClusRef & getITSABRef(GTrackID gid) const
void collectData(o2::framework::ProcessingContext &pc, const DataRequest &request)
void fillTrackMCLabels(const gsl::span< GTrackID > gids, std::vector< o2::MCCompLabel > &mcinfo) const
const o2::its::TrackITS & getITSTrack(GTrackID gid) const
uint16_t pattern
layers pattern
Definition TrkClusRef.h:29
static bool pathIsDirectory(const std::string_view p)
generic track with timestamp
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"