OPAL (Object Oriented Parallel Accelerator Library) 2022.1
OPAL
PyBindField.cpp
Go to the documentation of this file.
1#include <Python.h>
2#include <structmember.h>
3
4#include <exception>
5#include <iostream>
6#include <pybind11/pybind11.h>
7
10#include "AbsBeamline/Ring.h"
12
13//using namespace boost::python;
14
16std::string("Get the field value at a point in the field map.\n\n")+
17std::string(" x: x position [m]\n")+
18std::string(" y: y position [m]\n")+
19std::string(" z: z position [m]\n")+
20std::string(" t: time [ns]\n")+
21std::string("Returns a tuple containing 6 values:\n")+
22std::string(" out of bounds: 1 if the event was out of the field map\n")+
23std::string(" boundary, else 0.\n")+
24std::string(" Bx: x magnetic field [T]\n")+
25std::string(" By: y magnetic field [T]\n")+
26std::string(" Bz: z magnetic field [T]\n")+
27std::string(" Ex: x electric field\n")+
28std::string(" Ey: y electric field\n")+
29std::string(" Ez: z electric field\n");
30
31std::tuple<int, double, double, double, double, double, double>
32 get_field_value(double x, double y, double z, double t) {
33 Ring* ring = const_cast<Ring*>(Ring::getLastLockedRing());
34 if (ring == NULL) {
35 std::string err = "Could not find a ring object - maybe a "+
36 std::string("RingDefinition was not defined or KeepAlive was False");
37 // need to define proper exception handling
38 throw(err);
39 }
40 Vector_t R(x, y, z);
41 Vector_t P(0, 0, 0);
42 Vector_t E, B;
43 int outOfBounds = ring->apply(R, P, t, E, B);
44 if(outOfBounds);
45 std::tuple<int, double, double, double, double, double, double> value = {
46 outOfBounds, B[0]/10., B[1]/10., B[2]/10., E[0], E[1], E[2]
47 };
48 return value;
49}
50
51const char* module_docstring = "field module returns the field";
52
53PYBIND11_MODULE(bind_field, module) {
54 module.doc() = module_docstring;
55 module.def("get_field_value", &get_field_value,
57 pybind11::arg("x"),
58 pybind11::arg("y"),
59 pybind11::arg("z"),
60 pybind11::arg("t"));
61}
std::tuple< int, double, double, double, double, double, double > get_field_value(double x, double y, double z, double t)
Definition: PyBindField.cpp:32
const char * module_docstring
Definition: PyBindField.cpp:51
PYBIND11_MODULE(bind_field, module)
Definition: PyBindField.cpp:53
std::string get_field_value_docstring
Definition: PyBindField.cpp:15
arg(a))
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