WarpX
Public Member Functions | Public Attributes | List of all members
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 (int nlevs_max)
 
void ReadParameters ()
 
void AllocateMFs (int nlevs_max)
 
void AllocateLevelMFs (int lev, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm, int ncomps, const amrex::IntVect &ngJ, const amrex::IntVect &ngRho, const amrex::IntVect &jx_nodal_flag, const amrex::IntVect &jy_nodal_flag, const amrex::IntVect &jz_nodal_flag, const amrex::IntVect &rho_nodal_flag)
 
void ClearLevel (int lev)
 
void InitData ()
 
void GetCurrentExternal (amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &edge_lengths)
 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. More...
 
void GetCurrentExternal (std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &edge_lengths, int lev)
 
void CalculateCurrentAmpere (amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &Bfield, amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &edge_lengths)
 Function to calculate the total current based on Ampere's law while neglecting displacement current (J = curl x B). Used in the Ohm's law solver (kinetic-fluid hybrid model). More...
 
void CalculateCurrentAmpere (std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Bfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &edge_lengths, int lev)
 
void HybridPICSolveE (amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &Efield, amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &Jfield, amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &Bfield, amrex::Vector< std::unique_ptr< amrex::MultiFab >> const &rhofield, amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &edge_lengths, bool include_resistivity_term)
 Function to update the E-field using Ohm's law (hybrid-PIC model). More...
 
void HybridPICSolveE (std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Efield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Jfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Bfield, std::unique_ptr< amrex::MultiFab > const &rhofield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &edge_lengths, int lev, bool include_resistivity_term)
 
void HybridPICSolveE (std::array< std::unique_ptr< amrex::MultiFab >, 3 > &Efield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Jfield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &Bfield, std::unique_ptr< amrex::MultiFab > const &rhofield, std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &edge_lengths, int lev, PatchType patch_type, bool include_resistivity_term)
 
void BfieldEvolveRK (amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &Bfield, amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &Efield, amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &Jfield, amrex::Vector< std::unique_ptr< amrex::MultiFab >> const &rhofield, amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &edge_lengths, amrex::Real dt, DtType a_dt_type, amrex::IntVect ng, std::optional< bool > nodal_sync)
 
void BfieldEvolveRK (amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &Bfield, amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &Efield, amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &Jfield, amrex::Vector< std::unique_ptr< amrex::MultiFab >> const &rhofield, amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &edge_lengths, amrex::Real dt, int lev, DtType dt_type, amrex::IntVect ng, std::optional< bool > nodal_sync)
 
void FieldPush (amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &Bfield, amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &Efield, amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &Jfield, amrex::Vector< std::unique_ptr< amrex::MultiFab >> const &rhofield, amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &edge_lengths, amrex::Real dt, DtType dt_type, amrex::IntVect ng, std::optional< bool > nodal_sync)
 
void CalculateElectronPressure (DtType a_dt_type)
 Function to calculate the electron pressure at a given timestep type using the simulation charge density. Used in the Ohm's law solver (kinetic-fluid hybrid model). More...
 
void CalculateElectronPressure (int lev, DtType a_dt_type)
 
