Project
Loading...
Searching...
No Matches
PHOSTurnonCalibrator.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#include "PHOSBase/Geometry.h"
18#include "CCDB/CcdbApi.h"
19#include "CCDB/CcdbObjectInfo.h"
20
21#include "TF1.h"
22#include "TH1.h"
23#include "TGraphAsymmErrors.h"
24
25#include <fairlogger/Logger.h>
26#include <fstream> // std::ifstream
27
28using namespace o2::phos;
29
30PHOSTurnonSlot::PHOSTurnonSlot(bool useCCDB) : mUseCCDB(useCCDB)
31{
32 mFiredTiles.reset();
33 mNoisyTiles.reset();
34 mTurnOnHistos = std::make_unique<TurnOnHistos>();
35}
37{
38 mUseCCDB = other.mUseCCDB;
39 mRunStartTime = other.mUseCCDB;
40 mFiredTiles.reset();
41 mNoisyTiles.reset();
42 mTurnOnHistos = std::make_unique<TurnOnHistos>();
43}
44
46{
47 for (short ddl = 0; ddl < 14; ddl++) {
48 const std::array<float, TurnOnHistos::Npt>& all = mTurnOnHistos->getTotSpectrum(ddl);
49 const std::array<float, TurnOnHistos::Npt>& tr = mTurnOnHistos->getTrSpectrum(ddl);
50 float sumAll = 0, sumTr = 0.;
51 for (int i = 0; i < TurnOnHistos::Npt; i++) {
52 sumAll += all[i];
53 sumTr += tr[i];
54 }
55 LOG(info) << "DDL " << ddl << " total entries " << sumAll << " trigger clusters " << sumTr;
56 }
57}
58void PHOSTurnonSlot::fill(const gsl::span<const Cell>& cells, const gsl::span<const TriggerRecord>& cellTR,
59 const gsl::span<const Cluster>& clusters, const gsl::span<const TriggerRecord>& cluTR)
60{
61
62 auto ctr = cellTR.begin();
63 auto clutr = cluTR.begin();
64 while (ctr != cellTR.end() && clutr != cluTR.end()) {
65 //
66 // TODO! select NOT PHOS triggered events
67 // DataProcessingHeader::
68 //
69 if (ctr->getBCData() == clutr->getBCData()) {
70 scanClusters(cells, *ctr, clusters, *clutr);
71 } else {
72 LOG(error) << "Different TrigRecords for cells:" << ctr->getBCData() << " and clusters:" << clutr->getBCData();
73 }
74 ctr++;
75 clutr++;
76 }
77}
79{
80 mFiredTiles.reset();
81 mNoisyTiles.reset();
82 mTurnOnHistos.reset();
83}
84void PHOSTurnonSlot::scanClusters(const gsl::span<const Cell>& cells, const TriggerRecord& celltr,
85 const gsl::span<const Cluster>& clusters, const TriggerRecord& clutr)
86{
87 // First fill map of expected tiles from TRU cells
88 mNoisyTiles.reset();
89 int firstCellInEvent = celltr.getFirstEntry();
90 int lastCellInEvent = firstCellInEvent + celltr.getNumberOfObjects();
91 for (int i = firstCellInEvent; i < lastCellInEvent; i++) {
92 const Cell& c = cells[i];
93 if (c.getTRU()) {
94 mNoisyTiles.set(c.getTRUId() - Geometry::getTotalNCells() - 1);
95 }
96 }
97
98 // Copy to have good and noisy map
99 mFiredTiles.reset();
100 char mod;
101 float x, z;
102 short ddl;
103 int firstCluInEvent = clutr.getFirstEntry();
104 int lastCluInEvent = firstCluInEvent + clutr.getNumberOfObjects();
105 for (int i = firstCluInEvent; i < lastCluInEvent; i++) {
106 const Cluster& clu = clusters[i];
107 if (clu.getEnergy() < 1.e-4) {
108 continue;
109 }
110 mod = clu.module();
112 // TODO: do we need separate 2x2 and 4x4 spectra? Switch?
113 // short truId2x2 = Geometry::relPosToTruId(mod, x, z, 0);
114 // short truId4x4 = Geometry::relPosToTruId(mod, x, z, 1);
115 mTurnOnHistos->fillTotSp(ddl, clu.getEnergy());
116 // if (clu.firedTrigger() & 1) { //Bit 1: 2x2, bit 2 4x4
117 // mTurnOnHistos->fillFiredSp(ddl, clu.getEnergy());
118 // //Fill trigger map
119 // mFiredTiles.set(truId);
120 // }
121 // }
122 }
123 // Fill final good and noisy maps
124 mTurnOnHistos->fillFiredMap(mFiredTiles);
125 mNoisyTiles ^= mFiredTiles;
126 mTurnOnHistos->fillNoisyMap(mFiredTiles);
127}
128//==============================================
129
131{
132 // Extract results for the single slot
133 // if not ready yet, prepare containers
134 if (!mTurnOnHistos) {
135 mTurnOnHistos.reset(new TurnOnHistos());
136 }
137 PHOSTurnonSlot* c = slot.getContainer();
138 LOG(info) << "Finalize slot " << slot.getTFStart() << " <= TF <= " << slot.getTFEnd();
139 // Add histos
140 for (int mod = 0; mod < 8; mod++) {
141 mTurnOnHistos->merge(c->getCollectedHistos());
142 }
143 c->clear();
144}
146{
147
148 auto& cont = getSlots();
149 auto& slot = front ? cont.emplace_front(tstart, tend) : cont.emplace_back(tstart, tend);
150 slot.setContainer(std::make_unique<PHOSTurnonSlot>(mUseCCDB));
151 return slot;
152}
153
154bool PHOSTurnonCalibrator::process(uint64_t tf, const gsl::span<const Cell>& cells, const gsl::span<const TriggerRecord>& cellTR,
155 const gsl::span<const Cluster>& clusters, const gsl::span<const TriggerRecord>& cluTR)
156{
157 // process current TF
158 auto& slotTF = getSlotForTF(tf);
159 slotTF.getContainer()->setRunStartTime(tf);
160 slotTF.getContainer()->fill(cells, cellTR, clusters, cluTR);
161
162 return true;
163}
164
166{
167 // Use stored histos to calculate maps and turn-on curves
168 // return true of successful
169
170 // extract TOC
171 if (!mTriggerMap) {
172 mTriggerMap.reset(new TriggerMap());
173 }
174 TF1* th = new TF1("aTh", "[0]/(TMath::Exp(([1]-x)/[2])+1.)+(1.-[0])/(TMath::Exp(([3]-x)/[2])+1.)", 0., 40.);
175 std::array<std::array<float, 10>, TurnOnHistos::NDDL> params;
176 for (int ddl = 0; ddl < TurnOnHistos::NDDL; ddl++) {
177 TH1F hF("fired", "fired", 200, 0., 20.);
178 TH1F hA("all", "all", 200, 0., 20.);
179 const std::array<float, TurnOnHistos::Npt>& vf = mTurnOnHistos->getTrSpectrum(ddl);
180 const std::array<float, TurnOnHistos::Npt>& va = mTurnOnHistos->getTotSpectrum(ddl);
181 for (int i = 0; i < 200; i++) {
182 hF.SetBinContent(i + 1, vf[i]);
183 hA.SetBinContent(i + 1, va[i]);
184 }
185 hF.Sumw2();
186 hA.Sumw2();
187
188 TGraphAsymmErrors* gr = new TGraphAsymmErrors(&hF, &hA);
189 th->SetParameters(0.9, 3.5, 0.3, 7.5, 0.6);
190 gr->Fit(th, "Q", "", 2., 20.);
191 gr->SetName(Form("DDL_%d", ddl));
192 gr->SetTitle(Form("DDL %d", ddl));
193 // TODO!!! Add TGraph with fit to list of objects to send to QC
194 double* par = th->GetParameters();
195 for (int i = 0; i < 10; i++) {
196 params[ddl][i] = par[i];
197 }
198 }
199 std::string_view versionName{"default"};
200 mTriggerMap->addTurnOnCurvesParams(versionName, params);
201 // TODO: calculate bad map
202 // and fill object
203 // mTriggerMap->addBad2x2Channel(short cellID) ;
204}
Utils and constants for calibration and related workflows.
int32_t i
Device to calculate PHOS turn-on curve and bad map.
uint32_t c
Definition RawData.h:2
TFType getTFEnd() const
Definition TimeSlot.h:47
const Container * getContainer() const
Definition TimeSlot.h:53
TFType getTFStart() const
Definition TimeSlot.h:46
Contains PHOS cluster parameters.
Definition Cluster.h:39
float getEnergy() const
Definition Cluster.h:56
char module() const
Definition Cluster.h:92
void getLocalPosition(float &posX, float &posZ) const
Definition Cluster.h:76
static int getTotalNCells()
Definition Geometry.h:126
bool process(uint64_t tf, const gsl::span< const Cell > &cells, const gsl::span< const TriggerRecord > &trs, const gsl::span< const Cluster > &clusters, const gsl::span< const TriggerRecord > &cluTR)
Slot & emplaceNewSlot(bool front, TFType tstart, TFType tend) final
void fill(const gsl::span< const Cell > &cells, const gsl::span< const TriggerRecord > &trs, const gsl::span< const Cluster > &clusters, const gsl::span< const TriggerRecord > &cluTR)
Header for data corresponding to the same hardware trigger adapted from DataFormatsEMCAL/TriggerRecor...
int getNumberOfObjects() const
static constexpr short Npt
Number of bins in pt distribution.
static constexpr short NDDL
Number of DDLs.
GLint GLenum GLint x
Definition glcorearb.h:403
GLenum const GLfloat * params
Definition glcorearb.h:272
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
std::unique_ptr< GPUReconstructionTimeframe > tf
VectorOfTObjectPtrs other
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
Cluster clu
std::vector< Cluster > clusters
std::vector< Cell > cells