aboutsummaryrefslogtreecommitdiff
path: root/Source
diff options
context:
space:
mode:
Diffstat (limited to 'Source')
-rw-r--r--Source/FieldSolver/FiniteDifferenceSolver/EvolveB.cpp30
-rw-r--r--Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceSolver.H16
-rw-r--r--Source/FieldSolver/Make.package3
-rw-r--r--Source/WarpX.H6
-rw-r--r--Source/WarpX.cpp22
5 files changed, 36 insertions, 41 deletions
diff --git a/Source/FieldSolver/FiniteDifferenceSolver/EvolveB.cpp b/Source/FieldSolver/FiniteDifferenceSolver/EvolveB.cpp
index 9252237bd..647e27d8a 100644
--- a/Source/FieldSolver/FiniteDifferenceSolver/EvolveB.cpp
+++ b/Source/FieldSolver/FiniteDifferenceSolver/EvolveB.cpp
@@ -1,13 +1,13 @@
-#include <WarpXAlgorithmSelection.H>
-#include<FiniteDifferenceAlgorithms/YeeAlgorithm.H>
-#include<FiniteDifferenceSolver.H>
-#include<AMReX_Gpu.H>
+#include "WarpXAlgorithmSelection.H"
+#include "FiniteDifferenceAlgorithms/YeeAlgorithm.H"
+#include "FiniteDifferenceSolver.H"
+#include <AMReX_Gpu.H>
using namespace amrex;
void FiniteDifferenceSolver::EvolveB ( VectorField& Bfield,
- VectorField& Efield,
- amrex::Real dt ) {
+ VectorField const& Efield,
+ amrex::Real const dt ) {
// Select algorithm (The choice of algorithm is a runtime option,
// but we compile code for each algorithm, using templates)
if (m_fdtd_algo == MaxwellSolverAlgo::Yee){
@@ -19,10 +19,10 @@ void FiniteDifferenceSolver::EvolveB ( VectorField& Bfield,
}
}
-template<typename algo>
+template<typename T_Algo>
void FiniteDifferenceSolver::EvolveBwithAlgo ( VectorField& Bfield,
- VectorField& Efield,
- amrex::Real dt ) {
+ VectorField const& Efield,
+ amrex::Real const dt ) {
// Loop through the grids, and over the tiles within each grid
#ifdef _OPENMP
@@ -55,18 +55,18 @@ void FiniteDifferenceSolver::EvolveBwithAlgo ( VectorField& Bfield,
amrex::ParallelFor(tbx, tby, tbz,
[=] AMREX_GPU_DEVICE (int i, int j, int k){
- Bx(i, j, k) += dt * algo::UpwardDz(Ey, coefs_z, n_coefs_z, i, j, k)
- - dt * algo::UpwardDy(Ez, coefs_y, n_coefs_y, i, j, k);
+ Bx(i, j, k) += dt * T_Algo::UpwardDz(Ey, coefs_z, n_coefs_z, i, j, k)
+ - dt * T_Algo::UpwardDy(Ez, coefs_y, n_coefs_y, i, j, k);
},
[=] AMREX_GPU_DEVICE (int i, int j, int k){
- By(i, j, k) += dt * algo::UpwardDx(Ez, coefs_x, n_coefs_x, i, j, k)
- - dt * algo::UpwardDz(Ex, coefs_z, n_coefs_z, i, j, k);
+ By(i, j, k) += dt * T_Algo::UpwardDx(Ez, coefs_x, n_coefs_x, i, j, k)
+ - dt * T_Algo::UpwardDz(Ex, coefs_z, n_coefs_z, i, j, k);
},
[=] AMREX_GPU_DEVICE (int i, int j, int k){
- Bz(i, j, k) += dt * algo::UpwardDy(Ex, coefs_y, n_coefs_y, i, j, k)
- - dt * algo::UpwardDx(Ey, coefs_x, n_coefs_x, i, j, k);
+ Bz(i, j, k) += dt * T_Algo::UpwardDy(Ex, coefs_y, n_coefs_y, i, j, k)
+ - dt * T_Algo::UpwardDx(Ey, coefs_x, n_coefs_x, i, j, k);
}
);
diff --git a/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceSolver.H b/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceSolver.H
index 5496adff9..8653ce62b 100644
--- a/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceSolver.H
+++ b/Source/FieldSolver/FiniteDifferenceSolver/FiniteDifferenceSolver.H
@@ -1,9 +1,9 @@
#ifndef WARPX_FINITE_DIFFERENCE_SOLVER_H_
#define WARPX_FINITE_DIFFERENCE_SOLVER_H_
+#include "WarpXAlgorithmSelection.H"
+#include "FiniteDifferenceAlgorithms/YeeAlgorithm.H"
#include <AMReX_MultiFab.H>
-#include <WarpXAlgorithmSelection.H>
-#include<FiniteDifferenceAlgorithms/YeeAlgorithm.H>
/**
* \brief Top-level class for the electromagnetic finite-difference solver
@@ -17,7 +17,7 @@ class FiniteDifferenceSolver
using VectorField = std::array< std::unique_ptr<amrex::MultiFab>, 3 >;
// Constructor
- FiniteDifferenceSolver ( int fdtd_algo,
+ FiniteDifferenceSolver ( int const fdtd_algo,
std::array<amrex::Real,3> cell_size ) {
// Register the type of finite-difference algorithm
@@ -36,8 +36,8 @@ class FiniteDifferenceSolver
};
void EvolveB ( VectorField& Bfield,
- VectorField& Efield,
- amrex::Real dt );
+ VectorField const& Efield,
+ amrex::Real const dt );
private:
int m_fdtd_algo;
@@ -46,10 +46,10 @@ class FiniteDifferenceSolver
amrex::Gpu::ManagedVector<amrex::Real> stencil_coefs_y;
amrex::Gpu::ManagedVector<amrex::Real> stencil_coefs_z;
- template< typename fdtd_algo >
+ template< typename T_Algo >
void EvolveBwithAlgo ( VectorField& Bfield,
- VectorField& Efield,
- amrex::Real dt );
+ VectorField const& Efield,
+ amrex::Real const dt );
};
diff --git a/Source/FieldSolver/Make.package b/Source/FieldSolver/Make.package
index a7c8ce671..522c4c07a 100644
--- a/Source/FieldSolver/Make.package
+++ b/Source/FieldSolver/Make.package
@@ -6,9 +6,8 @@ ifeq ($(USE_PSATD),TRUE)
ifeq ($(USE_PSATD_PICSAR),TRUE)
include $(WARPX_HOME)/Source/FieldSolver/PicsarHybridSpectralSolver/Make.package
endif
-else
- include $(WARPX_HOME)/Source/FieldSolver/FiniteDifferenceSolver/Make.package
endif
+include $(WARPX_HOME)/Source/FieldSolver/FiniteDifferenceSolver/Make.package
INCLUDE_LOCATIONS += $(WARPX_HOME)/Source/FieldSolver
VPATH_LOCATIONS += $(WARPX_HOME)/Source/FieldSolver
diff --git a/Source/WarpX.H b/Source/WarpX.H
index 661b174b7..0b4aac0ce 100644
--- a/Source/WarpX.H
+++ b/Source/WarpX.H
@@ -8,10 +8,9 @@
#include <BilinearFilter.H>
#include <NCIGodfreyFilter.H>
+#include <FiniteDifferenceSolver.H>
#ifdef WARPX_USE_PSATD
# include <SpectralSolver.H>
-#else
-# include <FiniteDifferenceSolver.H>
#endif
#ifdef WARPX_USE_PSATD_HYBRID
# include <PicsarHybridFFTData.H>
@@ -696,10 +695,9 @@ private:
amrex::Vector<std::unique_ptr<SpectralSolver>> spectral_solver_fp;
amrex::Vector<std::unique_ptr<SpectralSolver>> spectral_solver_cp;
-#else
+#endif
amrex::Vector<std::unique_ptr<FiniteDifferenceSolver>> fdtd_solver_fp;
amrex::Vector<std::unique_ptr<FiniteDifferenceSolver>> fdtd_solver_cp;
-#endif
#ifdef WARPX_USE_PSATD_HYBRID
private:
diff --git a/Source/WarpX.cpp b/Source/WarpX.cpp
index c591a2856..c9084cd2a 100644
--- a/Source/WarpX.cpp
+++ b/Source/WarpX.cpp
@@ -893,11 +893,11 @@ WarpX::AllocLevelMFs (int lev, const BoxArray& ba, const DistributionMapping& dm
if (fft_hybrid_mpi_decomposition == false){
// Allocate and initialize the spectral solver
std::array<Real,3> dx = CellSize(lev);
- #if (AMREX_SPACEDIM == 3)
+#if (AMREX_SPACEDIM == 3)
RealVect dx_vect(dx[0], dx[1], dx[2]);
- #elif (AMREX_SPACEDIM == 2)
+#elif (AMREX_SPACEDIM == 2)
RealVect dx_vect(dx[0], dx[2]);
- #endif
+#endif
// Get the cell-centered box, with guard cells
BoxArray realspace_ba = ba; // Copy box
realspace_ba.enclosedCells().grow(ngE); // cell-centered + guard cells
@@ -905,11 +905,10 @@ WarpX::AllocLevelMFs (int lev, const BoxArray& ba, const DistributionMapping& dm
spectral_solver_fp[lev].reset( new SpectralSolver( realspace_ba, dm,
nox_fft, noy_fft, noz_fft, do_nodal, dx_vect, dt[lev] ) );
}
-#else
- std::array<Real,3> dx = CellSize(lev);
- fdtd_solver_fp[lev].reset(
- new FiniteDifferenceSolver(maxwell_fdtd_solver_id, dx) );
#endif
+ std::array<Real,3> const dx = CellSize(lev);
+ fdtd_solver_fp[lev].reset(
+ new FiniteDifferenceSolver(maxwell_fdtd_solver_id, dx) );
//
// The Aux patch (i.e., the full solution)
//
@@ -981,11 +980,11 @@ WarpX::AllocLevelMFs (int lev, const BoxArray& ba, const DistributionMapping& dm
if (fft_hybrid_mpi_decomposition == false){
// Allocate and initialize the spectral solver
std::array<Real,3> cdx = CellSize(lev-1);
- #if (AMREX_SPACEDIM == 3)
+#if (AMREX_SPACEDIM == 3)
RealVect cdx_vect(cdx[0], cdx[1], cdx[2]);
- #elif (AMREX_SPACEDIM == 2)
+#elif (AMREX_SPACEDIM == 2)
RealVect cdx_vect(cdx[0], cdx[2]);
- #endif
+#endif
// Get the cell-centered box, with guard cells
BoxArray realspace_ba = cba;// Copy box
realspace_ba.enclosedCells().grow(ngE);// cell-centered + guard cells
@@ -993,11 +992,10 @@ WarpX::AllocLevelMFs (int lev, const BoxArray& ba, const DistributionMapping& dm
spectral_solver_cp[lev].reset( new SpectralSolver( realspace_ba, dm,
nox_fft, noy_fft, noz_fft, do_nodal, cdx_vect, dt[lev] ) );
}
-#else
+#endif
std::array<Real,3> dx = CellSize(lev-1);
fdtd_solver_cp[lev].reset(
new FiniteDifferenceSolver( maxwell_fdtd_solver_id, dx ) );
-#endif
}
//