CbmRoot
Loading...
Searching...
No Matches
CbmLitFieldQa.cxx
Go to the documentation of this file.
1/* Copyright (C) 2009-2020 GSI/JINR-LIT, Darmstadt/Dubna
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Andrey Lebedev [committer] */
4
10#include "CbmLitFieldQa.h"
11
12#include "CbmHistManager.h"
13#include "CbmLitFieldQaReport.h"
14#include "CbmRichDetectorData.h" // for CbmRichPmtData, CbmRichPixelData
16#include "CbmRichGeoManager.h"
17#include "CbmUtils.h"
18#include "FairField.h"
19#include "FairRunAna.h"
20#include "FairRuntimeDb.h"
21#include "TGraph.h"
22#include "TGraph2D.h"
23#include "TVector3.h"
24
25#include <TFile.h>
26
27#include <boost/assign/list_of.hpp>
28
29#include <cmath>
30#include <sstream>
31
32using boost::assign::list_of;
33using Cbm::ToString;
34using namespace std;
35
37 : fField(NULL)
38 , fZSlicePosition()
39 , fXSlicePosition()
40 , fYSlicePosition()
41 , fNofSlices(0)
42 , fNofBinsX(100)
43 , fNofBinsY(100)
44 , fMinZFieldIntegral(171.)
45 , fMaxZFieldIntegral(330.)
46 , fAlongZAngles()
47 , fAlongZXY()
48 , fZMin(-10.)
49 , fZMax(300.)
50 , fZStep(1.)
51 //, fAcceptanceAngleX(35.)
52 //, fAcceptanceAngleY(35.)
53 , fHM(NULL)
54 , fOutputDir("./")
55{
56}
57
59
61{
63
64 // Calculate (X, Y) window for each slice
67 for (Int_t i = 0; i < fNofSlices; i++) {
68 // Double_t tanXangle = tan(fAcceptanceAngleX * 3.14159265 / 180); //
69 // Double_t tanYangle = tan(fAcceptanceAngleY * 3.14159265 / 180); //
70 // fXSlicePosition[i] = fZSlicePosition[i] * tanXangle;
71 // fYSlicePosition[i] = fZSlicePosition[i] * tanYangle;
72 fXSlicePosition[i] = 250.;
73 fYSlicePosition[i] = 250.;
74 }
75
76 vector<Double_t> tmp = list_of(0.)(10.)(20.)(30.);
77 fAlongZAngles = tmp;
78
79 fAlongZXY.push_back(make_pair(0., 0.));
80 fAlongZXY.push_back(make_pair(100., 0.));
81 fAlongZXY.push_back(make_pair(0., 100.));
82
83
84 fField = FairRunAna::Instance()->GetField();
85
86 fHM = new CbmHistManager();
87
90
92
94 report->Create(fHM, fOutputDir);
95 delete report;
96 TDirectory* oldir = gDirectory;
97 TFile* outFile = FairRootManager::Instance()->GetOutFile();
98 if (outFile != NULL) {
99 outFile->cd();
100 fHM->WriteToFile();
101 }
102 gDirectory->cd(oldir->GetPath());
103
104 // print some field values at RICH entrance
105 Double_t B[3];
106 vector<vector<Double_t>> pos = {{0., 0., 170.}, {0., 80., 170.}, {50., 0., 170.}, {0., 0., 250.}};
107 for (UInt_t i = 0; i < pos.size(); i++) {
108 fField->GetFieldValue(&pos[i][0], B);
109 Double_t magB = sqrt(B[0] * B[0] + B[1] * B[1] + B[2] * B[2]);
110 cout << "B at (" << pos[i][0] << ", " << pos[i][1] << ", " << pos[i][2] << ") = " << magB << " kGauss" << endl;
111 }
112
113 cout << "CbmLitFieldQa::Init() end" << endl;
114
115 return kSUCCESS;
116}
117
118void CbmLitFieldQa::Exec(Option_t*) {}
119
121
123{
124 string names[] = {"Bx", "By", "Bz", "Mod"};
125 string zTitle[] = {"B_{x} [kGauss]", "B_{y} [kGauss]", "B_{z} [kGauss]", "|B| [kGauss]"};
126 for (Int_t v = 0; v < 4; v++) {
127 for (Int_t i = 0; i < fNofSlices; i++) {
128 TGraph2D* graph = new TGraph2D();
129 graph->SetNpx(200);
130 graph->SetNpy(200);
131 string name = "hmf_" + names[v] + "_Graph2D_" + ToString<Int_t>(fZSlicePosition[i]);
132 string title = name + ";X [cm];Y [cm];" + zTitle[v];
133 graph->SetNameTitle(name.c_str(), title.c_str());
134 fHM->Add(name, graph);
135 }
136 }
137
138 for (Int_t v = 0; v < 4; v++) {
139 for (UInt_t i = 0; i < fAlongZAngles.size(); i++) {
140 TGraph* graph = new TGraph();
141 string name = "hmf_" + names[v] + "AlongZAngle_Graph_" + ToString<Int_t>(fAlongZAngles[i]);
142 string title = name + ";Z [cm];B [kGauss]";
143 graph->SetNameTitle(name.c_str(), title.c_str());
144 fHM->Add(name, graph);
145 }
146 for (UInt_t i = 0; i < fAlongZXY.size(); i++) {
147 TGraph* graph = new TGraph();
148 string name = "hmf_" + names[v] + "AlongZXY_Graph_" + ToString<Int_t>(fAlongZXY[i].first) + "_"
149 + ToString<Int_t>(fAlongZXY[i].second);
150 string title = name + ";Z [cm];B [kGauss]";
151 graph->SetNameTitle(name.c_str(), title.c_str());
152 fHM->Add(name, graph);
153 }
154 for (UInt_t i = 0; i < fAlongZXY.size(); i++) {
155 TGraph* graph = new TGraph();
156 string name = "hmf_" + names[v] + "AlongZXYIntegral_Graph_" + ToString<Int_t>(fAlongZXY[i].first) + "_"
157 + ToString<Int_t>(fAlongZXY[i].second);
158 string title = name + ";Z [cm];B_{Int_t} [kGauss*m]";
159 graph->SetNameTitle(name.c_str(), title.c_str());
160 fHM->Add(name, graph);
161 }
162 }
163
164 for (Int_t p = 0; p < 2; p++) {
165 for (Int_t v = 0; v < 4; v++) {
166 TGraph2D* graph = new TGraph2D();
167 graph->SetNpx(200);
168 graph->SetNpy(200);
169 string name = "hmf_" + names[v] + "_Rich_Pmt_" + ((p == 0) ? "up" : "down");
170 string title = name + ";X [cm];Y [cm];" + zTitle[v];
171 graph->SetNameTitle(name.c_str(), title.c_str());
172 fHM->Add(name, graph);
173 }
174 }
175}
176
178{
179 string names[] = {"Bx", "By", "Bz", "Mod"};
180 // Fill graphs for magnetic field for each (X, Y) slice
181 for (Int_t iSlice = 0; iSlice < fNofSlices; iSlice++) { // loop over slices
182 Double_t Z = fZSlicePosition[iSlice];
183
184 Int_t cnt = 0;
185 Double_t HX = 2 * fXSlicePosition[iSlice] / fNofBinsX; // step size for X position
186 Double_t HY = 2 * fYSlicePosition[iSlice] / fNofBinsY; // step size for Y position
187 for (Int_t iX = 0; iX < fNofBinsX; iX++) { // loop over x position
188 Double_t X = -fXSlicePosition[iSlice] + (iX + 0.5) * HX;
189 for (Int_t iY = 0; iY < fNofBinsY; iY++) { // loop over y position
190 Double_t Y = -fYSlicePosition[iSlice] + (iY + 0.5) * HY;
191
192 // Get field value
193 Double_t pos[3] = {X, Y, Z};
194 Double_t B[4];
195 fField->GetFieldValue(pos, B);
196
197 B[3] = sqrt(B[0] * B[0] + B[1] * B[1] + B[2] * B[2]);
198 for (Int_t v = 0; v < 4; v++) {
199 string name = "hmf_" + names[v] + "_Graph2D_" + ToString<Int_t>(fZSlicePosition[iSlice]);
200 fHM->G2(name)->SetPoint(cnt, X, Y, B[v]);
201 }
202 cnt++;
203 }
204 }
205 }
206
207 // Fill histograms for magnetic field along Z for different angles
208 for (UInt_t i = 0; i < fAlongZAngles.size(); i++) {
209 Int_t nofSteps = Int_t((fZMax - fZMin) / fZStep);
210 for (Int_t istep = 0; istep < nofSteps; istep++) {
211 Double_t Z = fZMin + istep * fZStep;
212 Double_t tanXangle = tan(fAlongZAngles[i] * 3.14159265 / 180); //
213 Double_t tanYangle = tan(fAlongZAngles[i] * 3.14159265 / 180); //
214 Double_t X = Z * tanXangle;
215 Double_t Y = Z * tanYangle;
216
217 // Get field value
218 Double_t pos[3] = {X, Y, Z};
219 Double_t B[4];
220 fField->GetFieldValue(pos, B);
221
222 B[3] = sqrt(B[0] * B[0] + B[1] * B[1] + B[2] * B[2]);
223 for (Int_t v = 0; v < 4; v++) {
224 string name = "hmf_" + names[v] + "AlongZAngle_Graph_" + ToString<Int_t>(fAlongZAngles[i]);
225 fHM->G1(name)->SetPoint(istep, Z, B[v]);
226 }
227 }
228 }
229 // Fill histograms for magnetic field along Z for different X position
230 for (UInt_t i = 0; i < fAlongZXY.size(); i++) {
231 Int_t nofSteps = Int_t((fZMax - fZMin) / fZStep);
232 Double_t integralB[4] = {0., 0., 0., 0.};
233 for (Int_t istep = 0; istep < nofSteps; istep++) {
234 Double_t Z = fZMin + istep * fZStep;
235 Double_t X = fAlongZXY[i].first;
236 Double_t Y = fAlongZXY[i].second;
237
238 // Get field value
239 Double_t pos[3] = {X, Y, Z};
240 Double_t B[4];
241 fField->GetFieldValue(pos, B);
242
243 B[3] = sqrt(B[0] * B[0] + B[1] * B[1] + B[2] * B[2]);
244
245 for (Int_t v = 0; v < 4; v++) {
246 string name = "hmf_" + names[v] + "AlongZXY_Graph_" + ToString<Int_t>(fAlongZXY[i].first) + "_"
247 + ToString<Int_t>(fAlongZXY[i].second);
248 fHM->G1(name)->SetPoint(istep, Z, B[v]);
249 // Calculate field integral in the RICH detector
250 if (Z >= fMinZFieldIntegral && Z <= fMaxZFieldIntegral) {
251 integralB[v] += 0.01 * fZStep * fabs(B[v]); // in kGauss * meter
252 name = "hmf_" + names[v] + "AlongZXYIntegral_Graph_" + ToString<Int_t>(fAlongZXY[i].first) + "_"
253 + ToString<Int_t>(fAlongZXY[i].second);
254 fHM->G1(name)->SetPoint(istep, Z, integralB[v]);
255 fHM->G1(name)->SetMaximum(1.1 * integralB[v]);
256 }
257 }
258 }
259 }
260}
261
263{
264 string names[] = {"Bx", "By", "Bz", "Mod"};
265 vector<Int_t> pixels = CbmRichDigiMapManager::GetInstance().GetPixelAddresses();
266 if (pixels.size() == 0) return;
267
268 Double_t maxModB = 0.;
269 //int iS = -1;
270 for (UInt_t i = 0; i < pixels.size(); i++) {
271 TVector3 inPos;
273 inPos.SetXYZ(pixelData->fX, pixelData->fY, pixelData->fZ);
274 string ud = "up";
275 if (inPos.Y() < 0) ud = "down";
276 TVector3 outPos;
278
279 // Get field value
280 // Take B-field values from real PMT position
281 Double_t pos[3] = {inPos.X(), inPos.Y(), inPos.Z()};
282
283 Double_t B[4];
284 fField->GetFieldValue(pos, B);
285
286 B[3] = sqrt(B[0] * B[0] + B[1] * B[1] + B[2] * B[2]);
287 if (maxModB < B[3]) maxModB = B[3];
288 for (Int_t v = 0; v < 4; v++) {
289 string name = "hmf_" + names[v] + "_Rich_Pmt_" + ud;
290 // Take B-field values from real PMT position, but display them for rotated X and Y
291 fHM->G2(name)->SetPoint(fHM->G2(name)->GetN(), outPos.X(), outPos.Y(), B[v]);
292 }
293 }
294
295 cout << "Maximum Bmod onto RICH PMT:" << maxModB << " kGauss" << endl;
296}
297
std::string ToString(ECbmModuleId modId)
Definition CbmDefs.cxx:70
Histogram manager.
Creates field QA report.
ClassImp(CbmLitFieldQa)
Field map QA.
friend fvec sqrt(const fvec &a)
fscal v[fmask::Size]
Definition KfSimdPseudo.h:4
bool first
Histogram manager.
TGraph2D * G2(const std::string &name) const
Return pointer to TGraph2D.
TGraph * G1(const std::string &name) const
Return pointer to TGraph.
void WriteToFile()
Write all objects to current opened file.
void Add(const std::string &name, TNamed *object)
Add new named object to manager.
Creates field QA report.
Field map QA.
vector< Double_t > fZSlicePosition
virtual void Finish()
Inherited from FairTask.
vector< Double_t > fAlongZAngles
vector< Double_t > fXSlicePosition
vector< std::pair< Double_t, Double_t > > fAlongZXY
vector< Double_t > fYSlicePosition
void CreateHistos()
Create histograms.
Double_t fMinZFieldIntegral
virtual void Exec(Option_t *opt)
Inherited from FairTask.
void FillBHistos()
Fill graphs and histos for field map for each field component (Bx, By, Bz, |B|).
Double_t fMaxZFieldIntegral
virtual InitStatus Init()
Inherited from FairTask.
void FillRichPmtPlaneBHistos()
Fill B-field histograms for RICH PMT plane.
CbmHistManager * fHM
virtual ~CbmLitFieldQa()
Destructor.
FairField * fField
CbmLitFieldQa()
Constructor.
static CbmRichDigiMapManager & GetInstance()
Return Instance of CbmRichGeoManager.
std::vector< Int_t > GetPixelAddresses()
Return addresses of all pixels.
CbmRichPixelData * GetPixelDataByAddress(Int_t address)
Return CbmRichDataPixel by digi address.
void RotatePoint(TVector3 *inPos, TVector3 *outPos, Bool_t noTilting=false)
static CbmRichGeoManager & GetInstance()
Base class for simulation reports.
void Create(CbmHistManager *histManager, const std::string &outputDir)
Main function which creates report data.
std::string ToString(const T &value)
Definition CbmUtils.h:26
Hash for CbmL1LinkKey.