Project
Loading...
Searching...
No Matches
PulseHeight.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
14
18#include <fairlogger/Logger.h>
19
20using namespace o2::trd;
21using namespace o2::trd::constants;
22
24{
25 mPHValues.clear();
26 mPHValuesHD.clear();
27 mDistances.clear();
28}
29
31{
32}
33
35{
36 mOutFile = std::make_unique<TFile>("trd_PH.root", "RECREATE");
37 if (mOutFile->IsZombie()) {
38 LOG(error) << "Failed to create output file!";
39 return;
40 }
41 mOutTree = std::make_unique<TTree>("ph", "Data points for PH histograms");
42 mOutTree->Branch("values", &mPHValuesPtr);
43 mOutTree->Branch("valuesHD", &mPHValuesHDPtr);
44 mOutTree->Branch("dist", &mDistancesPtr);
45 mWriteOutput = true;
46 LOG(info) << "Writing PH data points to local file trd_PH.root";
47}
48
50{
51 if (!mWriteOutput) {
52 return;
53 }
54 try {
55 mOutFile->cd();
56 mOutTree->Write();
57 mOutTree.reset();
58 mOutFile->Close();
59 mOutFile.reset();
60 } catch (std::exception const& e) {
61 LOG(error) << "Failed to write output file, reason: " << e.what();
62 }
63 mWriteOutput = false;
64}
65
66void PulseHeight::setInput(const o2::globaltracking::RecoContainer& input, gsl::span<const Digit>* digits)
67{
68 mTracksInITSTPCTRD = input.getITSTPCTRDTracks<TrackTRD>();
69 mTracksInTPCTRD = input.getTPCTRDTracks<TrackTRD>();
70 mTrackletsRaw = input.getTRDTracklets();
71 mTriggerRecords = input.getTRDTriggerRecords();
72 mTrackTriggerRecordsTPCTRD = input.getTPCTRDTriggers();
73 mTrackTriggerRecordsITSTPCTRD = input.getITSTPCTRDTriggers();
74 mDigits = digits;
75}
76
78{
79 LOGP(debug, "Processing {} tracklets and {} digits from {} triggers and {} ITS-TPC-TRD tracks and {} TPC-TRD tracks",
80 mTrackletsRaw.size(), mDigits->size(), mTriggerRecords.size(), mTracksInITSTPCTRD.size(), mTracksInTPCTRD.size());
81 if (mDigits->size() == 0) {
82 return;
83 }
84 size_t nTrkTrig[2] = {mTrackTriggerRecordsITSTPCTRD.size(), mTrackTriggerRecordsTPCTRD.size()};
85 int lastTrkTrig[2] = {0};
86 for (const auto& trig : mTriggerRecords) {
87 if (trig.getNumberOfDigits() == 0) {
88 continue;
89 }
90 for (int iTrackType = 0; iTrackType <= 1; ++iTrackType) {
91 const gsl::span<const TrackTRD>* tracks = (iTrackType == 0) ? &mTracksInITSTPCTRD : &mTracksInTPCTRD;
92 const gsl::span<const TrackTriggerRecord>* trackTriggers = (iTrackType == 0) ? &mTrackTriggerRecordsITSTPCTRD : &mTrackTriggerRecordsTPCTRD;
93 for (int iTrigTrack = lastTrkTrig[iTrackType]; iTrigTrack < nTrkTrig[iTrackType]; ++iTrigTrack) {
94 const auto& trigTrack = (*trackTriggers)[iTrigTrack];
95 if (trigTrack.getBCData().differenceInBC(trig.getBCData()) > 0) {
96 // aborting, since track trigger is later than digit trigger";
97 LOGP(debug, "Aborting, track trigger is too late");
98 break;
99 }
100 if (trigTrack.getBCData() != trig.getBCData()) {
101 // skipping, since track trigger earlier than digit trigger";
102 LOGP(debug, "Skipping, track trigger is too late");
103 ++lastTrkTrig[iTrackType];
104 continue;
105 }
106 if ((*tracks)[trigTrack.getFirstTrack()].hasPileUpInfo() && (*tracks)[trigTrack.getFirstTrack()].getPileUpTimeShiftMUS() < mParams.pileupCut) {
107 // rejecting triggers which are close to other collisions (avoid pile-up)
108 ++lastTrkTrig[iTrackType];
109 LOGP(debug, "Rececting trigger due to pileup of {}", (*tracks)[trigTrack.getFirstTrack()].getPileUpTimeShiftMUS());
110 break;
111 }
112 for (int iTrk = trigTrack.getFirstTrack(); iTrk < trigTrack.getFirstTrack() + trigTrack.getNumberOfTracks(); iTrk++) {
113 const auto& trk = (*tracks)[iTrk];
114 for (int iLayer = 0; iLayer < 6; iLayer++) {
115 int trkltIdx = trk.getTrackletIndex(iLayer);
116 if (trkltIdx < 0) {
117 continue;
118 }
119 findDigitsForTracklet(mTrackletsRaw[trkltIdx], trig, iTrackType);
120 }
121 }
122 }
123 }
124 }
125 if (mWriteOutput) {
126 mOutTree->Fill();
127 }
128}
129
131{
132 auto trkltDet = trklt.getDetector();
133 int iDigitFirst = trig.getFirstDigit();
134 int iDigitLast = trig.getFirstDigit() + trig.getNumberOfDigits();
135 for (int iDigit = iDigitFirst + 1; iDigit < iDigitLast - 1; ++iDigit) {
136 const auto& digit = (*mDigits)[iDigit];
137 if (digit.isSharedDigit()) {
138 continue; // avoid double-counting of the same digit
139 }
140 if (digit.getDetector() != trkltDet || digit.getPadRow() != trklt.getPadRow() || digit.getPadCol() != trklt.getPadCol(mApplyShift)) {
141 // for now we loose charge information from padrow-crossing tracklets (~15% of all tracklets)
142 continue;
143 }
144 int nNeighbours = 0;
145 bool left = false;
146 bool right = false;
147 const auto* digitLeft = &(*mDigits)[iDigit + 1];
148 const auto* digitRight = &(*mDigits)[iDigit - 1];
149 // due to shared digits the neighbouring element might still be in the same pad column
150 if (digitLeft->getPadCol() == digit.getPadCol() && iDigit < iDigitLast - 2) {
151 digitLeft = &(*mDigits)[iDigit + 2];
152 }
153 if (digitRight->getPadCol() == digit.getPadCol() && iDigit > iDigitFirst + 2) {
154 digitRight = &(*mDigits)[iDigit - 2];
155 }
156 LOG(debug) << "Central digit: " << digit;
157 LOG(debug) << "Left digit: " << *digitLeft;
158 LOG(debug) << "Right digit: " << *digitRight;
159 int direction = 0;
160 if (digitLeft->isNeighbour(digit)) {
161 ++nNeighbours;
162 left = true;
163 direction = digit.getPadCol() - digitLeft->getPadCol();
164 }
165 if (digitRight->isNeighbour(digit) && digit.getPadCol() - digitRight->getPadCol() != direction) {
166 ++nNeighbours;
167 right = true;
168 }
169 if (nNeighbours > 0) {
170 int digitTrackletDistance = 0;
171 auto adcSumMax = digit.getADCsum();
172 if (left && digitLeft->getADCsum() > adcSumMax) {
173 adcSumMax = digitLeft->getADCsum();
174 digitTrackletDistance = 1;
175 }
176 if (right && digitRight->getADCsum() > adcSumMax) {
177 adcSumMax = digitRight->getADCsum();
178 digitTrackletDistance = -1;
179 }
180 mDistances.push_back(digitTrackletDistance);
181 for (int iTb = 0; iTb < TIMEBINS; ++iTb) {
182 uint16_t phVal = digit.getADC()[iTb];
183 uint16_t phValHD = (digit.getADC()[iTb] << 2) + digit.getPreTrigPhase();
184 if (left) {
185 phVal += digitLeft->getADC()[iTb];
186 phValHD += (digitLeft->getADC()[iTb] << 2) + digit.getPreTrigPhase();
187 }
188 if (right) {
189 phVal += digitRight->getADC()[iTb];
190 phValHD += (digitRight->getADC()[iTb] << 2) + digit.getPreTrigPhase();
191 }
192 mPHValues.emplace_back(phVal, trkltDet, iTb, nNeighbours, type);
193 mPHValuesHD.emplace_back(phValHD, trkltDet, iTb, nNeighbours, type);
194 }
195 }
196 }
197}
Wrapper container for different reconstructed object types.
Global TRD definitions and constants.
std::ostringstream debug
Creates PH spectra from TRD digits found on tracks.
void findDigitsForTracklet(const Tracklet64 &trklt, const TriggerRecord &trig, int type)
Search for up to 3 neighbouring digits matching given tracklet and fill the PH values.
void setInput(const o2::globaltracking::RecoContainer &input, gsl::span< const Digit > *digits)
Initialize the input arrays.
void reset()
Reset the output.
void process()
Main processing function.
void init()
Initialize what is needed.
Header for data corresponding to the same hardware trigger adapted from DataFormatsITSMFT/ROFRecord.
int getNumberOfDigits() const
GLdouble GLdouble right
Definition glcorearb.h:4077
GLint GLint GLsizei GLint GLenum GLenum type
Definition glcorearb.h:275
constexpr int TIMEBINS
the number of time bins
Definition Constants.h:74
gsl::span< const o2::trd::TriggerRecord > getTRDTriggerRecords() const
gsl::span< const o2::trd::Tracklet64 > getTRDTracklets() const
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
std::vector< Digit > digits