Alamo
Expression.H
Go to the documentation of this file.
1 //
2 // Initialize a field using a mathematical expression.
3 // Expressions are imported as strings and are compiled real-time using the
4 // `AMReX Parser <https://amrex-codes.github.io/amrex/docs_html/Basics.html#parser>`_.
5 //
6 // Works for single or multiple-component fields.
7 // Use the :code:`regionN` (N=0,1,2, etc. up to number of components) to pass expression.
8 // For example:
9 //
10 // .. code-block:: bash
11 //
12 // ic.region0 = "sin(x*y*z)"
13 // ic.region1 = "3.0*(x > 0.5 and y > 0.5)"
14 //
15 // for a two-component field. It is up to you to make sure your expressions are parsed
16 // correctly; otherwise you will get undefined behavior.
17 //
18 // :bdg-primary-line:`Constants`
19 // You can add constants to your expressions using the :code:`constant` directive.
20 // For instance, in the following code
21 //
22 // .. code-block:: bash
23 //
24 // psi.ic.type=expression
25 // psi.ic.expression.constant.eps = 0.05
26 // psi.ic.expression.constant.R = 0.25
27 // psi.ic.expression.region0 = "0.5 + 0.5*tanh((x^2 + y^2 - R)/eps)"
28 //
29 // the constants :code:`eps` and :code:`R` are defined by the user and then used
30 // in the subsequent expression.
31 // The variables can have any name made up of characters that is not reserved.
32 // However, if multiple ICs are used, they must be defined each time for each IC.
33 //
34 
35 #ifndef IC_EXPRESSION_H_
36 #define IC_EXPRESSION_H_
37 #include "IC/IC.H"
38 #include "Util/Util.H"
39 #include "IO/ParmParse.H"
40 #include "AMReX_Parser.H"
41 
42 namespace IC
43 {
44 class Expression : public IC
45 {
46 private:
48  std::vector<amrex::Parser> parser;
49  std::vector<amrex::ParserExecutor<4>> f;
50  Expression::CoordSys coord = Expression::CoordSys::Cartesian;
51 public:
52  Expression(amrex::Vector<amrex::Geometry>& _geom) : IC(_geom) {}
53  Expression(amrex::Vector<amrex::Geometry>& _geom, IO::ParmParse& pp, std::string name) : IC(_geom)
54  {
55  pp_queryclass(name, *this);
56  }
57  void Add(const int& lev, Set::Field<Set::Scalar>& a_field, Set::Scalar a_time = 0.0)
58  {
59  Util::Assert(INFO, TEST(a_field[lev]->nComp() == (int)f.size()));
60  for (amrex::MFIter mfi(*a_field[lev], amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
61  {
62  amrex::Box bx;// = mfi.tilebox();
63  //bx.grow(a_field[lev]->nGrow());
64  amrex::IndexType type = a_field[lev]->ixType();
65  if (type == amrex::IndexType::TheCellType()) bx = mfi.growntilebox();
66  else if (type == amrex::IndexType::TheNodeType()) bx = mfi.grownnodaltilebox();
67  else Util::Abort(INFO, "Unkonwn index type");
68 
69  amrex::Array4<Set::Scalar> const& field = a_field[lev]->array(mfi);
70  for (unsigned int n = 0; n < f.size(); n++)
71  {
72  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
73  {
74  Set::Vector x = Set::Position(i, j, k, geom[lev], type);
75  if (coord == Expression::CoordSys::Cartesian)
76  {
77 #if AMREX_SPACEDIM == 1
78  field(i, j, k, n) = f[n](x(0), 0.0, 0.0, a_time);
79 #elif AMREX_SPACEDIM == 2
80  field(i, j, k, n) = f[n](x(0), x(1), 0.0, a_time);
81 #elif AMREX_SPACEDIM == 3
82  field(i, j, k, n) = f[n](x(0), x(1), x(2), a_time);
83 #endif
84  }
85 #if AMREX_SPACEDIM>1
86  else if (coord == Expression::CoordSys::Polar)
87  {
88  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);
89  }
90 #endif
91  });
92  }
93  }
94  a_field[lev]->FillBoundary();
95  };
96 
97  void Add(const int& lev, Set::Field<Set::Vector>& a_field, Set::Scalar a_time = 0.0)
98  {
99  Util::Assert(INFO, TEST(a_field[lev]->nComp() == 1));
100  Util::Assert(INFO, TEST(f.size() >= AMREX_SPACEDIM));
101  for (amrex::MFIter mfi(*a_field[lev], amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
102  {
103  amrex::Box bx;
104  amrex::IndexType type = a_field[lev]->ixType();
105  if (type == amrex::IndexType::TheCellType()) bx = mfi.growntilebox();
106  else if (type == amrex::IndexType::TheNodeType()) bx = mfi.grownnodaltilebox();
107  else Util::Abort(INFO, "Unkonwn index type");
108 
109  Set::Patch<Set::Vector> field = a_field.Patch(lev,mfi);
110  for (unsigned int n = 0; n < AMREX_SPACEDIM; n++)
111  {
112  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
113  {
114  Set::Vector x = Set::Position(i, j, k, geom[lev], type);
115  if (coord == Expression::CoordSys::Cartesian)
116  {
117 #if AMREX_SPACEDIM == 1
118  field(i, j, k)(n) = f[n](x(0), 0.0, 0.0, a_time);
119 #elif AMREX_SPACEDIM == 2
120  field(i, j, k)(n) = f[n](x(0), x(1), 0.0, a_time);
121 #elif AMREX_SPACEDIM == 3
122  field(i, j, k)(n) = f[n](x(0), x(1), x(2), a_time);
123 #endif
124  }
125 #if AMREX_SPACEDIM>1
126  else if (coord == Expression::CoordSys::Polar)
127  {
128  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);
129  }
130 #endif
131  });
132  }
133  }
134  a_field[lev]->FillBoundary();
135  };
136 
137  static void Parse(Expression& value, IO::ParmParse& pp)
138  {
139  for (int i = 0; true; i++)
140  {
141  std::string coordstr = "cartesian";
142  // coordinate system to use: "cartesian" (for x,y,z,t) and
143  // "polar" (for r, theta, z, t)
144  pp_query("coord", coordstr);
145  if (coordstr == "cartesian") value.coord = Expression::CoordSys::Cartesian;
146  else if (coordstr == "polar") value.coord = Expression::CoordSys::Polar;
147  else Util::Abort(INFO, "unsupported coordinates ", coordstr);
148 
149  std::string func = "0.0";
150  std::string name = "region" + std::to_string(i);
151 
152  if (!pp.contains(name.data())) break;
153  pp_query(name.data(), func);
154 
155  value.parser.push_back(amrex::Parser(func));
156 
157  //
158  // Read in user-defined constants and add them to the parser
159  //
160  std::string prefix = pp.getPrefix();
161  std::set<std::string> entries = pp.getEntries(prefix + ".constant");//"constant");
162  std::set<std::string>::iterator entry;
163  for (entry = entries.begin(); entry != entries.end(); entry++)
164  {
165  IO::ParmParse pp;
166  std::string fullname = *entry;
167  Set::Scalar val = NAN;
168  pp_query(fullname.data(),val);
169  std::string name = Util::String::Split(fullname,'.').back();
170  value.parser.back().setConstant(name,val);
171  }
172 
173  if (value.coord == Expression::CoordSys::Cartesian)
174  {
175  value.parser.back().registerVariables({ "x","y","z","t" });
176  value.f.push_back(value.parser.back().compile<4>());
177  }
178  else if (value.coord == Expression::CoordSys::Polar)
179  {
180  value.parser.back().registerVariables({ "r","theta","z","t" });
181  value.f.push_back(value.parser.back().compile<4>());
182  }
183  }
184 
185  };
186 };
187 }
188 
189 #endif
IC::Expression::parser
std::vector< amrex::Parser > parser
Definition: Expression.H:48
IC::IC::geom
amrex::Vector< amrex::Geometry > & geom
Definition: IC.H:45
Set::Position
AMREX_FORCE_INLINE Vector Position(const int &i, const int &j, const int &k, const amrex::Geometry &geom, const amrex::IndexType &ixType)
Definition: Base.H:121
Util.H
TEST
#define TEST(x)
Definition: Util.H:21
Util::String::Split
std::vector< std::string > Split(std::string &str, const char delim)
Definition: Util.cpp:291
IC::Expression::Polar
@ Polar
Definition: Expression.H:47
IC::Expression::Expression
Expression(amrex::Vector< amrex::Geometry > &_geom)
Definition: Expression.H:52
Set::Field< Set::Scalar >
Definition: Set.H:236
Set::Field::Patch
amrex::Array4< T > Patch(int lev, amrex::MFIter &mfi) const &
Definition: Set.H:76
Set::Vector
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition: Base.H:20
ParmParse.H
pp_query
#define pp_query(...)
Definition: ParmParse.H:104
Set::Scalar
amrex::Real Scalar
Definition: Base.H:19
IC::Expression::coord
Expression::CoordSys coord
Definition: Expression.H:50
pp_queryclass
#define pp_queryclass(...)
Definition: ParmParse.H:105
IC::Expression::Expression
Expression(amrex::Vector< amrex::Geometry > &_geom, IO::ParmParse &pp, std::string name)
Definition: Expression.H:53
IC::Expression::f
std::vector< amrex::ParserExecutor< 4 > > f
Definition: Expression.H:49
IC::Expression::Add
void Add(const int &lev, Set::Field< Set::Vector > &a_field, Set::Scalar a_time=0.0)
Definition: Expression.H:97
IO::ParmParse::contains
bool contains(std::string name)
Definition: ParmParse.H:151
Util::Assert
AMREX_FORCE_INLINE void Assert(std::string file, std::string func, int line, std::string smt, bool pass, Args const &... args)
Definition: Util.H:65
Util::Abort
void Abort(const char *msg)
Definition: Util.cpp:165
IC::Expression::Add
void Add(const int &lev, Set::Field< Set::Scalar > &a_field, Set::Scalar a_time=0.0)
Definition: Expression.H:57
IC::Expression
Definition: Expression.H:44
IC::Expression::Spherical
@ Spherical
Definition: Expression.H:47
IC::Expression::Cartesian
@ Cartesian
Definition: Expression.H:47
IC
Definition: Affine.H:18
IO::ParmParse::getPrefix
std::string getPrefix() const
Definition: ParmParse.H:131
IO::ParmParse
Definition: ParmParse.H:110
Set::Patch
amrex::Array4< T > const & Patch
Definition: Set.H:88
IC::Expression::Parse
static void Parse(Expression &value, IO::ParmParse &pp)
Definition: Expression.H:137
INFO
#define INFO
Definition: Util.H:20
IC.H
Set::Field< Set::Vector >
IC::Expression::CoordSys
CoordSys
Definition: Expression.H:47