WarpX
Loading...
Searching...
No Matches
FiniteDifferenceSolver Class Reference

Top-level class for the electromagnetic finite-difference solver. More...

#include <FiniteDifferenceSolver.H>

Public Member Functions

 FiniteDifferenceSolver (ElectromagneticSolverAlgo fdtd_algo, std::array< amrex::Real, 3 > cell_size, ablastr::utils::enums::GridType grid_type)
 Initialize the finite-difference Maxwell solver (for a given refinement level)
 
void EvolveB (ablastr::fields::MultiFabRegister &fields, int lev, PatchType patch_type, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > &flag_info_cell, std::array< std::unique_ptr< amrex::LayoutData< FaceInfoBox > >, 3 > &borrowing, amrex::Real dt)
 Update the B field, over one timestep.
 
void EvolveE (ablastr::fields::MultiFabRegister &fields, int lev, PatchType patch_type, ablastr::fields::VectorField const &Efield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const &eb_update_E, amrex::Real dt)
 Update the E field, over one timestep.
 
void EvolveF (amrex::MultiFab *Ffield, ablastr::fields::VectorField const &Efield, amrex::MultiFab *rhofield, int rho_comp, amrex::Real dt)
 Update the F field, over one timestep.
 
void EvolveG (amrex::MultiFab *Gfield, ablastr::fields::VectorField const &Bfield, amrex::Real dt)
 
void EvolveECTRho (ablastr::fields::VectorField const &Efield, ablastr::fields::VectorField const &edge_lengths, ablastr::fields::VectorField const &face_areas, ablastr::fields::VectorField const &ECTRhofield, int lev)
 Update the B field, over one timestep.
 
void ApplySilverMuellerBoundary (ablastr::fields::VectorField &Efield, ablastr::fields::VectorField &Bfield, amrex::Box domain_box, amrex::Real dt, amrex::Array< FieldBoundaryType, 3 > field_boundary_lo, amrex::Array< FieldBoundaryType, 3 > field_boundary_hi)
 Update the B field at the boundary, using the Silver-Mueller condition.
 
void ComputeDivE (ablastr::fields::VectorField const &Efield, amrex::MultiFab &divE)
 Update the F field, over one timestep.
 
void MacroscopicEvolveE (MacroscopicSolverAlgo macroscopic_solver_algo, ablastr::fields::VectorField const &Efield, ablastr::fields::VectorField const &Bfield, ablastr::fields::VectorField const &Jfield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const &eb_update_E, amrex::Real dt, std::unique_ptr< MacroscopicProperties > const &macroscopic_properties)
 Macroscopic E-update for non-vacuum medium using the user-selected finite-difference algorithm and macroscopic sigma-method defined in WarpXAlgorithmSelection.H.
 
void EvolveBPML (ablastr::fields::MultiFabRegister &fields, PatchType patch_type, int level, amrex::Real dt, bool dive_cleaning)
 Update the B field, over one timestep.
 
void EvolveEPML (ablastr::fields::MultiFabRegister &fields, PatchType patch_type, int level, MultiSigmaBox const &sigba, amrex::Real dt, bool pml_has_particles)
 Update the E field, over one timestep.
 
void EvolveFPML (amrex::MultiFab *Ffield, ablastr::fields::VectorField Efield, amrex::Real dt)
 Update the E field, over one timestep.
 
void HybridPICSolveE (ablastr::fields::VectorField const &Efield, ablastr::fields::VectorField &Jfield, ablastr::fields::VectorField const &Jifield, ablastr::fields::VectorField const &Bfield, amrex::MultiFab const &rhofield, amrex::MultiFab const &Pefield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const &eb_update_E, int lev, HybridPICModel const *hybrid_model, bool solve_for_Faraday)
 E-update in the hybrid PIC algorithm as described in Winske et al. (2003) Eq. 10. https://link.springer.com/chapter/10.1007/3-540-36530-3_8.
 
void CalculateCurrentAmpere (ablastr::fields::VectorField &Jfield, ablastr::fields::VectorField const &Bfield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const &eb_update_E, int lev)
 Calculation of total current using Ampere's law (without displacement current): J = (curl x B) / mu0.
 
void ComputeCurlA (ablastr::fields::VectorField &Bfield, ablastr::fields::VectorField const &Afield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const &eb_update_B, int lev)
 Calculation of B field from the vector potential A B = (curl x A) / mu0.
 
template<typename T_Algo>
void EvolveBCylindrical (ablastr::fields::VectorField const &Bfield, ablastr::fields::VectorField const &Efield, int lev, amrex::Real dt)
 
