WarpX
Loading...
Searching...
No Matches
ablastr::fields Namespace Reference

Namespaces

namespace  details
 

Classes

class  Direction
 
struct  MultiFabOwner
 
struct  MultiFabRegister
 

Typedefs

using ScalarField = amrex::MultiFab*
 
using ConstScalarField = amrex::MultiFab const *
 
using VectorField = std::array<amrex::MultiFab *, 3>
 
using ConstVectorField = std::array<amrex::MultiFab const *, 3>
 
using MultiLevelScalarField = amrex::Vector<ScalarField>
 
using ConstMultiLevelScalarField = amrex::Vector<ConstScalarField>
 
using MultiLevelVectorField = amrex::Vector<VectorField>
 
using ConstMultiLevelVectorField = amrex::Vector<ConstVectorField>
 

Functions

template<typename T_PostPhiCalculationFunctor = std::nullopt_t, typename T_BoundaryHandler = std::nullopt_t, typename T_FArrayBoxFactory = void>
void computeEffectivePotentialPhi (ablastr::fields::MultiLevelScalarField const &rho, ablastr::fields::MultiLevelScalarField const &phi, amrex::MultiFab const &sigma, amrex::Real relative_tolerance, amrex::Real absolute_tolerance, int max_iters, int verbosity, amrex::Vector< amrex::Geometry > const &geom, amrex::Vector< amrex::DistributionMapping > const &dmap, amrex::Vector< amrex::BoxArray > const &grids, utils::enums::GridType grid_type, bool is_solver_igf_on_lev0, bool eb_enabled=false, bool do_single_precision_comms=false, std::optional< amrex::Vector< amrex::IntVect > > rel_ref_ratio=std::nullopt, T_PostPhiCalculationFunctor post_phi_calculation=std::nullopt, T_BoundaryHandler const boundary_handler=std::nullopt, std::optional< amrex::Real const > current_time=std::nullopt, std::optional< amrex::Vector< T_FArrayBoxFactory const * > > eb_farray_box_factory=std::nullopt)
 
void computePhiIGF (amrex::MultiFab const &rho, amrex::MultiFab &phi, std::array< amrex::Real, 3 > const &cell_size, bool is_igf_2d_slices)
 Compute the electrostatic potential using the Integrated Green Function method as in http://dx.doi.org/10.1103/PhysRevSTAB.9.044204.
 
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real IntegratedPotential3D (amrex::Real x, amrex::Real y, amrex::Real z)
 Implements equation 2 in https://doi.org/10.1103/PhysRevSTAB.10.129901 with some modification to symmetrize the function.
 
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real IntegratedPotential2D (amrex::Real x, amrex::Real y)
 Implements equation 58 in https://doi.org/10.1016/j.jcp.2004.01.008.
 
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real SumOfIntegratedPotential3D (amrex::Real x, amrex::Real y, amrex::Real z, amrex::Real dx, amrex::Real dy, amrex::Real dz)
 add
 
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real SumOfIntegratedPotential2D (amrex::Real x, amrex::Real y, amrex::Real dx, amrex::Real dy)
 add
 
VectorField a2m (std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &old_vectorfield)
 
amrex::Real getMaxNormRho (ablastr::fields::ConstMultiLevelScalarField const &rho, int finest_level, amrex::Real &absolute_tolerance)
 
void interpolatePhiBetweenLevels (amrex::MultiFab const *phi_lev, amrex::MultiFab *phi_lev_plus_one, amrex::Geometry const &geom_lev, bool do_single_precision_comms, const amrex::IntVect &refratio, const int ncomp, const int ng)
 
template<typename T_PostPhiCalculationFunctor = std::nullopt_t, typename T_BoundaryHandler = std::nullopt_t, typename T_FArrayBoxFactory = void>
void computePhi (ablastr::fields::MultiLevelScalarField const &rho, ablastr::fields::MultiLevelScalarField const &phi, std::array< amrex::Real, 3 > const beta, amrex::Real relative_tolerance, amrex::Real absolute_tolerance, int max_iters, int verbosity, amrex::Vector< amrex::Geometry > const &geom, amrex::Vector< amrex::DistributionMapping > const &dmap, amrex::Vector< amrex::BoxArray > const &grids, utils::enums::GridType grid_type, bool is_solver_igf_on_lev0, bool const is_igf_2d, bool eb_enabled=false, bool do_single_precision_comms=false, std::optional< amrex::Vector< amrex::IntVect > > rel_ref_ratio=std::nullopt, T_PostPhiCalculationFunctor post_phi_calculation=std::nullopt, T_BoundaryHandler const boundary_handler=std::nullopt, std::optional< amrex::Real const > current_time=std::nullopt, std::optional< amrex::Vector< T_FArrayBoxFactory const * > > eb_farray_box_factory=std::nullopt)
 
