Project
Loading...
Searching...
No Matches
VertexTrackMatcherSpec.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
15
27#include "TStopwatch.h"
28
29using namespace o2::framework;
33
34namespace o2
35{
36namespace vertexing
37{
38
40{
41 public:
42 VertexTrackMatcherSpec(std::shared_ptr<DataRequest> dr, std::shared_ptr<o2::base::GRPGeomRequest> gr) : mDataRequest(dr), mGGCCDBRequest(gr){};
43 ~VertexTrackMatcherSpec() override = default;
44 void init(InitContext& ic) final;
45 void run(ProcessingContext& pc) final;
46 void endOfStream(EndOfStreamContext& ec) final;
47 void finaliseCCDB(ConcreteDataMatcher& matcher, void* obj) final;
48
49 private:
50 void updateTimeDependentParams(ProcessingContext& pc);
51 std::shared_ptr<DataRequest> mDataRequest;
52 std::shared_ptr<o2::base::GRPGeomRequest> mGGCCDBRequest;
53 o2::tpc::VDriftHelper mTPCVDriftHelper{};
55 TStopwatch mTimer;
56};
57
59{
60 //-------- init geometry and field --------//
61 mTimer.Stop();
62 mTimer.Reset();
63 mMatcher.setPrescaleLogs(ic.options().get<int>("prescale-logs"));
65}
66
68{
69 double timeCPU0 = mTimer.CpuTime(), timeReal0 = mTimer.RealTime();
70 mTimer.Start(false);
71
73 recoData.collectData(pc, *mDataRequest.get());
74 updateTimeDependentParams(pc); // make sure called after recoData.collectData as some objects might be fetched there
75
76 std::vector<o2::dataformats::VtxTrackIndex> trackIndex;
77 std::vector<o2::dataformats::VtxTrackRef> vtxRefs;
78
79 mMatcher.process(recoData, trackIndex, vtxRefs);
80
81 pc.outputs().snapshot(Output{"GLO", "PVTX_TRMTC", 0}, trackIndex);
82 pc.outputs().snapshot(Output{"GLO", "PVTX_TRMTCREFS", 0}, vtxRefs);
83
84 mTimer.Stop();
85 LOG(info) << "Made " << trackIndex.size() << " track associations for " << recoData.getPrimaryVertices().size()
86 << " vertices, timing: CPU: " << mTimer.CpuTime() - timeCPU0 << " Real: " << mTimer.RealTime() - timeReal0 << " s";
87}
88
89void VertexTrackMatcherSpec::updateTimeDependentParams(ProcessingContext& pc)
90{
92 mTPCVDriftHelper.extractCCDBInputs(pc);
93 static bool initOnceDone = false;
94 if (!initOnceDone) { // this params need to be queried only once
95 initOnceDone = true;
96 // put here init-once stuff
98 mMatcher.setITSROFrameLengthMUS(o2::base::GRPGeomHelper::instance().getGRPECS()->isDetContinuousReadOut(o2::detectors::DetID::ITS) ? alpParamsITS.roFrameLengthInBC * o2::constants::lhc::LHCBunchSpacingMUS : alpParamsITS.roFrameLengthTrig * 1.e-3);
100 mMatcher.setMFTROFrameLengthMUS(o2::base::GRPGeomHelper::instance().getGRPECS()->isDetContinuousReadOut(o2::detectors::DetID::MFT) ? alpParamsMFT.roFrameLengthInBC * o2::constants::lhc::LHCBunchSpacingMUS : alpParamsMFT.roFrameLengthTrig * 1.e-3);
101 mMatcher.init();
102 LOGP(info, "VertexTrackMatcher ITSROFrameLengthMUS:{} MFTROFrameLengthMUS:{}", mMatcher.getITSROFrameLengthMUS(), mMatcher.getMFTROFrameLengthMUS());
103 }
104 // we may have other params which need to be queried regularly
105 // VDrift may change from time to time
106 if (mTPCVDriftHelper.isUpdated()) {
108 auto& detParam = o2::tpc::ParameterDetector::Instance();
109 mMatcher.setTPCBin2MUS(elParam.ZbinWidth);
110 auto& vd = mTPCVDriftHelper.getVDriftObject();
111 mMatcher.setMaxTPCDriftTimeMUS(detParam.TPClength / (vd.refVDrift * vd.corrFact));
112 mMatcher.setTPCTDriftOffset(vd.getTimeOffset());
113 LOGP(info, "Updating TPC fast transform map with new VDrift factor of {} wrt reference {} and DriftTimeOffset correction {} wrt {} from source {}",
114 vd.corrFact, vd.refVDrift, vd.timeOffsetCorr, vd.refTimeOffset, mTPCVDriftHelper.getSourceName());
115 mTPCVDriftHelper.acknowledgeUpdate();
116 }
117}
118
120{
122 return;
123 }
124 if (mTPCVDriftHelper.accountCCDBInputs(matcher, obj)) {
125 return;
126 }
127 if (matcher == ConcreteDataMatcher("ITS", "ALPIDEPARAM", 0)) {
128 LOG(info) << "ITS Alpide param updated";
130 par.printKeyValues();
131 return;
132 }
133 if (matcher == ConcreteDataMatcher("MFT", "ALPIDEPARAM", 0)) {
134 LOG(info) << "MFT Alpide param updated";
136 par.printKeyValues();
137 return;
138 }
139}
140
142{
143 LOGF(info, "Primary vertex - track matching total timing: Cpu: %.3e Real: %.3e s in %d slots",
144 mTimer.CpuTime(), mTimer.RealTime(), mTimer.Counter() - 1);
145}
146
148{
149 std::vector<OutputSpec> outputs;
150 auto dataRequest = std::make_shared<DataRequest>();
151
152 dataRequest->requestTracks(src, false);
153 dataRequest->requestClusters(src & GTrackID::getSourcesMask("EMC,PHS,CPV"), false);
154 dataRequest->requestPrimaryVerticesTMP(false);
155
156 auto ggRequest = std::make_shared<o2::base::GRPGeomRequest>(false, // orbitResetTime
157 true, // GRPECS=true
158 true, // GRPLHCIF
159 false, // GRPMagField
160 false, // askMatLUT
162 dataRequest->inputs,
163 true);
164 o2::tpc::VDriftHelper::requestCCDBInputs(dataRequest->inputs);
165
166 outputs.emplace_back("GLO", "PVTX_TRMTC", 0, Lifetime::Timeframe);
167 outputs.emplace_back("GLO", "PVTX_TRMTCREFS", 0, Lifetime::Timeframe);
168
169 return DataProcessorSpec{
170 "pvertex-track-matching",
171 dataRequest->inputs,
172 outputs,
173 AlgorithmSpec{adaptFromTask<VertexTrackMatcherSpec>(dataRequest, ggRequest)},
174 Options{{"prescale-logs", VariantType::Int, 50, {"print vertex logs for each n-th TF"}}}};
175}
176
177} // namespace vertexing
178} // namespace o2
Wrapper container for different reconstructed object types.
Helper for geometry and GRP related CCDB requests.
Definition of the Names Generator class.
Definition of the parameter class for the detector.
Definition of the parameter class for the detector electronics.
Helper class to extract VDrift from different sources.
Specs for vertex track association device.
Class for vertex track association.
void checkUpdates(o2::framework::ProcessingContext &pc)
static GRPGeomHelper & instance()
void setRequest(std::shared_ptr< GRPGeomRequest > req)
static mask_t getSourcesMask(const std::string_view srcList)
Static class with identifiers, bitmasks and names for ALICE detectors.
Definition DetID.h:58
static constexpr ID ITS
Definition DetID.h:63
static constexpr ID MFT
Definition DetID.h:71
void snapshot(const Output &spec, T const &object)
ConfigParamRegistry const & options()
Definition InitContext.h:33
DataAllocator & outputs()
The data allocator is used to allocate memory for the output data.
static void requestCCDBInputs(std::vector< o2::framework::InputSpec > &inputs, bool laser=true, bool itstpcTgl=true)
void extractCCDBInputs(o2::framework::ProcessingContext &pc, bool laser=true, bool itstpcTgl=true)
const VDriftCorrFact & getVDriftObject() const
bool accountCCDBInputs(const o2::framework::ConcreteDataMatcher &matcher, void *obj)
static std::string_view getSourceName(Source s)
bool isUpdated() const
VertexTrackMatcherSpec(std::shared_ptr< DataRequest > dr, std::shared_ptr< o2::base::GRPGeomRequest > gr)
void endOfStream(EndOfStreamContext &ec) final
This is invoked whenever we have an EndOfStream event.
void finaliseCCDB(ConcreteDataMatcher &matcher, void *obj) final
void run(ProcessingContext &pc) final
~VertexTrackMatcherSpec() override=default
void process(const o2::globaltracking::RecoContainer &recoData, std::vector< VTIndex > &trackIndex, std::vector< VRef > &vtxRefs)
GLenum src
Definition glcorearb.h:1767
constexpr double LHCBunchSpacingMUS
Defining ITS Vertex explicitly as messageable.
Definition Cartesian.h:288
std::vector< ConfigParamSpec > Options
o2::framework::DataProcessorSpec getVertexTrackMatcherSpec(o2::dataformats::GlobalTrackID::mask_t src)
create a processor spec
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
void collectData(o2::framework::ProcessingContext &pc, const DataRequest &request)
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"