void FillElectronPressureMF (std::unique_ptr< amrex::MultiFab > const &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. More...
 
amrex::MultiFabget_pointer_current_fp_ampere (int lev, int direction) const
 
amrex::MultiFabget_pointer_current_fp_external (int lev, int direction) const
 
amrex::MultiFabget_pointer_electron_pressure_fp (int lev) const
 

Public Attributes

int m_substeps = 10
 
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
 
amrex::Real m_eta_h = 0.0
 
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_field_has_time_dependence = false
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > rho_fp_temp
 
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_fp_temp
 
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_fp_ampere
 
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > current_fp_external
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > electron_pressure_fp
 
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 ( int  nlevs_max)

Member Function Documentation

◆ AllocateLevelMFs()

void HybridPICModel::AllocateLevelMFs ( int  lev,
const amrex::BoxArray ba,
const amrex::DistributionMapping dm,
int  ncomps,
const amrex::IntVect ngJ,
const amrex::IntVect ngRho,
const amrex::IntVect jx_nodal_flag,
const amrex::IntVect jy_nodal_flag,
const amrex::IntVect jz_nodal_flag,
const amrex::IntVect rho_nodal_flag 
)

◆ AllocateMFs()

void HybridPICModel::AllocateMFs ( int  nlevs_max)

Allocate hybrid-PIC specific multifabs. Called in constructor.

◆ BfieldEvolveRK() [1/2]

void HybridPICModel::BfieldEvolveRK ( amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &  Bfield,
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &  Efield,
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &  Jfield,
amrex::Vector< std::unique_ptr< amrex::MultiFab >> const &  rhofield,
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &  edge_lengths,
amrex::Real  dt,
DtType  a_dt_type,
amrex::IntVect  ng,
std::optional< bool >  nodal_sync 
)

◆ BfieldEvolveRK() [2/2]

void HybridPICModel::BfieldEvolveRK ( amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &  Bfield,
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &  Efield,
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &  Jfield,
amrex::Vector< std::unique_ptr< amrex::MultiFab >> const &  rhofield,
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &  edge_lengths,
amrex::Real  dt,
int  lev,
DtType  dt_type,
amrex::IntVect  ng,
std::optional< bool >  nodal_sync 
)

◆ CalculateCurrentAmpere() [1/2]

void HybridPICModel::CalculateCurrentAmpere ( amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &  Bfield,
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &  edge_lengths 
)

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

Parameters
[in]BfieldMagnetic field from which the current is calculated.
[in]edge_lengthsLength of cell edges taking embedded boundaries into account

◆ CalculateCurrentAmpere() [2/2]

void HybridPICModel::CalculateCurrentAmpere ( std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &  Bfield,
std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &  edge_lengths,
int  lev 
)

◆ CalculateElectronPressure() [1/2]

void HybridPICModel::CalculateElectronPressure ( DtType  a_dt_type)

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

◆ CalculateElectronPressure() [2/2]

void HybridPICModel::CalculateElectronPressure ( int  lev,
DtType  a_dt_type 
)

◆ ClearLevel()

void HybridPICModel::ClearLevel ( int  lev)

Helper function to clear values from hybrid-PIC specific multifabs.

◆ FieldPush()

void HybridPICModel::FieldPush ( amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &  Bfield,
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &  Efield,
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &  Jfield,
amrex::Vector< std::unique_ptr< amrex::MultiFab >> const &  rhofield,
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &  edge_lengths,
amrex::Real  dt,
DtType  dt_type,
amrex::IntVect  ng,
std::optional< bool >  nodal_sync 
)

◆ FillElectronPressureMF()

void HybridPICModel::FillElectronPressureMF ( std::unique_ptr< amrex::MultiFab > const &  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_filedscalar electron pressure MultiFab at a given level
[in]rho_fieldscalar ion charge density Multifab at a given level

◆ get_pointer_current_fp_ampere()

amrex::MultiFab* HybridPICModel::get_pointer_current_fp_ampere ( int  lev,
int  direction 
) const
inline

◆ get_pointer_current_fp_external()

amrex::MultiFab* HybridPICModel::get_pointer_current_fp_external ( int  lev,
int  direction 
) const
inline

◆ get_pointer_electron_pressure_fp()

amrex::MultiFab* HybridPICModel::get_pointer_electron_pressure_fp ( int  lev) const
inline

◆ GetCurrentExternal() [1/2]

void HybridPICModel::GetCurrentExternal ( amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &  edge_lengths)

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.

◆ GetCurrentExternal() [2/2]

void HybridPICModel::GetCurrentExternal ( std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &  edge_lengths,
int  lev 
)

◆ HybridPICSolveE() [1/3]

void HybridPICModel::HybridPICSolveE ( amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> &  Efield,
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &  Jfield,
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &  Bfield,
amrex::Vector< std::unique_ptr< amrex::MultiFab >> const &  rhofield,
amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 >> const &  edge_lengths,
bool  include_resistivity_term 
)

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

◆ HybridPICSolveE() [2/3]

void HybridPICModel::HybridPICSolveE ( std::array< std::unique_ptr< amrex::MultiFab >, 3 > &  Efield,
std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &  Jfield,
std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &  Bfield,
std::unique_ptr< amrex::MultiFab > const &  rhofield,
std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &  edge_lengths,
int  lev,
bool  include_resistivity_term 
)

◆ HybridPICSolveE() [3/3]

void HybridPICModel::HybridPICSolveE ( std::array< std::unique_ptr< amrex::MultiFab >, 3 > &  Efield,
std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &  Jfield,
std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &  Bfield,
std::unique_ptr< amrex::MultiFab > const &  rhofield,
std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &  edge_lengths,
int  lev,
PatchType  patch_type,
bool  include_resistivity_term 
)

◆ InitData()

void HybridPICModel::InitData ( )

◆ 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

◆ current_fp_ampere

amrex::Vector<std::array< std::unique_ptr<amrex::MultiFab>, 3 > > HybridPICModel::current_fp_ampere

◆ current_fp_external

amrex::Vector<std::array< std::unique_ptr<amrex::MultiFab>, 3 > > HybridPICModel::current_fp_external

◆ current_fp_temp

amrex::Vector<std::array< std::unique_ptr<amrex::MultiFab>, 3 > > HybridPICModel::current_fp_temp

◆ electron_pressure_fp

amrex::Vector< std::unique_ptr<amrex::MultiFab> > HybridPICModel::electron_pressure_fp

◆ 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_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::Real HybridPICModel::m_eta_h = 0.0

Plasma hyper-resisitivity

◆ m_external_field_has_time_dependence

bool HybridPICModel::m_external_field_has_time_dependence = false

◆ m_gamma

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

Electron pressure scaling exponent

◆ 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

◆ rho_fp_temp

amrex::Vector< std::unique_ptr<amrex::MultiFab> > HybridPICModel::rho_fp_temp

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