WarpX
Loading...
Searching...
No Matches
RelativisticExplicitES Class Referencefinal

#include <RelativisticExplicitES.H>

Inheritance diagram for RelativisticExplicitES:
ElectrostaticSolver

Public Member Functions

 RelativisticExplicitES (int nlevs_max)
 
void InitData () override
 
void ComputeSpaceChargeField (ablastr::fields::MultiFabRegister &fields, MultiParticleContainer &mpc, MultiFluidContainer *mfl, int max_level) override
 Computes electrostatic fields for species that have initialize self fields turned on. A loop over all the species is performed and for each species (with self fields) the function AddSpaceChargeField is called. This function computes the electrostatic potential for species charge denisyt as source and then the electric and magnetic fields are updated to include the corresponding fields from the electrostatic potential. Then electric and magnetic fields are updated to include potential variation due to boundary conditions using the function AddBoundaryField
 
void AddSpaceChargeField (WarpXParticleContainer &pc, ablastr::fields::MultiLevelVectorField &Efield_fp, ablastr::fields::MultiLevelVectorField &Bfield_fp)
 
void AddBoundaryField (ablastr::fields::MultiLevelVectorField &Efield)
 
- Public Member Functions inherited from ElectrostaticSolver
 ElectrostaticSolver ()=default
 
 ElectrostaticSolver (int nlevs_max)
 
virtual ~ElectrostaticSolver ()
 
 ElectrostaticSolver (const ElectrostaticSolver &)=delete
 
ElectrostaticSolveroperator= (const ElectrostaticSolver &)=delete
 
 ElectrostaticSolver (ElectrostaticSolver &&)=delete
 
ElectrostaticSolveroperator= (ElectrostaticSolver &&)=delete
 
void ReadParameters ()
 
void setPhiBC (ablastr::fields::MultiLevelScalarField const &phi, amrex::Real t) const
 Set Dirichlet boundary conditions for the electrostatic solver. The given potential's values are fixed on the boundaries of the given dimension according to the desired values from the simulation input file, boundary.potential_lo and boundary.potential_hi.
 
void computePhi (ablastr::fields::MultiLevelScalarField const &rho, ablastr::fields::MultiLevelScalarField const &phi, std::array< amrex::Real, 3 > beta, amrex::Real required_precision, amrex::Real absolute_tolerance, int max_iters, int verbosity, bool is_igf_2d_slices, std::optional< ablastr::fields::MultiLevelVectorField > efield=std::nullopt) const
 
void computeE (ablastr::fields::MultiLevelVectorField const &E, ablastr::fields::MultiLevelScalarField const &phi, std::array< amrex::Real, 3 > beta) const
 Compute the electric field that corresponds to phi, and add it to the set of MultiFab E. The electric field is calculated by assuming that the source that produces the phi potential is moving with a constant speed $\vec{\beta}$:
 
void computeB (ablastr::fields::MultiLevelVectorField const &B, ablastr::fields::MultiLevelScalarField const &phi, std::array< amrex::Real, 3 > beta) const
 Compute the magnetic field that corresponds to phi, and add it to the set of MultiFab B. The magnetic field is calculated by assuming that the source that produces the phi potential is moving with a constant speed $\vec{\beta}$:
 

Additional Inherited Members

- Public Attributes inherited from ElectrostaticSolver
int num_levels
 
std::unique_ptr< PoissonBoundaryHandlerm_poisson_boundary_handler
 
amrex::Real self_fields_required_precision = 1e-11
 
amrex::Real self_fields_absolute_tolerance = 0.0
 
int self_fields_max_iters = 200
 
int self_fields_verbosity = 2
 
bool is_igf_2d_slices = false
 

Constructor & Destructor Documentation

◆ RelativisticExplicitES()

RelativisticExplicitES::RelativisticExplicitES ( int nlevs_max)
inline

Member Function Documentation

◆ AddBoundaryField()

void RelativisticExplicitES::AddBoundaryField ( ablastr::fields::MultiLevelVectorField & Efield)

Compute the potential phi by solving the Poisson equation with the simulation specific boundary conditions and boundary values, then add the E field due to that phi to Efield_fp.

Parameters
[in]EfieldEfield updated to include potential gradient from boundary condition

◆ AddSpaceChargeField()

void RelativisticExplicitES::AddSpaceChargeField ( WarpXParticleContainer & pc,
ablastr::fields::MultiLevelVectorField & Efield_fp,
ablastr::fields::MultiLevelVectorField & Bfield_fp )

Compute the charge density of the species paricle container, pc, and obtain the corresponding electrostatic potential to update the electric and magnetic fields.

Parameters
[in]pcparticle container for the selected species
[in]Efield_fpEfield updated to include potential computed for selected species charge density as source
[in]Bfield_fpBfield updated to include potential computed for selected species charge density as source

◆ ComputeSpaceChargeField()

void RelativisticExplicitES::ComputeSpaceChargeField ( ablastr::fields::MultiFabRegister & fields,
MultiParticleContainer & mpc,
MultiFluidContainer * mfl,
int max_level )
overridevirtual

Computes electrostatic fields for species that have initialize self fields turned on. A loop over all the species is performed and for each species (with self fields) the function AddSpaceChargeField is called. This function computes the electrostatic potential for species charge denisyt as source and then the electric and magnetic fields are updated to include the corresponding fields from the electrostatic potential. Then electric and magnetic fields are updated to include potential variation due to boundary conditions using the function AddBoundaryField

Parameters
fieldsMultiFab register of fields
mpcReference to multi-particle container
mflPointer to multi-fluid container
max_levelMaximum level of mesh refinement

Implements ElectrostaticSolver.

◆ InitData()

void RelativisticExplicitES::InitData ( )
overridevirtual

Reimplemented from ElectrostaticSolver.


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