1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
|
/* Copyright 2021-2022 The WarpX Community
*
* Authors: Axel Huebl
* License: BSD-3-Clause-LBNL
*/
#include "pyWarpX.H"
#include <WarpX.H>
// see WarpX.cpp - full includes for _fwd.H headers
#include <BoundaryConditions/PML.H>
#include <Diagnostics/MultiDiagnostics.H>
#include <Diagnostics/ReducedDiags/MultiReducedDiags.H>
#include <EmbeddedBoundary/WarpXFaceInfoBox.H>
#include <FieldSolver/FiniteDifferenceSolver/FiniteDifferenceSolver.H>
#include <FieldSolver/FiniteDifferenceSolver/MacroscopicProperties/MacroscopicProperties.H>
#include <FieldSolver/FiniteDifferenceSolver/HybridPICModel/HybridPICModel.H>
#ifdef WARPX_USE_PSATD
# include <FieldSolver/SpectralSolver/SpectralKSpace.H>
# ifdef WARPX_DIM_RZ
# include <FieldSolver/SpectralSolver/SpectralSolverRZ.H>
# include <BoundaryConditions/PML_RZ.H>
# else
# include <FieldSolver/SpectralSolver/SpectralSolver.H>
# endif // RZ ifdef
#endif // use PSATD ifdef
#include <FieldSolver/WarpX_FDTD.H>
#include <Filter/NCIGodfreyFilter.H>
#include <Particles/MultiParticleContainer.H>
#include <Particles/ParticleBoundaryBuffer.H>
#include <AcceleratorLattice/AcceleratorLattice.H>
#include <Utils/TextMsg.H>
#include <Utils/WarpXAlgorithmSelection.H>
#include <Utils/WarpXConst.H>
#include <Utils/WarpXProfilerWrapper.H>
#include <Utils/WarpXUtil.H>
#include <AMReX.H>
#include <AMReX_ParmParse.H>
#include <AMReX_ParallelDescriptor.H>
#if defined(AMREX_DEBUG) || defined(DEBUG)
# include <cstdio>
#endif
#include <string>
//using namespace warpx;
namespace warpx {
struct Config {};
}
void init_WarpX (py::module& m)
{
// Expose the WarpX instance
m.def("get_instance",
[] () { return &WarpX::GetInstance(); },
"Return a reference to the WarpX object.");
m.def("finalize", &WarpX::Finalize,
"Close out the WarpX related data");
py::class_<WarpX> warpx(m, "WarpX");
warpx
// WarpX is a Singleton Class with a private constructor
// https://github.com/ECP-WarpX/WarpX/pull/4104
// https://pybind11.readthedocs.io/en/stable/advanced/classes.html?highlight=singleton#custom-constructors
.def(py::init([]() {
return &WarpX::GetInstance();
}))
.def_static("get_instance",
[] () { return &WarpX::GetInstance(); },
"Return a reference to the WarpX object."
)
.def_static("finalize", &WarpX::Finalize,
"Close out the WarpX related data"
)
.def("initialize_data", &WarpX::InitData,
"Initializes the WarpX simulation"
)
.def("evolve", &WarpX::Evolve,
"Evolve the simulation the specified number of steps"
)
// from AmrCore->AmrMesh
.def("Geom",
//[](WarpX const & wx, int const lev) { return wx.Geom(lev); },
py::overload_cast< int >(&WarpX::Geom, py::const_),
py::arg("lev")
)
.def("DistributionMap",
[](WarpX const & wx, int const lev) { return wx.DistributionMap(lev); },
//py::overload_cast< int >(&WarpX::DistributionMap, py::const_),
py::arg("lev")
)
.def("boxArray",
[](WarpX const & wx, int const lev) { return wx.boxArray(lev); },
//py::overload_cast< int >(&WarpX::boxArray, py::const_),
py::arg("lev")
)
.def("multifab",
[](WarpX const & wx, std::string const multifab_name) {
if (wx.multifab_map.count(multifab_name) > 0) {
return wx.multifab_map.at(multifab_name);
} else {
throw std::runtime_error("The MultiFab '" + multifab_name + "' is unknown or is not allocated!");
}
},
py::arg("multifab_name"),
py::return_value_policy::reference_internal,
"Return MultiFabs by name, e.g., 'Efield_aux[x][l=0]', 'Efield_cp[x][l=0]', ..."
)
.def("multi_particle_container",
[](WarpX& wx){ return &wx.GetPartContainer(); },
py::return_value_policy::reference_internal
)
.def("get_particle_boundary_buffer",
[](WarpX& wx){ return &wx.GetParticleBoundaryBuffer(); },
py::return_value_policy::reference_internal
)
// Expose functions used to sync the charge density multifab
// accross tiles and apply appropriate boundary conditions
.def("sync_rho",
[](WarpX& wx){ wx.SyncRho(); }
)
#ifdef WARPX_DIM_RZ
.def("apply_inverse_volume_scaling_to_charge_density",
[](WarpX& wx, amrex::MultiFab* rho, int const lev) {
wx.ApplyInverseVolumeScalingToChargeDensity(rho, lev);
},
py::arg("rho"), py::arg("lev")
)
#endif
// Expose functions to get the current simulation step and time
.def("getistep",
[](WarpX const & wx, int lev){ return wx.getistep(lev); },
py::arg("lev")
)
.def("gett_new",
[](WarpX const & wx, int lev){ return wx.gett_new(lev); },
py::arg("lev")
)
.def("set_potential_on_eb",
[](WarpX& wx, std::string potential) {
wx.m_poisson_boundary_handler.setPotentialEB(potential);
},
py::arg("potential")
)
;
py::class_<warpx::Config>(m, "Config")
// .def_property_readonly_static(
// "warpx_version",
// [](py::object) { return Version(); },
// "WarpX version")
.def_property_readonly_static(
"have_mpi",
[](py::object){
#ifdef AMREX_USE_MPI
return true;
#else
return false;
#endif
})
.def_property_readonly_static(
"have_gpu",
[](py::object){
#ifdef AMREX_USE_GPU
return true;
#else
return false;
#endif
})
.def_property_readonly_static(
"have_omp",
[](py::object){
#ifdef AMREX_USE_OMP
return true;
#else
return false;
#endif
})
.def_property_readonly_static(
"gpu_backend",
[](py::object){
#ifdef AMREX_USE_CUDA
return "CUDA";
#elif defined(AMREX_USE_HIP)
return "HIP";
#elif defined(AMREX_USE_DPCPP)
return "SYCL";
#else
return py::none();
#endif
})
;
}
|