template<typename T_BoundaryHandler, typename T_PostACalculationFunctor = std::nullopt_t, typename T_FArrayBoxFactory = void>
void computeVectorPotential (amrex::Vector< amrex::Array< amrex::MultiFab *, 3 > > const &curr, amrex::Vector< amrex::Array< amrex::MultiFab *, 3 > > &A, amrex::Real relative_tolerance, amrex::Real absolute_tolerance, int max_iters, int verbosity, amrex::Vector< amrex::Geometry > const &geom, amrex::Vector< amrex::DistributionMapping > const &dmap, amrex::Vector< amrex::BoxArray > const &grids, T_BoundaryHandler const boundary_handler, bool eb_enabled=false, bool do_single_precision_comms=false, std::optional< amrex::Vector< amrex::IntVect > > rel_ref_ratio=std::nullopt, T_PostACalculationFunctor post_A_calculation=std::nullopt, std::optional< amrex::Real const > current_time=std::nullopt, std::optional< amrex::Vector< T_FArrayBoxFactory const * > > eb_farray_box_factory=std::nullopt)
 

Typedef Documentation

◆ ConstMultiLevelScalarField

A read-only multi-level scalar field

◆ ConstMultiLevelVectorField

A read-only multi-level vector field

◆ ConstScalarField

A read-only scalar field (a MultiFab)

Note: might still have components, e.g., for copies at different times.

◆ ConstVectorField

using ablastr::fields::ConstVectorField = std::array<amrex::MultiFab const *, 3>

A read-only vector field of three MultiFab

◆ MultiLevelScalarField

◆ MultiLevelVectorField

◆ ScalarField

A scalar field (a MultiFab)

Note: might still have components, e.g., for copies at different times.

◆ VectorField

A vector field of three MultiFab

Function Documentation

◆ a2m()

VectorField ablastr::fields::a2m ( std::array< std::unique_ptr< amrex::MultiFab >, 3 > const & old_vectorfield)

Little temporary helper function to pass temporary MultiFabs as VectorField.

Returns
pointers to externally managed vector field components (3 MultiFab)

◆ computeEffectivePotentialPhi()

template<typename T_PostPhiCalculationFunctor = std::nullopt_t, typename T_BoundaryHandler = std::nullopt_t, typename T_FArrayBoxFactory = void>
void ablastr::fields::computeEffectivePotentialPhi ( ablastr::fields::MultiLevelScalarField const & rho,
ablastr::fields::MultiLevelScalarField const & phi,
amrex::MultiFab const & sigma,
amrex::Real relative_tolerance,
amrex::Real absolute_tolerance,
int max_iters,
int verbosity,
amrex::Vector< amrex::Geometry > const & geom,
amrex::Vector< amrex::DistributionMapping > const & dmap,
amrex::Vector< amrex::BoxArray > const & grids,
utils::enums::GridType grid_type,
bool is_solver_igf_on_lev0,
bool eb_enabled = false,
bool do_single_precision_comms = false,
std::optional< amrex::Vector< amrex::IntVect > > rel_ref_ratio = std::nullopt,
T_PostPhiCalculationFunctor post_phi_calculation = std::nullopt,
T_BoundaryHandler const boundary_handler = std::nullopt,
std::optional< amrex::Real const > current_time = std::nullopt,
std::optional< amrex::Vector< T_FArrayBoxFactory const * > > eb_farray_box_factory = std::nullopt )

Compute the potential phi by solving the Poisson equation with a modifed dielectric function

Uses rho as a source. This uses the AMReX solver.

More specifically, this solves the equation

\[  \nabla \cdot \sigma \nabla \phi = - \rho/\epsilon_0
\]

