aboutsummaryrefslogtreecommitdiff
path: root/Source/Python/Particles/ParticleBoundaryBuffer.cpp
blob: 2a35faece9bef1f2d0840c2c70bb96bd3d9b9644 (plain) (blame)
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
/* Copyright 2021-2023 The WarpX Community
 *
 * Authors: Axel Huebl, Remi Lehe, Roelof Groenewald
 * License: BSD-3-Clause-LBNL
 */

#include "Python/pyWarpX.H"

#include <Particles/ParticleBoundaryBuffer.H>

namespace warpx {
    class BoundaryBufferParIter
        : public amrex::ParIter<0,0,PIdx::nattribs,0,amrex::PinnedArenaAllocator>
    {
    public:
        using amrex::ParIter<0,0,PIdx::nattribs,0,amrex::PinnedArenaAllocator>::ParIter;

        BoundaryBufferParIter(ContainerType& pc, int level) :
            amrex::ParIter<0,0,PIdx::nattribs,0,amrex::PinnedArenaAllocator>(pc, level) {}
    };
}

void init_BoundaryBufferParIter (py::module& m)
{
    py::class_<
        warpx::BoundaryBufferParIter,
        amrex::ParIter<0,0,PIdx::nattribs,0,amrex::PinnedArenaAllocator>
    >(m, "BoundaryBufferParIter")
        .def(py::init<amrex::ParIter<0,0,PIdx::nattribs,0,amrex::PinnedArenaAllocator>::ContainerType&, int>(),
            py::arg("particle_container"), py::arg("level")
        )
    ;
}

void init_ParticleBoundaryBuffer (py::module& m)
{
    py::class_<ParticleBoundaryBuffer>(m, "ParticleBoundaryBuffer")
        .def(py::init<>())
        .def("clear_particles",
            [](ParticleBoundaryBuffer& pbb) { pbb.clearParticles(); }
        )
        .def("get_particle_container",
            [](ParticleBoundaryBuffer& pbb,
            const std::string species_name, int boundary) {
                return &pbb.getParticleBuffer(species_name, boundary);
            },
            py::arg("species_name"), py::arg("boundary"),
            py::return_value_policy::reference_internal
        )
        .def("get_num_particles_in_container",
            [](ParticleBoundaryBuffer& pbb,
            const std::string species_name, int boundary, bool local)
            {
                return pbb.getNumParticlesInContainer(species_name, boundary, local);
            },
            py::arg("species_name"), py::arg("boundary"), py::arg("local")=true
        )
    ;
}