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

This class contains the parameters needed to evaluate hybrid field solutions (kinetic ions with fluid electrons). More...

#include <HybridPICModel.H>

Public Member Functions

 HybridPICModel ()
 
void ReadParameters ()
 
void AllocateLevelMFs (ablastr::fields::MultiFabRegister &fields, int lev, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm, int ncomps, const amrex::IntVect &ngJ, const amrex::IntVect &ngRho, const amrex::IntVect &ngEB, const amrex::IntVect &jx_nodal_flag, const amrex::IntVect &jy_nodal_flag, const amrex::IntVect &jz_nodal_flag, const amrex::IntVect &rho_nodal_flag, const amrex::IntVect &Ex_nodal_flag, const amrex::IntVect &Ey_nodal_flag, const amrex::IntVect &Ez_nodal_flag, const amrex::IntVect &Bx_nodal_flag, const amrex::IntVect &By_nodal_flag, const amrex::IntVect &Bz_nodal_flag) const
 
void InitData (const ablastr::fields::MultiFabRegister &fields)
 
void GetCurrentExternal ()
 Function to evaluate the external current expressions and populate the external current multifab. Note the external current can be a function of time and therefore this should be re-evaluated at every step.
 
void CalculatePlasmaCurrent (ablastr::fields::MultiLevelVectorField const &Bfield, amrex::Vector< std::array< std::unique_ptr< amrex::iMultiFab >, 3 > > &eb_update_E)
 Function to calculate the total plasma current based on Ampere's law while neglecting displacement current (J = curl x B). Any external current is subtracted as well. Used in the Ohm's law solver (kinetic-fluid hybrid model).
 
void CalculatePlasmaCurrent (ablastr::fields::VectorField const &Bfield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > &eb_update_E, int lev)
 
void HybridPICSolveE (ablastr::fields::MultiLevelVectorField const &Efield, ablastr::fields::MultiLevelVectorField const &Jfield, ablastr::fields::MultiLevelVectorField const &Bfield, ablastr::fields::MultiLevelScalarField const &rhofield, amrex::Vector< std::array< std::unique_ptr< amrex::iMultiFab >, 3 > > &eb_update_E, bool solve_for_Faraday) const
 Function to update the E-field using Ohm's law (hybrid-PIC model).
 
void HybridPICSolveE (ablastr::fields::VectorField const &Efield, ablastr::fields::VectorField const &Jfield, ablastr::fields::VectorField const &Bfield, amrex::MultiFab const &rhofield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > &eb_update_E, int lev, bool solve_for_Faraday) const
 
void HybridPICSolveE (ablastr::fields::VectorField const &Efield, ablastr::fields::VectorField const &Jfield, ablastr::fields::VectorField const &Bfield, amrex::MultiFab const &rhofield, std::array< std::unique_ptr< amrex::iMultiFab >, 3 > &eb_update_E, int lev, PatchType patch_type, bool solve_for_Faraday) const
 
void BfieldEvolveRK (ablastr::fields::MultiLevelVectorField const &Bfield, ablastr::fields::MultiLevelVectorField const &Efield, ablastr::fields::MultiLevelVectorField const &Jfield, ablastr::fields::MultiLevelScalarField const &rhofield, amrex::Vector< std::array< std::unique_ptr< amrex::iMultiFab >, 3 > > &eb_update_E, amrex::Real dt, SubcyclingHalf subcycling_half, amrex::IntVect ng, std::optional< bool > nodal_sync)
 
void BfieldEvolveRK (ablastr::fields::MultiLevelVectorField const &Bfield, ablastr::fields::MultiLevelVectorField const &Efield, ablastr::fields::MultiLevelVectorField const &Jfield, ablastr::fields::MultiLevelScalarField const &rhofield, amrex::Vector< std::array< std::unique_ptr< amrex::iMultiFab >, 3 > > &eb_update_E, amrex::Real dt, int lev, SubcyclingHalf subcycling_half, amrex::IntVect ng, std::optional< bool > nodal_sync)
 
