WarpX
Loading...
Searching...
No Matches
HybridPICModel.H
Go to the documentation of this file.
1/* Copyright 2023-2024 The WarpX Community
2 *
3 * This file is part of WarpX.
4 *
5 * Authors: Roelof Groenewald (TAE Technologies)
6 * S. Eric Clark (Helion Energy)
7 *
8 * License: BSD-3-Clause-LBNL
9 */
10
11#ifndef WARPX_HYBRIDPICMODEL_H_
12#define WARPX_HYBRIDPICMODEL_H_
13
14#include "HybridPICModel_fwd.H"
15
16#include "Fields.H"
17
20
23#include "Utils/WarpXConst.H"
24
27
28#include <AMReX_Array.H>
29#include <AMReX_REAL.H>
30#include <AMReX_BoxArray.H>
31#include <AMReX_IntVect.H>
33
34#include <optional>
35
41{
42public:
44
46 void ReadParameters ();
47
49 void AllocateLevelMFs (
51 int lev,
52 const amrex::BoxArray& ba,
54 int ncomps,
55 const amrex::IntVect& ngJ,
56 const amrex::IntVect& ngRho,
57 const amrex::IntVect& ngEB,
58 const amrex::IntVect& jx_nodal_flag,
59 const amrex::IntVect& jy_nodal_flag,
60 const amrex::IntVect& jz_nodal_flag,
61 const amrex::IntVect& rho_nodal_flag,
62 const amrex::IntVect& Ex_nodal_flag,
63 const amrex::IntVect& Ey_nodal_flag,
64 const amrex::IntVect& Ez_nodal_flag,
65 const amrex::IntVect& Bx_nodal_flag,
66 const amrex::IntVect& By_nodal_flag,
67 const amrex::IntVect& Bz_nodal_flag
68 ) const;
69
71
78 void GetCurrentExternal ();
79
91 amrex::Vector<std::array< std::unique_ptr<amrex::iMultiFab>,3 > >& eb_update_E
92 ) const;
94 ablastr::fields::VectorField const& Bfield,
95 std::array< std::unique_ptr<amrex::iMultiFab>,3 >& eb_update_E,
96 int lev
97 ) const;
98
103 void HybridPICSolveE (
108 amrex::Vector<std::array< std::unique_ptr<amrex::iMultiFab>,3 > >& eb_update_E,
109 bool solve_for_Faraday) const;
110
111 void HybridPICSolveE (
112 ablastr::fields::VectorField const& Efield,
113 ablastr::fields::VectorField const& Jfield,
114 ablastr::fields::VectorField const& Bfield,
115 amrex::MultiFab const& rhofield,
116 std::array< std::unique_ptr<amrex::iMultiFab>,3 >& eb_update_E,
117 int lev, bool solve_for_Faraday) const;
118
119 void HybridPICSolveE (
120 ablastr::fields::VectorField const& Efield,
121 ablastr::fields::VectorField const& Jfield,
122 ablastr::fields::VectorField const& Bfield,
123 amrex::MultiFab const& rhofield,
124 std::array< std::unique_ptr<amrex::iMultiFab>,3 >& eb_update_E,
125 int lev, PatchType patch_type, bool solve_for_Faraday) const;
126
127 void BfieldEvolveRK (
132 amrex::Vector<std::array< std::unique_ptr<amrex::iMultiFab>,3 > >& eb_update_E,
133 amrex::Real dt, SubcyclingHalf subcycling_half,
134 amrex::IntVect ng, std::optional<bool> nodal_sync);
135
136 void BfieldEvolveRK (
141 amrex::Vector<std::array< std::unique_ptr<amrex::iMultiFab>,3 > >& eb_update_E,
142 amrex::Real dt, int lev, SubcyclingHalf subcycling_half,
143 amrex::IntVect ng, std::optional<bool> nodal_sync);
144
145 void BfieldEvolveRKF45 (
150 amrex::Vector<std::array< std::unique_ptr<amrex::iMultiFab>,3 > >& eb_update_E,
151 amrex::Real dt_half, SubcyclingHalf subcycling_half,
152 amrex::IntVect ng, std::optional<bool> nodal_sync);
153
154 void BfieldEvolveRKF45 (
159 amrex::Vector<std::array< std::unique_ptr<amrex::iMultiFab>,3 > >& eb_update_E,
160 amrex::Real dt_half, int lev, SubcyclingHalf subcycling_half,
161 amrex::IntVect ng, std::optional<bool> nodal_sync);
162
163 void FieldPush (
168 amrex::Vector<std::array< std::unique_ptr<amrex::iMultiFab>,3 > >& eb_update_E,
169 amrex::Real dt, SubcyclingHalf subcycling_half,
170 amrex::IntVect ng, std::optional<bool> nodal_sync);
171
177 void CalculateElectronPressure () const;
178 void CalculateElectronPressure (int lev) const;
179
189 amrex::MultiFab& Pe_field,
190 amrex::MultiFab const& rho_field ) const;
191
192 // Declare variables to hold hybrid-PIC model parameters
195 int m_substeps = 10;
196
208 bool m_use_rkf45 = false;
209
211
218
221
223 std::string m_eta_expression = "0.0";
224 std::unique_ptr<amrex::Parser> m_resistivity_parser;
227
229 std::string m_eta_h_expression = "0.0";
230 std::unique_ptr<amrex::Parser> m_hyper_resistivity_parser;
234
236 std::string m_Jx_ext_grid_function = "0.0";
237 std::string m_Jy_ext_grid_function = "0.0";
238 std::string m_Jz_ext_grid_function = "0.0";
239 std::array< std::unique_ptr<amrex::Parser>, 3> m_J_external_parser;
240 std::array< amrex::ParserExecutor<4>, 3> m_J_external;
243
246 std::unique_ptr<ExternalVectorPotential> m_external_vector_potential;
247
266};
267
275
278 amrex::Real const T0,
279 amrex::Real const gamma,
280 amrex::Real const rho) {
281 return n0 * T0 * std::pow((rho/PhysConst::q_e)/n0, gamma);
282 }
283};
284
285#endif // WARPX_HYBRIDPICMODEL_H_
#define AMREX_FORCE_INLINE
#define AMREX_GPU_HOST_DEVICE
SubcyclingHalf
Subcycling half selector.
Definition WarpXAlgorithmSelection.H:166
amrex::ParserExecutor< 2 > m_eta_h
Definition HybridPICModel.H:231
void ReadParameters()
Definition HybridPICModel.cpp:31
amrex::GpuArray< int, 3 > Bz_IndexType
Definition HybridPICModel.H:259
amrex::GpuArray< int, 3 > Bx_IndexType
Definition HybridPICModel.H:255
amrex::Real m_n0_ref
Definition HybridPICModel.H:215
std::string m_Jz_ext_grid_function
Definition HybridPICModel.H:238
amrex::GpuArray< int, 3 > Ez_IndexType
Definition HybridPICModel.H:265
bool m_external_current_has_time_dependence
Definition HybridPICModel.H:242
amrex::Real m_substep_safety
Definition HybridPICModel.H:202
amrex::ParserExecutor< 3 > m_eta
Definition HybridPICModel.H:225
void GetCurrentExternal()
Function to evaluate the external current expressions and populate the external current multifab....
Definition HybridPICModel.cpp:290
std::array< amrex::ParserExecutor< 4 >, 3 > m_J_external
Definition HybridPICModel.H:240
void BfieldEvolveRKF45(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_half, SubcyclingHalf subcycling_half, amrex::IntVect ng, std::optional< bool > nodal_sync)
Definition HybridPICModel.cpp:664
amrex::Real m_elec_temp
Definition HybridPICModel.H:213
HybridPICModel()
Definition HybridPICModel.cpp:26
std::unique_ptr< ExternalVectorPotential > m_external_vector_potential
Definition HybridPICModel.H:246
void InitData(const ablastr::fields::MultiFabRegister &fields)
Definition HybridPICModel.cpp:184
amrex::GpuArray< int, 3 > Jx_IndexType
Definition HybridPICModel.H:249
bool m_add_external_fields
Definition HybridPICModel.H:245
amrex::GpuArray< int, 3 > Jz_IndexType
Definition HybridPICModel.H:253
std::string m_Jy_ext_grid_function
Definition HybridPICModel.H:237
bool m_use_rkf45
Definition HybridPICModel.H:208
std::unique_ptr< amrex::Parser > m_hyper_resistivity_parser
Definition HybridPICModel.H:230
bool m_hyper_resistivity_has_B_dependence
Definition HybridPICModel.H:233
std::unique_ptr< amrex::Parser > m_resistivity_parser
Definition HybridPICModel.H:224
std::string m_eta_expression
Definition HybridPICModel.H:223
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
Definition HybridPICModel.cpp:98
bool m_holmstrom_vacuum_region
Definition HybridPICModel.H:210
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 quas...
Definition HybridPICModel.cpp:436
amrex::Real m_substep_atol
Definition HybridPICModel.H:200
amrex::Real m_gamma
Definition HybridPICModel.H:217
std::array< std::unique_ptr< amrex::Parser >, 3 > m_J_external_parser
Definition HybridPICModel.H:239
amrex::Real m_substep_max_growth
Definition HybridPICModel.H:204
std::string m_Jx_ext_grid_function
Definition HybridPICModel.H:236
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)
Definition HybridPICModel.cpp:466
std::string m_eta_h_expression
Definition HybridPICModel.H:229
bool m_include_hyper_resistivity_term
Definition HybridPICModel.H:232
int m_max_substep_attempts
Definition HybridPICModel.H:206
int m_substeps
Definition HybridPICModel.H:195
bool m_resistivity_has_J_dependence
Definition HybridPICModel.H:226
amrex::GpuArray< int, 3 > Ex_IndexType
Definition HybridPICModel.H:261
amrex::Real m_n_floor
Definition HybridPICModel.H:220
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).
Definition HybridPICModel.cpp:342
amrex::Real m_substep_rtol
Definition HybridPICModel.H:198
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)
Definition HybridPICModel.cpp:1092
void CalculatePlasmaCurrent(ablastr::fields::MultiLevelVectorField const &Bfield, amrex::Vector< std::array< std::unique_ptr< amrex::iMultiFab >, 3 > > &eb_update_E) const
Function to calculate the total plasma current based on Ampere's law while neglecting displacement cu...
Definition HybridPICModel.cpp:307
bool m_has_external_current
Definition HybridPICModel.H:241
void CalculateElectronPressure() const
Function to calculate the electron pressure using the simulation charge density. Used in the Ohm's la...
Definition HybridPICModel.cpp:406
amrex::GpuArray< int, 3 > Jy_IndexType
Definition HybridPICModel.H:251
amrex::GpuArray< int, 3 > By_IndexType
Definition HybridPICModel.H:257
amrex::GpuArray< int, 3 > Ey_IndexType
Definition HybridPICModel.H:263
amrex_real Real
constexpr auto q_e
elementary charge [C]
Definition constant.H:169
Definition EffectivePotentialPoissonSolver.H:63
std::array< amrex::MultiFab *, 3 > VectorField
Definition MultiFabRegister.H:201
amrex::Vector< ScalarField > MultiLevelScalarField
Definition MultiFabRegister.H:210
amrex::Vector< VectorField > MultiLevelVectorField
Definition MultiFabRegister.H:218
PatchType
Definition Enums.H:30
IntVectND< 3 > IntVect
This struct contains only static functions to compute the electron pressure using the particle densit...
Definition HybridPICModel.H:274
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real get_pressure(amrex::Real const n0, amrex::Real const T0, amrex::Real const gamma, amrex::Real const rho)
Definition HybridPICModel.H:277
Definition MultiFabRegister.H:272