Template Parameters
T_PostPhiCalculationFunctora calculation per level directly after phi was calculated
T_BoundaryHandlerhandler for boundary conditions, for example
See also
ElectrostaticSolver::PoissonBoundaryHandler
Template Parameters
T_FArrayBoxFactoryusually nothing or an amrex::EBFArrayBoxFactory (EB ONLY)
Parameters
[in]rhoThe charge density a given species
[out]phiThe potential to be computed by this function
[in]sigmaThe matrix representing the mass operator used to lower the local plasma frequency
[in]relative_toleranceThe relative convergence threshold for the MLMG solver
[in]absolute_toleranceThe absolute convergence threshold for the MLMG solver
[in]max_itersThe maximum number of iterations allowed for the MLMG solver
[in]verbosityThe verbosity setting for the MLMG solver
[in]geomthe geometry per level (e.g., from AmrMesh)
[in]dmapthe distribution mapping per level (e.g., from AmrMesh)
[in]gridsthe grids per level (e.g., from AmrMesh)
[in]grid_typegrid type (e.g., collocated, staggered, hybrid)
[in]is_solver_igf_on_lev0boolean to select the Poisson solver: 1 for FFT on level 0 & Multigrid on other levels, 0 for Multigrid on all levels
[in]eb_enabledwhether embedded boundaries are enabled
[in]do_single_precision_commsperform communications in single precision
[in]rel_ref_ratiomesh refinement ratio between levels (default: 1)
[in]post_phi_calculationperform a calculation per level directly after phi was calculated; required for embedded boundaries (default: none)
[in]boundary_handlera handler for boundary conditions, for example
See also
ElectrostaticSolver::PoissonBoundaryHandler
Parameters
[in]current_timethe current time; required for embedded boundaries (default: none)
[in]eb_farray_box_factorya factory for field data,
See also
amrex::EBFArrayBoxFactory; required for embedded boundaries (default: none)

◆ computePhi()

template<typename T_PostPhiCalculationFunctor = std::nullopt_t, typename T_BoundaryHandler = std::nullopt_t, typename T_FArrayBoxFactory = void>
void ablastr::fields::computePhi ( ablastr::fields::MultiLevelScalarField const & rho,
ablastr::fields::MultiLevelScalarField const & phi,
std::array< amrex::Real, 3 > const beta,
amrex::Real relative_tolerance,
amrex::Real absolute_tolerance,
int max_iters,
int verbosity,
amrex::Vector< amrex::Geometry > const & geom,
amrex::Vector< amrex::DistributionMapping > const & dmap,
amrex::Vector< amrex::BoxArray > const & grids,
utils::enums::GridType grid_type,
bool is_solver_igf_on_lev0,
bool const is_igf_2d,
bool eb_enabled = false,
bool do_single_precision_comms = false,
std::optional< amrex::Vector< amrex::IntVect > > rel_ref_ratio = std::nullopt,
T_PostPhiCalculationFunctor post_phi_calculation = std::nullopt,
T_BoundaryHandler const boundary_handler = std::nullopt,
std::optional< amrex::Real const > current_time = std::nullopt,
std::optional< amrex::Vector< T_FArrayBoxFactory const * > > eb_farray_box_factory = std::nullopt )

Compute the potential phi by solving the Poisson equation

Uses rho as a source, assuming that the source moves at a constant speed $\vec{\beta}$. This uses the AMReX solver.

More specifically, this solves the equation

\[  \vec{\nabla}^2 r \phi - (\vec{\beta}\cdot\vec{\nabla})^2 r \phi = -\frac{r \rho}{\epsilon_0}
\]

