Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Mesh refinement: initialize the plasma on the correct level #554

Merged
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
3 changes: 3 additions & 0 deletions docs/source/run/parameters.rst
Original file line number Diff line number Diff line change
Expand Up @@ -156,6 +156,9 @@ plasma parameters for each plasma are specified via `<plasma name>.plasma_proper
Since in a quasi-static code, there is only a 2D plasma slice evolving along the longitudinal
coordinate, there is no need to specify a number of particles per cell in z.

* ``<plasma name>.level`` (`integer`) optional (default `0`)
Level of mesh refinement to which the plasma is assigned.

* ``<plasma name>.radius`` (`float`) optional (default `infinity`)
Radius of the plasma. Set a value to run simulations in a plasma column.

Expand Down
2 changes: 1 addition & 1 deletion src/Hipace.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -236,7 +236,7 @@ Hipace::InitData ()
AmrCore::InitFromScratch(0.0); // function argument is time
constexpr int lev = 0;
m_multi_beam.InitData(geom[lev]);
m_multi_plasma.InitData(lev, m_slice_ba[lev], m_slice_dm[lev], m_slice_geom[lev], geom[lev]);
m_multi_plasma.InitData(m_slice_ba, m_slice_dm, m_slice_geom, geom);
m_adaptive_time_step.Calculate(m_dt, m_multi_beam, m_multi_plasma.maxDensity());
#ifdef AMREX_USE_MPI
m_adaptive_time_step.WaitTimeStep(m_dt, m_comm_z);
Expand Down
6 changes: 3 additions & 3 deletions src/particles/MultiPlasma.H
Original file line number Diff line number Diff line change
Expand Up @@ -19,14 +19,14 @@ public:

/** \brief Loop over plasma species and initialize them.
*
* \param[in] lev MR level
* \param[in] slice_ba slice boxarray, on which plasma particles are defined
* \param[in] slice_dm DistributionMapping of the transverse slice domain
* \param[in] slice_gm slice geometry
* \param[in] gm Geometry of the simulation, to get the cell size
*/
void InitData (int lev, amrex::BoxArray slice_ba, amrex::DistributionMapping slice_dm,
amrex::Geometry slice_gm, amrex::Geometry gm);
void InitData (amrex::Vector<amrex::BoxArray> slice_ba,
amrex::Vector<amrex::DistributionMapping> slice_dm,
amrex::Vector<amrex::Geometry> slice_gm, amrex::Vector<amrex::Geometry> gm);


/** Loop over plasma species and depose their currents into the current 2D slice in fields
Expand Down
15 changes: 8 additions & 7 deletions src/particles/MultiPlasma.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -16,14 +16,15 @@ MultiPlasma::MultiPlasma (amrex::AmrCore* amr_core)
}

void
MultiPlasma::InitData (int lev, amrex::BoxArray slice_ba,
amrex::DistributionMapping slice_dm, amrex::Geometry slice_gm,
amrex::Geometry gm)
MultiPlasma::InitData (amrex::Vector<amrex::BoxArray> slice_ba,
amrex::Vector<amrex::DistributionMapping> slice_dm,
amrex::Vector<amrex::Geometry> slice_gm, amrex::Vector<amrex::Geometry> gm)
{
for (auto& plasma : m_all_plasmas) {
plasma.SetParticleBoxArray(lev, slice_ba);
plasma.SetParticleDistributionMap(lev, slice_dm);
plasma.SetParticleGeometry(lev, slice_gm);
const int lev = plasma.m_level;
plasma.SetParticleBoxArray(lev, slice_ba[lev]);
plasma.SetParticleDistributionMap(lev, slice_dm[lev]);
plasma.SetParticleGeometry(lev, slice_gm[lev]);
plasma.InitData();

if(plasma.m_can_ionize) {
Expand All @@ -35,7 +36,7 @@ MultiPlasma::InitData (int lev, amrex::BoxArray slice_ba,
}
AMREX_ALWAYS_ASSERT_WITH_MESSAGE(plasma_product != nullptr,
"Must specify a valid product plasma for Ionization using ionization_product");
plasma.InitIonizationModule(gm, plasma_product);
plasma.InitIonizationModule(gm[lev], plasma_product);
}
}
}
Expand Down
2 changes: 1 addition & 1 deletion src/particles/PlasmaParticleContainer.H
Original file line number Diff line number Diff line change
Expand Up @@ -90,7 +90,7 @@ public:
Fields& fields);

amrex::Real m_density {0}; /**< Density of the plasma */
int m_MR_level {0}; /**< mesh refinement level on which the plasma lives */
int m_level {0}; /**< mesh refinement level on which the plasma lives */
/** maximum weighting factor gamma/(Psi +1) before particle is regarded as violating
* the quasi-static approximation and is removed */
amrex::Real m_max_qsa_weighting_factor {35.};
Expand Down
2 changes: 1 addition & 1 deletion src/particles/PlasmaParticleContainer.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -70,7 +70,7 @@ PlasmaParticleContainer::ReadParameters ()
}
pp.query("ionization_product", m_product_name);
pp.query("density", m_density);
pp.query("MR_level", m_MR_level);
pp.query("level", m_level);
pp.query("radius", m_radius);
pp.query("hollow_core_radius", m_hollow_core_radius);
AMREX_ALWAYS_ASSERT_WITH_MESSAGE(m_hollow_core_radius < m_radius,
Expand Down
2 changes: 1 addition & 1 deletion src/particles/PlasmaParticleContainerInit.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -17,7 +17,7 @@ InitParticles (const amrex::IntVect& a_num_particles_per_cell,
{
HIPACE_PROFILE("PlasmaParticleContainer::InitParticles");

const int lev = 0;
const int lev = m_level;
const auto dx = ParticleGeom(lev).CellSizeArray();
const auto plo = ParticleGeom(lev).ProbLoArray();
const amrex::RealBox a_bounds = ParticleGeom(lev).ProbDomain();
Expand Down
2 changes: 1 addition & 1 deletion src/particles/deposition/PlasmaDepositCurrent.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -22,7 +22,7 @@ DepositCurrent (PlasmaParticleContainer& plasma, Fields & fields,
" (WhichSlice::Next) or for the ion charge deposition (WhichSLice::RhoIons)");

// only deposit plasma currents on their according MR level
if (plasma.m_MR_level != lev) return;
if (plasma.m_level != lev) return;

// Extract properties associated with physical size of the box
amrex::Real const * AMREX_RESTRICT dx = gm.CellSize();
Expand Down
2 changes: 1 addition & 1 deletion src/particles/pusher/PlasmaParticleAdvance.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -19,7 +19,7 @@ AdvancePlasmaParticles (PlasmaParticleContainer& plasma, Fields & fields,
using namespace amrex::literals;

// only push plasma particles on their according MR level
if (plasma.m_MR_level != lev) return;
if (plasma.m_level != lev) return;

// Extract properties associated with physical size of the box
amrex::Real const * AMREX_RESTRICT dx = gm.CellSize();
Expand Down