template<typename T_Algo>
void EvolveECylindrical (ablastr::fields::VectorField const &Efield, ablastr::fields::VectorField const &Bfield, ablastr::fields::VectorField const &Jfield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const &eb_update_E, amrex::MultiFab const *Ffield, int lev, amrex::Real dt)
 
template<typename T_Algo>
void EvolveFCylindrical (amrex::MultiFab *Ffield, ablastr::fields::VectorField const &Efield, amrex::MultiFab *rhofield, int rho_comp, amrex::Real dt)
 
template<typename T_Algo>
void ComputeDivECylindrical (ablastr::fields::VectorField const &Efield, amrex::MultiFab &divE)
 
template<typename T_Algo>
void HybridPICSolveECylindrical (ablastr::fields::VectorField const &Efield, ablastr::fields::VectorField const &Jfield, ablastr::fields::VectorField const &Jifield, ablastr::fields::VectorField const &Bfield, amrex::MultiFab const &rhofield, amrex::MultiFab const &Pefield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const &eb_update_E, int lev, HybridPICModel const *hybrid_model, bool solve_for_Faraday)
 
template<typename T_Algo>
void CalculateCurrentAmpereCylindrical (ablastr::fields::VectorField &Jfield, ablastr::fields::VectorField const &Bfield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const &eb_update_E, int lev)
 Calculate total current from Ampere's law without displacement current i.e. J = 1/mu_0 curl x B.
 
template<typename T_Algo>
void ComputeCurlACylindrical (ablastr::fields::VectorField &Bfield, ablastr::fields::VectorField const &Afield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const &eb_update_B, int lev)
 Calculate B from the curl of A i.e. B = curl(A) output field on B field mesh staggering.
 

Private Attributes

ElectromagneticSolverAlgo m_fdtd_algo
 
ablastr::utils::enums::GridType m_grid_type
 
amrex::Real m_dr
 
amrex::Real m_rmin
 
int m_nmodes
 
amrex::Vector< amrex::Real > m_h_stencil_coefs_r
 
amrex::Vector< amrex::Real > m_h_stencil_coefs_z
 
amrex::Gpu::DeviceVector< amrex::Real > m_stencil_coefs_r
 
amrex::Gpu::DeviceVector< amrex::Real > m_stencil_coefs_z
 

Detailed Description

Top-level class for the electromagnetic finite-difference solver.

Stores the coefficients of the finite-difference stencils, and has member functions to update fields over one time step.

Constructor & Destructor Documentation

◆ FiniteDifferenceSolver()

FiniteDifferenceSolver::FiniteDifferenceSolver ( ElectromagneticSolverAlgo fdtd_algo,
std::array< amrex::Real, 3 > cell_size,
ablastr::utils::enums::GridType grid_type )

Initialize the finite-difference Maxwell solver (for a given refinement level)

This function initializes the stencil coefficients for the chosen finite-difference algorithm

Parameters
fdtd_algoIdentifies the chosen algorithm, as defined in WarpXAlgorithmSelection.H
cell_sizeCell size along each dimension, for the chosen refinement level
grid_typeWhether the solver is applied to a collocated or staggered grid

Member Function Documentation

◆ ApplySilverMuellerBoundary()

void FiniteDifferenceSolver::ApplySilverMuellerBoundary ( ablastr::fields::VectorField & Efield,
ablastr::fields::VectorField & Bfield,
amrex::Box domain_box,
amrex::Real dt,
amrex::Array< FieldBoundaryType, 3 > field_boundary_lo,
amrex::Array< FieldBoundaryType, 3 > field_boundary_hi )

Update the B field at the boundary, using the Silver-Mueller condition.

◆ CalculateCurrentAmpere()

void FiniteDifferenceSolver::CalculateCurrentAmpere ( ablastr::fields::VectorField & Jfield,
ablastr::fields::VectorField const & Bfield,
std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const & eb_update_E,
int lev )

Calculation of total current using Ampere's law (without displacement current): J = (curl x B) / mu0.

Parameters
[out]Jfieldvector of current MultiFabs at a given level
[in]Bfieldvector of magnetic field MultiFabs at a given level
[in]eb_update_Eindicate in which cell E should be updated (related to embedded boundaries)
[in]levlevel number for the calculation

◆ CalculateCurrentAmpereCylindrical()

template<typename T_Algo>
void FiniteDifferenceSolver::CalculateCurrentAmpereCylindrical ( ablastr::fields::VectorField & Jfield,
ablastr::fields::VectorField const & Bfield,
std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const & eb_update_E,
int lev )

Calculate total current from Ampere's law without displacement current i.e. J = 1/mu_0 curl x B.

