aboutsummaryrefslogtreecommitdiff
path: root/Source/WarpXWrappers.cpp
diff options
context:
space:
mode:
Diffstat (limited to 'Source/WarpXWrappers.cpp')
-rw-r--r--Source/WarpXWrappers.cpp94
1 files changed, 93 insertions, 1 deletions
diff --git a/Source/WarpXWrappers.cpp b/Source/WarpXWrappers.cpp
index f8a839434..2df4cb549 100644
--- a/Source/WarpXWrappers.cpp
+++ b/Source/WarpXWrappers.cpp
@@ -5,6 +5,26 @@
#include <WarpXWrappers.h>
#include <WarpX.H>
+namespace
+{
+ double** getMultiFabPointers(const amrex::MultiFab& mf, int *num_boxes, int *ngrow, int **shapes)
+ {
+ *ngrow = mf.nGrow();
+ *num_boxes = mf.local_size();
+ *shapes = (int*) malloc(3*(*num_boxes)*sizeof(int));
+ double** data = (double**) malloc((*num_boxes)*sizeof(double*));
+
+ int i = 0;
+ for ( amrex::MFIter mfi(mf, false); mfi.isValid(); ++mfi, ++i ) {
+ data[i] = (double*) mf[mfi].dataPtr();
+ for (int j = 0; j < 3; ++j) {
+ (*shapes)[3*i+j] = mf[mfi].box().length(j);
+ }
+ }
+ return data;
+ }
+}
+
extern "C"
{
void amrex_init (int argc, char* argv[])
@@ -41,7 +61,10 @@ extern "C"
warpx.Evolve(numsteps);
}
- void addNParticles(int speciesnumber, int lenx, double* x, double* y, double* z, double* vx, double* vy, double* vz, int nattr, double* attr, int uniqueparticles)
+ void warpx_addNParticles(int speciesnumber, int lenx,
+ double* x, double* y, double* z,
+ double* vx, double* vy, double* vz,
+ int nattr, double* attr, int uniqueparticles)
{
auto & mypc = WarpX::GetInstance().GetPartContainer();
auto & myspc = mypc.GetParticleContainer(speciesnumber);
@@ -61,5 +84,74 @@ extern "C"
const amrex::Geometry& geom = warpx.Geom(0);
return geom.ProbHi(dir);
}
+
+ long warpx_getNumParticles(int speciesnumber) {
+ auto & mypc = WarpX::GetInstance().GetPartContainer();
+ auto & myspc = mypc.GetParticleContainer(speciesnumber);
+ return myspc.TotalNumberOfParticles();
+ }
+
+ double** warpx_getEfield(int lev, int direction,
+ int *return_size, int *ngrow, int **shapes) {
+ auto & mf = WarpX::GetInstance().getEfield(lev, direction);
+ return getMultiFabPointers(mf, return_size, ngrow, shapes);
+ }
+
+ double** warpx_getBfield(int lev, int direction,
+ int *return_size, int *ngrow, int **shapes) {
+
+ auto & mf = WarpX::GetInstance().getBfield(lev, direction);
+ return getMultiFabPointers(mf, return_size, ngrow, shapes);
+ }
+
+ double** warpx_getCurrentDensity(int lev, int direction,
+ int *return_size, int *ngrow, int **shapes) {
+
+ auto & mf = WarpX::GetInstance().getcurrent(lev, direction);
+ return getMultiFabPointers(mf, return_size, ngrow, shapes);
+ }
+
+ double** warpx_getParticleStructs(int speciesnumber,
+ int* num_tiles, int** particles_per_tile) {
+ auto & mypc = WarpX::GetInstance().GetPartContainer();
+ auto & myspc = mypc.GetParticleContainer(speciesnumber);
+
+ const int level = 0;
+
+ WarpXParIter pti(myspc, level);
+ *num_tiles = pti.numTiles();
+ *particles_per_tile = (int*) malloc(*num_tiles*sizeof(int));
+
+ double** data = (double**) malloc(*num_tiles*sizeof(typename WarpXParticleContainer::ParticleType*));
+ int i = 0;
+ for (WarpXParIter pti(myspc, level); pti.isValid(); ++pti, ++i) {
+ auto& aos = pti.GetArrayOfStructs();
+ data[i] = (double*) aos.data();
+ (*particles_per_tile)[i] = pti.numParticles();
+ }
+ return data;
+ }
+
+ double** warpx_getParticleArrays(int speciesnumber, int comp,
+ int* num_tiles, int** particles_per_tile) {
+ auto & mypc = WarpX::GetInstance().GetPartContainer();
+ auto & myspc = mypc.GetParticleContainer(speciesnumber);
+
+ const int level = 0;
+
+ WarpXParIter pti(myspc, level);
+ *num_tiles = pti.numTiles();
+ *particles_per_tile = (int*) malloc(*num_tiles*sizeof(int));
+
+ double** data = (double**) malloc(*num_tiles*sizeof(double*));
+ int i = 0;
+ for (WarpXParIter pti(myspc, level); pti.isValid(); ++pti, ++i) {
+ auto& soa = pti.GetStructOfArrays();
+ data[i] = (double*) soa[comp].dataPtr();
+ (*particles_per_tile)[i] = pti.numParticles();
+ }
+ return data;
+ }
+
}