void FieldPush (ablastr::fields::MultiLevelVectorField const &Bfield, ablastr::fields::MultiLevelVectorField const &Efield, ablastr::fields::MultiLevelVectorField const &Jfield, ablastr::fields::MultiLevelScalarField const &rhofield, amrex::Vector< std::array< std::unique_ptr< amrex::iMultiFab >, 3 > > &eb_update_E, amrex::Real dt, SubcyclingHalf subcycling_half, amrex::IntVect ng, std::optional< bool > nodal_sync)
 
void CalculateElectronPressure () const
 Function to calculate the electron pressure using the simulation charge density. Used in the Ohm's law solver (kinetic-fluid hybrid model).
 
void CalculateElectronPressure (int lev) const
 
void FillElectronPressureMF (amrex::MultiFab &Pe_field, amrex::MultiFab const &rho_field) const
 Fill the electron pressure multifab given the kinetic particle charge density (and assumption of quasi-neutrality) using the user specified electron equation of state.
 

Public Attributes

int m_substeps = 10
 
bool m_holmstrom_vacuum_region = false
 
amrex::Real m_elec_temp
 
amrex::Real m_n0_ref = 1.0
 
amrex::Real m_gamma = 5.0/3.0
 
amrex::Real m_n_floor = 1.0
 
std::string m_eta_expression = "0.0"
 
std::unique_ptr< amrex::Parserm_resistivity_parser
 
amrex::ParserExecutor< 2 > m_eta
 
bool m_resistivity_has_J_dependence = false
 
std::string m_eta_h_expression = "0.0"
 
std::unique_ptr< amrex::Parserm_hyper_resistivity_parser
 
amrex::ParserExecutor< 2 > m_eta_h
 
bool m_include_hyper_resistivity_term = false
 
bool m_hyper_resistivity_has_B_dependence = false
 
std::string m_Jx_ext_grid_function = "0.0"
 
std::string m_Jy_ext_grid_function = "0.0"
 
std::string m_Jz_ext_grid_function = "0.0"
 
std::array< std::unique_ptr< amrex::Parser >, 3 > m_J_external_parser
 
std::array< amrex::ParserExecutor< 4 >, 3 > m_J_external
 
bool m_external_current_has_time_dependence = false
 
bool m_add_external_fields = false
 
std::unique_ptr< ExternalVectorPotentialm_external_vector_potential
 
amrex::GpuArray< int, 3 > Jx_IndexType
 
amrex::GpuArray< int, 3 > Jy_IndexType
 
amrex::GpuArray< int, 3 > Jz_IndexType
 
amrex::GpuArray< int, 3 > Bx_IndexType
 
amrex::GpuArray< int, 3 > By_IndexType
 
amrex::GpuArray< int, 3 > Bz_IndexType
 
amrex::GpuArray< int, 3 > Ex_IndexType
 
amrex::GpuArray< int, 3 > Ey_IndexType
 
amrex::GpuArray< int, 3 > Ez_IndexType
 

Detailed Description

This class contains the parameters needed to evaluate hybrid field solutions (kinetic ions with fluid electrons).

Constructor & Destructor Documentation

◆ HybridPICModel()

HybridPICModel::HybridPICModel ( )

Member Function Documentation

◆ AllocateLevelMFs()

void HybridPICModel::AllocateLevelMFs ( ablastr::fields::MultiFabRegister & fields,
int lev,
const amrex::BoxArray & ba,
const amrex::DistributionMapping & dm,
int ncomps,
const amrex::IntVect & ngJ,
const amrex::IntVect & ngRho,
const amrex::IntVect & ngEB,
const amrex::IntVect & jx_nodal_flag,
const amrex::IntVect & jy_nodal_flag,
const amrex::IntVect & jz_nodal_flag,
const amrex::IntVect & rho_nodal_flag,
const amrex::IntVect & Ex_nodal_flag,
const amrex::IntVect & Ey_nodal_flag,
const amrex::IntVect & Ez_nodal_flag,
const amrex::IntVect & Bx_nodal_flag,
const amrex::IntVect & By_nodal_flag,
const amrex::IntVect & Bz_nodal_flag ) const

Allocate hybrid-PIC specific multifabs. Called in constructor.

◆ BfieldEvolveRK() [1/2]

