diff options
Diffstat (limited to 'Source/Parallelization')
-rw-r--r-- | Source/Parallelization/GuardCellManager.cpp | 9 | ||||
-rw-r--r-- | Source/Parallelization/WarpXComm_K.H | 6 | ||||
-rw-r--r-- | Source/Parallelization/WarpXRegrid.cpp | 4 |
3 files changed, 14 insertions, 5 deletions
diff --git a/Source/Parallelization/GuardCellManager.cpp b/Source/Parallelization/GuardCellManager.cpp index a4a3e78e2..05ea700d8 100644 --- a/Source/Parallelization/GuardCellManager.cpp +++ b/Source/Parallelization/GuardCellManager.cpp @@ -8,6 +8,7 @@ #include "Filter/NCIGodfreyFilter.H" #include <AMReX_Print.H> #include <AMReX_ParmParse.H> +#include <AMReX.H> using namespace amrex; @@ -111,10 +112,10 @@ guardCellManager::Init( pp.query("ny_guard", ngFFt_y); pp.query("nz_guard", ngFFt_z); -#if (AMREX_SPACEDIM == 3) - IntVect ngFFT = IntVect(ngFFt_x, ngFFt_y, ngFFt_z); +#if (AMREX_SPACEDIM == 3) + IntVect ngFFT = IntVect(ngFFt_x, ngFFt_y, ngFFt_z); #elif (AMREX_SPACEDIM == 2) - IntVect ngFFT = IntVect(ngFFt_x, ngFFt_z); + IntVect ngFFT = IntVect(ngFFt_x, ngFFt_z); #endif for (int i_dim=0; i_dim<AMREX_SPACEDIM; i_dim++ ){ @@ -132,6 +133,8 @@ guardCellManager::Init( ng_alloc_F_int = ng_required; } ng_alloc_F = IntVect(AMREX_D_DECL(ng_alloc_F_int, ng_alloc_F_int, ng_alloc_F_int)); +#else + ignore_unused(nox_fft, noy_fft, noz_fft); #endif ng_Extra = IntVect(static_cast<int>(aux_is_nodal and !do_nodal)); diff --git a/Source/Parallelization/WarpXComm_K.H b/Source/Parallelization/WarpXComm_K.H index 1b6eceb93..b490636ee 100644 --- a/Source/Parallelization/WarpXComm_K.H +++ b/Source/Parallelization/WarpXComm_K.H @@ -8,6 +8,7 @@ #define WARPX_COMM_K_H_ #include <AMReX_FArrayBox.H> +#include <AMReX.H> AMREX_GPU_DEVICE AMREX_FORCE_INLINE void warpx_interp_bfield_x (int j, int k, int l, @@ -384,6 +385,7 @@ void warpx_interp_nd_bfield_x (int j, int k, int l, { #if (AMREX_SPACEDIM == 2) Bxa(j,k,0) = 0.5*(Bxf(j,k-1,0) + Bxf(j,k,0)); + amrex::ignore_unused(l); #else Bxa(j,k,l) = 0.25*(Bxf(j,k-1,l-1) + Bxf(j,k,l-1) + Bxf(j,k-1,l) + Bxf(j,k,l)); #endif @@ -396,6 +398,7 @@ void warpx_interp_nd_bfield_y (int j, int k, int l, { #if (AMREX_SPACEDIM == 2) Bya(j,k,0) = 0.25*(Byf(j,k,0) + Byf(j-1,k,0) + Byf(j,k-1,0) + Byf(j-1,k-1,0)); + amrex::ignore_unused(l); #else Bya(j,k,l) = 0.25*(Byf(j-1,k,l-1) + Byf(j,k,l-1) + Byf(j-1,k,l) + Byf(j,k,l)); #endif @@ -408,6 +411,7 @@ void warpx_interp_nd_bfield_z (int j, int k, int l, { #if (AMREX_SPACEDIM == 2) Bza(j,k,0) = 0.5*(Bzf(j-1,k,0) + Bzf(j,k,0)); + amrex::ignore_unused(l); #else Bza(j,k,l) = 0.25*(Bzf(j-1,k-1,l) + Bzf(j,k-1,l) + Bzf(j-1,k,l) + Bzf(j,k,l)); #endif @@ -632,6 +636,7 @@ void warpx_interp_nd_efield_y (int j, int k, int l, { #if (AMREX_SPACEDIM == 2) Eya(j,k,0) = Eyf(j,k,0); + amrex::ignore_unused(l); #else Eya(j,k,l) = 0.5*(Eyf(j,k-1,l) + Eyf(j,k,l)); #endif @@ -644,6 +649,7 @@ void warpx_interp_nd_efield_z (int j, int k, int l, { #if (AMREX_SPACEDIM == 2) Eza(j,k,0) = 0.5*(Ezf(j,k-1,0) + Ezf(j,k,0)); + amrex::ignore_unused(l); #else Eza(j,k,l) = 0.5*(Ezf(j,k,l-1) + Ezf(j,k,l)); #endif diff --git a/Source/Parallelization/WarpXRegrid.cpp b/Source/Parallelization/WarpXRegrid.cpp index 6c6710121..9381bf60d 100644 --- a/Source/Parallelization/WarpXRegrid.cpp +++ b/Source/Parallelization/WarpXRegrid.cpp @@ -58,8 +58,8 @@ WarpX::LoadBalance () // As specified in the above calls to makeSFC and makeKnapSack, the new // distribution mapping is NOT communicated to all ranks; the loadbalanced // dm is up-to-date only on root, and we can decide whether to broadcast - if (load_balance_efficiency_ratio_threshold > 0.0 - & ParallelDescriptor::MyProc() == ParallelDescriptor::IOProcessorNumber()) + if ((load_balance_efficiency_ratio_threshold > 0.0) + && (ParallelDescriptor::MyProc() == ParallelDescriptor::IOProcessorNumber())) { doLoadBalance = (proposedEfficiency > load_balance_efficiency_ratio_threshold*currentEfficiency); } |