Skip to content

Commit 88f6194

Browse files
Pywrapper for ba convert, surroundingnodes, and mf.copy (#439)
1 parent dca9aec commit 88f6194

File tree

3 files changed

+37
-14
lines changed

3 files changed

+37
-14
lines changed

src/Base/BoxArray.cpp

Lines changed: 11 additions & 14 deletions
Original file line numberDiff line numberDiff line change
@@ -120,30 +120,27 @@ void init_BoxArray(py::module &m) {
120120
//! \brief Grow each Box in the BoxArray on the high end
121121
//! by n_cell cells in the idir direction.
122122
BoxArray& growHi (int idir, int n_cell);
123-
//! \brief Apply surroundingNodes(Box) to each Box in BoxArray.
124-
//! See the documentation of Box for details.
125-
BoxArray& surroundingNodes ();
126-
//!
123+
*/
124+
127125
//! \brief Apply surroundingNodes(Box,int) to each Box in
128126
//! BoxArray. See the documentation of Box for details.
129-
BoxArray& surroundingNodes (int dir);
130-
*/
127+
.def("surroundingNodes",
128+
py::overload_cast<>(&BoxArray::surroundingNodes))
129+
.def("surroundingNodes",
130+
py::overload_cast<int>(&BoxArray::surroundingNodes))
131131

132132
//! Apply Box::enclosedCells() to each Box in the BoxArray.
133133
.def("enclosed_cells",
134134
py::overload_cast<>(&BoxArray::enclosedCells))
135135
.def("enclosed_cells",
136136
py::overload_cast<int>(&BoxArray::enclosedCells))
137137

138+
//! Convert nodality of each box in the BoxArray
139+
.def("convert",
140+
py::overload_cast< IndexType >(&BoxArray::convert))
141+
.def("convert",
142+
py::overload_cast< IntVect const &>(&BoxArray::convert))
138143
/*
139-
//! Apply Box::convert(IndexType) to each Box in the BoxArray.
140-
BoxArray& convert (IndexType typ);
141-
142-
BoxArray& convert (const IntVect& typ);
143-
144-
//! Apply function (*fp)(Box) to each Box in the BoxArray.
145-
BoxArray& convert (Box (*fp)(const Box&));
146-
147144
//! Apply Box::shift(int,int) to each Box in the BoxArray.
148145
BoxArray& shift (int dir, int nzones);
149146

src/Base/MultiFab.H

Lines changed: 25 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -493,6 +493,31 @@ factory :
493493
;
494494
}
495495

496+
py_MultiFab
497+
.def("copymf",
498+
[](T &self, T const &src, int srccomp, int dstcomp, int numcomp, int nghost) {
499+
T::Copy(self, src, srccomp, dstcomp, numcomp, nghost);
500+
},
501+
py::arg("src"), py::arg("srccomp"), py::arg("dstcomp"), py::arg("numcomp"), py::arg("nghost"),
502+
"Copy from src to self including nghost ghost cells.\n"
503+
"The two MultiFabs MUST have the same underlying BoxArray. The copy is local"
504+
)
505+
;
506+
507+
// TODO: Missing in iMultiFab https://github.com/AMReX-Codes/amrex/issues/4317
508+
if constexpr (std::is_same_v<T, MultiFab>) {
509+
py_MultiFab
510+
.def("copymf",
511+
[](T &self, T const &src, int srccomp, int dstcomp, int numcomp, IntVect const &nghost) {
512+
T::Copy(self, src, srccomp, dstcomp, numcomp, nghost);
513+
},
514+
py::arg("src"), py::arg("srccomp"), py::arg("dstcomp"), py::arg("numcomp"), py::arg("nghost"),
515+
"Copy from src to self including nghost ghost cells.\n"
516+
"The two MultiFabs MUST have the same underlying BoxArray. The copy is local"
517+
)
518+
;
519+
}
520+
496521
py_MultiFab
497522
.def("subtract",
498523
[](T & self, T const & src, int srccomp, int comp, int numcomp, int nghost) {

src/Base/MultiFab.cpp

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -105,4 +105,5 @@ void init_MultiFab(py::module &m, py::class_< amrex::MFIter > & py_MFIter)
105105

106106
m.def("copy_mfab", py::overload_cast< MultiFab &, MultiFab const &, int, int, int, int >(&MultiFab::Copy), py::arg("dst"), py::arg("src"), py::arg("srccomp"), py::arg("dstcomp"), py::arg("numcomp"), py::arg("nghost"))
107107
.def("copy_mfab", py::overload_cast< MultiFab &, MultiFab const &, int, int, int, IntVect const & >(&MultiFab::Copy), py::arg("dst"), py::arg("src"), py::arg("srccomp"), py::arg("dstcomp"), py::arg("numcomp"), py::arg("nghost"));
108+
108109
}

0 commit comments

Comments
 (0)