void HybridPICModel::BfieldEvolveRK ( ablastr::fields::MultiLevelVectorField const & Bfield,
ablastr::fields::MultiLevelVectorField const & Efield,
ablastr::fields::MultiLevelVectorField const & Jfield,
ablastr::fields::MultiLevelScalarField const & rhofield,
amrex::Vector< std::array< std::unique_ptr< amrex::iMultiFab >, 3 > > & eb_update_E,
amrex::Real dt,
int lev,
SubcyclingHalf subcycling_half,
amrex::IntVect ng,
std::optional< bool > nodal_sync )

◆ BfieldEvolveRK() [2/2]

void HybridPICModel::BfieldEvolveRK ( ablastr::fields::MultiLevelVectorField const & Bfield,
ablastr::fields::MultiLevelVectorField const & Efield,
ablastr::fields::MultiLevelVectorField const & Jfield,
ablastr::fields::MultiLevelScalarField const & rhofield,
amrex::Vector< std::array< std::unique_ptr< amrex::iMultiFab >, 3 > > & eb_update_E,
amrex::Real dt,
SubcyclingHalf subcycling_half,
amrex::IntVect ng,
std::optional< bool > nodal_sync )

◆ CalculateElectronPressure() [1/2]

void HybridPICModel::CalculateElectronPressure ( ) const

Function to calculate the electron pressure using the simulation charge density. Used in the Ohm's law solver (kinetic-fluid hybrid model).

◆ CalculateElectronPressure() [2/2]

void HybridPICModel::CalculateElectronPressure ( int lev) const

◆ CalculatePlasmaCurrent() [1/2]

void HybridPICModel::CalculatePlasmaCurrent ( ablastr::fields::MultiLevelVectorField const & Bfield,
amrex::Vector< std::array< std::unique_ptr< amrex::iMultiFab >, 3 > > & eb_update_E )

Function to calculate the total plasma current based on Ampere's law while neglecting displacement current (J = curl x B). Any external current is subtracted as well. Used in the Ohm's law solver (kinetic-fluid hybrid model).

Parameters
[in]BfieldMagnetic field from which the current is calculated.
[in]eb_update_EIndicate in which cell J should be calculated (related to embedded boundaries).

◆ CalculatePlasmaCurrent() [2/2]

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

◆ FieldPush()

void HybridPICModel::FieldPush ( ablastr::fields::MultiLevelVectorField const & Bfield,
ablastr::fields::MultiLevelVectorField const & Efield,
ablastr::fields::MultiLevelVectorField const & Jfield,
ablastr::fields::MultiLevelScalarField const & rhofield,
amrex::Vector< std::array< std::unique_ptr< amrex::iMultiFab >, 3 > > & eb_update_E,
amrex::Real dt,
SubcyclingHalf subcycling_half,
amrex::IntVect ng,
std::optional< bool > nodal_sync )

◆ FillElectronPressureMF()

void HybridPICModel::FillElectronPressureMF ( amrex::MultiFab & Pe_field,
amrex::MultiFab const & rho_field ) const

Fill the electron pressure multifab given the kinetic particle charge density (and assumption of quasi-neutrality) using the user specified electron equation of state.

Parameters
[out]Pe_fieldscalar electron pressure MultiFab at a given level
[in]rho_fieldscalar ion charge density Multifab at a given level

◆ GetCurrentExternal()

void HybridPICModel::GetCurrentExternal ( )

Function to evaluate the external current expressions and populate the external current multifab. Note the external current can be a function of time and therefore this should be re-evaluated at every step.

◆ HybridPICSolveE() [1/3]

void HybridPICModel::HybridPICSolveE ( ablastr::fields::MultiLevelVectorField const & Efield,
ablastr::fields::MultiLevelVectorField const & Jfield,
ablastr::fields::MultiLevelVectorField const & Bfield,
ablastr::fields::MultiLevelScalarField const & rhofield,
amrex::Vector< std::array< std::unique_ptr< amrex::iMultiFab >, 3 > > & eb_update_E,
bool solve_for_Faraday ) const

Function to update the E-field using Ohm's law (hybrid-PIC model).

◆ HybridPICSolveE() [2/3]

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

◆ HybridPICSolveE() [3/3]

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

◆ InitData()

void HybridPICModel::InitData ( const ablastr::fields::MultiFabRegister & fields)

◆ ReadParameters()

void HybridPICModel::ReadParameters ( )