Template Parameters
T_PostPhiCalculationFunctora calculation per level directly after phi was calculated
T_BoundaryHandlerhandler for boundary conditions, for example
See also
ElectrostaticSolver::PoissonBoundaryHandler (EB ONLY)
Template Parameters
T_FArrayBoxFactoryusually nothing or an amrex::EBFArrayBoxFactory (EB ONLY)
Parameters
[in]rhoThe charge density a given species
[out]phiThe potential to be computed by this function
[in]betaRepresents the velocity of the source of phi
[in]relative_toleranceThe relative convergence threshold for the MLMG solver
[in]absolute_toleranceThe absolute convergence threshold for the MLMG solver
[in]max_itersThe maximum number of iterations allowed for the MLMG solver
[in]verbosityThe verbosity setting for the MLMG solver
[in]geomthe geometry per level (e.g., from AmrMesh)
[in]dmapthe distribution mapping per level (e.g., from AmrMesh)
[in]gridsthe grids per level (e.g., from AmrMesh)
[in]grid_typeInteger that corresponds to the type of grid used in the simulation (collocated, staggered, hybrid)
[in]is_solver_igf_on_lev0boolean to select the Poisson solver: 1 for FFT on level 0 & Multigrid on other levels, 0 for Multigrid on all levels
[in]is_igf_2dboolean to select between fully 3D Poisson solver and quasi-3D, i.e. one 2D Poisson solve on every z slice (default: false)
[in]eb_enabledsolve with embedded boundaries
[in]do_single_precision_commsperform communications in single precision
[in]rel_ref_ratiomesh refinement ratio between levels (default: 1)
[in]post_phi_calculationperform a calculation per level directly after phi was calculated; required for embedded boundaries (default: none)
[in]boundary_handlera handler for boundary conditions, for example
See also
ElectrostaticSolver::PoissonBoundaryHandler
Parameters
[in]current_timethe current time; required for embedded boundaries (default: none)
[in]eb_farray_box_factorya factory for field data,
See also
amrex::EBFArrayBoxFactory; required for embedded boundaries (default: none)

◆ computePhiIGF()

void ablastr::fields::computePhiIGF ( amrex::MultiFab const & rho,
amrex::MultiFab & phi,
std::array< amrex::Real, 3 > const & cell_size,
bool is_igf_2d_slices )

Compute the electrostatic potential using the Integrated Green Function method as in http://dx.doi.org/10.1103/PhysRevSTAB.9.044204.

Parameters
[in]rhothe charge density amrex::MultiFab
[out]phithe electrostatic potential amrex::MultiFab
[in]cell_sizean arreay of 3 reals dx dy dz
[in]is_igf_2d_slicesboolean to select between fully 3D Poisson solver and quasi-3D, i.e. one 2D Poisson solve on every z slice (default: false)

◆ computeVectorPotential()

template<typename T_BoundaryHandler, typename T_PostACalculationFunctor = std::nullopt_t, typename T_FArrayBoxFactory = void>
void ablastr::fields::computeVectorPotential ( amrex::Vector< amrex::Array< amrex::MultiFab *, 3 > > const & curr,
amrex::Vector< amrex::Array< amrex::MultiFab *, 3 > > & A,
amrex::Real relative_tolerance,
amrex::Real absolute_tolerance,
int max_iters,
int verbosity,
amrex::Vector< amrex::Geometry > const & geom,
amrex::Vector< amrex::DistributionMapping > const & dmap,
amrex::Vector< amrex::BoxArray > const & grids,
T_BoundaryHandler const boundary_handler,
bool eb_enabled = false,
bool do_single_precision_comms = false,
std::optional< amrex::Vector< amrex::IntVect > > rel_ref_ratio = std::nullopt,
T_PostACalculationFunctor post_A_calculation = std::nullopt,
std::optional< amrex::Real const > current_time = std::nullopt,
std::optional< amrex::Vector< T_FArrayBoxFactory const * > > eb_farray_box_factory = std::nullopt )

Compute the vector potential A by solving the Poisson equation

Uses J as a source, assuming that the source moves at a constant speed $\vec{\beta}$. This uses the AMReX solver.

More specifically, this solves the equation

\[  \vec{\nabla}^2 r \vec{A} - (\vec{\beta}\cdot\vec{\nabla})^2 r \vec{A} = - r \mu_0 \vec{J}
\]

