OPAL (Object Oriented Parallel Accelerator Library) 2024.2
OPAL
PyField.cpp
Go to the documentation of this file.
1//
2// Python API for Global (external) field access
3//
4// Copyright (c) 2023, Chris Rogers, STFC Rutherford Appleton Laboratory, Didcot, UK
5//
6// This file is part of OPAL.
7//
8// OPAL is free software: you can redistribute it and/or modify
9// it under the terms of the GNU General Public License as published by
10// the Free Software Foundation, either version 3 of the License, or
11// (at your option) any later version.
12//
13// You should have received a copy of the GNU General Public License
14// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
15//
16#include <boost/python.hpp>
17
18#include "AbsBeamline/Ring.h"
21#include "Track/TrackRun.h"
23
26
27namespace PyOpal {
28namespace Field {
29
30std::string field_docstring =
31 "field module enables user to get the field at a point";
32
34 "Get the field value at a point in the field map.\n"
35 "Only available in Cyclotron mode.\n"
36 "\n"
37 "The field lookup is performed against the last RINGDEFINITION that was\n"
38 "instantiated. This should be instantiated by calling\n"
39 "pyopal.parser.initialise_from_opal_file\n"
40 "\n"
41 "Parameters\n"
42 "----------\n"
43 "x : float\n"
44 " x position [m]\n"
45 "y : float\n"
46 " y position [m]\n"
47 "z : float\n"
48 " z position [m]\n"
49 "t: float\n"
50 " time [ns]\n"
51 "\n"
52 "Returns\n"
53 "-------\n"
54 "The function returns a tuple containing 7 values:\n"
55 "out of bounds : int\n"
56 " 1 if the event was out of the field map boundary, else 0.\n"
57 "Bx : float\n"
58 " x magnetic field [T]\n"
59 "By : float\n"
60 " y magnetic field [T]\n"
61 "Bz : float\n"
62 " z magnetic field [T]\n"
63 "Ex : float\n"
64 " x electric field\n"
65 "Ey : float\n"
66 " y electric field\n"
67 "Ez : float\n"
68 " z electric field\n";
69
70py::object get_field_value_cyclotron(double x,
71 double y,
72 double z,
73 double t,
74 ParallelCyclotronTracker* tracker) {
75 if (tracker == nullptr) {
76 throw(OpalException("PyField::get_field_value_cyclotron",
77 "ParallelCyclotronTracker was nullptr"));
78 }
79 Vector_t R(x, y, z);
80 Vector_t P, B, E;
81 int outOfBounds = tracker->computeExternalFields_m(R, P, t, E, B);
82 boost::python::tuple value = boost::python::make_tuple(outOfBounds,
83 B[0], B[1], B[2], E[0], E[1], E[2]);
84 return value;
85
86}
87
88py::object get_field_value(double x, double y, double z, double t) {
89 std::shared_ptr<Tracker> tracker = TrackRun::getTracker();
90 ParallelCyclotronTracker* trackerCycl =
91 dynamic_cast<ParallelCyclotronTracker*>(tracker.get());
92 if (trackerCycl != nullptr) {
93 return get_field_value_cyclotron(x, y, z, t, trackerCycl);
94 }
95 throw(OpalException("PyField::get_field_value",
96 "Could not find a ParallelCyclotronTracker - get_field_value only works in OPAL-CYCL mode"));
97}
98
99
100// returns a *borrowed* pointer
102 std::shared_ptr<Tracker> tracker = TrackRun::getTracker();
103 ParallelCyclotronTracker* trackerCycl =
104 dynamic_cast<ParallelCyclotronTracker*>(tracker.get());
105 Ring* ring = trackerCycl->getRing();
106 if (ring == nullptr) {
107 throw GeneralClassicException("PyRingDefinition::getSection",
108 "Internal PyOpal error - failed to cast to a Ring object");
109 }
110 return ring;
111}
112
113
115"Return a string holding the name of the i^th element [m].\n\n";
116std::string getElementName(int i) {
118 Component* component = sec->getComponent();
119 if (component == nullptr) {
120 throw GeneralClassicException("PyRingDefinition::getElementName",
121 "Internal PyOpal error - failed to cast to a Component");
122 }
123 return component->getName();
124}
125
127"Return a tuple holding the start position of the element (x, y, z) [m].\n\n";
128
129boost::python::object getElementStartPosition(int i) {
131 Vector_t pos = sec->getStartPosition();
132 return boost::python::make_tuple(pos[0], pos[1], pos[2]);
133}
134
135std::string end_pos_docstring =
136"Return a tuple holding the end position of the element (x, y, z) [m].\n\n";
137boost::python::object getElementEndPosition(int i) {
139 Vector_t pos = sec->getEndPosition();
140 return boost::python::make_tuple(pos[0], pos[1], pos[2]);
141}
142
144"Return a tuple holding the vector (x, y, z) normal to the face of the\n"
145"element start, pointing towards the element and having length 1.\n\n";
146boost::python::object getElementStartNormal(int i) {
148 Vector_t dir = sec->getStartNormal();
149 return boost::python::make_tuple(dir[0], dir[1], dir[2]);
150}
151
152std::string end_norm_docstring =
153"Return a tuple holding the vector (x, y, z) normal to the face of the\n"
154"element end, pointing towards the next element and having length 1.\n\n";
155boost::python::object getElementEndNormal(int i) {
157 Vector_t dir = sec->getEndNormal();
158 return boost::python::make_tuple(dir[0], dir[1], dir[2]);
159}
160
162"Return an integer corresponding to the number of elements stored in the Ring\n"
163"If this is 0, check that the track has been executed - the element\n"
164"placements are done during Track setup.\n\n";
167}
168
169
173 py::scope().attr("__doc__") = field_docstring.c_str();
174 py::def("get_field_value",
176 py::args("x", "y", "z", "t"),
178 );
179 py::def("get_number_of_elements",
181 num_elements_docstring.c_str());
182 py::def("get_element_start_position",
184 element_name_docstring.c_str());
185 py::def("get_element_start_normal", &getElementStartNormal);
186 py::def("get_element_end_position", &getElementEndPosition);
187 py::def("get_element_end_normal", &getElementEndNormal);
188 py::def("get_element_name", &getElementName);
189}
190
191}
192}
Tps< T > sec(const Tps< T > &x)
Secant.
Definition TpsMath.h:153
void Initialise()
Definition Globals.cpp:50
boost::python::object getElementStartNormal(int i)
Definition PyField.cpp:146
py::object get_field_value_cyclotron(double x, double y, double z, double t, ParallelCyclotronTracker *tracker)
Definition PyField.cpp:70
BOOST_PYTHON_MODULE(field)
Definition PyField.cpp:170
boost::python::object getElementEndPosition(int i)
Definition PyField.cpp:137
size_t getNumberOfElements()
Definition PyField.cpp:165
std::string end_pos_docstring
Definition PyField.cpp:135
boost::python::object getElementStartPosition(int i)
Definition PyField.cpp:129
std::string num_elements_docstring
Definition PyField.cpp:161
std::string start_pos_docstring
Definition PyField.cpp:126
std::string start_norm_docstring
Definition PyField.cpp:143
std::string field_docstring
Definition PyField.cpp:30
std::string element_name_docstring
Definition PyField.cpp:114
py::object get_field_value(double x, double y, double z, double t)
Definition PyField.cpp:88
Ring * getRing()
Definition PyField.cpp:101
std::string getElementName(int i)
Definition PyField.cpp:116
boost::python::object getElementEndNormal(int i)
Definition PyField.cpp:155
std::string get_field_value_docstring
Definition PyField.cpp:33
std::string end_norm_docstring
Definition PyField.cpp:152
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.
Interface for a single beam element.
Definition Component.h:50
virtual const std::string & getName() const
Get element name.
Ring describes a ring type geometry for tracking.
Definition Ring.h:53
size_t getNumberOfRingSections() const
Definition Ring.cpp:402
RingSection * getSection(int i) const
Definition Ring.cpp:387
Component placement handler in ring geometry.
Definition RingSection.h:58
static std::shared_ptr< Tracker > getTracker()
Definition TrackRun.cpp:733
The base class for all OPAL exceptions.
Vektor< double, 3 > Vector_t