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"
25
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 );
94 ablastr::fields::VectorField const& Bfield,
95 std::array< std::unique_ptr<amrex::iMultiFab>,3 >& eb_update_E,
96 int lev
97 );
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 FieldPush (
150 amrex::Vector<std::array< std::unique_ptr<amrex::iMultiFab>,3 > >& eb_update_E,
151 amrex::Real dt, SubcyclingHalf subcycling_half,
152 amrex::IntVect ng, std::optional<bool> nodal_sync);
153
159 void CalculateElectronPressure () const;
160 void CalculateElectronPressure (int lev) const;
161
171 amrex::MultiFab& Pe_field,
172 amrex::MultiFab const& rho_field ) const;
173
174 // Declare variables to hold hybrid-PIC model parameters
176 int m_substeps = 10;
177
179
181 amrex::Real m_elec_temp;
183 amrex::Real m_n0_ref = 1.0;
185 amrex::Real m_gamma = 5.0/3.0;
186
188 amrex::Real m_n_floor = 1.0;
189
191 std::string m_eta_expression = "0.0";
192 std::unique_ptr<amrex::Parser> m_resistivity_parser;
195
197 std::string m_eta_h_expression = "0.0";
198 std::unique_ptr<amrex::Parser> m_hyper_resistivity_parser;
202
204 std::string m_Jx_ext_grid_function = "0.0";
205 std::string m_Jy_ext_grid_function = "0.0";
206 std::string m_Jz_ext_grid_function = "0.0";
207 std::array< std::unique_ptr<amrex::Parser>, 3> m_J_external_parser;
208 std::array< amrex::ParserExecutor<4>, 3> m_J_external;
210
213 std::unique_ptr<ExternalVectorPotential> m_external_vector_potential;
214
233};
234
242
244 static amrex::Real get_pressure (amrex::Real const n0,
245 amrex::Real const T0,
246 amrex::Real const gamma,
247 amrex::Real const rho) {
248 return n0 * T0 * std::pow((rho/PhysConst::q_e)/n0, gamma);
249 }
250};
251
252#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:199
void ReadParameters()
Definition HybridPICModel.cpp:30
amrex::GpuArray< int, 3 > Bz_IndexType
Definition HybridPICModel.H:226
amrex::GpuArray< int, 3 > Bx_IndexType
Definition HybridPICModel.H:222
amrex::Real m_n0_ref
Definition HybridPICModel.H:183
std::string m_Jz_ext_grid_function
Definition HybridPICModel.H:206
amrex::GpuArray< int, 3 > Ez_IndexType
Definition HybridPICModel.H:232
bool m_external_current_has_time_dependence
Definition HybridPICModel.H:209
void GetCurrentExternal()
Function to evaluate the external current expressions and populate the external current multifab....
Definition HybridPICModel.cpp:259
std::array< amrex::ParserExecutor< 4 >, 3 > m_J_external
Definition HybridPICModel.H:208
amrex::Real m_elec_temp
Definition HybridPICModel.H:181
HybridPICModel()
Definition HybridPICModel.cpp:25
std::unique_ptr< ExternalVectorPotential > m_external_vector_potential
Definition HybridPICModel.H:213
void InitData(const ablastr::fields::MultiFabRegister &fields)
Definition HybridPICModel.cpp:157
amrex::GpuArray< int, 3 > Jx_IndexType
Definition HybridPICModel.H:216
bool m_add_external_fields
Definition HybridPICModel.H:212
amrex::GpuArray< int, 3 > Jz_IndexType
Definition HybridPICModel.H:220
std::string m_Jy_ext_grid_function
Definition HybridPICModel.H:205
std::unique_ptr< amrex::Parser > m_hyper_resistivity_parser
Definition HybridPICModel.H:198
bool m_hyper_resistivity_has_B_dependence
Definition HybridPICModel.H:201
std::unique_ptr< amrex::Parser > m_resistivity_parser
Definition HybridPICModel.H:192
std::string m_eta_expression
Definition HybridPICModel.H:191
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:73
bool m_holmstrom_vacuum_region
Definition HybridPICModel.H:178
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:410
amrex::Real m_gamma
Definition HybridPICModel.H:185
std::array< std::unique_ptr< amrex::Parser >, 3 > m_J_external_parser
Definition HybridPICModel.H:207
std::string m_Jx_ext_grid_function
Definition HybridPICModel.H:204
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:440
std::string m_eta_h_expression
Definition HybridPICModel.H:197
bool m_include_hyper_resistivity_term
Definition HybridPICModel.H:200
int m_substeps
Definition HybridPICModel.H:176
bool m_resistivity_has_J_dependence
Definition HybridPICModel.H:194
amrex::GpuArray< int, 3 > Ex_IndexType
Definition HybridPICModel.H:228
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 cu...
Definition HybridPICModel.cpp:276
amrex::Real m_n_floor
Definition HybridPICModel.H:188
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:316
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:573
amrex::ParserExecutor< 2 > m_eta
Definition HybridPICModel.H:193
void CalculateElectronPressure() const
Function to calculate the electron pressure using the simulation charge density. Used in the Ohm's la...
Definition HybridPICModel.cpp:380
amrex::GpuArray< int, 3 > Jy_IndexType
Definition HybridPICModel.H:218
amrex::GpuArray< int, 3 > By_IndexType
Definition HybridPICModel.H:224
amrex::GpuArray< int, 3 > Ey_IndexType
Definition HybridPICModel.H:230
constexpr auto q_e
elementary charge [C]
Definition constant.H:158
Definition EffectivePotentialPoissonSolver.H:63
std::array< amrex::MultiFab *, 3 > VectorField
Definition MultiFabRegister.H:191
amrex::Vector< ScalarField > MultiLevelScalarField
Definition MultiFabRegister.H:200
amrex::Vector< VectorField > MultiLevelVectorField
Definition MultiFabRegister.H:208
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:241
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:244
Definition MultiFabRegister.H:262