From 2a3fbdacbfb7edfd279c4c92543597bf0b475dd3 Mon Sep 17 00:00:00 2001 From: Luca Fedeli Date: Mon, 13 May 2024 13:15:14 +0200 Subject: [PATCH] ablastr::particles::compute_weights : implement 1D and use template parameter to specify if field is nodal (#4846) * use template parameter to specify if field is nodal in compute_weights and implement function also for 1D case * fix bug * remove unused variables * use amrex::IndexType::NODE as template parameter * Add Missing Include * Doc: Remove `(default)` Co-authored-by: Axel Huebl --- Source/Diagnostics/ParticleIO.cpp | 3 +- Source/EmbeddedBoundary/ParticleScraper.H | 7 +-- Source/Particles/ParticleBoundaryBuffer.cpp | 10 ++-- Source/ablastr/particles/NodalFieldGather.H | 58 +++++++++++++-------- 4 files changed, 48 insertions(+), 30 deletions(-) diff --git a/Source/Diagnostics/ParticleIO.cpp b/Source/Diagnostics/ParticleIO.cpp index c4164e383b1..f1f4d426a4a 100644 --- a/Source/Diagnostics/ParticleIO.cpp +++ b/Source/Diagnostics/ParticleIO.cpp @@ -283,7 +283,8 @@ storePhiOnParticles ( PinnedMemoryParticleContainer& tmp, getPosition(ip, xp, yp, zp); int i, j, k; amrex::Real W[AMREX_SPACEDIM][2]; - ablastr::particles::compute_weights(xp, yp, zp, plo, dxi, i, j, k, W); + ablastr::particles::compute_weights( + xp, yp, zp, plo, dxi, i, j, k, W); amrex::Real const phi_value = ablastr::particles::interp_field_nodal(i, j, k, W, phi_grid); phi_particle_arr[ip] = phi_value; } diff --git a/Source/EmbeddedBoundary/ParticleScraper.H b/Source/EmbeddedBoundary/ParticleScraper.H index ea4d0bc5c0e..1e915b39381 100644 --- a/Source/EmbeddedBoundary/ParticleScraper.H +++ b/Source/EmbeddedBoundary/ParticleScraper.H @@ -183,15 +183,16 @@ scrapeParticlesAtEB (PC& pc, const amrex::Vector& distan int i, j, k; amrex::Real W[AMREX_SPACEDIM][2]; - ablastr::particles::compute_weights(xp, yp, zp, plo, dxi, i, j, k, W); + ablastr::particles::compute_weights( + xp, yp, zp, plo, dxi, i, j, k, W); amrex::Real phi_value = ablastr::particles::interp_field_nodal(i, j, k, W, phi); if (phi_value < 0.0) { int ic, jc, kc; // Cell-centered indices - [[maybe_unused]] int nodal; amrex::Real Wc[AMREX_SPACEDIM][2]; // Cell-centered weights - ablastr::particles::compute_weights(xp, yp, zp, plo, dxi, ic, jc, kc, Wc, nodal=0); + ablastr::particles::compute_weights( + xp, yp, zp, plo, dxi, ic, jc, kc, Wc); amrex::RealVect normal = DistanceToEB::interp_normal(i, j, k, W, ic, jc, kc, Wc, phi, dxi); DistanceToEB::normalize(normal); amrex::RealVect pos; diff --git a/Source/Particles/ParticleBoundaryBuffer.cpp b/Source/Particles/ParticleBoundaryBuffer.cpp index c8c683f0abf..345a0bfa592 100644 --- a/Source/Particles/ParticleBoundaryBuffer.cpp +++ b/Source/Particles/ParticleBoundaryBuffer.cpp @@ -98,7 +98,8 @@ struct FindEmbeddedBoundaryIntersection { amrex::Real W[AMREX_SPACEDIM][2]; amrex::ParticleReal x_temp=xp, y_temp=yp, z_temp=zp; UpdatePosition(x_temp, y_temp, z_temp, ux, uy, uz, -dt_frac*dt); - ablastr::particles::compute_weights(x_temp, y_temp, z_temp, plo, dxi, i, j, k, W); + ablastr::particles::compute_weights( + x_temp, y_temp, z_temp, plo, dxi, i, j, k, W); amrex::Real phi_value = ablastr::particles::interp_field_nodal(i, j, k, W, phiarr); return phi_value; } ); @@ -115,11 +116,12 @@ struct FindEmbeddedBoundaryIntersection { // record the components of the normal on the destination int i, j, k; amrex::Real W[AMREX_SPACEDIM][2]; - ablastr::particles::compute_weights(x_temp, y_temp, z_temp, plo, dxi, i, j, k, W); + ablastr::particles::compute_weights( + x_temp, y_temp, z_temp, plo, dxi, i, j, k, W); int ic, jc, kc; // Cell-centered indices - int nodal; amrex::Real Wc[AMREX_SPACEDIM][2]; // Cell-centered weight - ablastr::particles::compute_weights(x_temp, y_temp, z_temp, plo, dxi, ic, jc, kc, Wc, nodal=0); // nodal=0 to calculate the weights with respect to the cell-centered nodes + ablastr::particles::compute_weights( + x_temp, y_temp, z_temp, plo, dxi, ic, jc, kc, Wc); amrex::RealVect normal = DistanceToEB::interp_normal(i, j, k, W, ic, jc, kc, Wc, phiarr, dxi); DistanceToEB::normalize(normal); diff --git a/Source/ablastr/particles/NodalFieldGather.H b/Source/ablastr/particles/NodalFieldGather.H index 5a10c73ae20..f324264db7f 100644 --- a/Source/ablastr/particles/NodalFieldGather.H +++ b/Source/ablastr/particles/NodalFieldGather.H @@ -12,16 +12,19 @@ #include #include #include +#include #include #include + namespace ablastr::particles { + /** * \brief Compute weight of each surrounding node (or cell-centered nodes) in interpolating a nodal ((or a cell-centered node) field - * to the given coordinates. If nodal=1, then the calculations will be done with respect to the nodes (default). If nodal=0, then the calculations will be done with respect to the cell-centered nodal) - * - * This currently only does linear order. + * to the given coordinates. If template parameter IdxType is amrex::IndexType::CellIndex::NODE, then the calculations will be done with respect to the nodes. + * If template parameter IdxType is amrex::IndexType::CellIndex::CELL, then the calculations will be done with respect to the cell-centered nodal. + * This currently only does linear order. * * \param xp,yp,zp Particle position coordinates * \param plo Index lower bounds of domain. @@ -30,24 +33,30 @@ namespace ablastr::particles * \param W 2D array of weights to store each neighbouring node (or cell-centered node) * \param nodal Int that tells if the weights are calculated in respect to the nodes (nodal=1) of the cell-centered nodes (nodal=0) */ +template AMREX_GPU_HOST_DEVICE AMREX_INLINE void compute_weights (const amrex::ParticleReal xp, const amrex::ParticleReal yp, const amrex::ParticleReal zp, amrex::GpuArray const& plo, amrex::GpuArray const& dxi, - int& i, int& j, int& k, amrex::Real W[AMREX_SPACEDIM][2], int nodal=1) noexcept + int& i, int& j, int& k, amrex::Real W[AMREX_SPACEDIM][2]) noexcept { using namespace amrex::literals; -#if !((nodal==0)||(nodal==1)) - ABLASTR_ABORT_WITH_MESSAGE("Error: 'nodal' has to be equal to 0 or 1"); -#endif + constexpr auto half_if_cell_centered = [](){ + if constexpr (IdxType == amrex::IndexType::CellIndex::NODE){ + return 0.0_rt; + } + else{ + return 0.5_rt; + } + }(); #if (defined WARPX_DIM_3D) - const amrex::Real x = (xp - plo[0]) * dxi[0] + static_cast(nodal-1)*0.5_rt; - const amrex::Real y = (yp - plo[1]) * dxi[1] + static_cast(nodal-1)*0.5_rt; - const amrex::Real z = (zp - plo[2]) * dxi[2] + static_cast(nodal-1)*0.5_rt; + const amrex::Real x = (xp - plo[0]) * dxi[0] - half_if_cell_centered; + const amrex::Real y = (yp - plo[1]) * dxi[1] - half_if_cell_centered; + const amrex::Real z = (zp - plo[2]) * dxi[2] - half_if_cell_centered; i = static_cast(amrex::Math::floor(x)); j = static_cast(amrex::Math::floor(y)); @@ -64,17 +73,17 @@ void compute_weights (const amrex::ParticleReal xp, #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ) # if (defined WARPX_DIM_XZ) - const amrex::Real x = (xp - plo[0]) * dxi[0] + static_cast(nodal-1)*0.5_rt; + const amrex::Real x = (xp - plo[0]) * dxi[0] - half_if_cell_centered; amrex::ignore_unused(yp); i = static_cast(amrex::Math::floor(x)); W[0][1] = x - i; # elif (defined WARPX_DIM_RZ) - const amrex::Real r = (std::sqrt(xp*xp+yp*yp) - plo[0]) * dxi[0] + static_cast(nodal-1)*0.5_rt; + const amrex::Real r = (std::sqrt(xp*xp+yp*yp) - plo[0]) * dxi[0] - half_if_cell_centered; i = static_cast(amrex::Math::floor(r)); W[0][1] = r - i; # endif - const amrex::Real z = (zp - plo[1]) * dxi[1] + static_cast(nodal-1)*0.5_rt; + const amrex::Real z = (zp - plo[1]) * dxi[1] - half_if_cell_centered; j = static_cast(amrex::Math::floor(z)); W[1][1] = z - j; @@ -83,8 +92,15 @@ void compute_weights (const amrex::ParticleReal xp, k = 0; #else - amrex::ignore_unused(xp, yp, zp, plo, dxi, i, j, k, W, nodal); - ABLASTR_ABORT_WITH_MESSAGE("Error: compute_weights not yet implemented in 1D"); + const amrex::Real z = (zp - plo[0]) * dxi[0] - half_if_cell_centered; + amrex::ignore_unused(xp, yp); + i = static_cast(amrex::Math::floor(z)); + + W[0][1] = z - i; + W[0][0] = 1.0_rt - W[0][1]; + + j = 0; + k = 0; #endif } @@ -100,8 +116,8 @@ amrex::Real interp_field_nodal (int i, int j, int k, const amrex::Real W[AMREX_SPACEDIM][2], amrex::Array4 const& scalar_field) noexcept { -#if (defined WARPX_DIM_3D) amrex::Real value = 0; +#if (defined WARPX_DIM_3D) value += scalar_field(i, j , k ) * W[0][0] * W[1][0] * W[2][0]; value += scalar_field(i+1, j , k ) * W[0][1] * W[1][0] * W[2][0]; value += scalar_field(i, j+1, k ) * W[0][0] * W[1][1] * W[2][0]; @@ -111,15 +127,13 @@ amrex::Real interp_field_nodal (int i, int j, int k, value += scalar_field(i , j+1, k+1) * W[0][0] * W[1][1] * W[2][1]; value += scalar_field(i+1, j+1, k+1) * W[0][1] * W[1][1] * W[2][1]; #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ) - amrex::Real value = 0; value += scalar_field(i, j , k) * W[0][0] * W[1][0]; value += scalar_field(i+1, j , k) * W[0][1] * W[1][0]; value += scalar_field(i, j+1, k) * W[0][0] * W[1][1]; value += scalar_field(i+1, j+1, k) * W[0][1] * W[1][1]; #else - const amrex::Real value = 0; - amrex::ignore_unused(i, j, k, W, scalar_field); - ABLASTR_ABORT_WITH_MESSAGE("Error: interp_field not yet implemented in 1D"); + value += scalar_field(i, j , k) * W[0][0]; + value += scalar_field(i+1, j , k) * W[0][1]; #endif return value; } @@ -144,7 +158,7 @@ amrex::Real doGatherScalarFieldNodal (const amrex::ParticleReal xp, // first find the weight of surrounding nodes to use during interpolation int ii, jj, kk; amrex::Real W[AMREX_SPACEDIM][2]; - compute_weights(xp, yp, zp, lo, dxi, ii, jj, kk, W); + compute_weights(xp, yp, zp, lo, dxi, ii, jj, kk, W); return interp_field_nodal(ii, jj, kk, W, scalar_field); } @@ -172,7 +186,7 @@ doGatherVectorFieldNodal (const amrex::ParticleReal xp, // first find the weight of surrounding nodes to use during interpolation int ii, jj, kk; amrex::Real W[AMREX_SPACEDIM][2]; - compute_weights(xp, yp, zp, lo, dxi, ii, jj, kk, W); + compute_weights(xp, yp, zp, lo, dxi, ii, jj, kk, W); amrex::GpuArray const field_interp = { interp_field_nodal(ii, jj, kk, W, vector_field_x),