OPAL (Object Oriented Parallel Accelerator Library) 2022.1
OPAL
PyField.cpp
Go to the documentation of this file.
1#include <Python.h>
2#include <structmember.h>
3
4#include <boost/python.hpp>
5
7#include "AbsBeamline/Ring.h"
8#include "Track/TrackRun.h"
12#include "PyOpal/PyField.h"
13
14
15namespace PyOpal {
16namespace Field {
17py::object get_field_value_parallelt(double x,
18 double y,
19 double z,
20 double t,
21 ParallelTTracker* tracker) {
22 if (tracker == NULL) {
23 throw(OpalException("PyField::get_field_value_parallelt",
24 "ParallelTTracker was NULL"));
25 }
26 int outOfBounds = 0;
27 Vector_t R(x, y, z);
28 Vector_t P(0, 0, 0);
29 Vector_t E, B;
30 outOfBounds = tracker->getFieldValue(R, t, B, E);
31 boost::python::tuple value = boost::python::make_tuple(outOfBounds,
32 B[0], B[1], B[2],
33 E[0], E[1], E[2]);
34 return value;
35}
36
37py::object get_field_value_cyclotron(double x,
38 double y,
39 double z,
40 double t,
41 ParallelCyclotronTracker* tracker) {
42 if (tracker == NULL) {
43 throw(OpalException("PyField::get_field_value_cyclotron",
44 "ParallelCyclotronTracker was NULL"));
45 }
46 boost::python::tuple value = boost::python::make_tuple(1,
47 -1, -1, -1, -1, -1, -1);
48 return value;
49
50}
51
52py::object get_field_value_ring(double x,
53 double y,
54 double z,
55 double t,
56 Ring* ring) {
57 if (ring == NULL) {
58 throw(OpalException("PyField::get_field_value_ring",
59 "Ring was NULL"));
60 }
61 if (ring == NULL) {
62 std::string err = "Could not find a ring object - maybe a "
63 "RingDefinition was not defined or KeepAlive was False";
64 throw(OpalException("PyField::get_field_value", err));
65 }
66 Vector_t R(x, y, z);
67 Vector_t P(0, 0, 0);
68 Vector_t E, B;
69 int outOfBounds = ring->apply(R, P, t, E, B);
70 boost::python::tuple value = boost::python::make_tuple(outOfBounds,
71 B[0]/10., B[1]/10., B[2]/10.,
72 E[0], E[1], E[2]);
73 return value;
74}
75
76py::object get_field_value(double x, double y, double z, double t) {
77 Ring* ring = const_cast<Ring*>(Ring::getLastLockedRing());
78 if (ring != NULL) {
79 return get_field_value_ring(x, y, z, t, ring);
80 }
81 Tracker* tracker = TrackRun::getTracker();
82 ParallelTTracker* trackerT = dynamic_cast<ParallelTTracker*>(tracker);
83 if (trackerT != NULL) {
84 return get_field_value_parallelt(x, y, z, t, trackerT);
85 }
86 ParallelCyclotronTracker* trackerCycl = dynamic_cast<ParallelCyclotronTracker*>(tracker);
87 if (trackerCycl != NULL) {
88 return get_field_value_cyclotron(x, y, z, t, trackerCycl);
89 }
90 throw(OpalException("PyField::get_field_value",
91 "Could not find a Ring, ParallelTTracker or ParallelCyclotronTracker"));
92
93}
94
97 py::scope().attr("__doc__") = field_docstring.c_str();
98 py::def("get_field_value",
100 py::args("x", "y", "z", "t"),
102 );
103}
104
105}
106}
107
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
py::object get_field_value_ring(double x, double y, double z, double t, Ring *ring)
Definition: PyField.cpp:52
py::object get_field_value_parallelt(double x, double y, double z, double t, ParallelTTracker *tracker)
Definition: PyField.cpp:17
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
Ring describes a ring type geometry for tracking.
Definition: Ring.h:64
virtual bool apply(const size_t &id, const double &t, Vector_t &E, Vector_t &B) override
Definition: Ring.cpp:99
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