forked from AliceO2Group/AliceO2
-
Notifications
You must be signed in to change notification settings - Fork 0
Expand file tree
/
Copy pathCheckChipResponseFile.C
More file actions
206 lines (184 loc) · 6.85 KB
/
CheckChipResponseFile.C
File metadata and controls
206 lines (184 loc) · 6.85 KB
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
// All rights not expressly granted are reserved.
//
// This software is distributed under the terms of the GNU General Public
// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
//
// In applying this license CERN does not waive the privileges and immunities
// granted to it by virtue of its status as an Intergovernmental Organization
// or submit itself to any jurisdiction.
/// \file CheckChipResponseFile.C
/// \brief Simple macro to check the chip response files
#if !defined(__CLING__) || defined(__ROOTCLING__)
#include <TFile.h>
#include <TGraph.h>
#include <TCanvas.h>
#include <TH1.h>
#include <TLegend.h>
#include <iostream>
#include <vector>
#include <string>
#define ENABLE_UPGRADES
#include "ITSMFTSimulation/AlpideSimResponse.h"
#include "ITS3Simulation/ChipSimResponse.h"
#include "ITS3Base/SegmentationMosaix.h"
#include "fairlogger/Logger.h"
#endif
using SegmentationMosaix = o2::its3::SegmentationMosaix;
double um2cm(double um) { return um * 1e-4; }
double cm2um(double cm) { return cm * 1e+4; }
std::unique_ptr<o2::its3::ChipSimResponse> mAlpSimResp0, mAlpSimResp1, mAptSimResp1;
std::unique_ptr<o2::its3::ChipSimResponse> loadResponse(const std::string& fileName, const std::string& respName)
{
TFile* f = TFile::Open(fileName.data());
if (!f) {
std::cerr << fileName << " not found" << std::endl;
return nullptr;
}
auto base = f->Get<o2::itsmft::AlpideSimResponse>(respName.c_str());
if (!base) {
std::cerr << respName << " not found in " << fileName << std::endl;
return nullptr;
}
return std::make_unique<o2::its3::ChipSimResponse>(base);
}
void LoadRespFunc()
{
std::string AptsFile = "$(O2_ROOT)/share/Detectors/Upgrades/ITS3/data/ITS3ChipResponseData/APTSResponseData.root";
std::string AlpideFile = "$(O2_ROOT)/share/Detectors/ITSMFT/data/AlpideResponseData/AlpideResponseData.root";
std::cout << "=====================\n";
LOGP(info, "ALPIDE Vbb=0V response");
mAlpSimResp0 = loadResponse(AlpideFile, "response0"); // Vbb=0V
mAlpSimResp0->computeCentreFromData();
mAlpSimResp0->print();
LOGP(info, "Response Centre {}", mAlpSimResp0->getRespCentreDep());
std::cout << "=====================\n";
LOGP(info, "ALPIDE Vbb=-3V response");
mAlpSimResp1 = loadResponse(AlpideFile, "response1"); // Vbb=-3V
mAlpSimResp1->computeCentreFromData();
mAlpSimResp1->print();
LOGP(info, "Response Centre {}", mAlpSimResp1->getRespCentreDep());
std::cout << "=====================\n";
LOGP(info, "APTS response");
mAptSimResp1 = loadResponse(AptsFile, "response1"); // APTS
mAptSimResp1->computeCentreFromData();
mAptSimResp1->print();
LOGP(info, "Response Centre {}", mAptSimResp1->getRespCentreDep());
std::cout << "=====================\n";
}
std::vector<float> getCollectionSeediciencies(o2::its3::ChipSimResponse* resp,
const std::vector<float>& depths)
{
std::vector<float> seed;
bool flipRow = false, flipCol = false;
for (auto depth : depths) {
auto rspmat = resp->getResponse(0.0, 0.0,
um2cm(depth) + 1.e-9,
flipRow, flipCol);
seed.push_back(rspmat ? rspmat->getValue(2, 2) : 0.f);
}
return seed;
}
std::vector<float> getShareValues(o2::its3::ChipSimResponse* resp,
const std::vector<float>& depths)
{
std::vector<float> share;
bool flipRow = false, flipCol = false;
for (auto depth : depths) {
auto rspmat = resp->getResponse(0.0, 0.0,
um2cm(depth) + 1.e-9,
flipRow, flipCol);
float s = 0;
int npix = resp->getNPix();
if (rspmat) {
for (int i = 0; i < npix; ++i)
for (int j = 0; j < npix; ++j)
if (!(i == npix / 2 && j == npix / 2))
s += rspmat->getValue(i, j);
}
share.push_back(s);
}
return share;
}
std::vector<float> getEffValues(o2::its3::ChipSimResponse* resp,
const std::vector<float>& depths)
{
std::vector<float> all;
bool flipRow = false, flipCol = false;
for (auto depth : depths) {
auto rspmat = resp->getResponse(0.0, 0.0,
um2cm(depth) + 1.e-9,
flipRow, flipCol);
float s = 0;
int npix = resp->getNPix();
if (rspmat) {
for (int i = 0; i < npix; ++i)
for (int j = 0; j < npix; ++j)
s += rspmat->getValue(i, j);
}
all.push_back(s);
}
return all;
}
void CheckChipResponseFile()
{
LoadRespFunc();
LOG(info) << "Response function loaded" << std::endl;
std::vector<float> vecDepth;
int numPoints = 100;
for (int i = 0; i < numPoints; ++i) {
float value = -50 + i * (100.0f / (numPoints - 1));
vecDepth.push_back(value);
}
int colors[] = {kOrange + 7, kRed + 1, kAzure + 4};
struct RespInfo {
std::unique_ptr<o2::its3::ChipSimResponse>& resp;
std::string title;
int color;
};
std::vector<RespInfo> responses = {
{mAptSimResp1, "APTS", colors[0]},
{mAlpSimResp0, "ALPIDE Vbb=0V", colors[1]},
{mAlpSimResp1, "ALPIDE Vbb=-3V", colors[2]}};
TCanvas* c1 = new TCanvas("c1", "c1", 800, 600);
TH1* frame = c1->DrawFrame(-50, -0.049, 50, 1.049);
frame->SetTitle(";Depth(um);Charge Collection Seed / Share / Eff");
TLegend* leg = new TLegend(0.15, 0.5, 0.4, 0.85);
leg->SetFillStyle(0);
leg->SetBorderSize(0);
for (auto& r : responses) {
if (!r.resp)
continue;
auto seed = getCollectionSeediciencies(r.resp.get(), vecDepth);
auto shr = getShareValues(r.resp.get(), vecDepth);
auto all = getEffValues(r.resp.get(), vecDepth);
auto grSeed = new TGraph(vecDepth.size(), vecDepth.data(), seed.data());
grSeed->SetTitle(Form("%s seed", r.title.c_str()));
grSeed->SetLineColor(r.color);
grSeed->SetLineWidth(2);
grSeed->SetMarkerColor(r.color);
grSeed->SetMarkerStyle(kFullCircle);
grSeed->SetMarkerSize(0.8);
grSeed->Draw("SAME LP");
leg->AddEntry(grSeed, Form("%s seed", r.title.c_str()), "lp");
auto grShare = new TGraph(vecDepth.size(), vecDepth.data(), shr.data());
grShare->SetLineColor(r.color);
grShare->SetLineWidth(2);
grShare->SetMarkerColor(r.color);
grShare->SetMarkerStyle(kOpenSquare);
grShare->SetMarkerSize(1);
grShare->Draw("SAME LP");
leg->AddEntry(grShare, Form("%s share", r.title.c_str()), "p");
auto grEff = new TGraph(vecDepth.size(), vecDepth.data(), all.data());
grEff->SetLineColor(r.color);
grEff->SetLineWidth(2);
grEff->SetMarkerColor(r.color);
grEff->SetMarkerStyle(kFullDiamond);
grEff->SetMarkerSize(1);
grEff->Draw("SAME LP");
leg->AddEntry(grEff, Form("%s eff", r.title.c_str()), "p");
}
leg->Draw();
c1->SaveAs("ChipResponse.pdf");
}