CbmRoot
CbmLitFieldApproximationQa.cxx
Go to the documentation of this file.
1 
7 #include "CbmDrawHist.h"
8 #include "CbmHistManager.h"
10 #include "CbmUtils.h"
11 #include "base/CbmLitFieldFitter.h"
13 #include "base/CbmLitFloat.h"
14 
15 #include "../../../parallel/LitFieldGrid.h"
16 #include "../../../parallel/LitFieldSlice.h"
17 #include "../../../parallel/LitFieldValue.h"
18 
19 #include "FairField.h"
20 #include "FairRunAna.h"
21 #include "FairRuntimeDb.h"
22 
23 #include "TCanvas.h"
24 #include "TF1.h"
25 #include "TF2.h"
26 #include "TGraph.h"
27 #include "TGraph2D.h"
28 #include "TH2D.h"
29 #include "TLegend.h"
30 #include "TPad.h"
31 #include "TPaveText.h"
32 #include "TStyle.h"
33 
34 #include <cmath>
35 #include <limits>
36 #include <sstream>
37 #include <string>
38 
39 #include <boost/assign/list_of.hpp>
40 
41 using boost::assign::list_of;
43 using Cbm::ToString;
47 using std::cout;
48 using std::endl;
49 using std::string;
50 
52  : fField(NULL)
53  , fNofSlices(0)
54  , fZSlicePosition()
55  , fXSlicePosition()
56  , fYSlicePosition()
57  , fFixedBounds(true)
58  , fAcceptanceAngleX(25.)
59  , fAcceptanceAngleY(25.)
60  , fNofBinsX(30)
61  , fNofBinsY(30)
62  , fUseEllipseAcc(true)
63  , fOutputDir("./field/")
64  , fFitter()
65  , fGridCreator()
66  , fPolynomDegreeIndex(1)
67  , fNofPolynoms(4)
68  , fPolynomDegrees()
69  , fHM(NULL) {}
70 
72 
74  fNofSlices = fZSlicePosition.size();
76 
77  // Calculate (X, Y) window for each slice
80  for (Int_t i = 0; i < fNofSlices; i++) {
81  Double_t tanXangle = tan(fAcceptanceAngleX * 3.14159265 / 180); //
82  Double_t tanYangle = tan(fAcceptanceAngleY * 3.14159265 / 180); //
83  fXSlicePosition[i] = fZSlicePosition[i] * tanXangle;
84  fYSlicePosition[i] = fZSlicePosition[i] * tanYangle;
85  }
86 
87  // Create field approximation tool for each polynom degree
88  fFitter.resize(fNofPolynoms);
89  for (UInt_t i = 0; i < fNofPolynoms; i++) {
91  fFitter[i]->SetXangle(fAcceptanceAngleX);
92  fFitter[i]->SetYangle(fAcceptanceAngleY);
93  fFitter[i]->SetNofBinsX(100); //fNofBinsX);
94  fFitter[i]->SetNofBinsY(100); //fNofBinsY);
95  fFitter[i]->SetUseEllipseAcc(fUseEllipseAcc);
96  }
97 
98  // Create grid creator tool
100 
101  fField = FairRunAna::Instance()->GetField();
102 
103  fHM = new CbmHistManager();
104 
105  // Fill and create graphs and histograms
106  CreateHistos();
107 
108  // This is always needed since field map is always drawn for comparison
109  FillBHistos();
110 
111  // Check and draw polynomial field approximation
113 
114  // Check and draw histograms for grid creator tool
116 
118  report->Create(fHM, fOutputDir);
119  delete report;
120 
121  fHM->WriteToFile();
122 
123  return kSUCCESS;
124 }
125 
126 void CbmLitFieldApproximationQa::Exec(Option_t* opt) {}
127 
129 
134  cout << fHM->ToString();
135 }
136 
138  string names[] = {"Bx", "By", "Bz", "Mod"};
139  string zTitle[] = {
140  "B_{x} [kGauss]", "B_{y} [kGauss]", "B_{z} [kGauss]", "|B| [kGauss]"};
141  for (Int_t v = 0; v < 4; v++) {
142  for (Int_t i = 0; i < fNofSlices; i++) {
143  TGraph2D* graph = new TGraph2D();
144  string name =
145  "hfa_" + names[v] + "_Graph2D_" + ToString<Int_t>(fZSlicePosition[i]);
146  string title = name + ";X [cm];Y [cm];" + zTitle[v];
147  graph->SetNameTitle(name.c_str(), title.c_str());
148  fHM->Add(name, graph);
149  }
150  }
151 }
152 
154  string names[] = {"Bx", "By", "Bz", "Mod"};
155 
156  Int_t nofBinsX = fNofBinsX;
157  Int_t nofBinsY = fNofBinsY;
158  Int_t nofBinsErrB = 100;
159  Int_t nofBinsRelErrB = 100;
160  Int_t nofBinsErrX = 100;
161  Int_t nofBinsErrY = 100;
162  Double_t minErrB = 0., maxErrB = 0., minRelErrB = 0., maxRelErrB = 0.;
163  if (fFixedBounds) {
164  minErrB = -0.5;
165  maxErrB = 0.5;
166  minRelErrB = -10.;
167  maxRelErrB = 10.;
168  }
169 
170  string zTitle[] = {
171  "B_{x} [kGauss]", "B_{y} [kGauss]", "B_{z} [kGauss]", "|B| [kGauss]"};
172  string errTitle[] = {"B_{x} error [kGauss]",
173  "B_{y} error [kGauss]",
174  "B_{z} error [kGauss]",
175  "|B| error [kGauss]"};
176  string relErrTitle[] = {"B_{x} relative error [%]",
177  "B_{y} relative error [%]",
178  "B_{z} relative error [%]",
179  "|B| relative error [%]"};
180 
181  // Create histograms
182  for (Int_t v = 0; v < 4; v++) {
183  for (Int_t i = 0; i < fNofSlices; i++) {
184  for (Int_t j = 0; j < fNofPolynoms; j++) {
185  TGraph2D* graph = new TGraph2D();
186  string name = "hfa_" + names[v] + "Apr_Graph2D" + "_"
187  + ToString<Int_t>(fZSlicePosition[i]) + "_"
188  + ToString<Int_t>(fPolynomDegrees[j]);
189  string title = name + ";X [cm];Y [cm];" + zTitle[v];
190  graph->SetNameTitle(name.c_str(), title.c_str());
191  fHM->Add(name, graph);
192 
193  name = "hfa_" + names[v] + "ErrApr_H1_"
194  + ToString<Int_t>(fZSlicePosition[i]) + "_"
195  + ToString<Int_t>(fPolynomDegrees[j]);
196  title = name + ";" + errTitle[v] + ";Counter";
197  fHM->Add(
198  name,
199  new TH1D(name.c_str(), title.c_str(), nofBinsErrB, minErrB, maxErrB));
200 
201  name = "hfa_" + names[v] + "ErrApr_H2_"
202  + ToString<Int_t>(fZSlicePosition[i]) + "_"
203  + ToString<Int_t>(fPolynomDegrees[j]);
204  title = name + ";X [cm];Y [cm];" + errTitle[v];
205  fHM->Add(name,
206  new TH2D(name.c_str(),
207  title.c_str(),
208  nofBinsErrX,
209  -fXSlicePosition[i],
211  nofBinsErrY,
212  -fYSlicePosition[i],
213  fYSlicePosition[i]));
214 
215  name = "hfa_" + names[v] + "RelErrApr_H1_"
216  + ToString<Int_t>(fZSlicePosition[i]) + "_"
217  + ToString<Int_t>(fPolynomDegrees[j]);
218  title = name + ";" + relErrTitle[v] + ";Counter";
219  fHM->Add(name,
220  new TH1D(name.c_str(),
221  title.c_str(),
222  nofBinsRelErrB,
223  minRelErrB,
224  maxRelErrB));
225 
226  name = "hfa_" + names[v] + "RelErrApr_H2_"
227  + ToString<Int_t>(fZSlicePosition[i]) + "_"
228  + ToString<Int_t>(fPolynomDegrees[j]);
229  title = name + ";X [cm];Y [cm];" + relErrTitle[v];
230  fHM->Add(name,
231  new TH2D(name.c_str(),
232  title.c_str(),
233  nofBinsErrX,
234  -fXSlicePosition[i],
236  nofBinsErrY,
237  -fYSlicePosition[i],
238  fYSlicePosition[i]));
239  }
240  }
241  }
242  cout << "-I- CbmLitFieldApproximationQa::CreateFitterErrHistos: Field fitter "
243  "error histograms created"
244  << endl;
245 }
246 
248  string names[] = {"Bx", "By", "Bz", "Mod"};
249 
250  Int_t nofBinsX = fNofBinsX;
251  Int_t nofBinsY = fNofBinsY;
252  Int_t nofBinsErrB = 100;
253  Int_t nofBinsRelErrB = 100;
254  Int_t nofBinsErrX = 100;
255  Int_t nofBinsErrY = 100;
256  Double_t minErrB = 0., maxErrB = 0., minRelErrB = 0., maxRelErrB = 0.;
257  if (fFixedBounds) {
258  minErrB = -0.5;
259  maxErrB = 0.5;
260  minRelErrB = -10.;
261  maxRelErrB = 10.;
262  }
263 
264  string zTitle[] = {
265  "B_{x} [kGauss]", "B_{y} [kGauss]", "B_{z} [kGauss]", "|B| [kGauss]"};
266  string errTitle[] = {"B_{x} error [kGauss]",
267  "B_{y} error [kGauss]",
268  "B_{z} error [kGauss]",
269  "|B| error [kGauss]"};
270  string relErrTitle[] = {"B_{x} relative error [%]",
271  "B_{y} relative error [%]",
272  "B_{z} relative error [%]",
273  "|B| relative error [%]"};
274 
275  // Create histograms
276  for (Int_t v = 0; v < 4; v++) {
277  for (Int_t i = 0; i < fNofSlices; i++) {
278  TGraph2D* graph = new TGraph2D();
279  string name = "hfa_" + names[v] + "Grid_Graph2D_"
280  + ToString<Int_t>(fZSlicePosition[i]);
281  string title = name + ";X [cm]; Y [cm];" + zTitle[v];
282  graph->SetNameTitle(name.c_str(), title.c_str());
283  fHM->Add(name, graph);
284 
285  name =
286  "hfa_" + names[v] + "ErrGrid_H1_" + ToString<Int_t>(fZSlicePosition[i]);
287  title = name + ";" + errTitle[v] + ";Counter";
288  fHM->Add(
289  name,
290  new TH1D(name.c_str(), title.c_str(), nofBinsErrB, minErrB, maxErrB));
291 
292  name =
293  "hfa_" + names[v] + "ErrGrid_H2_" + ToString<Int_t>(fZSlicePosition[i]);
294  title = name + ";X [cm];Y [cm];" + errTitle[v];
295  fHM->Add(name,
296  new TH2D(name.c_str(),
297  title.c_str(),
298  nofBinsErrX,
299  -fXSlicePosition[i],
301  nofBinsErrY,
302  -fYSlicePosition[i],
303  fYSlicePosition[i]));
304 
305  name = "hfa_" + names[v] + "RelErrGrid_H1_"
306  + ToString<Int_t>(fZSlicePosition[i]);
307  title = name + ";" + relErrTitle[v] + ";Counter";
308  fHM->Add(
309  name,
310  new TH1D(
311  name.c_str(), title.c_str(), nofBinsRelErrB, minRelErrB, maxRelErrB));
312 
313  name = "hfa_" + names[v] + "RelErrGrid_H2_"
314  + ToString<Int_t>(fZSlicePosition[i]);
315  title = name + ";X [cm];Y [cm];" + relErrTitle[v];
316  fHM->Add(name,
317  new TH2D(name.c_str(),
318  title.c_str(),
319  nofBinsErrX,
320  -fXSlicePosition[i],
322  nofBinsErrY,
323  -fYSlicePosition[i],
324  fYSlicePosition[i]));
325  }
326  }
327  cout << "-I- CbmLitFieldApproximationQa::CreateGridErrHistos(): Grid creator "
328  "error histograms created"
329  << endl;
330 }
331 
333  // Fill graphs for magnetic field for each (X, Y) slice
334  for (UInt_t iSlice = 0; iSlice < fNofSlices; iSlice++) { // loop over slices
335  Double_t Z = fZSlicePosition[iSlice];
336 
337  Int_t cnt = 0;
338  Double_t HX =
339  2 * fXSlicePosition[iSlice] / fNofBinsX; // step size for X position
340  Double_t HY =
341  2 * fYSlicePosition[iSlice] / fNofBinsY; // step size for Y position
342  for (Int_t iX = 0; iX < fNofBinsX; iX++) { // loop over x position
343  Double_t X = -fXSlicePosition[iSlice] + (iX + 0.5) * HX;
344  for (Int_t iY = 0; iY < fNofBinsY; iY++) { // loop over y position
345  Double_t Y = -fYSlicePosition[iSlice] + (iY + 0.5) * HY;
346 
347  // Check acceptance for ellipse
348  Double_t el =
349  (X * X) / (fXSlicePosition[iSlice] * fXSlicePosition[iSlice])
350  + (Y * Y) / (fYSlicePosition[iSlice] * fYSlicePosition[iSlice]);
351  if (fUseEllipseAcc && el > 1.) { continue; }
352 
353  // Get field value
354  Double_t pos[3] = {X, Y, Z};
355  Double_t B[3];
356  fField->GetFieldValue(pos, B);
357 
358  Double_t Bmod = sqrt(B[0] * B[0] + B[1] * B[1] + B[2] * B[2]);
359 
360  string s = ToString<Int_t>(fZSlicePosition[iSlice]);
361  fHM->G2(string("hfa_Bx_Graph2D_") + s)->SetPoint(cnt, X, Y, B[0]);
362  fHM->G2(string("hfa_By_Graph2D_") + s)->SetPoint(cnt, X, Y, B[1]);
363  fHM->G2(string("hfa_Bz_Graph2D_") + s)->SetPoint(cnt, X, Y, B[2]);
364  fHM->G2(string("hfa_Mod_Graph2D_") + s)->SetPoint(cnt, X, Y, Bmod);
365  cnt++;
366  }
367  }
368  }
369 }
370 
372  vector<vector<LitFieldSliceScal>> slices;
373  slices.resize(fNofPolynoms);
374  for (UInt_t i = 0; i < fNofPolynoms; i++) {
375  slices[i].resize(fNofSlices);
376  }
377 
378  // Approximate field in each slice for each polynomial degree
379  for (UInt_t i = 0; i < fNofPolynoms; i++) {
380  for (Int_t j = 0; j < fNofSlices; j++) {
381  fFitter[i]->FitSlice(fZSlicePosition[j], slices[i][j]);
382  cout << "-I- CbmLitFieldApproximationQa::FillFieldApproximationHistos: "
383  << " field approximation (degree=" << fPolynomDegrees[i]
384  << ", Z=" << fZSlicePosition[j] << ")" << endl;
385  }
386  }
387 
388  // Fill graph for approximated field map
389  for (UInt_t iSlice = 0; iSlice < fNofSlices; iSlice++) { // Loop over slices
390  Double_t Z = fZSlicePosition[iSlice];
391  Int_t cnt = 0;
392 
393  Double_t HX =
394  2 * fXSlicePosition[iSlice] / fNofBinsX; // Step size for X position
395  Double_t HY =
396  2 * fYSlicePosition[iSlice] / fNofBinsY; // Step size for Y position
397  for (Int_t iX = 0; iX < fNofBinsX; iX++) { // Loop over x position
398  Double_t X = -fXSlicePosition[iSlice] + (iX + 0.5) * HX;
399  for (Int_t iY = 0; iY < fNofBinsY; iY++) { // Loop over y position
400  Double_t Y = -fYSlicePosition[iSlice] + (iY + 0.5) * HY;
401 
402  // Check acceptance for ellipse
403  Double_t el =
404  (X * X) / (fXSlicePosition[iSlice] * fXSlicePosition[iSlice])
405  + (Y * Y) / (fYSlicePosition[iSlice] * fYSlicePosition[iSlice]);
406  if (fUseEllipseAcc && el > 1.) { continue; }
407 
408  for (Int_t p = 0; p < fNofPolynoms; p++) {
410  slices[p][iSlice].GetFieldValue(X, Y, v);
411  Double_t mod = sqrt(v.Bx * v.Bx + v.By * v.By + v.Bz * v.Bz);
412  string s = ToString<Int_t>(fZSlicePosition[iSlice]) + "_"
413  + ToString<Int_t>(fPolynomDegrees[p]);
414  fHM->G2(string("hfa_BxApr_Graph2D_") + s)->SetPoint(cnt, X, Y, v.Bx);
415  fHM->G2(string("hfa_ByApr_Graph2D_") + s)->SetPoint(cnt, X, Y, v.By);
416  fHM->G2(string("hfa_BzApr_Graph2D_") + s)->SetPoint(cnt, X, Y, v.Bz);
417  fHM->G2(string("hfa_ModApr_Graph2D_") + s)->SetPoint(cnt, X, Y, mod);
418  }
419  cnt++;
420  } // End loop over y position
421  } // End loop over x position
422  } // End loop over slices
423 
424  // Fill error histograms
425  Int_t nofBinsX = 100;
426  Int_t nofBinsY = 100;
427  for (Int_t iSlice = 0; iSlice < fNofSlices; iSlice++) {
428  Double_t Z = fZSlicePosition[iSlice];
429  Double_t HX =
430  2 * fXSlicePosition[iSlice] / nofBinsX; // step size for X position
431  Double_t HY =
432  2 * fYSlicePosition[iSlice] / nofBinsY; // step size for Y position
433  for (Int_t iX = 0; iX < nofBinsX; iX++) { // loop over x position
434  Double_t X = -fXSlicePosition[iSlice] + (iX + 0.5) * HX;
435  for (Int_t iY = 0; iY < nofBinsY; iY++) { // loop over y position
436  Double_t Y = -fYSlicePosition[iSlice] + (iY + 0.5) * HY;
437 
438  // Check acceptance for ellipse
439  Double_t el =
440  (X * X) / (fXSlicePosition[iSlice] * fXSlicePosition[iSlice])
441  + (Y * Y) / (fYSlicePosition[iSlice] * fYSlicePosition[iSlice]);
442  if (fUseEllipseAcc && el > 1.) { continue; }
443 
444  // Get field value
445  Double_t pos[3] = {X, Y, Z};
446  Double_t B[3];
447  fField->GetFieldValue(pos, B);
448 
449  Double_t Bmod = sqrt(B[0] * B[0] + B[1] * B[1] + B[2] * B[2]);
450 
451  for (Int_t p = 0; p < fNofPolynoms; p++) {
453  slices[p][iSlice].GetFieldValue(X, Y, v);
454  Double_t mod = sqrt(v.Bx * v.Bx + v.By * v.By + v.Bz * v.Bz);
455 
456  Double_t errBx = B[0] - v.Bx;
457  Double_t errBy = B[1] - v.By;
458  Double_t errBz = B[2] - v.Bz;
459  Double_t errMod = Bmod - mod;
460  Double_t relErrBx = (B[0] != 0.) ? (errBx / B[0]) * 100. : 0.;
461  Double_t relErrBy = (B[1] != 0.) ? (errBy / B[1]) * 100. : 0.;
462  Double_t relErrBz = (B[2] != 0.) ? (errBz / B[2]) * 100. : 0.;
463  Double_t relErrMod = (Bmod != 0.) ? (errMod / Bmod) * 100. : 0;
464 
465  string s = ToString<Int_t>(fZSlicePosition[iSlice]) + "_"
466  + ToString<Int_t>(fPolynomDegrees[p]);
467  fHM->H1(string("hfa_BxErrApr_H1_") + s)->Fill(errBx);
468  fHM->H1(string("hfa_BxRelErrApr_H1_") + s)->Fill(relErrBx);
469  fHM->H2(string("hfa_BxErrApr_H2_") + s)->Fill(X, Y, errBx);
470  fHM->H2(string("hfa_BxRelErrApr_H2_") + s)->Fill(X, Y, relErrBx);
471  fHM->H1(string("hfa_ByErrApr_H1_") + s)->Fill(errBy);
472  fHM->H1(string("hfa_ByRelErrApr_H1_") + s)->Fill(relErrBy);
473  fHM->H2(string("hfa_ByErrApr_H2_") + s)->Fill(X, Y, errBy);
474  fHM->H2(string("hfa_ByRelErrApr_H2_") + s)->Fill(X, Y, relErrBy);
475  fHM->H1(string("hfa_BzErrApr_H1_") + s)->Fill(errBz);
476  fHM->H1(string("hfa_BzRelErrApr_H1_") + s)->Fill(relErrBz);
477  fHM->H2(string("hfa_BzErrApr_H2_") + s)->Fill(X, Y, errBz);
478  fHM->H2(string("hfa_BzRelErrApr_H2_") + s)->Fill(X, Y, relErrBz);
479  fHM->H1(string("hfa_ModErrApr_H1_") + s)->Fill(errMod);
480  fHM->H1(string("hfa_ModRelErrApr_H1_") + s)->Fill(relErrMod);
481  fHM->H2(string("hfa_ModErrApr_H2_") + s)->Fill(X, Y, errMod);
482  fHM->H2(string("hfa_ModRelErrApr_H2_") + s)->Fill(X, Y, relErrMod);
483  }
484  }
485  }
486  }
487 }
488 
490  vector<LitFieldGrid> grids;
491  grids.resize(fNofSlices);
492  for (Int_t iSlice = 0; iSlice < fNofSlices; iSlice++) {
493  fGridCreator->CreateGrid(fZSlicePosition[iSlice], grids[iSlice]);
494  }
495 
496  // Fill graph
497  for (Int_t iSlice = 0; iSlice < fNofSlices; iSlice++) {
498  Int_t cnt = 0;
499  Double_t Z = fZSlicePosition[iSlice];
500  Double_t HX =
501  2 * fXSlicePosition[iSlice] / fNofBinsX; // step size for X position
502  Double_t HY =
503  2 * fYSlicePosition[iSlice] / fNofBinsY; // step size for Y position
504  for (Int_t iX = 0; iX < fNofBinsX; iX++) { // loop over x position
505  Double_t X = -fXSlicePosition[iSlice] + (iX + 0.5) * HX;
506  for (Int_t iY = 0; iY < fNofBinsY; iY++) { // loop over y position
507  Double_t Y = -fYSlicePosition[iSlice] + (iY + 0.5) * HY;
509  grids[iSlice].GetFieldValue(X, Y, v);
510  Double_t mod = sqrt(v.Bx * v.Bx + v.By * v.By + v.Bz * v.Bz);
511  string s = ToString<Int_t>(fZSlicePosition[iSlice]);
512  fHM->G2(string("hfa_BxGrid_Graph2D_") + s)->SetPoint(cnt, X, Y, v.Bx);
513  fHM->G2(string("hfa_ByGrid_Graph2D_") + s)->SetPoint(cnt, X, Y, v.By);
514  fHM->G2(string("hfa_BzGrid_Graph2D_") + s)->SetPoint(cnt, X, Y, v.Bz);
515  fHM->G2(string("hfa_ModGrid_Graph2D_") + s)->SetPoint(cnt, X, Y, mod);
516  cnt++;
517  }
518  }
519  }
520 
521  // Fill error histograms
522  Int_t nofBinsX = 100;
523  Int_t nofBinsY = 100;
524  for (Int_t iSlice = 0; iSlice < fNofSlices; iSlice++) {
525  Double_t Z = fZSlicePosition[iSlice];
526  Double_t HX =
527  2 * fXSlicePosition[iSlice] / nofBinsX; // step size for X position
528  Double_t HY =
529  2 * fYSlicePosition[iSlice] / nofBinsY; // step size for Y position
530  for (Int_t iX = 0; iX < nofBinsX; iX++) { // loop over x position
531  Double_t X = -fXSlicePosition[iSlice] + (iX + 0.5) * HX;
532  for (Int_t iY = 0; iY < nofBinsY; iY++) { // loop over y position
533  Double_t Y = -fYSlicePosition[iSlice] + (iY + 0.5) * HY;
534 
535  // Get field value
536  Double_t pos[3] = {X, Y, Z};
537  Double_t B[3];
538  fField->GetFieldValue(pos, B);
539 
540  Double_t Bmod = sqrt(B[0] * B[0] + B[1] * B[1] + B[2] * B[2]);
541 
543  grids[iSlice].GetFieldValue(X, Y, v);
544  Double_t mod = sqrt(v.Bx * v.Bx + v.By * v.By + v.Bz * v.Bz);
545 
546  Double_t errBx = B[0] - v.Bx;
547  Double_t errBy = B[1] - v.By;
548  Double_t errBz = B[2] - v.Bz;
549  Double_t errMod = Bmod - mod;
550  Double_t relErrBx = (B[0] != 0.) ? (errBx / B[0]) * 100. : 0.;
551  Double_t relErrBy = (B[1] != 0.) ? (errBy / B[1]) * 100. : 0.;
552  Double_t relErrBz = (B[2] != 0.) ? (errBz / B[2]) * 100. : 0.;
553  Double_t relErrMod = (Bmod != 0.) ? (errMod / Bmod) * 100. : 0;
554 
555  string s = ToString<Int_t>(fZSlicePosition[iSlice]);
556  fHM->H1(string("hfa_BxErrGrid_H1_") + s)->Fill(errBx);
557  fHM->H1(string("hfa_BxRelErrGrid_H1_") + s)->Fill(relErrBx);
558  fHM->H2(string("hfa_BxErrGrid_H2_") + s)->Fill(X, Y, errBx);
559  fHM->H2(string("hfa_BxRelErrGrid_H2_") + s)->Fill(X, Y, relErrBx);
560  fHM->H1(string("hfa_ByErrGrid_H1_") + s)->Fill(errBy);
561  fHM->H1(string("hfa_ByRelErrGrid_H1_") + s)->Fill(relErrBy);
562  fHM->H2(string("hfa_ByErrGrid_H2_") + s)->Fill(X, Y, errBy);
563  fHM->H2(string("hfa_ByRelErrGrid_H2_") + s)->Fill(X, Y, relErrBy);
564  fHM->H1(string("hfa_BzErrGrid_H1_") + s)->Fill(errBz);
565  fHM->H1(string("hfa_BzRelErrGrid_H1_") + s)->Fill(relErrBz);
566  fHM->H2(string("hfa_BzErrGrid_H2_") + s)->Fill(X, Y, errBz);
567  fHM->H2(string("hfa_BzRelErrGrid_H2_") + s)->Fill(X, Y, relErrBz);
568  fHM->H1(string("hfa_ModErrGrid_H1_") + s)->Fill(errMod);
569  fHM->H1(string("hfa_ModRelErrGrid_H1_") + s)->Fill(relErrMod);
570  fHM->H2(string("hfa_ModErrGrid_H2_") + s)->Fill(X, Y, errMod);
571  fHM->H2(string("hfa_ModRelErrGrid_H2_") + s)->Fill(X, Y, relErrMod);
572  }
573  }
574  }
575 }
576 
CbmLitFieldApproximationQa::fXSlicePosition
vector< Double_t > fXSlicePosition
Definition: CbmLitFieldApproximationQa.h:126
CbmLitFieldApproximationQa::CreateFieldHistos
void CreateFieldHistos()
Create field histograms.
Definition: CbmLitFieldApproximationQa.cxx:137
CbmLitFieldApproximationQa::Init
virtual InitStatus Init()
Inherited from FairTask.
Definition: CbmLitFieldApproximationQa.cxx:73
CbmLitFieldApproximationQa::Finish
virtual void Finish()
Inherited from FairTask.
Definition: CbmLitFieldApproximationQa.cxx:128
CbmLitFieldApproximationQa::fUseEllipseAcc
Bool_t fUseEllipseAcc
Definition: CbmLitFieldApproximationQa.h:144
CbmLitFieldApproximationQa::fPolynomDegrees
vector< UInt_t > fPolynomDegrees
Definition: CbmLitFieldApproximationQa.h:150
CbmHistManager::ToString
std::string ToString() const
Return string representation of class.
Definition: core/base/CbmHistManager.cxx:258
CbmHistManager::WriteToFile
void WriteToFile()
Write all histograms to current opened file.
Definition: core/base/CbmHistManager.cxx:103
CbmLitFieldApproximationQa::FillBHistos
void FillBHistos()
Fill graphs and histos for field map for each field component (Bx, By, Bz, |B|).
Definition: CbmLitFieldApproximationQa.cxx:332
CbmLitFieldApproximationQaReport
Creates field QA report.
Definition: CbmLitFieldApproximationQaReport.h:20
CbmLitFieldApproximationQa::fHM
CbmHistManager * fHM
Definition: CbmLitFieldApproximationQa.h:160
ClassImp
ClassImp(CbmLitFieldApproximationQa)
sqrt
friend F32vec4 sqrt(const F32vec4 &a)
Definition: L1/vectors/P4_F32vec4.h:41
i
int i
Definition: L1/vectors/P4_F32vec4.h:25
CbmLitFieldApproximationQa::fYSlicePosition
vector< Double_t > fYSlicePosition
Definition: CbmLitFieldApproximationQa.h:127
CbmHistManager::H2
TH2 * H2(const std::string &name) const
Return pointer to TH2 histogram.
Definition: CbmHistManager.h:190
CbmLitFieldApproximationQa::fNofBinsX
Int_t fNofBinsX
Definition: CbmLitFieldApproximationQa.h:140
CbmLitFieldApproximationQa::fNofSlices
Int_t fNofSlices
Definition: CbmLitFieldApproximationQa.h:119
CbmLitFieldGridCreator
Definition: CbmLitFieldGridCreator.h:17
CbmDrawHist.h
Helper functions for drawing 1D and 2D histograms and graphs.
CbmLitFieldApproximationQa::fNofBinsY
Int_t fNofBinsY
Definition: CbmLitFieldApproximationQa.h:142
CbmHistManager.h
Histogram manager.
CbmLitFieldApproximationQa::FillFieldApproximationHistos
void FillFieldApproximationHistos()
Fill histograms for polynomial field approximation.
Definition: CbmLitFieldApproximationQa.cxx:371
CbmLitFieldApproximationQa::Exec
virtual void Exec(Option_t *opt)
Inherited from FairTask.
Definition: CbmLitFieldApproximationQa.cxx:126
lit::parallel::LitFieldGrid
Class stores a grid of magnetic field values in XY slice at Z position.
Definition: LitFieldGrid.h:44
CbmLitFieldFitter.h
Implementation of the polynomial field approximation.
CbmLitFieldApproximationQa::fNofPolynoms
UInt_t fNofPolynoms
Definition: CbmLitFieldApproximationQa.h:148
CbmHistManager
Histogram manager.
Definition: CbmHistManager.h:41
lit::parallel::LitFieldSliceScal
LitFieldSlice< fscal > LitFieldSliceScal
Scalar version of LitFieldSlice.
Definition: LitFieldSlice.h:389
CbmLitFieldApproximationQa::CbmLitFieldApproximationQa
CbmLitFieldApproximationQa()
Constructor.
Definition: CbmLitFieldApproximationQa.cxx:51
CbmLitFieldApproximationQa::~CbmLitFieldApproximationQa
virtual ~CbmLitFieldApproximationQa()
Destructor.
Definition: CbmLitFieldApproximationQa.cxx:71
CbmHistManager::G2
TGraph2D * G2(const std::string &name) const
Return pointer to TGraph2D.
Definition: CbmHistManager.h:243
CbmLitFieldApproximationQa::fField
FairField * fField
Definition: CbmLitFieldApproximationQa.h:116
CbmSimulationReport::Create
void Create(CbmHistManager *histManager, const std::string &outputDir)
Main function which creates report data.
CbmLitFieldApproximationQa::fOutputDir
string fOutputDir
Definition: CbmLitFieldApproximationQa.h:130
CbmLitFieldApproximationQa
Field map approximation QA.
Definition: CbmLitFieldApproximationQa.h:34
CbmLitFieldApproximationQa::fGridCreator
CbmLitFieldGridCreator * fGridCreator
Definition: CbmLitFieldApproximationQa.h:158
CbmLitFieldApproximationQa::CreateGridHistos
void CreateGridHistos()
Create histograms for grid creator.
Definition: CbmLitFieldApproximationQa.cxx:247
CbmLitFloat.h
Define floating point number type litfloat.
CbmHistManager::H1
TH1 * H1(const std::string &name) const
Return pointer to TH1 histogram.
Definition: CbmHistManager.h:170
CbmLitFieldGridCreator.h
Class creates grid with magnetic field values at a certain Z position.
CbmLitFieldApproximationQaReport.h
Creates field QA report.
CbmUtils.h
CbmLitFieldApproximationQa::CreateFitterHistos
void CreateFitterHistos()
Create histograms for field approximation.
Definition: CbmLitFieldApproximationQa.cxx:153
CbmLitFieldApproximationQa::fAcceptanceAngleY
Double_t fAcceptanceAngleY
Definition: CbmLitFieldApproximationQa.h:138
CbmLitFieldApproximationQa::fFitter
vector< CbmLitFieldFitter * > fFitter
Definition: CbmLitFieldApproximationQa.h:152
v
__m128 v
Definition: L1/vectors/P4_F32vec4.h:1
CbmLitFieldApproximationQa::FillGridCreatorHistos
void FillGridCreatorHistos()
fill histograms for grid creator tool.
Definition: CbmLitFieldApproximationQa.cxx:489
CbmLitFieldFitter
Implementation of the polynomial field approximation.
Definition: CbmLitFieldFitter.h:60
pos
TVector3 pos
Definition: CbmMvdSensorDigiToHitTask.cxx:60
CbmSimulationReport
Base class for simulation reports.
Definition: CbmSimulationReport.h:28
Cbm::SaveCanvasAsImage
void SaveCanvasAsImage(TCanvas *c, const std::string &dir, const std::string &option)
Definition: CbmUtils.cxx:20
CbmLitFieldApproximationQa::fFixedBounds
Bool_t fFixedBounds
Definition: CbmLitFieldApproximationQa.h:162
lit::parallel::LitFieldValueScal
LitFieldValue< fscal > LitFieldValueScal
Scalar version of LitFieldValue.
Definition: LitFieldValue.h:64
CbmLitFieldApproximationQa::fZSlicePosition
vector< Double_t > fZSlicePosition
Definition: CbmLitFieldApproximationQa.h:125
CbmLitFieldApproximationQa::CreateHistos
void CreateHistos()
Create histograms.
Definition: CbmLitFieldApproximationQa.cxx:130
Cbm::ToString
std::string ToString(const T &value)
Definition: CbmUtils.h:16
CbmLitFieldApproximationQa.h
Field map approximation QA.
CbmLitFieldApproximationQa::fAcceptanceAngleX
Double_t fAcceptanceAngleX
Definition: CbmLitFieldApproximationQa.h:136
CbmLitFieldGridCreator::CreateGrid
void CreateGrid(fscal Z, lit::parallel::LitFieldGrid &grid)
Main function which creates grid with magnetic field values in (X, Y) slice.
Definition: CbmLitFieldGridCreator.cxx:23
CbmHistManager::Add
void Add(const std::string &name, TNamed *object)
Add new named object to manager.
Definition: CbmHistManager.h:58