OPAL (Object Oriented Parallel Accelerator Library) 2022.1
OPAL
PyField.cpp
Go to the documentation of this file.
1#include <boost/python.hpp>
2
4#include "AbsBeamline/Ring.h"
5#include "Track/TrackRun.h"
8
11
12namespace PyOpal {
13namespace Field {
14
15std::string field_docstring =
16 "field module enables user to get the field at a point";
17
18std::string get_field_value_docstring =
19 "Get the field value at a point in the field map.\n"
20 "\n"
21 "The field lookup is performed against the last RINGDEFINITION that was\n"
22 "instantiated. This should be instantiated by calling\n"
23 "pyopal.parser.initialise_from_opal_file\n"
24 "\n"
25 "Parameters\n"
26 "----------\n"
27 "x : float\n"
28 " x position [m]\n"
29 "y : float\n"
30 " y position [m]\n"
31 "z : float\n"
32 " z position [m]\n"
33 "t: float\n"
34 " time [ns]\n"
35 "\n"
36 "Returns\n"
37 "-------\n"
38 "The function returns a tuple containing 7 values:\n"
39 "out of bounds : int\n"
40 " 1 if the event was out of the field map boundary, else 0.\n"
41 "Bx : float\n"
42 " x magnetic field [T]\n"
43 "By : float\n"
44 " y magnetic field [T]\n"
45 "Bz : float\n"
46 " z magnetic field [T]\n"
47 "Ex : float\n"
48 " x electric field\n"
49 "Ey : float\n"
50 " y electric field\n"
51 "Ez : float\n"
52 " z electric field\n";
53
54py::object get_field_value_cyclotron(double x,
55 double y,
56 double z,
57 double t,
58 ParallelCyclotronTracker* tracker) {
59 if (tracker == NULL) {
60 throw(OpalException("PyField::get_field_value_cyclotron",
61 "ParallelCyclotronTracker was NULL"));
62 }
63 Vector_t R(x, y, z);
64 Vector_t P, B, E;
65 int outOfBounds = tracker->computeExternalFields_m(R, P, t, E, B);
66 boost::python::tuple value = boost::python::make_tuple(outOfBounds,
67 B[0], B[1], B[2], E[0], E[1], E[2]);
68 return value;
69
70}
71
72py::object get_field_value(double x, double y, double z, double t) {
73 /*
74 Ring* ring = const_cast<Ring*>(Ring::getLastLockedRing());
75 if (ring != NULL) {
76 return get_field_value_ring(x, y, z, t, ring);
77 }*/
78 std::shared_ptr<Tracker> tracker = TrackRun::getTracker();
79 ParallelCyclotronTracker* trackerCycl =
80 dynamic_cast<ParallelCyclotronTracker*>(tracker.get());
81 if (trackerCycl != nullptr) {
82 return get_field_value_cyclotron(x, y, z, t, trackerCycl);
83 }
84 throw(OpalException("PyField::get_field_value",
85 "Could not find a ParallelCyclotronTracker"));
86}
87
90 py::scope().attr("__doc__") = field_docstring.c_str();
91 py::def("get_field_value",
93 py::args("x", "y", "z", "t"),
95 );
96}
97
98}
99}
100
py::object get_field_value_cyclotron(double x, double y, double z, double t, ParallelCyclotronTracker *tracker)
Definition: PyField.cpp:37
BOOST_PYTHON_MODULE(field)
Definition: PyField.cpp:95
std::string field_docstring
Definition: PyField.h:7
py::object get_field_value(double x, double y, double z, double t)
Definition: PyField.cpp:76
std::string get_field_value_docstring
Definition: PyField.h:10
bool computeExternalFields_m(const Vector_t &R, const Vector_t &P, const double &t, Vector_t &Efield, Vector_t &Bfield)
Calculate the field map at an arbitrary point.
static std::shared_ptr< Tracker > getTracker()
Definition: TrackRun.cpp:737
The base class for all OPAL exceptions.
Definition: OpalException.h:28
Definition: Field.h:33