|
| 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) |
| |
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
- Template Parameters
-
| T_PostPhiCalculationFunctor | a calculation per level directly after phi was calculated |
| T_BoundaryHandler | handler for boundary conditions, for example |
- See also
- ElectrostaticSolver::PoissonBoundaryHandler
- Template Parameters
-
- Parameters
-
| [in] | rho | The charge density a given species |
| [out] | phi | The potential to be computed by this function |
| [in] | sigma | The matrix representing the mass operator used to lower the local plasma frequency |
| [in] | relative_tolerance | The relative convergence threshold for the MLMG solver |
| [in] | absolute_tolerance | The absolute convergence threshold for the MLMG solver |
| [in] | max_iters | The maximum number of iterations allowed for the MLMG solver |
| [in] | verbosity | The verbosity setting for the MLMG solver |
| [in] | geom | the geometry per level (e.g., from AmrMesh) |
| [in] | dmap | the distribution mapping per level (e.g., from AmrMesh) |
| [in] | grids | the grids per level (e.g., from AmrMesh) |
| [in] | grid_type | grid type (e.g., collocated, staggered, hybrid) |
| [in] | is_solver_igf_on_lev0 | boolean to select the Poisson solver: 1 for FFT on level 0 & Multigrid on other levels, 0 for Multigrid on all levels |
| [in] | eb_enabled | whether embedded boundaries are enabled |
| [in] | do_single_precision_comms | perform communications in single precision |
| [in] | rel_ref_ratio | mesh refinement ratio between levels (default: 1) |
| [in] | post_phi_calculation | perform a calculation per level directly after phi was calculated; required for embedded boundaries (default: none) |
| [in] | boundary_handler | a handler for boundary conditions, for example |
- See also
- ElectrostaticSolver::PoissonBoundaryHandler
- Parameters
-
| [in] | current_time | the current time; required for embedded boundaries (default: none) |
| [in] | eb_farray_box_factory | a factory for field data, |
- See also
- amrex::EBFArrayBoxFactory; required for embedded boundaries (default: none)
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
. This uses the AMReX solver.
More specifically, this solves the equation
- Template Parameters
-
| T_PostPhiCalculationFunctor | a calculation per level directly after phi was calculated |
| T_BoundaryHandler | handler for boundary conditions, for example |
- See also
- ElectrostaticSolver::PoissonBoundaryHandler (EB ONLY)
- Template Parameters
-
- Parameters
-
| [in] | rho | The charge density a given species |
| [out] | phi | The potential to be computed by this function |
| [in] | beta | Represents the velocity of the source of phi |
| [in] | relative_tolerance | The relative convergence threshold for the MLMG solver |
| [in] | absolute_tolerance | The absolute convergence threshold for the MLMG solver |
| [in] | max_iters | The maximum number of iterations allowed for the MLMG solver |
| [in] | verbosity | The verbosity setting for the MLMG solver |
| [in] | geom | the geometry per level (e.g., from AmrMesh) |
| [in] | dmap | the distribution mapping per level (e.g., from AmrMesh) |
| [in] | grids | the grids per level (e.g., from AmrMesh) |
| [in] | grid_type | Integer that corresponds to the type of grid used in the simulation (collocated, staggered, hybrid) |
| [in] | is_solver_igf_on_lev0 | boolean to select the Poisson solver: 1 for FFT on level 0 & Multigrid on other levels, 0 for Multigrid on all levels |
| [in] | is_igf_2d | boolean to select between fully 3D Poisson solver and quasi-3D, i.e. one 2D Poisson solve on every z slice (default: false) |
| [in] | eb_enabled | solve with embedded boundaries |
| [in] | do_single_precision_comms | perform communications in single precision |
| [in] | rel_ref_ratio | mesh refinement ratio between levels (default: 1) |
| [in] | post_phi_calculation | perform a calculation per level directly after phi was calculated; required for embedded boundaries (default: none) |
| [in] | boundary_handler | a handler for boundary conditions, for example |
- See also
- ElectrostaticSolver::PoissonBoundaryHandler
- Parameters
-
| [in] | current_time | the current time; required for embedded boundaries (default: none) |
| [in] | eb_farray_box_factory | a factory for field data, |
- See also
- amrex::EBFArrayBoxFactory; required for embedded boundaries (default: none)
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
. This uses the AMReX solver.
More specifically, this solves the equation
- Template Parameters
-
| T_BoundaryHandler | handler for boundary conditions, for example |
- See also
- MagnetostaticSolver::MultiPoissonBoundaryHandler
- Template Parameters
-
| T_PostACalculationFunctor | a calculation per level directly after A was calculated |
| T_FArrayBoxFactory | usually nothing or an amrex::EBFArrayBoxFactory (EB ONLY) |
- Parameters
-
| [in] | curr | The current density a given species |
| [out] | A | The vector potential to be computed by this function |
| [in] | relative_tolerance | The relative convergence threshold for the MLMG solver |
| [in] | absolute_tolerance | The absolute convergence threshold for the MLMG solver |
| [in] | max_iters | The maximum number of iterations allowed for the MLMG solver |
| [in] | verbosity | The verbosity setting for the MLMG solver |
| [in] | geom | the geometry per level (e.g., from AmrMesh) |
| [in] | dmap | the distribution mapping per level (e.g., from AmrMesh) |
| [in] | grids | the grids per level (e.g., from AmrMesh) |
| [in] | boundary_handler | a handler for boundary conditions, for example |
- See also
- MagnetostaticSolver::VectorPoissonBoundaryHandler
- Parameters
-
| [in] | eb_enabled | solve with embedded boundaries |
| [in] | do_single_precision_comms | perform communications in single precision |
| [in] | rel_ref_ratio | mesh refinement ratio between levels (default: 1) |
| [in] | post_A_calculation | perform a calculation per level directly after A was calculated; required for embedded boundaries (default: none) |
| [in] | current_time | the current time; required for embedded boundaries (default: none) |
| [in] | eb_farray_box_factory | a factory for field data, |
- See also
- amrex::EBFArrayBoxFactory; required for embedded boundaries (default: none)