Belle II Software  release-06-00-14
CreateFieldMapModule.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #include <geometry/modules/CreateFieldMapModule.h>
10 #include <framework/geometry/BFieldManager.h>
11 #include <framework/utilities/Utils.h>
12 #include <framework/gearbox/Unit.h>
13 #include <framework/geometry/B2Vector3.h>
14 
15 #include <TFile.h>
16 #include <TH2D.h>
17 #include <TTree.h>
18 
19 #include <boost/algorithm/string.hpp>
20 #include <boost/format.hpp>
21 
22 using namespace Belle2;
23 
24 //-----------------------------------------------------------------
25 // Register the Module
26 //-----------------------------------------------------------------
27 REG_MODULE(CreateFieldMap)
28 
29 //-----------------------------------------------------------------
30 // Implementation
31 //-----------------------------------------------------------------
33 {
34  // Set module properties
35  setDescription("Create Field maps of the Belle II magnetic field used in the simulation");
36 
37  // Parameter definitions
38  addParam("filename", m_filename, "ROOT filename for the map", std::string("FieldMap.root"));
39  addParam("type", m_type, "type of the fieldmap (xy, zx, zy, zr)", std::string("zx"));
40  addParam("nU", m_nU, "number of steps along first coordinate", 1800);
41  addParam("minU", m_minU, "minimum value for the first coordinate", -400.);
42  addParam("maxU", m_maxU, "maximum value for the first coordinate", 500.);
43  addParam("nV", m_nV, "number of steps along second coordinate", 1600);
44  addParam("minV", m_minV, "minimum value for the second coordinate", -400.);
45  addParam("maxV", m_maxV, "maximum value for the second coordinate", 400.);
46  addParam("phi", m_phi, "phi angle for the scan in radians", 0.);
47  addParam("wOffset", m_wOffset, "value of third coordinate", 0.);
48  addParam("nPhi", m_nPhi, "number of phi steps for zr averaging", m_nPhi);
49  addParam("saveAllPoints", m_createTree, "save all sampled points in a TTree, "
50  "WARNING: output can be huge", false);
51 }
52 
54 {
55  if (m_filename.empty()) {
56  B2WARNING("CreateFieldMap: No Filename given, just sampling for fun");
57  }
58  boost::to_lower(m_type);
59  boost::trim(m_type);
60  if (m_type != "xy" && m_type != "zx" && m_type != "zy" && m_type != "zr") {
61  B2ERROR("CreateFieldMap type '" << m_type << "' not valid, use one of 'xy', 'zx', 'zy' or 'zr'");
62  }
63  if (m_nU < 0 || m_nV < 0) {
64  B2ERROR("Number of steps has to be positive");
65  }
66  if (m_maxU == m_minU) {
67  B2ERROR("Range for first coordinate is zero");
68  } else if (m_maxU < m_minU) {
69  std::swap(m_maxU, m_minU);
70  }
71  if (m_maxV == m_minV) {
72  B2ERROR("Range for second coordinate is zero");
73  } else if (m_maxV < m_minV) {
74  std::swap(m_maxV, m_minV);
75  }
76  if (m_type == "zr") {
77  if (m_minV < 0) {
78  B2ERROR("R values (minV, maxV) must be positive");
79  }
80  if (m_nPhi <= 0) {
81  B2ERROR("Number of steps in phi must be positive");
82  }
83  }
84 }
85 
87 {
88  const bool save{!m_filename.empty()};
89  TFile* outfile{nullptr};
90  if (save) {
91  //Create histograms
92  outfile = new TFile(m_filename.c_str(), "RECREATE");
93  outfile->cd();
94  }
95 
96  //Determine type of scan
97  EFieldTypes type = c_XY;
98  if (m_type == "zx")
99  type = c_ZX;
100  else if (m_type == "zy")
101  type = c_ZY;
102  else if (m_type == "zr")
103  type = c_ZR;
104 
105  //some values needed for output
106  int lastPercent(-1);
107  uint64_t nSteps = m_nU * m_nV;
108  uint64_t curStep{0};
109  double startTime{0};
110  auto showProgress = [&]() {
111  if (curStep == 0) startTime = Utils::getClock();
112  const int64_t donePercent = 100 * ++curStep / nSteps;
113  if (donePercent > lastPercent) {
114  const double totalTime = Utils::getClock() - startTime;
115  const double perStep = totalTime / curStep;
116  if (donePercent == 100) {
117  B2INFO(boost::format("BField %s Scan: %d samples, %.3f us per sample, total: %.2f seconds")
118  % m_type % curStep
119  % (perStep / Unit::us) % (totalTime / Unit::s));
120  } else {
121  B2INFO(boost::format("BField %s Scan: %3d%%, %.3f us per sample")
122  % m_type % donePercent
123  % (perStep / Unit::us));
124  }
125  lastPercent = donePercent;
126  }
127  };
128 
129  struct { float x{0}, y{0}, z{0}, bx{0}, by{0}, bz{0}; } field_point;
130  TTree* all_values{nullptr};
131  if (save && m_createTree) {
132  all_values = new TTree("bfield_values", "All B field values");
133  all_values->Branch("x", &field_point.x, "x/F");
134  all_values->Branch("y", &field_point.y, "y/F");
135  all_values->Branch("z", &field_point.z, "z/F");
136  all_values->Branch("bx", &field_point.bx, "bx/F");
137  all_values->Branch("by", &field_point.by, "by/F");
138  all_values->Branch("bz", &field_point.bz, "bz/F");
139  }
140  auto fillTree = [&](const B2Vector3D & p, const B2Vector3D & b) {
141  if (!all_values) return;
142  field_point.x = p.X();
143  field_point.y = p.Y();
144  field_point.z = p.Z();
145  field_point.bx = b.X();
146  field_point.by = b.Y();
147  field_point.bz = b.Z();
148  all_values->Fill();
149  };
150 
151  if (type == c_ZR) {
152  nSteps *= m_nPhi;
153  TH2D* h_b = new TH2D("B", "$B$ average;$z$/cm;$r$/cm", m_nU, m_minU, m_maxU, m_nV, m_minV, m_maxV);
154  TH2D* h_br = new TH2D("Br", "$B_r$ average;$z$/cm;$r$/cm", m_nU, m_minU, m_maxU, m_nV, m_minV, m_maxV);
155  TH2D* h_bz = new TH2D("Bz", "$B_z$ average;$z$/cm;$r$/cm", m_nU, m_minU, m_maxU, m_nV, m_minV, m_maxV);
156  for (int iU = 0; iU < m_nU; ++iU) {
157  for (int iV = 0; iV < m_nV; ++iV) {
158  //find value for first and second coordinate
159  const double u = h_b->GetXaxis()->GetBinCenter(iU + 1);
160  const double v = h_b->GetYaxis()->GetBinCenter(iV + 1);
161  //Determine global coordinates
162  for (int iPhi = 0; iPhi < m_nPhi; ++iPhi) {
163  B2Vector3D pos(v, 0, u);
164  pos.RotateZ(2 * M_PI * iPhi / m_nPhi);
165  //Obtain magnetic field
167  //And fill histograms
168  if (save) {
169  h_br->Fill(u, v, bfield.Perp());
170  h_bz->Fill(u, v, bfield.Z());
171  h_b->Fill(u, v, bfield.Mag());
172  fillTree(pos, bfield);
173  }
174  showProgress();
175  }
176  }
177  }
178  //Write histograms and close file.
179  if (save) {
180  for (TH2D* h : {h_br, h_bz, h_b}) {
181  h->Scale(1. / m_nPhi);
182  h->Write();
183  }
184  }
185  } else {
186  const std::string nu = m_type.substr(0, 1);
187  const std::string nv = m_type.substr(1, 1);
188  TH2D* h_b = new TH2D("B", ("$B$;$" + nu + "$/cm;$" + nv + "$/cm").c_str(), m_nU, m_minU, m_maxU, m_nV, m_minV, m_maxV);
189  TH2D* h_bx = new TH2D("Bx", ("$B_x$;$" + nu + "$/cm;$" + nv + "$/cm").c_str(), m_nU, m_minU, m_maxU, m_nV, m_minV, m_maxV);
190  TH2D* h_by = new TH2D("By", ("$B_y$;$" + nu + "$/cm;$" + nv + "$/cm").c_str(), m_nU, m_minU, m_maxU, m_nV, m_minV, m_maxV);
191  TH2D* h_bz = new TH2D("Bz", ("$B_z$;$" + nu + "$/cm;$" + nv + "$/cm").c_str(), m_nU, m_minU, m_maxU, m_nV, m_minV, m_maxV);
192 
193  //Loop over all bins
194  for (int iU = 0; iU < m_nU; ++iU) {
195  for (int iV = 0; iV < m_nV; ++iV) {
196  B2Vector3D pos(0, 0, 0);
197  //find value for first and second coordinate
198  const double u = h_bx->GetXaxis()->GetBinCenter(iU + 1);
199  const double v = h_bx->GetYaxis()->GetBinCenter(iV + 1);
200  //Determine global coordinates
201  switch (type) {
202  case c_XY:
203  pos.SetXYZ(u, v , m_wOffset);
204  break;
205  case c_ZX:
206  pos.SetXYZ(v, m_wOffset, u);
207  break;
208  case c_ZY:
209  pos.SetXYZ(m_wOffset, v, u);
210  break;
211  default:
212  break;
213  }
214  pos.RotateZ(m_phi);
215  //Obtain magnetic field
217  //And fill histograms
218  if (save) {
219  h_bx->Fill(u, v, bfield.X());
220  h_by->Fill(u, v, bfield.Y());
221  h_bz->Fill(u, v, bfield.Z());
222  h_b->Fill(u, v, bfield.Mag());
223  fillTree(pos, bfield);
224  }
225  showProgress();
226  }
227  }
228  //Write histograms.
229  if (save) {
230  for (TH2D* h : {h_bx, h_by, h_bz, h_b}) {
231  h->Write();
232  }
233  }
234  }
235  if (save) {
236  if (all_values) all_values->Write();
237  outfile->Close();
238  delete outfile;
239  }
240  //histograms seem to be deleted when file is closed
241 }
DataType Z() const
access variable Z (= .at(2) without boundary check)
Definition: B2Vector3.h:420
DataType X() const
access variable X (= .at(0) without boundary check)
Definition: B2Vector3.h:416
DataType Y() const
access variable Y (= .at(1) without boundary check)
Definition: B2Vector3.h:418
DataType Mag() const
The magnitude (rho in spherical coordinate system).
Definition: B2Vector3.h:144
DataType Perp() const
The transverse component (R in cylindrical coordinate system).
Definition: B2Vector3.h:185
static B2Vector3D getFieldInTesla(const B2Vector3D &pos)
return the magnetic field at a given position in Tesla.
Definition: BFieldManager.h:70
Create Field maps of the Belle II magnetic field used in the simulation.
int m_nPhi
number of steps in phi
virtual void initialize() override
Check input parameters.
double m_minV
start value for the first coordinate
double m_wOffset
offset on the third coordinate
double m_minU
start value for the first coordinate
double m_phi
phi rotation when sampling magnetic field (value of pi will convert ZX to ZY scan)
int m_nV
number of steps along the second coordinate
int m_nU
number of steps along the first coordinate
std::string m_type
type of the fieldmap (zx, zy, zr)
double m_maxV
end value for the first coordinate
virtual void beginRun() override
Create the fieldmap.
std::string m_filename
output filename for the fieldmap
EFieldTypes
Types of Fieldmap to be created.
@ c_ZR
scan along Z and R, averaging over phi
double m_maxU
end value for the first coordinate
bool m_createTree
if true create a TTree with all sampled points
Base class for Modules.
Definition: Module.h:72
static const double us
[microsecond]
Definition: Unit.h:97
static const double s
[second]
Definition: Unit.h:95
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
double getClock()
Return current value of the real-time clock.
Definition: Utils.cc:65
Abstract base class for different kinds of events.