35#ifndef IC_EXPRESSION_H_
36#define IC_EXPRESSION_H_
40#include "AMReX_Parser.H"
49 std::vector<amrex::ParserExecutor<4>>
f;
52 static constexpr const char*
name =
"expression";
54 IC<
Set::Scalar>(_geom),
IC<
Set::Vector>(_geom) {}
63 for (amrex::MFIter mfi(*a_field[lev], amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
67 amrex::IndexType type = a_field[lev]->ixType();
68 if (type == amrex::IndexType::TheCellType()) bx = mfi.growntilebox();
69 else if (type == amrex::IndexType::TheNodeType()) bx = mfi.grownnodaltilebox();
72 amrex::Array4<Set::Scalar>
const& field = a_field[lev]->array(mfi);
73 for (
unsigned int n = 0; n <
f.size(); n++)
75 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
80#if AMREX_SPACEDIM == 1
81 field(i, j, k, n) =
f[n](x(0), 0.0, 0.0, a_time);
82#elif AMREX_SPACEDIM == 2
83 field(i, j, k, n) =
f[n](x(0), x(1), 0.0, a_time);
84#elif AMREX_SPACEDIM == 3
85 field(i, j, k, n) =
f[n](x(0), x(1), x(2), a_time);
91 field(i, j, k, n) =
f[n](sqrt(x(0)* x(0) + x(1) * x(1)), std::atan2(x(1), x(0)), x(2), a_time);
97 a_field[lev]->FillBoundary();
104 for (amrex::MFIter mfi(*a_field[lev], amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
107 amrex::IndexType type = a_field[lev]->ixType();
108 if (type == amrex::IndexType::TheCellType()) bx = mfi.growntilebox();
109 else if (type == amrex::IndexType::TheNodeType()) bx = mfi.grownnodaltilebox();
113 for (
unsigned int n = 0; n < AMREX_SPACEDIM; n++)
115 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
120#if AMREX_SPACEDIM == 1
121 field(i, j, k)(n) =
f[n](x(0), 0.0, 0.0, a_time);
122#elif AMREX_SPACEDIM == 2
123 field(i, j, k)(n) =
f[n](x(0), x(1), 0.0, a_time);
124#elif AMREX_SPACEDIM == 3
125 field(i, j, k)(n) =
f[n](x(0), x(1), x(2), a_time);
131 field(i, j, k)(n) =
f[n](sqrt(x(0)* x(0) + x(1) * x(1)), std::atan2(x(1), x(0)), x(2), a_time);
137 a_field[lev]->FillBoundary();
143 std::string coordstr =
"";
150 std::vector<std::string> expression_strs;
155 for (
unsigned int i = 0; i < expression_strs.size(); i++)
157 value.
parser.push_back(amrex::Parser(expression_strs[i]));
163 std::set<std::string> entries = pp.getEntries(prefix +
".constant");
164 std::set<std::string>::iterator entry;
165 for (entry = entries.begin(); entry != entries.end(); entry++)
168 std::string fullname = *entry;
177 value.
parser.back().registerVariables({
"x",
"y",
"z",
"t" });
178 value.
f.push_back(value.
parser.back().compile<4>());
182 value.
parser.back().registerVariables({
"r",
"theta",
"z",
"t" });
183 value.
f.push_back(value.
parser.back().compile<4>());
#define pp_query_validate(...)
#define pp_queryclass(...)
static constexpr const char * name
virtual void Add(const int &lev, Set::Field< Set::Vector > &a_field, Set::Scalar a_time=0.0) override
virtual void Add(const int &lev, Set::Field< Set::Scalar > &a_field, Set::Scalar a_time=0.0) override
Expression::CoordSys coord
Expression(amrex::Vector< amrex::Geometry > &_geom, IO::ParmParse &pp, std::string name)
std::vector< amrex::ParserExecutor< 4 > > f
std::vector< amrex::Parser > parser
Expression(amrex::Vector< amrex::Geometry > &_geom)
static void Parse(Expression &value, IO::ParmParse &pp)
Pure abstract IC object from which all other IC objects inherit.
std::string getPrefix() const
int query_enumerate(std::string a_name, std::vector< T > &value, int number=1, std::string file="", std::string func="", int line=__LINE__)
amrex::Array4< T > Patch(int lev, amrex::MFIter &mfi) const &
Initialize a spherical inclusion.
A collection of data types and symmetry-reduced data structures.
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
AMREX_FORCE_INLINE Vector Position(const int &i, const int &j, const int &k, const amrex::Geometry &geom, const amrex::IndexType &ixType)
std::vector< std::string > Split(std::string &str, const char delim)
void Abort(const char *msg)
AMREX_FORCE_INLINE void Assert(std::string file, std::string func, int line, std::string smt, bool pass, Args const &... args)