Read user-defined model parameters. Called in constructor.

Member Data Documentation

◆ Bx_IndexType

amrex::GpuArray<int, 3> HybridPICModel::Bx_IndexType

Gpu Vector with index type of the Bx multifab

◆ By_IndexType

amrex::GpuArray<int, 3> HybridPICModel::By_IndexType

Gpu Vector with index type of the By multifab

◆ Bz_IndexType

amrex::GpuArray<int, 3> HybridPICModel::Bz_IndexType

Gpu Vector with index type of the Bz multifab

◆ Ex_IndexType

amrex::GpuArray<int, 3> HybridPICModel::Ex_IndexType

Gpu Vector with index type of the Ex multifab

◆ Ey_IndexType

amrex::GpuArray<int, 3> HybridPICModel::Ey_IndexType

Gpu Vector with index type of the Ey multifab

◆ Ez_IndexType

amrex::GpuArray<int, 3> HybridPICModel::Ez_IndexType

Gpu Vector with index type of the Ez multifab

◆ Jx_IndexType

amrex::GpuArray<int, 3> HybridPICModel::Jx_IndexType

Gpu Vector with index type of the Jx multifab

◆ Jy_IndexType

amrex::GpuArray<int, 3> HybridPICModel::Jy_IndexType

Gpu Vector with index type of the Jy multifab

◆ Jz_IndexType

amrex::GpuArray<int, 3> HybridPICModel::Jz_IndexType

Gpu Vector with index type of the Jz multifab

◆ m_add_external_fields

bool HybridPICModel::m_add_external_fields = false

External E/B fields

◆ m_elec_temp

amrex::Real HybridPICModel::m_elec_temp

Electron temperature in eV

◆ m_eta

amrex::ParserExecutor<2> HybridPICModel::m_eta

◆ m_eta_expression

std::string HybridPICModel::m_eta_expression = "0.0"

Plasma resistivity

◆ m_eta_h

amrex::ParserExecutor<2> HybridPICModel::m_eta_h

◆ m_eta_h_expression

std::string HybridPICModel::m_eta_h_expression = "0.0"

Plasma hyper-resisitivity

◆ m_external_current_has_time_dependence

bool HybridPICModel::m_external_current_has_time_dependence = false

◆ m_external_vector_potential

std::unique_ptr<ExternalVectorPotential> HybridPICModel::m_external_vector_potential

◆ m_gamma

amrex::Real HybridPICModel::m_gamma = 5.0/3.0

Electron pressure scaling exponent

◆ m_holmstrom_vacuum_region

bool HybridPICModel::m_holmstrom_vacuum_region = false

◆ m_hyper_resistivity_has_B_dependence

bool HybridPICModel::m_hyper_resistivity_has_B_dependence = false

◆ m_hyper_resistivity_parser

std::unique_ptr<amrex::Parser> HybridPICModel::m_hyper_resistivity_parser

◆ m_include_hyper_resistivity_term

bool HybridPICModel::m_include_hyper_resistivity_term = false

◆ m_J_external

std::array< amrex::ParserExecutor<4>, 3> HybridPICModel::m_J_external

◆ m_J_external_parser

std::array< std::unique_ptr<amrex::Parser>, 3> HybridPICModel::m_J_external_parser

◆ m_Jx_ext_grid_function

std::string HybridPICModel::m_Jx_ext_grid_function = "0.0"

External current

◆ m_Jy_ext_grid_function

std::string HybridPICModel::m_Jy_ext_grid_function = "0.0"

◆ m_Jz_ext_grid_function

std::string HybridPICModel::m_Jz_ext_grid_function = "0.0"

◆ m_n0_ref

amrex::Real HybridPICModel::m_n0_ref = 1.0

Reference electron density

◆ m_n_floor

amrex::Real HybridPICModel::m_n_floor = 1.0

Plasma density floor - if n < n_floor it will be set to n_floor

◆ m_resistivity_has_J_dependence

bool HybridPICModel::m_resistivity_has_J_dependence = false

◆ m_resistivity_parser

std::unique_ptr<amrex::Parser> HybridPICModel::m_resistivity_parser

◆ m_substeps

int HybridPICModel::m_substeps = 10

Number of substeps to take when evolving B


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