diff --git a/Docs/source/usage/parameters.rst b/Docs/source/usage/parameters.rst index b9d82d5014a..54b3493c978 100644 --- a/Docs/source/usage/parameters.rst +++ b/Docs/source/usage/parameters.rst @@ -971,16 +971,21 @@ Particle initialization The ``external_file`` option is currently implemented for 2D, 3D and RZ geometries, with record components in the cartesian coordinates ``(x,y,z)`` for 3D and RZ, and ``(x,z)`` for 2D. For more information on the `openPMD format `__ and how to build WarpX with it, please visit :ref:`the install section `. - * ``NFluxPerCell``: Continuously inject a flux of macroparticles from a planar surface. + * ``NFluxPerCell``: Continuously inject a flux of macroparticles from a surface. The emitting surface can be chosen to be either a plane + defined by the user (using some of the parameters listed below), or the embedded boundary (see :ref:`Embedded Boundary Conditions `). This requires the additional parameters: * ``.flux_profile`` (see the description of this parameter further below) - * ``.surface_flux_pos`` (`double`, location of the injection plane [meter]) + * ``.inject_from_embedded_boundary`` (`0` or `1`, default `0` ; whether to inject from the embedded boundary or from a user-specified plane. + When injecting from the embedded boundary, the momentum distribution specified by the user along ``z`` (see e.g. ``uz_m``, ``uz_th`` below) is interpreted + as the momentum distribution along the local normal to the embedded boundary.) - * ``.flux_normal_axis`` (`x`, `y`, or `z` for 3D, `x` or `z` for 2D, or `r`, `t`, or `z` for RZ. When `flux_normal_axis` is `r` or `t`, the `x` and `y` components of the user-specified momentum distribution are interpreted as the `r` and `t` components respectively) + * ``.surface_flux_pos`` (only used when injecting from a plane, `double`, location of the injection plane [meter]) - * ``.flux_direction`` (`-1` or `+1`, direction of flux relative to the plane) + * ``.flux_normal_axis`` (only used when injecting from a plane, `x`, `y`, or `z` for 3D, `x` or `z` for 2D, or `r`, `t`, or `z` for RZ. When `flux_normal_axis` is `r` or `t`, the `x` and `y` components of the user-specified momentum distribution are interpreted as the `r` and `t` components respectively) + + * ``.flux_direction`` (only used when injecting from a plane, `-1` or `+1`, direction of flux relative to the plane) * ``.num_particles_per_cell`` (`double`) diff --git a/Examples/Tests/flux_injection/CMakeLists.txt b/Examples/Tests/flux_injection/CMakeLists.txt index d09b83d7618..0929fc3d4c4 100644 --- a/Examples/Tests/flux_injection/CMakeLists.txt +++ b/Examples/Tests/flux_injection/CMakeLists.txt @@ -20,3 +20,33 @@ add_warpx_test( diags/diag1000120 # output OFF # dependency ) + +add_warpx_test( + test_3d_flux_injection_from_eb # name + 3 # dims + 2 # nprocs + inputs_test_3d_flux_injection_from_eb # inputs + analysis_flux_injection_from_eb.py # analysis + diags/diag1000010 # output + OFF # dependency +) + +add_warpx_test( + test_rz_flux_injection_from_eb # name + RZ # dims + 2 # nprocs + inputs_test_rz_flux_injection_from_eb # inputs + analysis_flux_injection_from_eb.py # analysis + diags/diag1000010 # output + OFF # dependency +) + +add_warpx_test( + test_2d_flux_injection_from_eb # name + 2 # dims + 2 # nprocs + inputs_test_2d_flux_injection_from_eb # inputs + analysis_flux_injection_from_eb.py # analysis + diags/diag1000010 # output + OFF # dependency +) diff --git a/Examples/Tests/flux_injection/analysis_flux_injection_from_eb.py b/Examples/Tests/flux_injection/analysis_flux_injection_from_eb.py new file mode 100755 index 00000000000..36ff50bea06 --- /dev/null +++ b/Examples/Tests/flux_injection/analysis_flux_injection_from_eb.py @@ -0,0 +1,161 @@ +#!/usr/bin/env python3 +# +# Copyright 2024 Remi Lehe +# +# This file is part of WarpX. +# +# License: BSD-3-Clause-LBNL + +""" +This script tests the emission of particles from the embedded boundary. +(In this case, the embedded boundary is a sphere in 3D and RZ, a cylinder in 2D.) +We check that the embedded boundary emits the correct number of particles, and that +the particle distributions are consistent with the expected distributions. +""" + +import os +import re +import sys + +import matplotlib.pyplot as plt +import numpy as np +import yt +from scipy.constants import c, m_e +from scipy.special import erf + +sys.path.insert(1, "../../../../warpx/Regression/Checksum/") +import checksumAPI + +yt.funcs.mylog.setLevel(0) + +# Open plotfile specified in command line +fn = sys.argv[1] +ds = yt.load(fn) +ad = ds.all_data() +t_max = ds.current_time.item() # time of simulation + +# Extract the dimensionality of the simulation +with open("./warpx_used_inputs", "r") as f: + warpx_used_inputs = f.read() +if re.search("geometry.dims = 2", warpx_used_inputs): + dims = "2D" +elif re.search("geometry.dims = RZ", warpx_used_inputs): + dims = "RZ" +elif re.search("geometry.dims = 3", warpx_used_inputs): + dims = "3D" + +# Total number of electrons expected: +# Simulation parameters determine the total number of particles emitted (Ntot) +flux = 1.0 # in m^-2.s^-1, from the input script +R = 2.0 # in m, radius of the sphere +if dims == "3D" or dims == "RZ": + emission_surface = 4 * np.pi * R**2 # in m^2 +elif dims == "2D": + emission_surface = 2 * np.pi * R # in m +Ntot = flux * emission_surface * t_max + +# Parameters of the histogram +hist_bins = 50 +hist_range = [-0.5, 0.5] + + +# Define function that histograms and checks the data +def gaussian_dist(u, u_th): + return 1.0 / ((2 * np.pi) ** 0.5 * u_th) * np.exp(-(u**2) / (2 * u_th**2)) + + +def gaussian_flux_dist(u, u_th, u_m): + normalization_factor = u_th**2 * np.exp(-(u_m**2) / (2 * u_th**2)) + ( + np.pi / 2 + ) ** 0.5 * u_m * u_th * (1 + erf(u_m / (2**0.5 * u_th))) + result = ( + 1.0 + / normalization_factor + * np.where(u > 0, u * np.exp(-((u - u_m) ** 2) / (2 * u_th**2)), 0) + ) + return result + + +def compare_gaussian(u, w, u_th, label=""): + du = (hist_range[1] - hist_range[0]) / hist_bins + w_hist, u_hist = np.histogram(u, bins=hist_bins, weights=w / du, range=hist_range) + u_hist = 0.5 * (u_hist[1:] + u_hist[:-1]) + w_th = Ntot * gaussian_dist(u_hist, u_th) + plt.plot(u_hist, w_hist, label=label + ": simulation") + plt.plot(u_hist, w_th, "--", label=label + ": theory") + assert np.allclose(w_hist, w_th, atol=0.07 * w_th.max()) + + +def compare_gaussian_flux(u, w, u_th, u_m, label=""): + du = (hist_range[1] - hist_range[0]) / hist_bins + w_hist, u_hist = np.histogram(u, bins=hist_bins, weights=w / du, range=hist_range) + u_hist = 0.5 * (u_hist[1:] + u_hist[:-1]) + w_th = Ntot * gaussian_flux_dist(u_hist, u_th, u_m) + plt.plot(u_hist, w_hist, label=label + ": simulation") + plt.plot(u_hist, w_th, "--", label=label + ": theory") + assert np.allclose(w_hist, w_th, atol=0.05 * w_th.max()) + + +# Load data and perform check + +plt.figure() + +if dims == "3D": + x = ad["electron", "particle_position_x"].to_ndarray() + y = ad["electron", "particle_position_y"].to_ndarray() + z = ad["electron", "particle_position_z"].to_ndarray() +elif dims == "2D": + x = ad["electron", "particle_position_x"].to_ndarray() + y = np.zeros_like(x) + z = ad["electron", "particle_position_y"].to_ndarray() +elif dims == "RZ": + theta = ad["electron", "particle_theta"].to_ndarray() + r = ad["electron", "particle_position_x"].to_ndarray() + x = r * np.cos(theta) + y = r * np.sin(theta) + z = ad["electron", "particle_position_y"].to_ndarray() +ux = ad["electron", "particle_momentum_x"].to_ndarray() / (m_e * c) +uy = ad["electron", "particle_momentum_y"].to_ndarray() / (m_e * c) +uz = ad["electron", "particle_momentum_z"].to_ndarray() / (m_e * c) +w = ad["electron", "particle_weight"].to_ndarray() + +# Check that the total number of particles emitted is correct +Ntot_sim = np.sum(w) +print("Ntot_sim = ", Ntot_sim) +print("Ntot = ", Ntot) +assert np.isclose(Ntot_sim, Ntot, rtol=0.01) + +# Check that none of the particles are inside the EB +# A factor 0.98 is applied to accomodate +# the cut-cell approximation of the sphere +assert np.all(x**2 + y**2 + z**2 > (0.98 * R) ** 2) + +# Check that the normal component of the velocity is consistent with the expected distribution +r = np.sqrt(x**2 + y**2 + z**2) +nx = x / r +ny = y / r +nz = z / r +u_n = ux * nx + uy * ny + uz * nz # normal component +compare_gaussian_flux(u_n, w, u_th=0.1, u_m=0.07, label="u_n") + +# Pick a direction that is orthogonal to the normal direction, and check the distribution +vx = ny / np.sqrt(nx**2 + ny**2) +vy = -nx / np.sqrt(nx**2 + ny**2) +vz = 0 +u_perp = ux * vx + uy * vy + uz * vz +compare_gaussian(u_perp, w, u_th=0.01, label="u_perp") + +# Pick the other perpendicular direction, and check the distribution +# The third direction is obtained by the cross product (n x v) +wx = ny * vz - nz * vy +wy = nz * vx - nx * vz +wz = nx * vy - ny * vx +u_perp2 = ux * wx + uy * wy + uz * wz +compare_gaussian(u_perp2, w, u_th=0.01, label="u_perp") + +plt.tight_layout() +plt.savefig("Distribution.png") + +# Verify checksum +test_name = os.path.split(os.getcwd())[1] +checksumAPI.evaluate_checksum(test_name, fn) diff --git a/Examples/Tests/flux_injection/inputs_base_from_eb b/Examples/Tests/flux_injection/inputs_base_from_eb new file mode 100644 index 00000000000..3e32d8799b6 --- /dev/null +++ b/Examples/Tests/flux_injection/inputs_base_from_eb @@ -0,0 +1,42 @@ +# Maximum number of time steps +max_step = 10 + +# The lo and hi ends of grids are multipliers of blocking factor +amr.blocking_factor = 8 + +# Maximum allowable size of each subdomain in the problem domain; +# this is used to decompose the domain for parallel calculations. +amr.max_grid_size = 8 + +# Maximum level in hierarchy (for now must be 0, i.e., one level in total) +amr.max_level = 0 + +# Deactivate Maxwell solver +algo.maxwell_solver = none +warpx.const_dt = 1e-9 + +# Embedded boundary +warpx.eb_implicit_function = "-(x**2+y**2+z**2-2**2)" + +# particles +particles.species_names = electron +algo.particle_shape = 3 + +electron.charge = -q_e +electron.mass = m_e +electron.injection_style = NFluxPerCell +electron.inject_from_embedded_boundary = 1 +electron.num_particles_per_cell = 100 +electron.flux_profile = parse_flux_function +electron.flux_function(x,y,z,t) = "1." +electron.momentum_distribution_type = gaussianflux +electron.ux_th = 0.01 +electron.uy_th = 0.01 +electron.uz_th = 0.1 +electron.uz_m = 0.07 + +# Diagnostics +diagnostics.diags_names = diag1 +diag1.intervals = 10 +diag1.diag_type = Full +diag1.fields_to_plot = none diff --git a/Examples/Tests/flux_injection/inputs_test_2d_flux_injection_from_eb b/Examples/Tests/flux_injection/inputs_test_2d_flux_injection_from_eb new file mode 100644 index 00000000000..f2e6f177887 --- /dev/null +++ b/Examples/Tests/flux_injection/inputs_test_2d_flux_injection_from_eb @@ -0,0 +1,13 @@ +FILE = inputs_base_from_eb + +# number of grid points +amr.n_cell = 16 16 + +# Geometry +geometry.dims = 2 +geometry.prob_lo = -4 -4 +geometry.prob_hi = 4 4 + +# Boundary condition +boundary.field_lo = periodic periodic +boundary.field_hi = periodic periodic diff --git a/Examples/Tests/flux_injection/inputs_test_3d_flux_injection_from_eb b/Examples/Tests/flux_injection/inputs_test_3d_flux_injection_from_eb new file mode 100644 index 00000000000..81ddc039977 --- /dev/null +++ b/Examples/Tests/flux_injection/inputs_test_3d_flux_injection_from_eb @@ -0,0 +1,13 @@ +FILE = inputs_base_from_eb + +# number of grid points +amr.n_cell = 16 16 16 + +# Geometry +geometry.dims = 3 +geometry.prob_lo = -4 -4 -4 +geometry.prob_hi = 4 4 4 + +# Boundary condition +boundary.field_lo = periodic periodic periodic +boundary.field_hi = periodic periodic periodic diff --git a/Examples/Tests/flux_injection/inputs_test_rz_flux_injection_from_eb b/Examples/Tests/flux_injection/inputs_test_rz_flux_injection_from_eb new file mode 100644 index 00000000000..4c970257f57 --- /dev/null +++ b/Examples/Tests/flux_injection/inputs_test_rz_flux_injection_from_eb @@ -0,0 +1,15 @@ +FILE = inputs_base_from_eb + +# number of grid points +amr.n_cell = 8 16 + +# Geometry +geometry.dims = RZ +geometry.prob_lo = 0 -4 +geometry.prob_hi = 4 4 + +# Boundary condition +boundary.field_lo = none periodic +boundary.field_hi = pec periodic + +electron.num_particles_per_cell = 300 diff --git a/Regression/Checksum/benchmarks_json/test_2d_flux_injection_from_eb.json b/Regression/Checksum/benchmarks_json/test_2d_flux_injection_from_eb.json new file mode 100644 index 00000000000..dd489f16e05 --- /dev/null +++ b/Regression/Checksum/benchmarks_json/test_2d_flux_injection_from_eb.json @@ -0,0 +1,11 @@ +{ + "lev=0": {}, + "electron": { + "particle_momentum_x": 6.990772711451971e-19, + "particle_momentum_y": 5.4131306169803364e-20, + "particle_momentum_z": 6.997294931789925e-19, + "particle_position_x": 35518.95120597846, + "particle_position_y": 35517.855675902414, + "particle_weight": 1.25355e-07 + } +} diff --git a/Regression/Checksum/benchmarks_json/test_3d_flux_injection_from_eb.json b/Regression/Checksum/benchmarks_json/test_3d_flux_injection_from_eb.json new file mode 100644 index 00000000000..e947a8af07b --- /dev/null +++ b/Regression/Checksum/benchmarks_json/test_3d_flux_injection_from_eb.json @@ -0,0 +1,12 @@ +{ + "lev=0": {}, + "electron": { + "particle_momentum_x": 4.371688233196277e-18, + "particle_momentum_y": 4.368885079657374e-18, + "particle_momentum_z": 4.367429424105371e-18, + "particle_position_x": 219746.94401890738, + "particle_position_y": 219690.7015248918, + "particle_position_z": 219689.45580938633, + "particle_weight": 4.954974999999999e-07 + } +} \ No newline at end of file diff --git a/Regression/Checksum/benchmarks_json/test_rz_flux_injection_from_eb.json b/Regression/Checksum/benchmarks_json/test_rz_flux_injection_from_eb.json new file mode 100644 index 00000000000..23884de9725 --- /dev/null +++ b/Regression/Checksum/benchmarks_json/test_rz_flux_injection_from_eb.json @@ -0,0 +1,12 @@ +{ + "lev=0": {}, + "electron": { + "particle_momentum_x": 6.734984863106283e-19, + "particle_momentum_y": 6.786279785869023e-19, + "particle_momentum_z": 1.0527983828124758e-18, + "particle_position_x": 53309.270966506396, + "particle_position_y": 53302.3776094842, + "particle_theta": 58707.74469425615, + "particle_weight": 4.991396867417661e-07 + } +} \ No newline at end of file diff --git a/Source/Initialization/PlasmaInjector.H b/Source/Initialization/PlasmaInjector.H index b9fe2323290..f14720d271c 100644 --- a/Source/Initialization/PlasmaInjector.H +++ b/Source/Initialization/PlasmaInjector.H @@ -131,6 +131,8 @@ public: int flux_normal_axis; int flux_direction; // -1 for left, +1 for right + bool m_inject_from_eb = false; // whether to inject from the embedded boundary + bool radially_weighted = true; std::string str_flux_function; diff --git a/Source/Initialization/PlasmaInjector.cpp b/Source/Initialization/PlasmaInjector.cpp index 3d846375a99..76bb7a5be42 100644 --- a/Source/Initialization/PlasmaInjector.cpp +++ b/Source/Initialization/PlasmaInjector.cpp @@ -9,6 +9,7 @@ */ #include "PlasmaInjector.H" +#include "EmbeddedBoundary/Enabled.H" #include "Initialization/GetTemperature.H" #include "Initialization/GetVelocity.H" #include "Initialization/InjectorDensity.H" @@ -303,50 +304,65 @@ void PlasmaInjector::setupNFluxPerCell (amrex::ParmParse const& pp_species) "(Please visit PR#765 for more information.)"); } #endif - utils::parser::getWithParser(pp_species, source_name, "surface_flux_pos", surface_flux_pos); - utils::parser::queryWithParser(pp_species, source_name, "flux_tmin", flux_tmin); - utils::parser::queryWithParser(pp_species, source_name, "flux_tmax", flux_tmax); - std::string flux_normal_axis_string; - utils::parser::get(pp_species, source_name, "flux_normal_axis", flux_normal_axis_string); - flux_normal_axis = -1; + + // Check whether injection from the embedded boundary is requested + utils::parser::queryWithParser(pp_species, source_name, "inject_from_embedded_boundary", m_inject_from_eb); + if (m_inject_from_eb) { + AMREX_ALWAYS_ASSERT_WITH_MESSAGE( EB::enabled(), + "Error: Embedded boundary injection is only available when " + "embedded boundaries are enabled."); + flux_normal_axis = 2; // Interpret z as the normal direction to the EB + flux_direction = 1; + } else { + // Injection is through a plane in this case. + // Parse the parameters of the plane (position, normal direction, etc.) + + utils::parser::getWithParser(pp_species, source_name, "surface_flux_pos", surface_flux_pos); + utils::parser::queryWithParser(pp_species, source_name, "flux_tmin", flux_tmin); + utils::parser::queryWithParser(pp_species, source_name, "flux_tmax", flux_tmax); + std::string flux_normal_axis_string; + utils::parser::get(pp_species, source_name, "flux_normal_axis", flux_normal_axis_string); + flux_normal_axis = -1; #ifdef WARPX_DIM_RZ - if (flux_normal_axis_string == "r" || flux_normal_axis_string == "R") { - flux_normal_axis = 0; - } - if (flux_normal_axis_string == "t" || flux_normal_axis_string == "T") { - flux_normal_axis = 1; - } + if (flux_normal_axis_string == "r" || flux_normal_axis_string == "R") { + flux_normal_axis = 0; + } + if (flux_normal_axis_string == "t" || flux_normal_axis_string == "T") { + flux_normal_axis = 1; + } #else # ifndef WARPX_DIM_1D_Z - if (flux_normal_axis_string == "x" || flux_normal_axis_string == "X") { - flux_normal_axis = 0; - } + if (flux_normal_axis_string == "x" || flux_normal_axis_string == "X") { + flux_normal_axis = 0; + } # endif #endif #ifdef WARPX_DIM_3D - if (flux_normal_axis_string == "y" || flux_normal_axis_string == "Y") { - flux_normal_axis = 1; - } + if (flux_normal_axis_string == "y" || flux_normal_axis_string == "Y") { + flux_normal_axis = 1; + } #endif - if (flux_normal_axis_string == "z" || flux_normal_axis_string == "Z") { - flux_normal_axis = 2; - } + if (flux_normal_axis_string == "z" || flux_normal_axis_string == "Z") { + flux_normal_axis = 2; + } #ifdef WARPX_DIM_3D - const std::string flux_normal_axis_help = "'x', 'y', or 'z'."; + const std::string flux_normal_axis_help = "'x', 'y', or 'z'."; #else # ifdef WARPX_DIM_RZ - const std::string flux_normal_axis_help = "'r' or 'z'."; + const std::string flux_normal_axis_help = "'r' or 'z'."; # elif WARPX_DIM_XZ - const std::string flux_normal_axis_help = "'x' or 'z'."; + const std::string flux_normal_axis_help = "'x' or 'z'."; # else - const std::string flux_normal_axis_help = "'z'."; + const std::string flux_normal_axis_help = "'z'."; # endif -#endif - WARPX_ALWAYS_ASSERT_WITH_MESSAGE(flux_normal_axis >= 0, - "Error: Invalid value for flux_normal_axis. It must be " + flux_normal_axis_help); - utils::parser::getWithParser(pp_species, source_name, "flux_direction", flux_direction); - WARPX_ALWAYS_ASSERT_WITH_MESSAGE(flux_direction == +1 || flux_direction == -1, - "Error: flux_direction must be -1 or +1."); + #endif + WARPX_ALWAYS_ASSERT_WITH_MESSAGE(flux_normal_axis >= 0, + "Error: Invalid value for flux_normal_axis. It must be " + flux_normal_axis_help); + utils::parser::getWithParser(pp_species, source_name, "flux_direction", flux_direction); + WARPX_ALWAYS_ASSERT_WITH_MESSAGE(flux_direction == +1 || flux_direction == -1, + "Error: flux_direction must be -1 or +1."); + } + // Construct InjectorPosition with InjectorPositionRandom. h_flux_pos = std::make_unique( (InjectorPositionRandomPlane*)nullptr, diff --git a/Source/Particles/AddPlasmaUtilities.H b/Source/Particles/AddPlasmaUtilities.H index 8f0489e3921..bb05d7be3c8 100644 --- a/Source/Particles/AddPlasmaUtilities.H +++ b/Source/Particles/AddPlasmaUtilities.H @@ -22,6 +22,29 @@ #include #include +struct PDim3 { + amrex::ParticleReal x, y, z; + + AMREX_GPU_HOST_DEVICE + PDim3(const amrex::XDim3& a): + x{static_cast(a.x)}, + y{static_cast(a.y)}, + z{static_cast(a.z)} + {} + + AMREX_GPU_HOST_DEVICE + ~PDim3() = default; + + AMREX_GPU_HOST_DEVICE + PDim3(PDim3 const &) = default; + AMREX_GPU_HOST_DEVICE + PDim3& operator=(PDim3 const &) = default; + AMREX_GPU_HOST_DEVICE + PDim3(PDim3&&) = default; + AMREX_GPU_HOST_DEVICE + PDim3& operator=(PDim3&&) = default; +}; + /* Finds the overlap region between the given tile_realbox and part_realbox, returning true if an overlap exists and false if otherwise. This also sets the parameters overlap_realbox, @@ -71,12 +94,124 @@ int compute_area_weights (const amrex::IntVect& iv, const int normal_axis) { return r; } + +#ifdef AMREX_USE_EB +/* + * \brief This computes the scale_fac (used for setting the particle weights) on a on area basis + * (used for flux injection from the embedded boundary). + * + * \param[in] dx: cell size in each direction + * \param[in] num_ppc_real: number of particles per cell + * \param[in] eb_bnd_normal_arr: array containing the normal to the embedded boundary + * \param[in] i, j, k: indices of the cell + * + * \return scale_fac: the scaling factor to be applied to the weight of the particles + */ +AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE +amrex::Real compute_scale_fac_area_eb ( + const amrex::GpuArray& dx, + const amrex::Real num_ppc_real, + amrex::Array4 const& eb_bnd_normal_arr, + int i, int j, int k ) { + using namespace amrex::literals; + // Scale particle weight by the area of the emitting surface, within one cell + // By definition, eb_bnd_area_arr is normalized (unitless). + // Here we undo the normalization (i.e. multiply by the surface used for normalization in amrex: + // see https://amrex-codes.github.io/amrex/docs_html/EB.html#embedded-boundary-data-structures) +#if defined(WARPX_DIM_3D) + const amrex::Real nx = eb_bnd_normal_arr(i,j,k,0); + const amrex::Real ny = eb_bnd_normal_arr(i,j,k,1); + const amrex::Real nz = eb_bnd_normal_arr(i,j,k,2); + amrex::Real scale_fac = std::sqrt(amrex::Math::powi<2>(nx*dx[1]*dx[2]) + + amrex::Math::powi<2>(ny*dx[0]*dx[2]) + + amrex::Math::powi<2>(nz*dx[0]*dx[1])); + +#elif defined(WARPX_DIM_RZ) || defined(WARPX_DIM_XZ) + const amrex::Real nx = eb_bnd_normal_arr(i,j,k,0); + const amrex::Real nz = eb_bnd_normal_arr(i,j,k,1); + amrex::Real scale_fac = std::sqrt(amrex::Math::powi<2>(nx*dx[1]) + + amrex::Math::powi<2>(nz*dx[0])); +#else + amrex::ignore_unused(dx, eb_bnd_normal_arr, i, j, k); + amrex::Real scale_fac = 0.0_rt; +#endif + // Do not multiply by eb_bnd_area_arr(i,j,k) here because this + // already taken into account by emitting a number of macroparticles + // that is proportional to the area of eb_bnd_area_arr(i,j,k). + scale_fac /= num_ppc_real; + return scale_fac; +} + +/* \brief Rotate the momentum of the particle so that the z direction + * transforms to the normal of the embedded boundary. + * + * More specifically, before calling this function, `pu.z` has the + * momentum distribution that is meant for the direction normal to the + * embedded boundary, and `pu.x`/`pu.y` have the momentum distribution that + * is meant for the tangentional direction. After calling this function, + * `pu.x`, `pu.y`, `pu.z` will have the correct momentum distribution, + * consistent with the local normal to the embedded boundary. + * + * \param[inout] pu momentum of the particle + * \param[in] eb_bnd_normal_arr: array containing the normal to the embedded boundary + * \param[in] i, j, k: indices of the cell + * */ +AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE +void rotate_momentum_eb ( + PDim3 & pu, + amrex::Array4 const& eb_bnd_normal_arr, + int i, int j, int k ) +{ + using namespace amrex::literals; + +#if defined(WARPX_DIM_3D) + // The minus sign below takes into account the fact that eb_bnd_normal_arr + // points towards the covered region, while particles are to be emitted + // *away* from the covered region. + amrex::Real const nx = -eb_bnd_normal_arr(i,j,k,0); + amrex::Real const ny = -eb_bnd_normal_arr(i,j,k,1); + amrex::Real const nz = -eb_bnd_normal_arr(i,j,k,2); + + // Rotate the momentum in theta and phi + amrex::Real const cos_theta = nz; + amrex::Real const sin_theta = std::sqrt(1._rt-nz*nz); + amrex::Real const nperp = std::sqrt(nx*nx + ny*ny); + amrex::Real cos_phi = 1; + amrex::Real sin_phi = 0; + if ( nperp > 0.0 ) { + cos_phi = nx/nperp; + sin_phi = ny/nperp; + } + // Apply rotation matrix + amrex::Real const ux = pu.x*cos_theta*cos_phi - pu.y*sin_phi + pu.z*sin_theta*cos_phi; + amrex::Real const uy = pu.x*cos_theta*sin_phi + pu.y*cos_phi + pu.z*sin_theta*sin_phi; + amrex::Real const uz = -pu.x*sin_theta + pu.z*cos_theta; + pu.x = ux; + pu.y = uy; + pu.z = uz; + +#elif defined(WARPX_DIM_RZ) || defined(WARPX_DIM_XZ) + // The minus sign below takes into account the fact that eb_bnd_normal_arr + // points towards the covered region, while particles are to be emitted + // *away* from the covered region. + amrex::Real const sin_theta = -eb_bnd_normal_arr(i,j,k,0); + amrex::Real const cos_theta = -eb_bnd_normal_arr(i,j,k,1); + amrex::Real const uz = pu.z*cos_theta - pu.x*sin_theta; + amrex::Real const ux = pu.x*cos_theta + pu.z*sin_theta; + pu.x = ux; + pu.z = uz; +#else + amrex::ignore_unused(pu, eb_bnd_normal_arr, i, j, k); +#endif +} +#endif //AMREX_USE_EB + /* This computes the scale_fac (used for setting the particle weights) on a on area basis - (used for flux injection). + (used for flux injection from a plane). */ AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE -amrex::Real compute_scale_fac_area (const amrex::GpuArray& dx, +amrex::Real compute_scale_fac_area_plane (const amrex::GpuArray& dx, const amrex::Real num_ppc_real, const int flux_normal_axis) { using namespace amrex::literals; amrex::Real scale_fac = AMREX_D_TERM(dx[0],*dx[1],*dx[2])/num_ppc_real; diff --git a/Source/Particles/PhysicalParticleContainer.cpp b/Source/Particles/PhysicalParticleContainer.cpp index 26f9fee38d3..7c70c9a35c4 100644 --- a/Source/Particles/PhysicalParticleContainer.cpp +++ b/Source/Particles/PhysicalParticleContainer.cpp @@ -147,29 +147,6 @@ namespace return z0; } - struct PDim3 { - ParticleReal x, y, z; - - AMREX_GPU_HOST_DEVICE - PDim3(const amrex::XDim3& a): - x{static_cast(a.x)}, - y{static_cast(a.y)}, - z{static_cast(a.z)} - {} - - AMREX_GPU_HOST_DEVICE - ~PDim3() = default; - - AMREX_GPU_HOST_DEVICE - PDim3(PDim3 const &) = default; - AMREX_GPU_HOST_DEVICE - PDim3& operator=(PDim3 const &) = default; - AMREX_GPU_HOST_DEVICE - PDim3(PDim3&&) = default; - AMREX_GPU_HOST_DEVICE - PDim3& operator=(PDim3&&) = default; - }; - AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE XDim3 getCellCoords (const GpuArray& lo_corner, const GpuArray& dx, @@ -1371,6 +1348,22 @@ PhysicalParticleContainer::AddPlasmaFlux (PlasmaInjector const& plasma_injector, const auto dx = geom.CellSizeArray(); const auto problo = geom.ProbLoArray(); +#ifdef AMREX_USE_EB + bool const inject_from_eb = plasma_injector.m_inject_from_eb; // whether to inject from EB or from a plane + // Extract data structures for embedded boundaries + amrex::FabArray const* eb_flag = nullptr; + amrex::MultiCutFab const* eb_bnd_area = nullptr; + amrex::MultiCutFab const* eb_bnd_normal = nullptr; + amrex::MultiCutFab const* eb_bnd_cent = nullptr; + if (inject_from_eb) { + amrex::EBFArrayBoxFactory const& eb_box_factory = WarpX::GetInstance().fieldEBFactory(0); + eb_flag = &eb_box_factory.getMultiEBCellFlagFab(); + eb_bnd_area = &eb_box_factory.getBndryArea(); + eb_bnd_normal = &eb_box_factory.getBndryNormal(); + eb_bnd_cent = &eb_box_factory.getBndryCent(); + } +#endif + amrex::LayoutData* cost = WarpX::getCosts(0); // Create temporary particle container to which particles will be added; @@ -1428,9 +1421,20 @@ PhysicalParticleContainer::AddPlasmaFlux (PlasmaInjector const& plasma_injector, RealBox overlap_realbox; Box overlap_box; IntVect shifted; - const bool no_overlap = find_overlap_flux(tile_realbox, part_realbox, dx, problo, plasma_injector, overlap_realbox, overlap_box, shifted); - if (no_overlap) { - continue; // Go to the next tile +#ifdef AMREX_USE_EB + if (inject_from_eb) { + // Injection from EB + const amrex::FabType fab_type = (*eb_flag)[mfi].getType(tile_box); + if (fab_type == amrex::FabType::regular) { continue; } // Go to the next tile + if (fab_type == amrex::FabType::covered) { continue; } // Go to the next tile + overlap_box = tile_box; + overlap_realbox = part_realbox; + } else +#endif + { + // Injection from a plane + const bool no_overlap = find_overlap_flux(tile_realbox, part_realbox, dx, problo, plasma_injector, overlap_realbox, overlap_box, shifted); + if (no_overlap) { continue; } // Go to the next tile } const int grid_id = mfi.index(); @@ -1450,24 +1454,57 @@ PhysicalParticleContainer::AddPlasmaFlux (PlasmaInjector const& plasma_injector, if (refine_injection) { fine_overlap_box = overlap_box & amrex::shift(fine_injection_box, -shifted); } + +#ifdef AMREX_USE_EB + // Extract data structures for embedded boundaries + amrex::Array4::value_type> eb_flag_arr; + amrex::Array4 eb_bnd_area_arr; + amrex::Array4 eb_bnd_normal_arr; + amrex::Array4 eb_bnd_cent_arr; + if (inject_from_eb) { + eb_flag_arr = eb_flag->array(mfi); + eb_bnd_area_arr = eb_bnd_area->array(mfi); + eb_bnd_normal_arr = eb_bnd_normal->array(mfi); + eb_bnd_cent_arr = eb_bnd_cent->array(mfi); + } +#endif + amrex::ParallelForRNG(overlap_box, [=] AMREX_GPU_DEVICE (int i, int j, int k, amrex::RandomEngine const& engine) noexcept { const IntVect iv(AMREX_D_DECL(i, j, k)); amrex::ignore_unused(j,k); - auto lo = getCellCoords(overlap_corner, dx, {0._rt, 0._rt, 0._rt}, iv); - auto hi = getCellCoords(overlap_corner, dx, {1._rt, 1._rt, 1._rt}, iv); - - if (flux_pos->overlapsWith(lo, hi)) + // Determine the number of macroparticles to inject in this cell (num_ppc_int) +#ifdef AMREX_USE_EB + amrex::Real num_ppc_real_in_this_cell = num_ppc_real; // user input: number of macroparticles per cell + if (inject_from_eb) { + // Injection from EB + // Skip cells that are not partially covered by the EB + if (eb_flag_arr(i,j,k).isRegular() || eb_flag_arr(i,j,k).isCovered()) { return; } + // Scale by the (normalized) area of the EB surface in this cell + num_ppc_real_in_this_cell *= eb_bnd_area_arr(i,j,k); + } else +#else + amrex::Real const num_ppc_real_in_this_cell = num_ppc_real; // user input: number of macroparticles per cell +#endif { - auto index = overlap_box.index(iv); - int r = 1; - if (fine_overlap_box.ok() && fine_overlap_box.contains(iv)) { - r = compute_area_weights(rrfac, flux_normal_axis); - } - const int num_ppc_int = static_cast(num_ppc_real*r + amrex::Random(engine)); - pcounts[index] = num_ppc_int; + // Injection from a plane + auto lo = getCellCoords(overlap_corner, dx, {0._rt, 0._rt, 0._rt}, iv); + auto hi = getCellCoords(overlap_corner, dx, {1._rt, 1._rt, 1._rt}, iv); + // Skip cells that do not overlap with the plane + if (!flux_pos->overlapsWith(lo, hi)) { return; } } + + auto index = overlap_box.index(iv); + // Take into account refined injection region + int r = 1; + if (fine_overlap_box.ok() && fine_overlap_box.contains(iv)) { + r = compute_area_weights(rrfac, flux_normal_axis); + } + const int num_ppc_int = static_cast(num_ppc_real_in_this_cell*r + amrex::Random(engine)); + pcounts[index] = num_ppc_int; + + amrex::ignore_unused(j,k); }); // Max number of new particles. All of them are created, @@ -1537,7 +1574,15 @@ PhysicalParticleContainer::AddPlasmaFlux (PlasmaInjector const& plasma_injector, amrex::ignore_unused(j,k); const auto index = overlap_box.index(iv); - Real scale_fac = compute_scale_fac_area(dx, num_ppc_real, flux_normal_axis); + Real scale_fac; +#ifdef AMREX_USE_EB + if (inject_from_eb) { + scale_fac = compute_scale_fac_area_eb(dx, num_ppc_real, eb_bnd_normal_arr, i, j, k ); + } else +#endif + { + scale_fac = compute_scale_fac_area_plane(dx, num_ppc_real, flux_normal_axis); + } if (fine_overlap_box.ok() && fine_overlap_box.contains(iv)) { scale_fac /= compute_area_weights(rrfac, flux_normal_axis); @@ -1548,13 +1593,32 @@ PhysicalParticleContainer::AddPlasmaFlux (PlasmaInjector const& plasma_injector, const long ip = poffset[index] + i_part; pa_idcpu[ip] = amrex::SetParticleIDandCPU(pid+ip, cpuid); - // This assumes the flux_pos is of type InjectorPositionRandomPlane - const XDim3 r = (fine_overlap_box.ok() && fine_overlap_box.contains(iv)) ? - // In the refined injection region: use refinement ratio `rrfac` - flux_pos->getPositionUnitBox(i_part, rrfac, engine) : - // Otherwise: use 1 as the refinement ratio - flux_pos->getPositionUnitBox(i_part, amrex::IntVect::TheUnitVector(), engine); - auto pos = getCellCoords(overlap_corner, dx, r, iv); + // Determine the position of the particle within the cell + XDim3 pos; + XDim3 r; +#ifdef AMREX_USE_EB + if (inject_from_eb) { +#if defined(WARPX_DIM_3D) + pos.x = overlap_corner[0] + (iv[0] + 0.5_rt + eb_bnd_cent_arr(i,j,k,0))*dx[0]; + pos.y = overlap_corner[1] + (iv[1] + 0.5_rt + eb_bnd_cent_arr(i,j,k,1))*dx[1]; + pos.z = overlap_corner[2] + (iv[2] + 0.5_rt + eb_bnd_cent_arr(i,j,k,2))*dx[2]; +#elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ) + pos.x = overlap_corner[0] + (iv[0] + 0.5_rt + eb_bnd_cent_arr(i,j,k,0))*dx[0]; + pos.y = 0.0_rt; + pos.z = overlap_corner[1] + (iv[1] + 0.5_rt + eb_bnd_cent_arr(i,j,k,1))*dx[1]; +#endif + } else +#endif + { + // Injection from a plane + // This assumes the flux_pos is of type InjectorPositionRandomPlane + r = (fine_overlap_box.ok() && fine_overlap_box.contains(iv)) ? + // In the refined injection region: use refinement ratio `rrfac` + flux_pos->getPositionUnitBox(i_part, rrfac, engine) : + // Otherwise: use 1 as the refinement ratio + flux_pos->getPositionUnitBox(i_part, amrex::IntVect::TheUnitVector(), engine); + pos = getCellCoords(overlap_corner, dx, r, iv); + } auto ppos = PDim3(pos); // inj_mom would typically be InjectorMomentumGaussianFlux @@ -1595,6 +1659,15 @@ PhysicalParticleContainer::AddPlasmaFlux (PlasmaInjector const& plasma_injector, continue; } +#ifdef AMREX_USE_EB + if (inject_from_eb) { + // Injection from EB: rotate momentum according to the normal of the EB surface + // (The above code initialized the momentum by assuming that z is the direction + // normal to the EB surface. Thus we need to rotate from z to the normal.) + rotate_momentum_eb(pu, eb_bnd_normal_arr, i, j , k); + } +#endif + #ifdef WARPX_DIM_RZ // Conversion from cylindrical to Cartesian coordinates // Replace the x and y, setting an angle theta. @@ -1610,7 +1683,11 @@ PhysicalParticleContainer::AddPlasmaFlux (PlasmaInjector const& plasma_injector, const amrex::Real radial_position = ppos.x; ppos.x = radial_position*cos_theta; ppos.y = radial_position*sin_theta; - if (loc_flux_normal_axis != 2) { + if ((loc_flux_normal_axis != 2) +#ifdef AMREX_USE_EB + || (inject_from_eb) +#endif + ) { // Rotate the momentum // This because, when the flux direction is e.g. "r" // the `inj_mom` objects generates a v*Gaussian distribution