diff options
-rwxr-xr-x | Python/pywarpx/_libwarpx.py | 91 | ||||
-rw-r--r-- | Python/pywarpx/fields.py | 4 | ||||
-rw-r--r-- | Source/Python/WarpXWrappers.cpp | 73 | ||||
-rw-r--r-- | Source/Python/WarpXWrappers.h | 50 |
4 files changed, 127 insertions, 91 deletions
diff --git a/Python/pywarpx/_libwarpx.py b/Python/pywarpx/_libwarpx.py index 2b5688933..6c3e54619 100755 --- a/Python/pywarpx/_libwarpx.py +++ b/Python/pywarpx/_libwarpx.py @@ -57,14 +57,35 @@ try: except OSError: raise Exception('libwarpx%s.so was not installed. It can be installed by running "make" in the Python directory of WarpX'%geometry_dim) +# WarpX can be compiled using either double or float +libwarpx.warpx_Real_size.restype = ctypes.c_int +libwarpx.warpx_ParticleReal_size.restype = ctypes.c_int + +_Real_size = libwarpx.warpx_Real_size() +_ParticleReal_size = libwarpx.warpx_ParticleReal_size() + +if _Real_size == 8: + c_real = ctypes.c_double + _numpy_real_dtype = 'f8' +else: + c_real = ctypes.c_float + _numpy_real_dtype = 'f4' + +if _ParticleReal_size == 8: + c_particlereal = ctypes.c_double + _numpy_particlereal_dtype = 'f8' +else: + c_particlereal = ctypes.c_float + _numpy_particlereal_dtype = 'f4' + dim = libwarpx.warpx_SpaceDim() -# our particle data type -_p_struct = [(d, 'f8') for d in 'xyz'[:dim]] + [('id', 'i4'), ('cpu', 'i4')] +# our particle data type, depends on _ParticleReal_size +_p_struct = [(d, _numpy_particlereal_dtype) for d in 'xyz'[:dim]] + [('id', 'i4'), ('cpu', 'i4')] _p_dtype = np.dtype(_p_struct, align=True) _numpy_to_ctypes = {} -_numpy_to_ctypes['f8'] = ctypes.c_double +_numpy_to_ctypes[_numpy_particlereal_dtype] = c_particlereal _numpy_to_ctypes['i4'] = ctypes.c_int class Particle(ctypes.Structure): @@ -75,8 +96,10 @@ class Particle(ctypes.Structure): _LP_particle_p = ctypes.POINTER(ctypes.POINTER(Particle)) _LP_c_int = ctypes.POINTER(ctypes.c_int) _LP_c_void_p = ctypes.POINTER(ctypes.c_void_p) -_LP_c_double = ctypes.POINTER(ctypes.c_double) -_LP_LP_c_double = ctypes.POINTER(_LP_c_double) +_LP_c_real = ctypes.POINTER(c_real) +_LP_LP_c_real = ctypes.POINTER(_LP_c_real) +_LP_c_particlereal = ctypes.POINTER(c_particlereal) +_LP_LP_c_particlereal = ctypes.POINTER(_LP_c_particlereal) _LP_c_char = ctypes.POINTER(ctypes.c_char) _LP_LP_c_char = ctypes.POINTER(_LP_c_char) @@ -100,63 +123,63 @@ def _array1d_from_pointer(pointer, dtype, size): # set the arg and return types of the wrapped functions libwarpx.amrex_init.argtypes = (ctypes.c_int, _LP_LP_c_char) libwarpx.warpx_getParticleStructs.restype = _LP_particle_p -libwarpx.warpx_getParticleArrays.restype = _LP_LP_c_double -libwarpx.warpx_getEfield.restype = _LP_LP_c_double +libwarpx.warpx_getParticleArrays.restype = _LP_LP_c_particlereal +libwarpx.warpx_getEfield.restype = _LP_LP_c_real libwarpx.warpx_getEfieldLoVects.restype = _LP_c_int -libwarpx.warpx_getEfieldCP.restype = _LP_LP_c_double +libwarpx.warpx_getEfieldCP.restype = _LP_LP_c_real libwarpx.warpx_getEfieldCPLoVects.restype = _LP_c_int -libwarpx.warpx_getEfieldFP.restype = _LP_LP_c_double +libwarpx.warpx_getEfieldFP.restype = _LP_LP_c_real libwarpx.warpx_getEfieldFPLoVects.restype = _LP_c_int -libwarpx.warpx_getBfield.restype = _LP_LP_c_double +libwarpx.warpx_getBfield.restype = _LP_LP_c_real libwarpx.warpx_getBfieldLoVects.restype = _LP_c_int -libwarpx.warpx_getBfieldCP.restype = _LP_LP_c_double +libwarpx.warpx_getBfieldCP.restype = _LP_LP_c_real libwarpx.warpx_getBfieldCPLoVects.restype = _LP_c_int -libwarpx.warpx_getBfieldFP.restype = _LP_LP_c_double +libwarpx.warpx_getBfieldFP.restype = _LP_LP_c_real libwarpx.warpx_getBfieldFPLoVects.restype = _LP_c_int -libwarpx.warpx_getCurrentDensity.restype = _LP_LP_c_double +libwarpx.warpx_getCurrentDensity.restype = _LP_LP_c_real libwarpx.warpx_getCurrentDensityLoVects.restype = _LP_c_int -libwarpx.warpx_getCurrentDensityCP.restype = _LP_LP_c_double +libwarpx.warpx_getCurrentDensityCP.restype = _LP_LP_c_real libwarpx.warpx_getCurrentDensityCPLoVects.restype = _LP_c_int -libwarpx.warpx_getCurrentDensityFP.restype = _LP_LP_c_double +libwarpx.warpx_getCurrentDensityFP.restype = _LP_LP_c_real libwarpx.warpx_getCurrentDensityFPLoVects.restype = _LP_c_int -#libwarpx.warpx_getPMLSigma.restype = _LP_c_double -#libwarpx.warpx_getPMLSigmaStar.restype = _LP_c_double -#libwarpx.warpx_ComputePMLFactors.argtypes = (ctypes.c_int, ctypes.c_double) +#libwarpx.warpx_getPMLSigma.restype = _LP_c_real +#libwarpx.warpx_getPMLSigmaStar.restype = _LP_c_real +#libwarpx.warpx_ComputePMLFactors.argtypes = (ctypes.c_int, c_real) libwarpx.warpx_addNParticles.argtypes = (ctypes.c_int, ctypes.c_int, - _ndpointer(ctypes.c_double, flags="C_CONTIGUOUS"), - _ndpointer(ctypes.c_double, flags="C_CONTIGUOUS"), - _ndpointer(ctypes.c_double, flags="C_CONTIGUOUS"), - _ndpointer(ctypes.c_double, flags="C_CONTIGUOUS"), - _ndpointer(ctypes.c_double, flags="C_CONTIGUOUS"), - _ndpointer(ctypes.c_double, flags="C_CONTIGUOUS"), + _ndpointer(c_particlereal, flags="C_CONTIGUOUS"), + _ndpointer(c_particlereal, flags="C_CONTIGUOUS"), + _ndpointer(c_particlereal, flags="C_CONTIGUOUS"), + _ndpointer(c_particlereal, flags="C_CONTIGUOUS"), + _ndpointer(c_particlereal, flags="C_CONTIGUOUS"), + _ndpointer(c_particlereal, flags="C_CONTIGUOUS"), ctypes.c_int, - _ndpointer(ctypes.c_double, flags="C_CONTIGUOUS"), + _ndpointer(c_particlereal, flags="C_CONTIGUOUS"), ctypes.c_int) -libwarpx.warpx_getProbLo.restype = ctypes.c_double -libwarpx.warpx_getProbHi.restype = ctypes.c_double +libwarpx.warpx_getProbLo.restype = c_real +libwarpx.warpx_getProbHi.restype = c_real libwarpx.warpx_getistep.restype = ctypes.c_int -libwarpx.warpx_gett_new.restype = ctypes.c_double -libwarpx.warpx_getdt.restype = ctypes.c_double +libwarpx.warpx_gett_new.restype = c_real +libwarpx.warpx_getdt.restype = c_real libwarpx.warpx_maxStep.restype = ctypes.c_int -libwarpx.warpx_stopTime.restype = ctypes.c_double +libwarpx.warpx_stopTime.restype = c_real libwarpx.warpx_checkInt.restype = ctypes.c_int libwarpx.warpx_plotInt.restype = ctypes.c_int libwarpx.warpx_finestLevel.restype = ctypes.c_int -libwarpx.warpx_EvolveE.argtypes = [ctypes.c_double] -libwarpx.warpx_EvolveB.argtypes = [ctypes.c_double] +libwarpx.warpx_EvolveE.argtypes = [c_real] +libwarpx.warpx_EvolveB.argtypes = [c_real] libwarpx.warpx_FillBoundaryE.argtypes = [] libwarpx.warpx_FillBoundaryB.argtypes = [] libwarpx.warpx_UpdateAuxilaryData.argtypes = [] libwarpx.warpx_SyncCurrent.argtypes = [] -libwarpx.warpx_PushParticlesandDepose.argtypes = [ctypes.c_double] +libwarpx.warpx_PushParticlesandDepose.argtypes = [c_real] libwarpx.warpx_getistep.argtypes = [ctypes.c_int] libwarpx.warpx_setistep.argtypes = [ctypes.c_int, ctypes.c_int] libwarpx.warpx_gett_new.argtypes = [ctypes.c_int] -libwarpx.warpx_sett_new.argtypes = [ctypes.c_int, ctypes.c_double] +libwarpx.warpx_sett_new.argtypes = [ctypes.c_int, c_real] libwarpx.warpx_getdt.argtypes = [ctypes.c_int] def get_nattr(): diff --git a/Python/pywarpx/fields.py b/Python/pywarpx/fields.py index 8b7283a46..fa247db07 100644 --- a/Python/pywarpx/fields.py +++ b/Python/pywarpx/fields.py @@ -127,7 +127,7 @@ class _MultiFABWrapper(object): sss = (max(0, ixstop - ixstart), max(0, iystop - iystart), max(0, izstop - izstart)) - resultglobal = np.zeros(sss) + resultglobal = np.zeros(sss, dtype=_libwarpx._numpy_real_dtype) datalist = [] for i in range(len(fields)): @@ -222,7 +222,7 @@ class _MultiFABWrapper(object): max(0, izstop - izstart)) if ncomps > 1 and ic is None: sss = tuple(list(sss) + [ncomps]) - resultglobal = np.zeros(sss) + resultglobal = np.zeros(sss, dtype=_libwarpx._numpy_real_dtype) datalist = [] for i in range(len(fields)): diff --git a/Source/Python/WarpXWrappers.cpp b/Source/Python/WarpXWrappers.cpp index c6d7dfdb8..be9ec9519 100644 --- a/Source/Python/WarpXWrappers.cpp +++ b/Source/Python/WarpXWrappers.cpp @@ -1,7 +1,4 @@ -#include <AMReX.H> -#include <AMReX_BLProfiler.H> - #include <WarpXWrappers.h> #include <WarpXParticleContainer.H> #include <WarpX.H> @@ -10,7 +7,7 @@ namespace { - double** getMultiFabPointers(const amrex::MultiFab& mf, int *num_boxes, int *ncomps, int *ngrow, int **shapes) + amrex::Real** getMultiFabPointers(const amrex::MultiFab& mf, int *num_boxes, int *ncomps, int *ngrow, int **shapes) { *ncomps = mf.nComp(); *ngrow = mf.nGrow(); @@ -18,14 +15,14 @@ namespace int shapesize = AMREX_SPACEDIM; if (mf.nComp() > 1) shapesize += 1; *shapes = (int*) malloc(shapesize * (*num_boxes) * sizeof(int)); - double** data = (double**) malloc((*num_boxes) * sizeof(double*)); + amrex::Real** data = (amrex::Real**) malloc((*num_boxes) * sizeof(amrex::Real*)); #ifdef _OPENMP #pragma omp parallel #endif for ( amrex::MFIter mfi(mf, false); mfi.isValid(); ++mfi ) { int i = mfi.LocalIndex(); - data[i] = (double*) mf[mfi].dataPtr(); + data[i] = (amrex::Real*) mf[mfi].dataPtr(); for (int j = 0; j < AMREX_SPACEDIM; ++j) { (*shapes)[shapesize*i+j] = mf[mfi].box().length(j); } @@ -53,6 +50,16 @@ namespace extern "C" { + int warpx_Real_size() + { + return (int)sizeof(amrex::Real); + } + + int warpx_ParticleReal_size() + { + return (int)sizeof(amrex::ParticleReal); + } + int warpx_nSpecies() { auto & mypc = WarpX::GetInstance().GetPartContainer(); @@ -165,9 +172,9 @@ extern "C" } 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) + amrex::ParticleReal* x, amrex::ParticleReal* y, amrex::ParticleReal* z, + amrex::ParticleReal* vx, amrex::ParticleReal* vy, amrex::ParticleReal* vz, + int nattr, amrex::ParticleReal* attr, int uniqueparticles) { auto & mypc = WarpX::GetInstance().GetPartContainer(); auto & myspc = mypc.GetParticleContainer(speciesnumber); @@ -180,14 +187,14 @@ extern "C" ConvertLabParamsToBoost(); } - double warpx_getProbLo(int dir) + amrex::Real warpx_getProbLo(int dir) { WarpX& warpx = WarpX::GetInstance(); const amrex::Geometry& geom = warpx.Geom(0); return geom.ProbLo(dir); } - double warpx_getProbHi(int dir) + amrex::Real warpx_getProbHi(int dir) { WarpX& warpx = WarpX::GetInstance(); const amrex::Geometry& geom = warpx.Geom(0); @@ -200,7 +207,7 @@ extern "C" return myspc.TotalNumberOfParticles(); } - double** warpx_getEfield(int lev, int direction, + amrex::Real** warpx_getEfield(int lev, int direction, int *return_size, int *ncomps, int *ngrow, int **shapes) { auto & mf = WarpX::GetInstance().getEfield(lev, direction); return getMultiFabPointers(mf, return_size, ncomps, ngrow, shapes); @@ -212,7 +219,7 @@ extern "C" return getMultiFabLoVects(mf, return_size, ngrow); } - double** warpx_getEfieldCP(int lev, int direction, + amrex::Real** warpx_getEfieldCP(int lev, int direction, int *return_size, int *ncomps, int *ngrow, int **shapes) { auto & mf = WarpX::GetInstance().getEfield_cp(lev, direction); return getMultiFabPointers(mf, return_size, ncomps, ngrow, shapes); @@ -224,7 +231,7 @@ extern "C" return getMultiFabLoVects(mf, return_size, ngrow); } - double** warpx_getEfieldFP(int lev, int direction, + amrex::Real** warpx_getEfieldFP(int lev, int direction, int *return_size, int *ncomps, int *ngrow, int **shapes) { auto & mf = WarpX::GetInstance().getEfield_fp(lev, direction); return getMultiFabPointers(mf, return_size, ncomps, ngrow, shapes); @@ -236,7 +243,7 @@ extern "C" return getMultiFabLoVects(mf, return_size, ngrow); } - double** warpx_getBfield(int lev, int direction, + amrex::Real** warpx_getBfield(int lev, int direction, int *return_size, int *ncomps, int *ngrow, int **shapes) { auto & mf = WarpX::GetInstance().getBfield(lev, direction); return getMultiFabPointers(mf, return_size, ncomps, ngrow, shapes); @@ -248,7 +255,7 @@ extern "C" return getMultiFabLoVects(mf, return_size, ngrow); } - double** warpx_getBfieldCP(int lev, int direction, + amrex::Real** warpx_getBfieldCP(int lev, int direction, int *return_size, int *ncomps, int *ngrow, int **shapes) { auto & mf = WarpX::GetInstance().getBfield_cp(lev, direction); return getMultiFabPointers(mf, return_size, ncomps, ngrow, shapes); @@ -260,7 +267,7 @@ extern "C" return getMultiFabLoVects(mf, return_size, ngrow); } - double** warpx_getBfieldFP(int lev, int direction, + amrex::Real** warpx_getBfieldFP(int lev, int direction, int *return_size, int *ncomps, int *ngrow, int **shapes) { auto & mf = WarpX::GetInstance().getBfield_fp(lev, direction); return getMultiFabPointers(mf, return_size, ncomps, ngrow, shapes); @@ -272,7 +279,7 @@ extern "C" return getMultiFabLoVects(mf, return_size, ngrow); } - double** warpx_getCurrentDensity(int lev, int direction, + amrex::Real** warpx_getCurrentDensity(int lev, int direction, int *return_size, int *ncomps, int *ngrow, int **shapes) { auto & mf = WarpX::GetInstance().getcurrent(lev, direction); return getMultiFabPointers(mf, return_size, ncomps, ngrow, shapes); @@ -284,7 +291,7 @@ extern "C" return getMultiFabLoVects(mf, return_size, ngrow); } - double** warpx_getCurrentDensityCP(int lev, int direction, + amrex::Real** warpx_getCurrentDensityCP(int lev, int direction, int *return_size, int *ncomps, int *ngrow, int **shapes) { auto & mf = WarpX::GetInstance().getcurrent_cp(lev, direction); return getMultiFabPointers(mf, return_size, ncomps, ngrow, shapes); @@ -296,7 +303,7 @@ extern "C" return getMultiFabLoVects(mf, return_size, ngrow); } - double** warpx_getCurrentDensityFP(int lev, int direction, + amrex::Real** warpx_getCurrentDensityFP(int lev, int direction, int *return_size, int *ncomps, int *ngrow, int **shapes) { auto & mf = WarpX::GetInstance().getcurrent_fp(lev, direction); return getMultiFabPointers(mf, return_size, ncomps, ngrow, shapes); @@ -308,7 +315,7 @@ extern "C" return getMultiFabLoVects(mf, return_size, ngrow); } - double** warpx_getParticleStructs(int speciesnumber, int lev, + amrex::ParticleReal** warpx_getParticleStructs(int speciesnumber, int lev, int* num_tiles, int** particles_per_tile) { auto & mypc = WarpX::GetInstance().GetPartContainer(); auto & myspc = mypc.GetParticleContainer(speciesnumber); @@ -320,17 +327,17 @@ extern "C" *num_tiles = i; *particles_per_tile = (int*) malloc(*num_tiles*sizeof(int)); - double** data = (double**) malloc(*num_tiles*sizeof(typename WarpXParticleContainer::ParticleType*)); + amrex::ParticleReal** data = (amrex::ParticleReal**) malloc(*num_tiles*sizeof(typename WarpXParticleContainer::ParticleType*)); i = 0; for (WarpXParIter pti(myspc, lev); pti.isValid(); ++pti, ++i) { auto& aos = pti.GetArrayOfStructs(); - data[i] = (double*) aos.data(); + data[i] = (amrex::ParticleReal*) aos.data(); (*particles_per_tile)[i] = pti.numParticles(); } return data; } - double** warpx_getParticleArrays(int speciesnumber, int comp, int lev, + amrex::ParticleReal** warpx_getParticleArrays(int speciesnumber, int comp, int lev, int* num_tiles, int** particles_per_tile) { auto & mypc = WarpX::GetInstance().GetPartContainer(); auto & myspc = mypc.GetParticleContainer(speciesnumber); @@ -342,11 +349,11 @@ extern "C" *num_tiles = i; *particles_per_tile = (int*) malloc(*num_tiles*sizeof(int)); - double** data = (double**) malloc(*num_tiles*sizeof(double*)); + amrex::ParticleReal** data = (amrex::ParticleReal**) malloc(*num_tiles*sizeof(amrex::ParticleReal*)); i = 0; for (WarpXParIter pti(myspc, lev); pti.isValid(); ++pti, ++i) { auto& soa = pti.GetStructOfArrays(); - data[i] = (double*) soa.GetRealData(comp).dataPtr(); + data[i] = (amrex::ParticleReal*) soa.GetRealData(comp).dataPtr(); (*particles_per_tile)[i] = pti.numParticles(); } return data; @@ -361,11 +368,11 @@ extern "C" warpx.MoveWindow (true); } - void warpx_EvolveE (double dt) { + void warpx_EvolveE (amrex::Real dt) { WarpX& warpx = WarpX::GetInstance(); warpx.EvolveE (dt); } - void warpx_EvolveB (double dt) { + void warpx_EvolveB (amrex::Real dt) { WarpX& warpx = WarpX::GetInstance(); warpx.EvolveB (dt); } @@ -385,7 +392,7 @@ extern "C" WarpX& warpx = WarpX::GetInstance(); warpx.UpdateAuxilaryData (); } - void warpx_PushParticlesandDepose (double cur_time) { + void warpx_PushParticlesandDepose (amrex::Real cur_time) { WarpX& warpx = WarpX::GetInstance(); warpx.PushParticlesandDepose (cur_time); } @@ -398,15 +405,15 @@ extern "C" WarpX& warpx = WarpX::GetInstance(); warpx.setistep (lev, ii); } - double warpx_gett_new (int lev) { + amrex::Real warpx_gett_new (int lev) { WarpX& warpx = WarpX::GetInstance(); return warpx.gett_new (lev); } - void warpx_sett_new (int lev, double time) { + void warpx_sett_new (int lev, amrex::Real time) { WarpX& warpx = WarpX::GetInstance(); warpx.sett_new (lev, time); } - double warpx_getdt (int lev) { + amrex::Real warpx_getdt (int lev) { WarpX& warpx = WarpX::GetInstance(); return warpx.getdt (lev); } @@ -415,7 +422,7 @@ extern "C" WarpX& warpx = WarpX::GetInstance(); return warpx.maxStep (); } - double warpx_stopTime () { + amrex::Real warpx_stopTime () { WarpX& warpx = WarpX::GetInstance(); return warpx.stopTime (); } diff --git a/Source/Python/WarpXWrappers.h b/Source/Python/WarpXWrappers.h index 233a0b930..a272b9250 100644 --- a/Source/Python/WarpXWrappers.h +++ b/Source/Python/WarpXWrappers.h @@ -1,6 +1,9 @@ #ifndef WARPX_WRAPPERS_H_ #define WARPX_WRAPPERS_H_ +#include <AMReX.H> +#include <AMReX_BLProfiler.H> + #ifdef BL_USE_MPI #include <mpi.h> #endif @@ -9,6 +12,9 @@ extern "C" { #endif + int warpx_Real_size(); + int warpx_ParticleReal_size(); + int warpx_nSpecies(); bool warpx_use_fdtd_nci_corr(); @@ -49,61 +55,61 @@ extern "C" { void warpx_evolve (int numsteps); // -1 means the inputs parameter will be used. 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); + amrex::ParticleReal* x, amrex::ParticleReal* y, amrex::ParticleReal* z, + amrex::ParticleReal* vx, amrex::ParticleReal* vy, amrex::ParticleReal* vz, + int nattr, amrex::ParticleReal* attr, int uniqueparticles); void warpx_ConvertLabParamsToBoost(); - double warpx_getProbLo(int dir); + amrex::Real warpx_getProbLo(int dir); - double warpx_getProbHi(int dir); + amrex::Real warpx_getProbHi(int dir); long warpx_getNumParticles(int speciesnumber); - double** warpx_getEfield(int lev, int direction, - int *return_size, int* ncomps, int* ngrow, int **shapes); + amrex::Real** warpx_getEfield(int lev, int direction, + int *return_size, int* ncomps, int* ngrow, int **shapes); int* warpx_getEfieldLoVects(int lev, int direction, int *return_size, int* ngrow); - double** warpx_getBfield(int lev, int direction, - int *return_size, int* ncomps, int* ngrow, int **shapes); + amrex::Real** warpx_getBfield(int lev, int direction, + int *return_size, int* ncomps, int* ngrow, int **shapes); int* warpx_getBfieldLoVects(int lev, int direction, int *return_size, int* ngrow); - double** warpx_getCurrentDensity(int lev, int direction, - int *return_size, int* ncomps, int* ngrow, int **shapes); + amrex::Real** warpx_getCurrentDensity(int lev, int direction, + int *return_size, int* ncomps, int* ngrow, int **shapes); int* warpx_getCurrentDensityLoVects(int lev, int direction, int *return_size, int* ngrow); - double** warpx_getParticleStructs(int speciesnumber, int lev, - int* num_tiles, int** particles_per_tile); + amrex::ParticleReal** warpx_getParticleStructs(int speciesnumber, int lev, + int* num_tiles, int** particles_per_tile); - double** warpx_getParticleArrays(int speciesnumber, int comp, int lev, - int* num_tiles, int** particles_per_tile); + amrex::ParticleReal** warpx_getParticleArrays(int speciesnumber, int comp, int lev, + int* num_tiles, int** particles_per_tile); void warpx_ComputeDt (); void warpx_MoveWindow (); - void warpx_EvolveE (double dt); - void warpx_EvolveB (double dt); + void warpx_EvolveE (amrex::Real dt); + void warpx_EvolveB (amrex::Real dt); void warpx_FillBoundaryE (); void warpx_FillBoundaryB (); void warpx_SyncCurrent (); void warpx_UpdateAuxilaryData (); - void warpx_PushParticlesandDepose (double cur_time); + void warpx_PushParticlesandDepose (amrex::Real cur_time); int warpx_getistep (int lev); void warpx_setistep (int lev, int ii); - double warpx_gett_new (int lev); - void warpx_sett_new (int lev, double time); - double warpx_getdt (int lev); + amrex::Real warpx_gett_new (int lev); + void warpx_sett_new (int lev, amrex::Real time); + amrex::Real warpx_getdt (int lev); int warpx_maxStep (); - double warpx_stopTime (); + amrex::Real warpx_stopTime (); int warpx_checkInt (); int warpx_plotInt (); |