Parameters
[out]Jfieldvector of total current MultiFabs at a given level
[in]Bfieldvector of magnetic field MultiFabs at a given level
[in]eb_update_Especifies where the plasma current should be calculated.
[in]levrefinement level

◆ ComputeCurlA()

void FiniteDifferenceSolver::ComputeCurlA ( ablastr::fields::VectorField & Bfield,
ablastr::fields::VectorField const & Afield,
std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const & eb_update_B,
int lev )

Calculation of B field from the vector potential A B = (curl x A) / mu0.

Parameters
[out]Bfieldvector of current MultiFabs at a given level
[in]Afieldvector of magnetic field MultiFabs at a given level
[in]eb_update_Barray indicating where the field should be updated with respect to the position of the embedded boundary
[in]levlevel number for the calculation

◆ ComputeCurlACylindrical()

template<typename T_Algo>
void FiniteDifferenceSolver::ComputeCurlACylindrical ( ablastr::fields::VectorField & Bfield,
ablastr::fields::VectorField const & Afield,
std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const & eb_update_B,
int lev )

Calculate B from the curl of A i.e. B = curl(A) output field on B field mesh staggering.

Parameters
[out]Bfieldoutput of curl operation
[in]Afieldinput staggered field, should be on E/J/A mesh staggering
[in]eb_update_Bspecifies where the plasma current should be calculated.
[in]levrefinement level

◆ ComputeDivE()

void FiniteDifferenceSolver::ComputeDivE ( ablastr::fields::VectorField const & Efield,
amrex::MultiFab & divE )

Update the F field, over one timestep.

◆ ComputeDivECylindrical()

template<typename T_Algo>
void FiniteDifferenceSolver::ComputeDivECylindrical ( ablastr::fields::VectorField const & Efield,
amrex::MultiFab & divE )

◆ EvolveB()

void FiniteDifferenceSolver::EvolveB ( ablastr::fields::MultiFabRegister & fields,
int lev,
PatchType patch_type,
std::array< std::unique_ptr< amrex::iMultiFab >, 3 > & flag_info_cell,
std::array< std::unique_ptr< amrex::LayoutData< FaceInfoBox > >, 3 > & borrowing,
amrex::Real dt )

Update the B field, over one timestep.

◆ EvolveBCylindrical()

template<typename T_Algo>
void FiniteDifferenceSolver::EvolveBCylindrical ( ablastr::fields::VectorField const & Bfield,
ablastr::fields::VectorField const & Efield,
int lev,
amrex::Real dt )

◆ EvolveBPML()

void FiniteDifferenceSolver::EvolveBPML ( ablastr::fields::MultiFabRegister & fields,
PatchType patch_type,
int level,
amrex::Real dt,
bool dive_cleaning )

Update the B field, over one timestep.

◆ EvolveE()

void FiniteDifferenceSolver::EvolveE ( ablastr::fields::MultiFabRegister & fields,
int lev,
PatchType patch_type,
ablastr::fields::VectorField const & Efield,
std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const & eb_update_E,
amrex::Real dt )

Update the E field, over one timestep.

◆ EvolveECTRho()

void FiniteDifferenceSolver::EvolveECTRho ( ablastr::fields::VectorField const & Efield,
ablastr::fields::VectorField const & edge_lengths,
ablastr::fields::VectorField const & face_areas,
ablastr::fields::VectorField const & ECTRhofield,
int lev )

Update the B field, over one timestep.

◆ EvolveECylindrical()

template<typename T_Algo>
void FiniteDifferenceSolver::EvolveECylindrical ( ablastr::fields::VectorField const & Efield,
ablastr::fields::VectorField const & Bfield,
ablastr::fields::VectorField const & Jfield,
std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const & eb_update_E,
amrex::MultiFab const * Ffield,
int lev,
amrex::Real dt )

◆ EvolveEPML()

void FiniteDifferenceSolver::EvolveEPML ( ablastr::fields::MultiFabRegister & fields,
PatchType patch_type,
int level,
MultiSigmaBox const & sigba,
amrex::Real dt,
bool pml_has_particles )

Update the E field, over one timestep.

◆ EvolveF()

void FiniteDifferenceSolver::EvolveF ( amrex::MultiFab * Ffield,
ablastr::fields::VectorField const & Efield,
amrex::MultiFab * rhofield,
int rho_comp,
amrex::Real dt )

Update the F field, over one timestep.

◆ EvolveFCylindrical()

template<typename T_Algo>
void FiniteDifferenceSolver::EvolveFCylindrical ( amrex::MultiFab * Ffield,
ablastr::fields::VectorField const & Efield,
amrex::MultiFab * rhofield,
int rho_comp,
amrex::Real dt )

◆ EvolveFPML()

