10 #include <dqm/analysis/modules/DQMHistAnalysisCDCMonObj.h>
13 #include <cdc/geometry/CDCGeometryPar.h>
32 DQMHistAnalysisCDCMonObjModule::DQMHistAnalysisCDCMonObjModule()
36 setDescription(
"Modify and analyze the data quality histograms of CDCMonObj");
38 for (
int i = 0; i < 300; i++) {
42 for (
int i = 0; i < 56; i++)
m_hHits[i] =
nullptr;
53 if (!(*m_channelMapFromDB).isValid()) {
54 B2FATAL(
"Channel map is not valid");
59 B2FATAL(
"CDCGeometryp is not valid");
66 gStyle->SetOptStat(0);
67 gStyle->SetPalette(kViridis);
68 gStyle->SetPadTopMargin(0.1);
69 gStyle->SetPadRightMargin(0.05);
70 gStyle->SetPadBottomMargin(0.1);
71 gStyle->SetPadLeftMargin(0.15);
73 m_cMain =
new TCanvas(
"cdc_main",
"cdc_main", 1500, 1000);
76 m_cADC =
new TCanvas(
"cdc_adc",
"cdc_adc", 2000, 10000);
78 m_cTDC =
new TCanvas(
"cdc_tdc",
"cdc_tdc", 2000, 10000);
80 m_cHit =
new TCanvas(
"cdc_hit",
"cdc_hit", 1500, 6000);
83 B2DEBUG(20,
"DQMHistAnalysisCDCMonObj: initialized.");
90 const int isl = cm.getISuperLayer();
91 const int il = cm.getILayer();
92 const int iw = cm.getIWire();
93 const int iBoard = cm.getBoardID();
94 const int iCh = cm.getBoardChannel();
95 const WireID wireId(isl, il, iw);
96 m_chMap.insert(std::make_pair(wireId, std::make_pair(iBoard, iCh)));
102 for (
const auto& sense : geom.getSenseLayers()) {
103 int i = sense.getId();
104 if (i < 0 || i > 55) {
105 B2FATAL(
"no such sense layer");
112 for (
const auto& field : geom.getFieldLayers()) {
113 int i = field.getId();
114 if (i < 0 || i > 54) {
115 B2FATAL(
"no such sense layer");
119 m_fieldR[56] = geom.getOuterWall(0).getRmin();
120 m_fieldR[0] = geom.getInnerWall(0).getRmax();
127 for (
int ilayer = 0; ilayer < 56; ++ilayer) {
131 for (
int iwire = 0; iwire <
m_nSenseWires[ilayer]; ++iwire) {
132 float phi = dPhi * (iwire +
m_offset[ilayer]);
133 float phi1 = phi - dPhi * 0.5;
134 float phi2 = phi + dPhi * 0.5;
135 Double_t x_pos[] = {r1* (sin(phi)*
tan(phi - phi1) + cos(phi)),
138 r1* (sin(phi)*
tan(phi - phi2) + cos(phi))
140 Double_t y_pos[] = {r1* (-cos(phi)*
tan(phi - phi1) + sin(phi)),
143 r1* (-cos(phi)*
tan(phi - phi2) + sin(phi))
145 h->AddBin(4, x_pos, y_pos);
154 for (
int il = 0; il < 56; ++il) {
156 const int y =
m_hHits[il]->GetBinContent(iw + 1);
167 TH1D* hist = (TH1D*)h->Clone();
168 hist->SetBinContent(1, 0.0);
169 float m = hist->GetMean();
175 TH1D* hist = (TH1D*)h->Clone();
176 hist->SetBinContent(1, 0.0);
177 if (hist->GetMean() == 0) {
return 0.0;}
178 double quantiles[1] = {0.0};
179 double probSums[1] = {0.5};
180 hist->GetQuantiles(1, quantiles, probSums);
181 float median = quantiles[0];
187 const WireID w(layer, wire);
192 B2ERROR(
"no corresponding board/channel found layer " << layer <<
" wire " << wire);
193 return std::make_pair(-1, -1);
200 B2DEBUG(20,
"end run");
207 B2INFO(
"Histogram named m_hADC is not found.");
211 TF1* fitFunc[300] = {};
212 for (
int i = 0; i < 300; ++i) {
213 fitFunc[i] =
new TF1(Form(
"f%d", i),
"[0]+[6]*x+[1]*(exp([2]*(x-[3]))/(1+exp(-([4]-x)/[5])))",
214 4921 - 100, 4921 + 100);
215 fitFunc[i]->SetParLimits(6, 0, 0.1);
216 fitFunc[i]->SetParLimits(4, 4850., 5000.0);
217 fitFunc[i]->SetParLimits(5, 0, 50.0);
220 int neve =
m_hTDC->GetEntries();
221 if (neve == 0)neve = 1;
223 B2DEBUG(20,
"adc related");
226 TH1F* hADCMean =
new TH1F(
"hADCMean",
"ADC mean;board;adc mean", 300, 0, 300);
227 TH1F* hADC1000 =
new TH1F(
"ADC1000",
"ADC1000", 300, 0, 300);
228 TH1F* hADC0 =
new TH1F(
"ADC0",
"ADC0", 300, 0, 300);
231 std::vector<float> means = {};
232 std::vector<float> medians = {};
234 for (
int i = 0; i < 300; ++i) {
235 m_hADCs[i] =
m_hADC->ProjectionY(Form(
"hADC%d", i), i + 1, i + 1,
"");
236 m_hADCs[i]->SetTitle(Form(
"hADC%d", i));
237 float n =
static_cast<float>(
m_hADCs[i]->GetEntries());
240 hADC0->SetBinContent(i + 1, -0.1);
242 float n0 =
static_cast<float>(
m_hADCs[i]->GetBinContent(1));
244 B2DEBUG(21,
"bad adc bid " << i <<
" " << n0 <<
" " << n);
247 float bin1 =
m_hADCs[i]->GetBinContent(1);
251 medians.push_back(md);
252 hADCMean->SetBinContent(i + 1, m);
253 hADCMean->SetBinError(i + 1, 0);
254 double overflow =
m_hADCs[i]->GetBinContent(
m_hADCs[i]->GetNbinsX() + 1);
255 hADC1000->SetBinContent(i + 1, overflow / (overflow + n));
256 hADC0->SetBinContent(i + 1, bin1 / (overflow + n));
261 B2DEBUG(20,
"tdc related");
263 TH1F* hTDCEdge =
new TH1F(
"hTDCEdge",
"TDC edge;board;tdc edge [nsec]", 300, 0, 300);
264 TH1F* hTDCSlope =
new TH1F(
"hTDCSlope",
"TDC slope;board;tdc slope [nsec]", 300, 0, 300);
265 std::vector<float> tdcEdges = {};
266 std::vector<float> tdcSlopes = {};
267 for (
int i = 0; i < 300; ++i) {
268 m_hTDCs[i] =
m_hTDC->ProjectionY(Form(
"hTDC%d", i), i + 1, i + 1);
269 m_hTDCs[i]->SetTitle(Form(
"hTDC%d", i));
272 tdcEdges.push_back(0);
273 tdcSlopes.push_back(0);
275 double init_p0 =
m_hTDCs[i]->GetBinContent(700 + 60);
276 fitFunc[i]->SetParameters(init_p0, 100, 0.01, 4700, 4900, 2, 0.01);
277 fitFunc[i]->SetParameter(6, 0.02);
278 fitFunc[i]->SetParLimits(0, init_p0 - 200, init_p0 + 200);
281 xxx =
m_hTDCs[i]->Fit(fitFunc[i],
"qM0",
"", 4850, 5000);
283 xxx =
m_hTDCs[i]->Fit(fitFunc[i],
"qM0",
"", 4800, 5000);
285 float p4 = fitFunc[i]->GetParameter(4);
286 float p5 = fitFunc[i]->GetParameter(5);
288 if (xxx != -1 && 4850 < p4 && p4 < 5000) {
289 hTDCEdge->SetBinContent(i + 1, p4);
290 hTDCEdge->SetBinError(i + 1, 0);
291 hTDCSlope->SetBinContent(i + 1, p5);
292 hTDCSlope->SetBinError(i + 1, 0);
294 tdcEdges.push_back(p4);
295 tdcSlopes.push_back(p5);
301 B2DEBUG(20,
"hit related");
302 TH1F* hHitPerLayer =
new TH1F(
"hHitPerLayer",
"hit/Layer;layer", 56, 0, 56);
304 for (
int i = 0; i < 56; ++i) {
305 m_hHits[i] =
m_hHit->ProjectionY(Form(
"hHit%d", i), i + 1, i + 1);
306 m_hHits[i]->SetTitle(Form(
"hHit%d", i));
310 for (
int j = 0; j < nBins; ++j) {
311 nhitSumL +=
m_hHits[i]->GetBinContent(j + 1);
314 hHitPerLayer->SetBinContent(i + 1,
static_cast<float>(1.0 * nhitSumL / neve));
315 }
else hHitPerLayer->SetBinContent(i + 1,
static_cast<float>(nhitSumL));
316 hHitPerLayer->SetBinError(i + 1, 0);
322 B2DEBUG(20,
"bad wire related");
323 hBadChannel =
new TH2F(
"hbadch",
"bad channel map;wire;layer", 400, 0, 400, 56, 0, 56);
324 for (
int i = 0; i < 400; ++i) {
325 for (
int j = 0; j < 56; ++j) {
330 hBadChannelBC =
new TH2F(
"hbadchBC",
"bad channel map per board/channel;board;channel", 300, 0, 300, 48, 0, 48);
331 for (
int i = 0; i < 300; ++i) {
332 for (
int j = 0; j < 48; ++j) {
339 h2p->SetTitle(
"bad wires in xy view");
340 h2p->GetXaxis()->SetTitle(
"X [cm]");
341 h2p->GetYaxis()->SetTitle(
"Y [cm]");
344 const int l = lw.first;
345 const int w = lw.second;
346 B2DEBUG(21,
"l " << l <<
" w " << w);
351 float dPhi =
static_cast<float>(2.0 * M_PI /
m_nSenseWires[l]);
352 float phi = dPhi * (w +
m_offset[l]);
353 float x = r * cos(phi);
354 float y = r * sin(phi);
355 h2p->Fill(x, y, 1.1);
358 B2DEBUG(20,
"writing");
362 hADCMean->SetMinimum(0);
363 hADCMean->SetMaximum(300);
364 hADCMean->DrawCopy();
367 hTDCEdge->SetMinimum(4800);
368 hTDCEdge->SetMaximum(5000);
369 hTDCEdge->DrawCopy();
372 hTDCSlope->SetMinimum(0);
373 hTDCSlope->SetMaximum(50);
374 hTDCSlope->DrawCopy();
381 hHitPerLayer->DrawCopy();
383 hADC1000->DrawCopy();
388 for (
int i = 0; i < 56; i++) {
394 m_cADC->Divide(6, 50, 0.0002, 0.0002);
395 m_cTDC->Divide(6, 50, 0.0002, 0.0002);
397 for (
int i = 0; i < 300; i++) {
399 Double_t max =
m_hADCs[i]->GetMaximum();
400 m_hADCs[i]->GetYaxis()->SetRangeUser(0, 3 * max);
405 fitFunc[i]->SetLineColor(kRed);
406 fitFunc[i]->Draw(
"same");
407 max =
m_hTDCs[i]->GetMaximum();
408 TLine* l1 =
new TLine(tdcEdges[i], 0, tdcEdges[i], max * 1.05);
409 l1->SetLineColor(kRed);
410 TLine* l0 =
new TLine(4910, 0, 4910, max * 1.05);
416 h2p->DrawCopy(
"col");
417 float superLayerR[10] = {16.3, 24.3, 35.66, 46.63, 57.55, 68.47,
418 79.39, 90.31, 101.23, 112.05
422 for (
int i = 0; i < 10; ++i) {
423 circs[i] =
new TEllipse(0, 0, superLayerR[i], superLayerR[i]);
424 circs[i]->SetFillStyle(4000);
425 circs[i]->SetLineStyle(kDashed);
426 circs[i]->SetLineColor(0);
427 circs[i]->Draw(
"same");
433 m_monObj->
setVariable(
"adcMean", std::accumulate(means.begin(), means.end(), 0.0) / means.size());
434 m_monObj->
setVariable(
"adcMeanMedianBoard", std::accumulate(medians.begin(), medians.end(), 0.0) / medians.size());
437 m_monObj->
setVariable(
"tdcEdge", std::accumulate(tdcEdges.begin(), tdcEdges.end(), 0.0) / (tdcEdges.size() - 1 - nDeadTDC));
439 m_monObj->
setVariable(
"tdcSlope", std::accumulate(tdcSlopes.begin(), tdcSlopes.end(), 0.0) / (tdcSlopes.size() - 1 - nDeadTDC));
453 B2DEBUG(20,
"terminate called");
The Class for CDC geometry.
Class for accessing arrays of objects in the database.
Class for accessing objects in the database.
~DQMHistAnalysisCDCMonObjModule()
Destructor.
std::map< WireID, std::pair< int, int > > m_chMap
Channel map retrieved
TCanvas * m_cTDC
TDC panel.
void initialize() override final
Initialize the Module.
DBArray< CDCChannelMap > * m_channelMapFromDB
Channel map retrieved from DB.
std::vector< std::pair< int, int > > m_badChannels
bad wires list
int m_nSenseWires[56]
number of wires for each layer.
TCanvas * m_cHit
Hit panel.
TH1D * m_hTDCs[300]
TDC histograms with track associated hits for each board (0-299)
MonitoringObject * m_monObj
monitoring object
void terminate() override final
Termination action.
DBObjPtr< CDCGeometry > * m_cdcGeo
Geometry of CDC.
TH2F * m_hTDC
Summary of TDC histograms with track associated hits.
TH1D * m_hADCs[300]
ADC histograms with track associated hits for each board (0-299)
void makeBadChannelList()
make bad channel list.
TH2F * hBadChannel
bad channel map;wire;layer
std::pair< int, int > getBoardChannel(unsigned short layer, unsigned short wire)
Get board/channel from layer/wire.
TH2Poly * h2p
bad wires in xy view
TH1D * m_hHits[56]
hit histograms for each layer (0-55)
void endRun() override final
End-of-run action.
void configureBins(TH2Poly *h)
Configure bins of TH2Poly.
float getHistMean(TH1D *h) const
Get mean of ADC histogram excluding 0-th bin.
float m_offset[56]
Offset of sense layer
void beginRun() override final
Called when entering a new run.
TH2F * m_hADC
Summary of ADC histograms with track associated hits.
TH2F * hBadChannelBC
bad channel map per board/channel;board;channel
float m_fieldR[57]
Radius of field layer.
float getHistMedian(TH1D *h) const
Get median of ADC histogram excluding 0-th bin.
float m_senseR[56]
Radius of sense (+field) layer.
TCanvas * m_cADC
ADC panel.
TCanvas * m_cMain
main panel
TH2F * m_hHit
Summary of hit histograms.
The base class for the histogram analysis module.
static TH1 * findHist(const std::string &histname, bool onlyIfUpdated=false)
Get histogram from list (no other search).
static MonitoringObject * getMonitoringObject(const std::string &histname)
Get MonitoringObject with given name (new object is created if non-existing)
void setDescription(const std::string &description)
Sets the description of the module.
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
void setVariable(const std::string &var, float val, float upErr=-1., float dwErr=-1)
set value to float variable (new variable is made if not yet existing)
void addCanvas(TCanvas *canv)
Add Canvas to monitoring object.
Class to identify a wire inside the CDC.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
double tan(double a)
tan for double
Abstract base class for different kinds of events.