OPALX (Object Oriented Parallel Accelerator Library for Exascale) MINIorX
OPALX
DumpEMFields.cpp
Go to the documentation of this file.
1//
2// Class DumpEMFields
3// DumpEMFields dumps the dynamically changing fields of a Ring in a user-
4// defined grid.
5//
6// Copyright (c) 2017, Chris Rogers
7// All rights reserved
8//
9// This file is part of OPAL.
10//
11// OPAL is free software: you can redistribute it and/or modify
12// it under the terms of the GNU General Public License as published by
13// the Free Software Foundation, either version 3 of the License, or
14// (at your option) any later version.
15//
16// You should have received a copy of the GNU General Public License
17// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
18//
20
25#include "Physics/Units.h"
27#include "Utilities/Util.h"
28
29#include <boost/filesystem.hpp>
30
31#include <cmath>
32#include <fstream>
33#include <map>
34
35extern Inform* gmsg;
36
37std::unordered_set<DumpEMFields*> DumpEMFields::dumpsSet_m;
38
40 Action(SIZE, "DUMPEMFIELDS",
41 "The \"DUMPEMFIELDS\" statement dumps a field map to a user-defined "
42 "field file, for checking that fields are generated correctly. "
43 "The fields are written out on a grid in space and time."),
44 grid_m(nullptr),
45 filename_m("") {
46
47 // would be nice if "steps" could be integer
49 ("FILE_NAME", "Name of the file to which field data is dumped");
50
52 ("COORDINATE_SYSTEM", "Choose to use CARTESIAN or CYLINDRICAL coordinates", {"CARTESIAN", "CYLINDRICAL"}, "CARTESIAN");
53
55 ("X_START", "(Cartesian) Start point in the grid in x [m]");
56
58 ("DX", "(Cartesian) Grid step size in x [m]");
59
61 ("X_STEPS", "(Cartesian) Number of steps in x");
62
64 ("Y_START", "(Cartesian) Start point in the grid in y [m]");
65
67 ("DY", "(Cartesian) Grid step size in y [m]");
68
70 ("Y_STEPS", "(Cartesian) Number of steps in y");
71
73 ("Z_START", "Start point in the grid in z [m]");
74
76 ("DZ", "Grid step size in z [m]");
77
79 ("Z_STEPS", "Number of steps in z");
80
82 ("T_START", "Start point in the grid in time [ns]");
83
85 ("DT", "Grid step size in time [ns]");
86
88 ("T_STEPS", "Number of steps in time");
89
91 ("R_START", "(Cylindrical) Start point in the grid in radius [m]");
92
94 ("DR", "(Cylindrical) Grid step size in radius [m]");
95
97 ("R_STEPS", "(Cylindrical) Number of steps in radius");
98
100 ("PHI_START", "(Cylindrical) Start point in the grid in phi [rad]");
101
103 ("DPHI", "(Cylindrical) Grid step size in phi [rad]");
104
106 ("PHI_STEPS", "(Cylindrical) Number of steps in phi");
107
109}
110
111DumpEMFields::DumpEMFields(const std::string& name, DumpEMFields* parent):
112 Action(name, parent), grid_m(nullptr)
113{}
114
116 delete grid_m;
117 dumpsSet_m.erase(this);
118}
119
120DumpEMFields* DumpEMFields::clone(const std::string& name) {
121 DumpEMFields* dumper = new DumpEMFields(name, this);
122 if (grid_m != nullptr) {
123 dumper->grid_m = grid_m->clone();
124 }
125 dumper->filename_m = filename_m;
127 if (dumpsSet_m.find(this) != dumpsSet_m.end()) {
128 dumpsSet_m.insert(dumper);
129 }
130 return dumper;
131}
132
134 static const std::map<std::string, CoordinateSystem> stringCoordinateSystem_s = {
135 {"CARTESIAN", CoordinateSystem::CARTESIAN},
136 {"CYLINDRICAL", CoordinateSystem::CYLINDRICAL}
137 };
138 coordinates_m = stringCoordinateSystem_s.at(Attributes::getString(itsAttr[COORDINATE_SYSTEM]));
139}
140
142 buildGrid();
143 // the routine for action (OpalParser/OpalParser) calls execute and then
144 // deletes 'this'; so we must build a copy that lasts until the field maps
145 // are constructed and we are ready for tracking (which is when the field
146 // maps are written). Hence the clone call below.
147 dumpsSet_m.insert(this->clone(""));
148}
149
151 std::vector<double> spacing(4);
152 std::vector<double> origin(4);
153 std::vector<int> gridSize(4);
155
156 switch (coordinates_m) {
158 origin[0] = Attributes::getReal(itsAttr[X_START]);
159 spacing[0] = Attributes::getReal(itsAttr[DX]);
160 double nx = Attributes::getReal(itsAttr[X_STEPS]);
161 checkInt(nx, "X_STEPS");
162 gridSize[0] = nx;
163
164 origin[1] = Attributes::getReal(itsAttr[Y_START]);
165 spacing[1] = Attributes::getReal(itsAttr[DY]);
166 double ny = Attributes::getReal(itsAttr[Y_STEPS]);
167 checkInt(ny, "Y_STEPS");
168 gridSize[1] = ny;
169
170 break;
171 }
173 origin[0] = Attributes::getReal(itsAttr[R_START]);
174 spacing[0] = Attributes::getReal(itsAttr[DR]);
176 checkInt(nr, "R_STEPS");
177 gridSize[0] = nr;
178
180 spacing[1] = Attributes::getReal(itsAttr[DPHI]);
181 double nphi = Attributes::getReal(itsAttr[PHI_STEPS]);
182 checkInt(nphi, "PHI_STEPS");
183 gridSize[1] = nphi;
184
185 break;
186 }
187 }
188
189 origin[2] = Attributes::getReal(itsAttr[Z_START]);
190 spacing[2] = Attributes::getReal(itsAttr[DZ]);
191 double nz = Attributes::getReal(itsAttr[Z_STEPS]);
192 checkInt(nz, "Z_STEPS");
193 gridSize[2] = nz;
194
195 origin[3] = Attributes::getReal(itsAttr[T_START]);
196 spacing[3] = Attributes::getReal(itsAttr[DT]);
197 double nt = Attributes::getReal(itsAttr[T_STEPS]);
198 checkInt(nt, "T_STEPS");
199 gridSize[3] = nt;
200
201 if (grid_m != nullptr) {
202 delete grid_m;
203 }
204
205 grid_m = new interpolation::NDGrid(4, &gridSize[0], &spacing[0], &origin[0]);
206
208}
209
212 for (dump_iter it = dumpsSet_m.begin(); it != dumpsSet_m.end(); ++it) {
213 (*it)->writeFieldThis(field);
214 }
215}
216
217void DumpEMFields::checkInt(double real, std::string name, double tolerance) {
218 real += tolerance; // prevent rounding error
219 if (std::abs(std::floor(real) - real) > 2*tolerance) {
220 throw OpalException("DumpEMFields::checkInt",
221 "Value for " + name +
222 " should be an integer but a real value was found");
223 }
224 if (std::floor(real) < 0.5) {
225 throw OpalException("DumpEMFields::checkInt",
226 "Value for " + name + " should be 1 or more");
227 }
228}
229
230void DumpEMFields::writeHeader(std::ofstream& fout) const {
231 fout << grid_m->end().toInteger() << "\n";
232 switch (coordinates_m) {
234 fout << 1 << " x [m]\n";
235 fout << 2 << " y [m]\n";
236 fout << 3 << " z [m]\n";
237 fout << 4 << " t [ns]\n";
238 fout << 5 << " Bx [kGauss]\n";
239 fout << 6 << " By [kGauss]\n";
240 fout << 7 << " Bz [kGauss]\n";
241 fout << 8 << " Ex [MV/m]\n";
242 fout << 9 << " Ey [MV/m]\n";
243 fout << 10 << " Ez [MV/m]\n";
244 break;
245 }
247 fout << 1 << " r [m]\n";
248 fout << 2 << " phi [deg]\n";
249 fout << 3 << " z [m]\n";
250 fout << 4 << " t [ns]\n";
251 fout << 5 << " Br [kGauss]\n";
252 fout << 6 << " Bphi [kGauss]\n";
253 fout << 7 << " Bz [kGauss]\n";
254 fout << 8 << " Er [MV/m]\n";
255 fout << 9 << " Ephi [MV/m]\n";
256 fout << 10 << " Ez [MV/m]\n";
257 break;
258 }
259 }
260 fout << 0 << std::endl;
261}
262
264 const Vector_t<double, 3>& pointIn,
265 const double& time,
266 std::ofstream& fout) const {
267 Vector_t<double, 3> centroid(0., 0., 0.);
268 Vector_t<double, 3> E(0., 0., 0.);
269 Vector_t<double, 3> B(0., 0., 0.);
270 Vector_t<double, 3> point = pointIn;
272 // pointIn is r, phi, z
273 point[0] = std::cos(pointIn[1])*pointIn[0];
274 point[1] = std::sin(pointIn[1])*pointIn[0];
275 }
276
277 field->apply(point, centroid, time, E, B);
278 Vector_t<double, 3> Bout = B;
279 Vector_t<double, 3> Eout = E;
281 // pointIn is r, phi, z
282 Bout[0] = B[0]*std::cos(pointIn[1])+B[1]*std::sin(pointIn[1]);
283 Bout[1] = -B[0]*std::sin(pointIn[1])+B[1]*std::cos(pointIn[1]);
284 Eout[0] = E[0]*std::cos(pointIn[1])+E[1]*std::sin(pointIn[1]);
285 Eout[1] = -E[0]*std::sin(pointIn[1])+E[1]*std::cos(pointIn[1]);
286 fout << pointIn[0] << " " << pointIn[1]*Units::rad2deg << " " << pointIn[2] << " " << time << " ";
287 } else {
288 fout << pointIn[0] << " " << pointIn[1] << " " << pointIn[2] << " " << time << " ";
289 }
290
291 fout << Bout[0] << " " << Bout[1] << " " << Bout[2] << " ";
292 fout << Eout[0] << " " << Eout[1] << " " << Eout[2] << "\n";
293}
294
296 if (grid_m == nullptr) {
297 throw OpalException("DumpEMFields::writeFieldThis",
298 "The grid was nullptr; there was a problem with the DumpEMFields initialisation.");
299 }
300 if (field == nullptr) {
301 throw OpalException("DumpEMFields::writeFieldThis",
302 "The field to be written was nullptr.");
303 }
304
305 *gmsg << *this << endl;
306
307 std::string fname;
308 if (boost::filesystem::path(filename_m).is_absolute() == true) {
309 fname = filename_m;
310 } else {
311 fname = Util::combineFilePath({
314 });
315 }
316
317 std::vector<double> point_std(4);
318 Vector_t<double, 3> point(0., 0., 0.);
319 std::ofstream fout;
320 try {
321 fout.open(fname.c_str(), std::ofstream::out);
322 } catch (std::exception& exc) {
323 throw OpalException("DumpEMFields::writeFieldThis",
324 "Failed to open DumpEMFields file " + filename_m);
325 }
326 if (!fout.good()) {
327 throw OpalException("DumpEMFields::writeFieldThis",
328 "Failed to open DumpEMFields file " + filename_m);
329 }
330 // set precision
331 writeHeader(fout);
333 it < grid_m->end();
334 ++it) {
335 it.getPosition(&point_std[0]);
336 for (size_t i = 0; i < 3; ++i) {
337 point[i] = point_std[i];
338 }
339 double time = point_std[3];
340 writeFieldLine(field, point, time, fout);
341 }
342 if (!fout.good()) {
343 throw OpalException("DumpEMFields::writeFieldThis",
344 "Something went wrong during writing " + filename_m);
345 }
346 fout.close();
347}
348
349void DumpEMFields::print(std::ostream& os) const {
350 os << "* ************* D U M P E M F I E L D S ****************************************** " << std::endl;
351 os << "* File name: '" << filename_m << "'\n";
353 os << "* Coordinate system: " << Attributes::getString(itsAttr[COORDINATE_SYSTEM]) << '\n'
354 << "* X_START = " << Attributes::getReal(itsAttr[X_START]) << " [m]\n"
355 << "* DX = " << Attributes::getReal(itsAttr[DX]) << " [m]\n"
356 << "* X_STEPS = " << Attributes::getReal(itsAttr[X_STEPS]) << '\n'
357 << "* Y_START = " << Attributes::getReal(itsAttr[Y_START]) << " [m]\n"
358 << "* DY = " << Attributes::getReal(itsAttr[DY]) << " [m]\n"
359 << "* Y_STEPS = " << Attributes::getReal(itsAttr[Y_STEPS]) << '\n';
361 os << "* Coordinate system: " << Attributes::getString(itsAttr[COORDINATE_SYSTEM]) << '\n'
362 << "* R_START = " << Attributes::getReal(itsAttr[R_START]) << " [m]\n"
363 << "* DR = " << Attributes::getReal(itsAttr[DR]) << " [m]\n"
364 << "* R_STEPS = " << Attributes::getReal(itsAttr[R_STEPS]) << '\n'
365 << "* PHI_START = " << Attributes::getReal(itsAttr[PHI_START]) << " [rad]\n"
366 << "* DPHI = " << Attributes::getReal(itsAttr[DPHI]) << " [rad]\n"
367 << "* PHI_STEPS = " << Attributes::getReal(itsAttr[PHI_STEPS]) << '\n';
368 }
369 os << "* Z_START = " << Attributes::getReal(itsAttr[Z_START]) << " [m]\n"
370 << "* DZ = " << Attributes::getReal(itsAttr[DZ]) << " [m]\n"
371 << "* Z_STEPS = " << Attributes::getReal(itsAttr[Z_STEPS]) << '\n'
372 << "* T_START = " << Attributes::getReal(itsAttr[T_START]) << " [ns]\n"
373 << "* DT = " << Attributes::getReal(itsAttr[DT]) << " [ns]\n"
374 << "* T_STEPS = " << Attributes::getReal(itsAttr[T_STEPS]) << '\n';
375 os << "* ********************************************************************************** " << std::endl;
376}
Inform * gmsg
Definition: changes.cpp:7
@ SIZE
Definition: IndexMap.cpp:173
PartBunch< T, Dim >::ConstIterator end(PartBunch< T, Dim > const &bunch)
const int nr
Definition: ClassicRandom.h:24
Inform & endl(Inform &inf)
Definition: Inform.cpp:42
double getReal(const Attribute &attr)
Return real value.
Definition: Attributes.cpp:252
Attribute makePredefinedString(const std::string &name, const std::string &help, const std::initializer_list< std::string > &predefinedStrings)
Make predefined string attribute.
Definition: Attributes.cpp:409
Attribute makeReal(const std::string &name, const std::string &help)
Make real attribute.
Definition: Attributes.cpp:240
std::string getString(const Attribute &attr)
Get string value.
Definition: Attributes.cpp:343
Attribute makeString(const std::string &name, const std::string &help)
Make string attribute.
Definition: Attributes.cpp:332
constexpr double rad2deg
Definition: Units.h:146
std::string::iterator iterator
Definition: MSLang.h:15
std::string combineFilePath(std::initializer_list< std::string > ilist)
Definition: Util.cpp:197
Interface for a single beam element.
Definition: Component.h:50
virtual bool apply(const size_t &i, const double &t, Vector_t< double, 3 > &E, Vector_t< double, 3 > &B)
Definition: Component.cpp:80
The base class for all OPAL actions.
Definition: Action.h:30
void registerOwnership(const AttributeHandler::OwnerType &itsClass) const
Definition: Object.cpp:189
std::vector< Attribute > itsAttr
The object attributes.
Definition: Object.h:216
static OpalData * getInstance()
Definition: OpalData.cpp:195
std::string getAuxiliaryOutputDirectory() const
get the name of the the additional data directory
Definition: OpalData.cpp:677
DumpEMFields dumps the dynamically changing fields of a Ring in a user- defined grid.
Definition: DumpEMFields.h:54
virtual void writeFieldThis(Component *field)
virtual void buildGrid()
void print(std::ostream &os) const
Print the attributes of DumpEMFields to standard out.
static std::unordered_set< DumpEMFields * > dumpsSet_m
Definition: DumpEMFields.h:148
void parseCoordinateSystem()
std::string filename_m
Definition: DumpEMFields.h:144
static void writeFields(Component *field)
Write the fields for all defined DumpEMFields objects.
interpolation::NDGrid * grid_m
Definition: DumpEMFields.h:143
CoordinateSystem coordinates_m
Definition: DumpEMFields.h:146
void writeHeader(std::ofstream &fout) const
virtual void execute()
Builds the grid but does not write the field map.
static void checkInt(double value, std::string name, double tolerance=1e-9)
virtual DumpEMFields * clone(const std::string &name)
Make a clone (overloadable copy-constructor).
DumpEMFields()
Constructor.
virtual ~DumpEMFields()
Destructor deletes grid_m and if in the dumps set, take it out.
void writeFieldLine(Component *field, const Vector_t< double, 3 > &point, const double &time, std::ofstream &fout) const
Used to loop over some, or all, points in the mesh, as in stl Enables e.g.
int toInteger() const
Return an integer representation of the iterator.
NDGrid holds grid information for a rectilinear grid in some arbitrary dimensional space.
Definition: NDGrid.h:57
Mesh::Iterator begin() const
Returns the origin of the mesh (lowest possible index in all dimensions)
Definition: NDGrid.h:358
Mesh::Iterator end() const
Returns the end of the mesh (highest possible index in all dimensions)
Definition: NDGrid.h:362
NDGrid * clone()
Inheritable copy constructor.
Definition: NDGrid.h:350
The base class for all OPAL exceptions.
Definition: OpalException.h:28
Definition: Inform.h:40