void FiniteDifferenceSolver::EvolveFPML ( amrex::MultiFab * Ffield,
ablastr::fields::VectorField Efield,
amrex::Real dt )

Update the E field, over one timestep.

◆ EvolveG()

void FiniteDifferenceSolver::EvolveG ( amrex::MultiFab * Gfield,
ablastr::fields::VectorField const & Bfield,
amrex::Real dt )

◆ HybridPICSolveE()

void FiniteDifferenceSolver::HybridPICSolveE ( ablastr::fields::VectorField const & Efield,
ablastr::fields::VectorField & Jfield,
ablastr::fields::VectorField const & Jifield,
ablastr::fields::VectorField const & Bfield,
amrex::MultiFab const & rhofield,
amrex::MultiFab const & Pefield,
std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const & eb_update_E,
int lev,
HybridPICModel const * hybrid_model,
bool solve_for_Faraday )

E-update in the hybrid PIC algorithm as described in Winske et al. (2003) Eq. 10. https://link.springer.com/chapter/10.1007/3-540-36530-3_8.

Parameters
[out]Efieldvector of electric field MultiFabs updated at a given level
[in]Jfieldvector of total plasma current MultiFabs at a given level
[in]Jifieldvector of ion current density MultiFabs at a given level
[in]Bfieldvector of magnetic field MultiFabs at a given level
[in]rhofieldscalar ion charge density Multifab at a given level
[in]Pefieldscalar electron pressure MultiFab at a given level
[in]eb_update_Eindicate in which cell E should be updated (related to embedded boundaries)
[in]levlevel number for the calculation
[in]hybrid_modelinstance of the hybrid-PIC model
[in]solve_for_Faradayboolean flag for whether the E-field is solved to be used in Faraday's equation

◆ HybridPICSolveECylindrical()

template<typename T_Algo>
void FiniteDifferenceSolver::HybridPICSolveECylindrical ( ablastr::fields::VectorField const & Efield,
ablastr::fields::VectorField const & Jfield,
ablastr::fields::VectorField const & Jifield,
ablastr::fields::VectorField const & Bfield,
amrex::MultiFab const & rhofield,
amrex::MultiFab const & Pefield,
std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const & eb_update_E,
int lev,
HybridPICModel const * hybrid_model,
bool solve_for_Faraday )

◆ MacroscopicEvolveE()

void FiniteDifferenceSolver::MacroscopicEvolveE ( MacroscopicSolverAlgo macroscopic_solver_algo,
ablastr::fields::VectorField const & Efield,
ablastr::fields::VectorField const & Bfield,
ablastr::fields::VectorField const & Jfield,
std::array< std::unique_ptr< amrex::iMultiFab >, 3 > const & eb_update_E,
amrex::Real dt,
std::unique_ptr< MacroscopicProperties > const & macroscopic_properties )

Macroscopic E-update for non-vacuum medium using the user-selected finite-difference algorithm and macroscopic sigma-method defined in WarpXAlgorithmSelection.H.

Parameters
[in]macroscopic_solver_algomacroscopic Maxwell solver algorithm (BackwardEuler or Lax-Wendroff)
[out]Efieldvector of electric field MultiFabs updated at a given level
[in]Bfieldvector of magnetic field MultiFabs at a given level
[in]Jfieldvector of current density MultiFabs at a given level
[in]eb_update_Eindicate in which cell E should be updated (related to embedded boundaries)
[in]dttimestep of the simulation
[in]macroscopic_propertiescontains user-defined properties of the medium.

Member Data Documentation

◆ m_dr

amrex::Real FiniteDifferenceSolver::m_dr
private

◆ m_fdtd_algo

ElectromagneticSolverAlgo FiniteDifferenceSolver::m_fdtd_algo
private

◆ m_grid_type

ablastr::utils::enums::GridType FiniteDifferenceSolver::m_grid_type
private

◆ m_h_stencil_coefs_r

amrex::Vector<amrex::Real> FiniteDifferenceSolver::m_h_stencil_coefs_r
private

◆ m_h_stencil_coefs_z

amrex::Vector<amrex::Real> FiniteDifferenceSolver::m_h_stencil_coefs_z
private

◆ m_nmodes

int FiniteDifferenceSolver::m_nmodes
private

◆ m_rmin

amrex::Real FiniteDifferenceSolver::m_rmin
private

◆ m_stencil_coefs_r

amrex::Gpu::DeviceVector<amrex::Real> FiniteDifferenceSolver::m_stencil_coefs_r
private

◆ m_stencil_coefs_z

amrex::Gpu::DeviceVector<amrex::Real> FiniteDifferenceSolver::m_stencil_coefs_z
private

The documentation for this class was generated from the following files: