aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
-rw-r--r--Docs/source/usage/parameters.rst13
-rw-r--r--Examples/Modules/boosted_diags/inputs_3d_slice15
-rw-r--r--Source/Diagnostics/Diagnostics.cpp7
-rw-r--r--Source/Initialization/CustomDensityProb.H4
-rw-r--r--Source/Initialization/InjectorDensity.cpp2
-rw-r--r--Source/Initialization/PlasmaInjector.cpp18
-rw-r--r--Source/Initialization/WarpXInitData.cpp4
-rw-r--r--Source/Particles/LaserParticleContainer.cpp10
-rw-r--r--Source/Particles/MultiParticleContainer.cpp4
-rw-r--r--Source/Utils/WarpXUtil.H6
-rw-r--r--Source/Utils/WarpXUtil.cpp82
-rw-r--r--Source/WarpX.cpp12
-rw-r--r--Source/main.cpp2
13 files changed, 132 insertions, 47 deletions
diff --git a/Docs/source/usage/parameters.rst b/Docs/source/usage/parameters.rst
index c7ede31c8..8b57ed75e 100644
--- a/Docs/source/usage/parameters.rst
+++ b/Docs/source/usage/parameters.rst
@@ -4,7 +4,7 @@ Input Parameters
================
.. note::
- The WarpXParser (see :ref:`running-cpp-parameters-parser`) is used for the right-hand-side of all input parameters that consist in a single real number, so expressions like ``<species_name>.density_max = "2.+1."`` and/or using user-defined constants are accepted. See below for more detail.
+ The WarpXParser (see :ref:`running-cpp-parameters-parser`) is used for the right-hand-side of all input parameters that consist of one or more floats, so expressions like ``<species_name>.density_max = "2.+1."`` and/or using user-defined constants are accepted. See below for more detail.
.. _running-cpp-parameters-overall:
@@ -130,7 +130,7 @@ Setting up the field mesh
* ``geometry.coord_sys`` (`integer`) optional (default `0`)
Coordinate system used by the simulation. 0 for Cartesian, 1 for cylindrical.
-* ``geometry.prob_lo`` and ``geometry.prob_hi`` (`2 floats in 2D`, `3 integers in 3D`; in meters)
+* ``geometry.prob_lo`` and ``geometry.prob_hi`` (`2 floats in 2D`, `3 floats in 3D`; in meters)
The extent of the full simulation box. This box is rectangular, and thus its
extent is given here by the coordinates of the lower corner (``geometry.prob_lo``) and
upper corner (``geometry.prob_hi``). The first axis of the coordinates is x
@@ -146,7 +146,7 @@ Setting up the field mesh
The speed of moving window, in units of the speed of light
(i.e. use ``1.0`` for a moving window that moves exactly at the speed of light)
-* ``warpx.fine_tag_lo`` and ``warpx.fine_tag_hi`` (`2 floats in 2D`, `3 integers in 3D`; in meters) optional
+* ``warpx.fine_tag_lo`` and ``warpx.fine_tag_hi`` (`2 floats in 2D`, `3 floats in 3D`; in meters) optional
**When using static mesh refinement with 1 level**, the extent of the refined patch.
This patch is rectangular, and thus its extent is given here by the coordinates
of the lower corner (``warpx.fine_tag_lo``) and upper corner (``warpx.fine_tag_hi``).
@@ -330,12 +330,13 @@ Math parser and user-defined constants
--------------------------------------
WarpX provides a math parser that reads expressions in the input file.
-It can be used in all input parameters that consist in one real number.
+It can be used in all input parameters that consist of one or more floats.
+Note that when multiple floats are expected, the expressions are space delimited.
WarpX constants
^^^^^^^^^^^^^^^
-WarpX provides a few pre-defined constants, that can be used for any parameter that consists in one real number.
+WarpX provides a few pre-defined constants, that can be used for any parameter that consists of one or more floats.
======== ===================
q_e elementary charge
@@ -352,7 +353,7 @@ User-defined constants
^^^^^^^^^^^^^^^^^^^^^^
Users can define their own constants in the input file.
-These constants can be used for any parameter that consists in one real number.
+These constants can be used for any parameter that consists of one or more floats.
User-defined constant names can contain only letters, numbers and the character ``_``.
The name of each constant has to begin with a letter. The following names are used
by WarpX, and cannot be used as user-defined constants: ``x``, ``y``, ``z``, ``X``, ``Y``, ``t``.
diff --git a/Examples/Modules/boosted_diags/inputs_3d_slice b/Examples/Modules/boosted_diags/inputs_3d_slice
index d2663eece..196c15641 100644
--- a/Examples/Modules/boosted_diags/inputs_3d_slice
+++ b/Examples/Modules/boosted_diags/inputs_3d_slice
@@ -5,10 +5,17 @@ amr.max_grid_size = 64
amr.blocking_factor = 32
amr.max_level = 0
+my_constants.xmin = -128.e-6
+my_constants.ymin = -128.e-6
+my_constants.zmin = -40.e-6
+my_constants.xmax = +128.e-6
+my_constants.ymax = +128.e-6
+my_constants.zmax = 0.96e-6
+
geometry.coord_sys = 0 # 0: Cartesian
geometry.is_periodic = 1 1 0 # Is periodic?
-geometry.prob_lo = -128.e-6 -128.e-6 -40.e-6
-geometry.prob_hi = 128.e-6 128.e-6 0.96e-6
+geometry.prob_lo = xmin ymin zmin
+geometry.prob_hi = xmax ymax zmax
algo.current_deposition = esirkepov
algo.charge_deposition = standard
@@ -100,8 +107,8 @@ laser1.profile_t_peak = 40.e-15 # The time at which the laser reaches its pea
laser1.profile_focal_distance = 0.5e-3 # Focal distance from the antenna (in meters)
laser1.wavelength = 0.81e-6 # The wavelength of the laser (in meters)
-slice.dom_lo = -128.e-6 0.0 -40.e-6
-slice.dom_hi = 128.e-6 0.0 0.96e-6
+slice.dom_lo = xmin 0.0 zmin
+slice.dom_hi = xmax 0.0 zmax
slice.coarsening_ratio = 1 1 1
slice.plot_int = -1
slice.num_slice_snapshots_lab = 4
diff --git a/Source/Diagnostics/Diagnostics.cpp b/Source/Diagnostics/Diagnostics.cpp
index 1e305d0fc..d58457a70 100644
--- a/Source/Diagnostics/Diagnostics.cpp
+++ b/Source/Diagnostics/Diagnostics.cpp
@@ -81,14 +81,14 @@ Diagnostics::BaseReadParameters ()
m_lo.resize(AMREX_SPACEDIM);
m_hi.resize(AMREX_SPACEDIM);
- bool lo_specified = pp_diag_name.queryarr("diag_lo", m_lo);
+ bool lo_specified = queryArrWithParser(pp_diag_name, "diag_lo", m_lo, 0, AMREX_SPACEDIM);
if (!lo_specified) {
for (int idim=0; idim < AMREX_SPACEDIM; ++idim) {
m_lo[idim] = warpx.Geom(0).ProbLo(idim);
}
}
- bool hi_specified = pp_diag_name.queryarr("diag_hi", m_hi);
+ bool hi_specified = queryArrWithParser(pp_diag_name, "diag_hi", m_hi, 0, AMREX_SPACEDIM);
if (!hi_specified) {
for (int idim =0; idim < AMREX_SPACEDIM; ++idim) {
m_hi[idim] = warpx.Geom(0).ProbHi(idim);
@@ -193,7 +193,8 @@ Diagnostics::InitData ()
amrex::ParmParse pp_diag_name(m_diag_name);
amrex::Vector <amrex::Real> dummy_val(AMREX_SPACEDIM);
- if ( pp_diag_name.queryarr("diag_lo", dummy_val) || pp_diag_name.queryarr("diag_hi", dummy_val) ) {
+ if ( queryArrWithParser(pp_diag_name, "diag_lo", dummy_val, 0, AMREX_SPACEDIM) ||
+ queryArrWithParser(pp_diag_name, "diag_hi", dummy_val, 0, AMREX_SPACEDIM) ) {
// set geometry filter for particle-diags to true when the diagnostic domain-extent
// is specified by the user
for (int i = 0; i < m_output_species.size(); ++i) {
diff --git a/Source/Initialization/CustomDensityProb.H b/Source/Initialization/CustomDensityProb.H
index 3151d4904..935345bd9 100644
--- a/Source/Initialization/CustomDensityProb.H
+++ b/Source/Initialization/CustomDensityProb.H
@@ -12,6 +12,8 @@
#include <AMReX_Gpu.H>
#include <AMReX_Dim3.H>
+#include "Utils/WarpXUtil.H"
+
// An example of Custom Density Profile
// struct whose getDensity returns density at a given position computed from
@@ -25,7 +27,7 @@ struct InjectorDensityCustom
std::vector<amrex::Real> v;
AMREX_ALWAYS_ASSERT_WITH_MESSAGE(v.size() <= 6,
"Too many parameters for InjectorDensityCustom");
- pp_species_name.getarr("custom_profile_params", v);
+ getArrWithParser(pp_species_name, "custom_profile_params", v);
for (int i = 0; i < static_cast<int>(v.size()); ++i) {
p[i] = v[i];
}
diff --git a/Source/Initialization/InjectorDensity.cpp b/Source/Initialization/InjectorDensity.cpp
index a6df6ae13..c2b2d3605 100644
--- a/Source/Initialization/InjectorDensity.cpp
+++ b/Source/Initialization/InjectorDensity.cpp
@@ -43,7 +43,7 @@ InjectorDensityPredefined::InjectorDensityPredefined (
std::vector<amrex::Real> v;
// Read parameters for the predefined plasma profile.
- pp_species_name.getarr("predefined_profile_params", v);
+ getArrWithParser(pp_species_name, "predefined_profile_params", v);
AMREX_ALWAYS_ASSERT_WITH_MESSAGE(v.size() <= 6,
"Too many parameters for InjectorDensityPredefined");
for (int i = 0; i < static_cast<int>(v.size()); ++i) {
diff --git a/Source/Initialization/PlasmaInjector.cpp b/Source/Initialization/PlasmaInjector.cpp
index e838cfcb4..5c99bb03a 100644
--- a/Source/Initialization/PlasmaInjector.cpp
+++ b/Source/Initialization/PlasmaInjector.cpp
@@ -143,8 +143,8 @@ PlasmaInjector::PlasmaInjector (int ispecies, const std::string& name)
if (part_pos_s == "python") {
return;
} else if (part_pos_s == "singleparticle") {
- pp_species_name.getarr("single_particle_pos", single_particle_pos, 0, 3);
- pp_species_name.getarr("single_particle_vel", single_particle_vel, 0, 3);
+ getArrWithParser(pp_species_name, "single_particle_pos", single_particle_pos, 0, 3);
+ getArrWithParser(pp_species_name, "single_particle_vel", single_particle_vel, 0, 3);
for (auto& x : single_particle_vel) {
x *= PhysConst::c;
}
@@ -152,13 +152,13 @@ PlasmaInjector::PlasmaInjector (int ispecies, const std::string& name)
add_single_particle = true;
return;
} else if (part_pos_s == "multipleparticles") {
- pp_species_name.getarr("multiple_particles_pos_x", multiple_particles_pos_x);
- pp_species_name.getarr("multiple_particles_pos_y", multiple_particles_pos_y);
- pp_species_name.getarr("multiple_particles_pos_z", multiple_particles_pos_z);
- pp_species_name.getarr("multiple_particles_vel_x", multiple_particles_vel_x);
- pp_species_name.getarr("multiple_particles_vel_y", multiple_particles_vel_y);
- pp_species_name.getarr("multiple_particles_vel_z", multiple_particles_vel_z);
- pp_species_name.getarr("multiple_particles_weight", multiple_particles_weight);
+ getArrWithParser(pp_species_name, "multiple_particles_pos_x", multiple_particles_pos_x);
+ getArrWithParser(pp_species_name, "multiple_particles_pos_y", multiple_particles_pos_y);
+ getArrWithParser(pp_species_name, "multiple_particles_pos_z", multiple_particles_pos_z);
+ getArrWithParser(pp_species_name, "multiple_particles_vel_x", multiple_particles_vel_x);
+ getArrWithParser(pp_species_name, "multiple_particles_vel_y", multiple_particles_vel_y);
+ getArrWithParser(pp_species_name, "multiple_particles_vel_z", multiple_particles_vel_z);
+ getArrWithParser(pp_species_name, "multiple_particles_weight", multiple_particles_weight);
AMREX_ALWAYS_ASSERT_WITH_MESSAGE(
((multiple_particles_pos_x.size() == multiple_particles_pos_y.size()) &&
(multiple_particles_pos_x.size() == multiple_particles_pos_z.size()) &&
diff --git a/Source/Initialization/WarpXInitData.cpp b/Source/Initialization/WarpXInitData.cpp
index b6182b1ff..892172745 100644
--- a/Source/Initialization/WarpXInitData.cpp
+++ b/Source/Initialization/WarpXInitData.cpp
@@ -340,12 +340,12 @@ WarpX::InitLevelData (int lev, Real /*time*/)
// if the input string is "constant", the values for the
// external grid must be provided in the input.
if (B_ext_grid_s == "constant")
- pp_warpx.getarr("B_external_grid", B_external_grid);
+ getArrWithParser(pp_warpx, "B_external_grid", B_external_grid);
// if the input string is "constant", the values for the
// external grid must be provided in the input.
if (E_ext_grid_s == "constant")
- pp_warpx.getarr("E_external_grid", E_external_grid);
+ getArrWithParser(pp_warpx, "E_external_grid", E_external_grid);
// initialize the averaged fields only if the averaged algorithm
// is activated ('psatd.do_time_averaging=1')
diff --git a/Source/Particles/LaserParticleContainer.cpp b/Source/Particles/LaserParticleContainer.cpp
index c5fe3e4fb..3e573f33e 100644
--- a/Source/Particles/LaserParticleContainer.cpp
+++ b/Source/Particles/LaserParticleContainer.cpp
@@ -46,9 +46,9 @@ LaserParticleContainer::LaserParticleContainer (AmrCore* amr_core, int ispecies,
std::transform(laser_type_s.begin(), laser_type_s.end(), laser_type_s.begin(), ::tolower);
// Parse the properties of the antenna
- pp_laser_name.getarr("position", m_position);
- pp_laser_name.getarr("direction", m_nvec);
- pp_laser_name.getarr("polarization", m_p_X);
+ getArrWithParser(pp_laser_name, "position", m_position);
+ getArrWithParser(pp_laser_name, "direction", m_nvec);
+ getArrWithParser(pp_laser_name, "polarization", m_p_X);
pp_laser_name.query("pusher_algo", m_pusher_algo);
getWithParser(pp_laser_name, "wavelength", m_wavelength);
@@ -127,10 +127,10 @@ LaserParticleContainer::LaserParticleContainer (AmrCore* amr_core, int ispecies,
m_laser_injection_box= Geom(0).ProbDomain();
{
Vector<Real> lo, hi;
- if (pp_laser_name.queryarr("prob_lo", lo)) {
+ if (queryArrWithParser(pp_laser_name, "prob_lo", lo, 0, AMREX_SPACEDIM)) {
m_laser_injection_box.setLo(lo);
}
- if (pp_laser_name.queryarr("prob_hi", hi)) {
+ if (queryArrWithParser(pp_laser_name, "prob_hi", hi, 0, AMREX_SPACEDIM)) {
m_laser_injection_box.setHi(hi);
}
}
diff --git a/Source/Particles/MultiParticleContainer.cpp b/Source/Particles/MultiParticleContainer.cpp
index 4f6bdf05c..eb5293bee 100644
--- a/Source/Particles/MultiParticleContainer.cpp
+++ b/Source/Particles/MultiParticleContainer.cpp
@@ -117,13 +117,13 @@ MultiParticleContainer::ReadParameters ()
// then the values for the external B on particles must
// be provided in the input file.
if (m_B_ext_particle_s == "constant")
- pp_particles.getarr("B_external_particle", m_B_external_particle);
+ getArrWithParser(pp_particles, "B_external_particle", m_B_external_particle);
// if the input string for E_external on particles is "constant"
// then the values for the external E on particles must
// be provided in the input file.
if (m_E_ext_particle_s == "constant")
- pp_particles.getarr("E_external_particle", m_E_external_particle);
+ getArrWithParser(pp_particles, "E_external_particle", m_E_external_particle);
// if the input string for B_ext_particle_s is
// "parse_b_ext_particle_function" then the mathematical expression
diff --git a/Source/Utils/WarpXUtil.H b/Source/Utils/WarpXUtil.H
index 057b3193c..45b6de45f 100644
--- a/Source/Utils/WarpXUtil.H
+++ b/Source/Utils/WarpXUtil.H
@@ -19,6 +19,7 @@
#include <cstdint>
#include <string>
+void ParseGeometryInput();
void ReadBoostedFrameParameters(amrex::Real& gamma_boost, amrex::Real& beta_boost,
amrex::Vector<int>& boost_direction);
@@ -179,6 +180,8 @@ WarpXParser makeParser (std::string const& parse_function, std::vector<std::stri
* \param[out] val where the value queried and parsed is stored
*/
int queryWithParser (const amrex::ParmParse& a_pp, char const * const str, amrex::Real& val);
+int queryArrWithParser (const amrex::ParmParse& a_pp, char const * const str, std::vector<amrex::Real>& val,
+ const int start_ix, const int num_val);
/**
* \brief Similar to amrex::ParmParse::get, but also supports math expressions for the value.
@@ -193,6 +196,9 @@ int queryWithParser (const amrex::ParmParse& a_pp, char const * const str, amrex
* \param[out] val where the value queried and parsed is stored
*/
void getWithParser (const amrex::ParmParse& a_pp, char const * const str, amrex::Real& val);
+void getArrWithParser (const amrex::ParmParse& a_pp, char const * const str, std::vector<amrex::Real>& val);
+void getArrWithParser (const amrex::ParmParse& a_pp, char const * const str, std::vector<amrex::Real>& val,
+ const int start_ix, const int num_val);
namespace WarpXUtilMsg{
diff --git a/Source/Utils/WarpXUtil.cpp b/Source/Utils/WarpXUtil.cpp
index ac76b8ce6..d85cd17a4 100644
--- a/Source/Utils/WarpXUtil.cpp
+++ b/Source/Utils/WarpXUtil.cpp
@@ -20,6 +20,22 @@
using namespace amrex;
+void ParseGeometryInput()
+{
+ ParmParse pp_geometry("geometry");
+
+ Vector<Real> prob_lo(AMREX_SPACEDIM);
+ Vector<Real> prob_hi(AMREX_SPACEDIM);
+
+ getArrWithParser(pp_geometry, "prob_lo", prob_lo, 0, AMREX_SPACEDIM);
+ AMREX_ALWAYS_ASSERT(prob_lo.size() == AMREX_SPACEDIM);
+ getArrWithParser(pp_geometry, "prob_hi", prob_hi, 0, AMREX_SPACEDIM);
+ AMREX_ALWAYS_ASSERT(prob_hi.size() == AMREX_SPACEDIM);
+
+ pp_geometry.addarr("prob_lo", prob_lo);
+ pp_geometry.addarr("prob_hi", prob_hi);
+}
+
void ReadBoostedFrameParameters(Real& gamma_boost, Real& beta_boost,
Vector<int>& boost_direction)
{
@@ -72,21 +88,19 @@ void ConvertLabParamsToBoost()
ParmParse pp_amr("amr");
ParmParse pp_slice("slice");
- pp_geometry.getarr("prob_lo",prob_lo,0,AMREX_SPACEDIM);
- AMREX_ALWAYS_ASSERT(prob_lo.size() == AMREX_SPACEDIM);
- pp_geometry.getarr("prob_hi",prob_hi,0,AMREX_SPACEDIM);
- AMREX_ALWAYS_ASSERT(prob_hi.size() == AMREX_SPACEDIM);
+ getArrWithParser(pp_geometry, "prob_lo", prob_lo, 0, AMREX_SPACEDIM);
+ getArrWithParser(pp_geometry, "prob_hi", prob_hi, 0, AMREX_SPACEDIM);
- pp_slice.queryarr("dom_lo",slice_lo,0,AMREX_SPACEDIM);
+ queryArrWithParser(pp_slice, "dom_lo", slice_lo, 0, AMREX_SPACEDIM);
AMREX_ALWAYS_ASSERT(slice_lo.size() == AMREX_SPACEDIM);
- pp_slice.queryarr("dom_hi",slice_hi,0,AMREX_SPACEDIM);
+ queryArrWithParser(pp_slice, "dom_hi", slice_hi, 0, AMREX_SPACEDIM);
AMREX_ALWAYS_ASSERT(slice_hi.size() == AMREX_SPACEDIM);
pp_amr.query("max_level", max_level);
if (max_level > 0){
- pp_warpx.getarr("fine_tag_lo", fine_tag_lo);
- pp_warpx.getarr("fine_tag_hi", fine_tag_hi);
+ getArrWithParser(pp_warpx, "fine_tag_lo", fine_tag_lo);
+ getArrWithParser(pp_warpx, "fine_tag_hi", fine_tag_hi);
}
@@ -273,6 +287,58 @@ getWithParser (const amrex::ParmParse& a_pp, char const * const str, amrex::Real
val = parser.eval();
}
+int
+queryArrWithParser (const amrex::ParmParse& a_pp, char const * const str, std::vector<amrex::Real>& val,
+ const int start_ix, const int num_val)
+{
+ // call amrex::ParmParse::query, check if the user specified str.
+ std::vector<std::string> tmp_str_arr;
+ int is_specified = a_pp.queryarr(str, tmp_str_arr, start_ix, num_val);
+ if (is_specified)
+ {
+ // If so, create parser objects and apply them to the values provided by the user.
+ int const n = static_cast<int>(tmp_str_arr.size());
+ val.resize(n);
+ for (int i=0 ; i < n ; i++) {
+ auto parser = makeParser(tmp_str_arr[i], {});
+ val[i] = parser.eval();
+ }
+ }
+ // return the same output as amrex::ParmParse::query
+ return is_specified;
+}
+
+void
+getArrWithParser (const amrex::ParmParse& a_pp, char const * const str, std::vector<amrex::Real>& val)
+{
+ // Create parser objects and apply them to the values provided by the user.
+ std::vector<std::string> tmp_str_arr;
+ a_pp.getarr(str, tmp_str_arr);
+
+ int const n = static_cast<int>(tmp_str_arr.size());
+ val.resize(n);
+ for (int i=0 ; i < n ; i++) {
+ auto parser = makeParser(tmp_str_arr[i], {});
+ val[i] = parser.eval();
+ }
+}
+
+void
+getArrWithParser (const amrex::ParmParse& a_pp, char const * const str, std::vector<amrex::Real>& val,
+ const int start_ix, const int num_val)
+{
+ // Create parser objects and apply them to the values provided by the user.
+ std::vector<std::string> tmp_str_arr;
+ a_pp.getarr(str, tmp_str_arr, start_ix, num_val);
+
+ int const n = static_cast<int>(tmp_str_arr.size());
+ val.resize(n);
+ for (int i=0 ; i < n ; i++) {
+ auto parser = makeParser(tmp_str_arr[i], {});
+ val[i] = parser.eval();
+ }
+}
+
/**
* \brief Ensures that the blocks are setup correctly for the RZ spectral solver
* When using the RZ spectral solver, the Hankel transform cannot be
diff --git a/Source/WarpX.cpp b/Source/WarpX.cpp
index 5bd68be51..fc066472b 100644
--- a/Source/WarpX.cpp
+++ b/Source/WarpX.cpp
@@ -536,9 +536,9 @@ WarpX::ReadParameters ()
pp_warpx.query("num_mirrors", num_mirrors);
if (num_mirrors>0){
mirror_z.resize(num_mirrors);
- pp_warpx.getarr("mirror_z", mirror_z, 0, num_mirrors);
+ getArrWithParser(pp_warpx, "mirror_z", mirror_z, 0, num_mirrors);
mirror_z_width.resize(num_mirrors);
- pp_warpx.getarr("mirror_z_width", mirror_z_width, 0, num_mirrors);
+ getArrWithParser(pp_warpx, "mirror_z_width", mirror_z_width, 0, num_mirrors);
mirror_z_npoints.resize(num_mirrors);
pp_warpx.getarr("mirror_z_npoints", mirror_z_npoints, 0, num_mirrors);
}
@@ -703,8 +703,8 @@ WarpX::ReadParameters ()
if (maxLevel() > 0) {
Vector<Real> lo, hi;
- pp_warpx.getarr("fine_tag_lo", lo);
- pp_warpx.getarr("fine_tag_hi", hi);
+ getArrWithParser(pp_warpx, "fine_tag_lo", lo);
+ getArrWithParser(pp_warpx, "fine_tag_hi", hi);
fine_tag_lo = RealVect{lo};
fine_tag_hi = RealVect{hi};
}
@@ -1005,8 +1005,8 @@ WarpX::ReadParameters ()
{
slice_crse_ratio[idim] = 1;
}
- pp_slice.queryarr("dom_lo",slice_lo,0,AMREX_SPACEDIM);
- pp_slice.queryarr("dom_hi",slice_hi,0,AMREX_SPACEDIM);
+ queryArrWithParser(pp_slice, "dom_lo", slice_lo, 0, AMREX_SPACEDIM);
+ queryArrWithParser(pp_slice, "dom_hi", slice_hi, 0, AMREX_SPACEDIM);
pp_slice.queryarr("coarsening_ratio",slice_crse_ratio,0,AMREX_SPACEDIM);
pp_slice.query("plot_int",slice_plot_int);
slice_realbox.setLo(slice_lo);
diff --git a/Source/main.cpp b/Source/main.cpp
index 477a0d495..276b7c31b 100644
--- a/Source/main.cpp
+++ b/Source/main.cpp
@@ -34,6 +34,8 @@ int main(int argc, char* argv[])
rocfft_setup();
#endif
+ ParseGeometryInput();
+
ConvertLabParamsToBoost();
ReadBCParams();