Template Parameters
T_BoundaryHandlerhandler for boundary conditions, for example
See also
MagnetostaticSolver::MultiPoissonBoundaryHandler
Template Parameters
T_PostACalculationFunctora calculation per level directly after A was calculated
T_FArrayBoxFactoryusually nothing or an amrex::EBFArrayBoxFactory (EB ONLY)
Parameters
[in]currThe current density a given species
[out]AThe vector potential to be computed by this function
[in]relative_toleranceThe relative convergence threshold for the MLMG solver
[in]absolute_toleranceThe absolute convergence threshold for the MLMG solver
[in]max_itersThe maximum number of iterations allowed for the MLMG solver
[in]verbosityThe verbosity setting for the MLMG solver
[in]geomthe geometry per level (e.g., from AmrMesh)
[in]dmapthe distribution mapping per level (e.g., from AmrMesh)
[in]gridsthe grids per level (e.g., from AmrMesh)
[in]boundary_handlera handler for boundary conditions, for example
See also
MagnetostaticSolver::VectorPoissonBoundaryHandler
Parameters
[in]eb_enabledsolve with embedded boundaries
[in]do_single_precision_commsperform communications in single precision
[in]rel_ref_ratiomesh refinement ratio between levels (default: 1)
[in]post_A_calculationperform a calculation per level directly after A was calculated; required for embedded boundaries (default: none)
[in]current_timethe current time; required for embedded boundaries (default: none)
[in]eb_farray_box_factorya factory for field data,
See also
amrex::EBFArrayBoxFactory; required for embedded boundaries (default: none)

◆ getMaxNormRho()

amrex::Real ablastr::fields::getMaxNormRho ( ablastr::fields::ConstMultiLevelScalarField const & rho,
int finest_level,
amrex::Real & absolute_tolerance )
inline

Compute the L-infinity norm of the charge density rho across all MR levels to determine if rho is zero everywhere

Parameters
[in]rhoThe charge density a given species
[in]finest_levelThe most refined mesh refinement level
[in]absolute_toleranceThe absolute convergence threshold for the MLMG solver

◆ IntegratedPotential2D()

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real ablastr::fields::IntegratedPotential2D ( amrex::Real x,
amrex::Real y )

Implements equation 58 in https://doi.org/10.1016/j.jcp.2004.01.008.

Parameters
[in]xx-coordinate of given location
[in]yy-coordinate of given location
Returns
the integrated Green function G in 2D

◆ IntegratedPotential3D()

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real ablastr::fields::IntegratedPotential3D ( amrex::Real x,
amrex::Real y,
amrex::Real z )

Implements equation 2 in https://doi.org/10.1103/PhysRevSTAB.10.129901 with some modification to symmetrize the function.

Parameters
[in]xx-coordinate of given location
[in]yy-coordinate of given location
[in]zz-coordinate of given location
Returns
the integrated Green function G in 3D

◆ interpolatePhiBetweenLevels()

void ablastr::fields::interpolatePhiBetweenLevels ( amrex::MultiFab const * phi_lev,
amrex::MultiFab * phi_lev_plus_one,
amrex::Geometry const & geom_lev,
bool do_single_precision_comms,
const amrex::IntVect & refratio,
const int ncomp,
const int ng )
inline

Interpolate the potential phi from level lev to lev+1

Needed to solve Poisson equation on lev+1 The coarser level lev provides both the boundary values and initial guess for phi on the finer level lev+1

Parameters
[in]phi_levThe potential on a given mesh refinement level lev
[in]phi_lev_plus_oneThe potential on the next level lev+1
[in]geom_levThe geometry of level lev
[in]do_single_precision_commsperform communications in single precision
[in]refratiomesh refinement ratio between level lev and lev+1
[in]ncompNumber of components of the multifab (1)
[in]ngNumber of ghost cells (1 if collocated, 0 otherwise)

◆ SumOfIntegratedPotential2D()

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real ablastr::fields::SumOfIntegratedPotential2D ( amrex::Real x,
amrex::Real y,
amrex::Real dx,
amrex::Real dy )

add

Parameters
[in]xx-coordinate of given location
[in]yy-coordinate of given location
[in]dxcell size along x
[in]dycell size along y
Returns
the sum of integrated Green function G in 2D

◆ SumOfIntegratedPotential3D()

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real ablastr::fields::SumOfIntegratedPotential3D ( amrex::Real x,
amrex::Real y,
amrex::Real z,
amrex::Real dx,
amrex::Real dy,
amrex::Real dz )

add

Parameters
[in]xx-coordinate of given location
[in]yy-coordinate of given location
[in]zz-coordinate of given location
[in]dxcell size along x
[in]dycell size along y
[in]dzcell size along z
Returns
the sum of integrated Green function G in 3D