CbmRoot
Loading...
Searching...
No Matches
CbmCaInputQaSts.cxx
Go to the documentation of this file.
1/* Copyright (C) 2023 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Sergei Zharko [committer] */
4
9
10#include "CbmCaInputQaSts.h"
11
12#include "CbmAddress.h"
13#include "CbmMCDataArray.h"
14#include "CbmMCEventList.h"
15#include "CbmMCTrack.h"
16#include "CbmMatch.h"
17#include "CbmQaCanvas.h"
18#include "CbmQaTable.h"
19#include "CbmQaUtil.h"
20#include "CbmStsCluster.h"
21#include "CbmStsHit.h"
22#include "CbmStsPoint.h"
24#include "CbmTimeSlice.h"
25#include "FairRootManager.h"
26#include "Logger.h"
27#include "TBox.h"
28#include "TClonesArray.h"
29#include "TEllipse.h"
30#include "TF1.h"
31#include "TFormula.h"
32#include "TGraphAsymmErrors.h"
33#include "TH1F.h"
34#include "TH2F.h"
35#include "TMath.h"
36#include "TProfile.h"
37#include "TProfile2D.h"
38#include "TStyle.h"
39
40#include <algorithm>
41#include <fstream>
42#include <iomanip>
43#include <numeric>
44#include <tuple>
45
47
48// ---------------------------------------------------------------------------------------------------------------------
49//
50CbmCaInputQaSts::CbmCaInputQaSts(int verbose, bool isMCUsed) : CbmCaInputQaBase("CbmCaInputQaSts", verbose, isMCUsed)
51{
52 // Default parameters of task
54}
55
56// ---------------------------------------------------------------------------------------------------------------------
57//
59{
61 if (IsMCUsed()) {
62 for (int idig = 0; idig <= fkMaxDigisInClusterForPulls; idig++) {
64 auto [msgU, resU] = CheckRangePull(fvph_pull_u_Ndig[idig]);
65 StoreCheckResult(Form("pull_u_%d_digis", idig), resU, msgU);
66
68 auto [msgV, resV] = CheckRangePull(fvph_pull_v_Ndig[idig]);
69 StoreCheckResult(Form("pull_v_%d_digis", idig), resV, msgV);
70 }
71 } // McUsed
72}
73
74
75// ---------------------------------------------------------------------------------------------------------------------
76//
78{
79 gStyle->SetOptFit(1);
80
81 if (IsMCUsed()) {
82 { // u-coordinate vs N digis
83 auto* canv = MakeQaObject<CbmQaCanvas>("vs N digi/pull_u", "Pulls for u coordinate vs N digis", 1600, 800);
84 canv->Divide2D(fkMaxDigisInClusterForPulls + 1);
85 for (int i = 0; i <= fkMaxDigisInClusterForPulls; ++i) {
86 canv->cd((i > 0) ? i : fkMaxDigisInClusterForPulls + 1);
87 fvph_pull_u_Ndig[i]->DrawCopy("", "");
88 }
89 }
90
91 { // v-coordinate vs N digis
92 auto* canv = MakeQaObject<CbmQaCanvas>("vs N digi/pull_v", "Pulls for v coordinate vs N digis", 1600, 800);
93 canv->Divide2D(fkMaxDigisInClusterForPulls + 1);
94 for (int i = 0; i <= fkMaxDigisInClusterForPulls; ++i) {
95 canv->cd((i > 0) ? i : fkMaxDigisInClusterForPulls + 1);
96 fvph_pull_v_Ndig[i]->DrawCopy("", "");
97 }
98 }
99 }
100
102}
103
104// ---------------------------------------------------------------------------------------------------------------------
105//
107{
109
110 // Vectors with pointers to histograms
111 fvph_pull_u_Ndig.clear();
112 fvph_pull_v_Ndig.clear();
113}
114
115// ---------------------------------------------------------------------------------------------------------------------
116//
118{
119 auto SetRange = [](std::array<double, 2>& range, double min, double max) {
120 range[0] = min;
121 range[1] = max;
122 };
123 // Hit errors
124 SetRange(fRHitDx, 0.0000, 0.0050); // [cm]
125 SetRange(fRHitDy, 0.0000, 0.0200); // [cm]
126 SetRange(fRHitDu, 0.0000, 0.0050); // [cm]
127 SetRange(fRHitDv, 0.0000, 0.0050); // [cm]
128 SetRange(fRHitDt, 0.0000, 10.000); // [ns]
129 // Residuals
130 SetRange(fRResX, -0.02, 0.02);
131 SetRange(fRResY, -0.10, 0.10);
132 SetRange(fRResU, -0.02, 0.02);
133 SetRange(fRResV, -0.02, 0.02);
134 SetRange(fRResT, -25.0, 25.0);
135}
136
137// ---------------------------------------------------------------------------------------------------------------------
138//
140{
141 const auto* pHit = dynamic_cast<const CbmStsHit*>(fpHits->At(fHitQaData.GetHitIndex()));
142
143 if (IsMCUsed()) { // u-coordinate residual per number of digis
144 const auto* pCluster = dynamic_cast<const CbmStsCluster*>(fpClusters->At(pHit->GetFrontClusterId()));
145 assert(pCluster);
146 int nDigis = pCluster->GetNofDigis();
147 if (nDigis > fkMaxDigisInClusterForPulls) {
148 nDigis = 0;
149 }
150 fvph_pull_u_Ndig[nDigis]->Fill(fHitQaData.GetPullU());
151 }
152
153 if (IsMCUsed()) { // v-coordinate residual per number of digis
154 const auto* pCluster = dynamic_cast<const CbmStsCluster*>(fpClusters->At(pHit->GetBackClusterId()));
155 assert(pCluster);
156 int nDigis = pCluster->GetNofDigis();
157 if (nDigis > fkMaxDigisInClusterForPulls) {
158 nDigis = 0;
159 }
160 fvph_pull_v_Ndig[nDigis]->Fill(fHitQaData.GetPullV());
161 }
162}
163
164// ---------------------------------------------------------------------------------------------------------------------
165//
167
168// ---------------------------------------------------------------------------------------------------------------------
169//
171{
172 // ----- Detector intereface
173 //
175
176 auto initStatus = CbmCaInputQaBase::InitQa();
177 if (initStatus != kSUCCESS) {
178 return initStatus;
179 }
180
181 // ----- Input data initialization
182 //
183 fpClusters = dynamic_cast<TClonesArray*>(FairRootManager::Instance()->GetObject("StsCluster"));
184 LOG_IF(fatal, !fpClusters) << "\033[1;31m" << fName << ": container of hit clusters in STS is not found\033[0m";
185
186 // ----- Histogram initialization
187 //
189 std::string detName = fpDetInterface->GetDetectorName();
190
191 // rename some histogramms with respect to STS specifics
192 for (int iSt = 0; iSt <= nSt; ++iSt) {
193 TString tsuff = (iSt == nSt) ? "" : Form(" in %s station %d", detName.c_str(), iSt); // Histogram title suffix
194 fvph_hit_du[iSt]->SetTitle((TString) "Hit position error across front strips" + tsuff);
195 fvph_hit_dv[iSt]->SetTitle((TString) "Hit position error across back strips" + tsuff);
196 fvph_hit_kuv[iSt]->SetTitle((TString) "Hit error correlation between front and back strips" + tsuff);
197 if (IsMCUsed()) {
198 fvph_res_u[iSt]->SetTitle((TString) "Residuals for Front strip coordinate" + tsuff);
199 fvph_res_v[iSt]->SetTitle((TString) "Residuals for Back strip coordinate" + tsuff);
200 fvph_pull_u[iSt]->SetTitle((TString) "Pulls for Front strip coordinate" + tsuff);
201 fvph_pull_v[iSt]->SetTitle((TString) "Pulls for Back strip coordinate" + tsuff);
202 fvph_res_u_vs_u[iSt]->SetTitle((TString) "Residuals for Front strip coordinate" + tsuff);
203 fvph_res_v_vs_v[iSt]->SetTitle((TString) "Residuals for Back strip coordinate" + tsuff);
204 fvph_pull_u_vs_u[iSt]->SetTitle((TString) "Pulls for Front strip coordinate" + tsuff);
205 fvph_pull_v_vs_v[iSt]->SetTitle((TString) "Pulls for Back strip coordinate" + tsuff);
206 }
207 }
208
209 if (IsMCUsed()) {
212
213 for (int idig = 0; idig <= fkMaxDigisInClusterForPulls; idig++) {
214 TString sN = "All stations/pull/Ndigi/pull_u_";
215 TString sT = "Pulls for Front strip coordinate, hits with ";
216 if (idig == 0) {
217 sN += Form("%d+_digi", fkMaxDigisInClusterForPulls + 1);
218 sT += Form(">=%d digi", fkMaxDigisInClusterForPulls + 1);
219 }
220 else {
221 sN += Form("%d_digi", idig);
222 sT += Form("%d digi", idig);
223 }
224 sT += "; (u_{reco} - u_{MC}) / #sigma_{u}^{reco}";
226 }
227
228 for (int idig = 0; idig <= fkMaxDigisInClusterForPulls; idig++) {
229 TString sN = "All stations/pull/Ndigi/pull_v_";
230 TString sT = "Pulls for Back strip coordinate, hits with ";
231 if (idig == 0) {
232 sN += Form("%d+_digi", fkMaxDigisInClusterForPulls + 1);
233 sT += Form(">=%d digi", fkMaxDigisInClusterForPulls + 1);
234 }
235 else {
236 sN += Form("%d_digi", idig);
237 sT += Form("%d digi", idig);
238 }
239 sT += "; (v_{reco} - v_{MC}) / #sigma_{v}^{reco}";
241 }
242 }
243
244 return kSUCCESS;
245}
ClassImp(CbmCaInputQaSts)
QA-task for CA tracking input from MuCh detector (header)
Definition of the CbmQaCanvas class.
Definition of CbmQaTable class.
Useful utilities for CBM QA tasks.
Data class for STS clusters.
Data class for a reconstructed hit in the STS.
friend fscal max(fscal x, fscal y)
friend fscal min(fscal x, fscal y)
A QA-task class, which provides assurance of MuCh hits and geometry.
CbmTrackingDetectorInterfaceBase * fpDetInterface
void CreateSummary() override
Creates summary cavases, tables etc.
std::pair< std::string, bool > CheckRangePull(TH1 *h) const
void Check() override
Method to check, if the QA results are acceptable.
InitStatus InitQa() override
Initializes QA task.
void DeInit() override
De-initializes histograms.
A QA-task class, which provides assurance of MuCh hits and geometry.
void FillHistogramsPerHit() override
Fills histograms per hit.
CbmCaInputQaSts(int verbose, bool isMCUsed)
void FillHistogramsPerPoint() override
Fills histograms per MC point.
std::vector< TH1F * > fvph_pull_u_Ndig
pull for u coordinate, depending on N digis in the cluster
InitStatus InitQa() override
Initializes QA.
void DefineParameters() override
Defines parameters of the task.
const int fkMaxDigisInClusterForPulls
max digis in cluster for separate histogramming of puls
void CreateSummary() override
Creates summary cavases, tables etc.
std::vector< TH1F * > fvph_pull_v_Ndig
pull for v coordinate, depending on N digis in the cluster
void DeInit() override
De-initializes histograms.
void Check() override
Method to check, if the QA results are acceptable.
TClonesArray * fpClusters
Array of hit clusters.
T * MakeQaObject(TString sName, TString sTitle, Args... args)
Definition CbmQaIO.h:261
void StoreCheckResult(const std::string &tag, bool result, const std::string &msg="")
Stores check flag to the check-list.
bool IsMCUsed() const
Returns flag, whether MC information is used or not in the task.
Definition CbmQaTask.h:126
Data class for STS clusters.
data class for a reconstructed 3-d hit in the STS
Definition CbmStsHit.h:35
static CbmStsTrackingInterface * Instance()
Gets pointer to the instance of the CbmStsTrackingInterface class.
int GetNtrackingStations() const
Gets actual number of stations, provided by the current geometry setup.
virtual std::string GetDetectorName() const =0
Returns the name of the detector subsystem.
double GetPullU() const
Gets pull of u-coordinate.
double GetPullV() const
Gets pull of v-coordinate.
int GetHitIndex() const
Gets hit index.
void SetLargeStats(TH1 *pHist)
Set large stat. window.