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