Alamo
Expression.H
Go to the documentation of this file.
1 #ifndef BC_OPERATOR_ELASTIC_EXPRESSION_H
2 #define BC_OPERATOR_ELASTIC_EXPRESSION_H
3 
4 // #include "Operator/Elastic.H"
5 #include "IO/ParmParse.H"
8 
9 #include "AMReX_Parser.H"
10 
11 namespace BC
12 {
13 namespace Operator
14 {
15 namespace Elastic
16 {
17 class Expression : public Elastic
18 {
19 public:
20  //static constexpr const char* const Elastic::strings[];
21  #if AMREX_SPACEDIM==2
22  static const constexpr char * const facestr[] = {
23  "xlo","ylo","xhi","yhi","xloylo","xloyhi","xhiylo","xhiyhi","int"
24  };
25  #elif AMREX_SPACEDIM==3
26  static const constexpr char * const facestr[] = {
27  "xlo","ylo","zlo","xhi","yhi","zhi",
28  "ylozlo","ylozhi","yhizlo","yhizhi",
29  "zloxlo","zloxhi","zhixlo","zhixhi",
30  "xloylo","xloyhi","xhiylo","xhiyhi",
31  "xloylozlo","xloylozhi","xloyhizlo","xloyhizhi",
32  "xhiylozlo","xhiylozhi","xhiyhizlo","xhiyhizhi",
33  "int"
34  };
35  #endif
36 
38  {
39  // By default, all boundary conditions are displacement
40  for (int face = 0; face < m_nfaces; face++)
41  for (int direction = 0; direction < AMREX_SPACEDIM; direction++)
42  {
43  m_bc_type[face][direction] = Type::Displacement;
44  }
45  };
47 
48 
49  using Elastic::Init;
50 
51  virtual void
52  Init(amrex::FabArray<amrex::BaseFab<Set::Vector>> *a_rhs,
53  const amrex::Geometry &a_geom,
54  bool a_homogeneous = false) const override
55  {
56  amrex::Box domain(a_geom.Domain());
57  domain.convert(amrex::IntVect::TheNodeVector());
58  const amrex::Dim3 lo= amrex::lbound(domain), hi = amrex::ubound(domain);
59  for (amrex::MFIter mfi(*a_rhs, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
60  {
61  amrex::Box bx = mfi.tilebox();
62  bx.grow(2);
63  bx = bx & domain;
64  amrex::Array4<Set::Vector> const& rhs = a_rhs->array(mfi);
65  amrex::ParallelFor (bx,[=] AMREX_GPU_DEVICE(int i, int j, int k) {
66 
67  for (int dir = 0; dir < AMREX_SPACEDIM; dir++)
68  {
69  Face face = Face::INT;
70 
71  #if AMREX_SPACEDIM == 2
72 
73  if (i==lo.x && j==lo.y) face = Face::XLO_YLO;
74  else if (i==lo.x && j==hi.y) face = Face::XLO_YHI;
75  else if (i==hi.x && j==lo.y) face = Face::XHI_YLO;
76  else if (i==hi.x && j==hi.y) face = Face::XHI_YHI;
77 
78  else if (i==lo.x) face = Face::XLO;
79  else if (i==hi.x) face = Face::XHI;
80  else if (j==lo.y) face = Face::YLO;
81  else if (j==hi.y) face = Face::YHI;
82 
83  #elif AMREX_SPACEDIM == 3
84 
85  if (i==lo.x && j==lo.y && k==lo.z) face = Face::XLO_YLO_ZLO;
86  else if (i==lo.x && j==lo.y && k==hi.z) face = Face::XLO_YLO_ZHI;
87  else if (i==lo.x && j==hi.y && k==lo.z) face = Face::XLO_YHI_ZLO;
88  else if (i==lo.x && j==hi.y && k==hi.z) face = Face::XLO_YHI_ZHI;
89  else if (i==hi.x && j==lo.y && k==lo.z) face = Face::XHI_YLO_ZLO;
90  else if (i==hi.x && j==lo.y && k==hi.z) face = Face::XHI_YLO_ZHI;
91  else if (i==hi.x && j==hi.y && k==lo.z) face = Face::XHI_YHI_ZLO;
92  else if (i==hi.x && j==hi.y && k==hi.z) face = Face::XHI_YHI_ZHI;
93 
94  else if (j==lo.y && k==lo.z) face = Face::YLO_ZLO;
95  else if (j==lo.y && k==hi.z) face = Face::YLO_ZHI;
96  else if (j==hi.y && k==lo.z) face = Face::YHI_ZLO;
97  else if (j==hi.y && k==hi.z) face = Face::YHI_ZHI;
98  else if (k==lo.z && i==lo.x) face = Face::ZLO_XLO;
99  else if (k==lo.z && i==hi.x) face = Face::ZLO_XHI;
100  else if (k==hi.z && i==lo.x) face = Face::ZHI_XLO;
101  else if (k==hi.z && i==hi.x) face = Face::ZHI_XHI;
102  else if (i==lo.x && j==lo.y) face = Face::XLO_YLO;
103  else if (i==lo.x && j==hi.y) face = Face::XLO_YHI;
104  else if (i==hi.x && j==lo.y) face = Face::XHI_YLO;
105  else if (i==hi.x && j==hi.y) face = Face::XHI_YHI;
106 
107  else if (i==lo.x) face = Face::XLO;
108  else if (i==hi.x) face = Face::XHI;
109  else if (j==lo.y) face = Face::YLO;
110  else if (j==hi.y) face = Face::YHI;
111  else if (k==lo.z) face = Face::ZLO;
112  else if (k==hi.z) face = Face::ZHI;
113 
114  #endif
115 
116  amrex::IndexType type = a_rhs->ixType();
117  if (!(face == Face::INT))
118  {
119  if (a_homogeneous && m_bc_type[face][dir] == Type::Displacement)
120  rhs(i,j,k)(dir) = 0.0;
121  else
122  {
123  Set::Vector x = Set::Position(i, j, k, a_geom, type);
124  #if AMREX_SPACEDIM==1
125  rhs(i,j,k)(dir) = m_bc_func[face][dir](x(0),0.0,0.0,m_time);
126  #elif AMREX_SPACEDIM==2
127  rhs(i,j,k)(dir) = m_bc_func[face][dir](x(0),x(1),0.0,m_time);
128  #elif AMREX_SPACEDIM==3
129  rhs(i,j,k)(dir) = m_bc_func[face][dir](x(0),x(1),x(2),m_time);
130  #endif
131  }
132  }
133  }
134  });
135  }
136  }
137 
138 
139  virtual void
140  Init(amrex::MultiFab * a_rhs,
141  const amrex::Geometry &a_geom,
142  bool a_homogeneous = false) const override
143  {
144  amrex::Box domain(a_geom.Domain());
145  domain.convert(amrex::IntVect::TheNodeVector());
146  const amrex::Dim3 lo= amrex::lbound(domain), hi = amrex::ubound(domain);
147  for (amrex::MFIter mfi(*a_rhs, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
148  {
149  amrex::Box bx = mfi.tilebox();
150  bx.grow(2);
151  bx = bx & domain;
152  amrex::Array4<amrex::Real> const& rhs = a_rhs->array(mfi);
153  amrex::ParallelFor (bx,[=] AMREX_GPU_DEVICE(int i, int j, int k) {
154 
155  for (int dir = 0; dir < AMREX_SPACEDIM; dir++)
156  {
157  Face face = Face::INT;
158 
159  #if AMREX_SPACEDIM == 2
160 
161  if (i==lo.x && j==lo.y) face = Face::XLO_YLO;
162  else if (i==lo.x && j==hi.y) face = Face::XLO_YHI;
163  else if (i==hi.x && j==lo.y) face = Face::XHI_YLO;
164  else if (i==hi.x && j==hi.y) face = Face::XHI_YHI;
165 
166  else if (i==lo.x) face = Face::XLO;
167  else if (i==hi.x) face = Face::XHI;
168  else if (j==lo.y) face = Face::YLO;
169  else if (j==hi.y) face = Face::YHI;
170 
171  #elif AMREX_SPACEDIM == 3
172 
173  if (i==lo.x && j==lo.y && k==lo.z) face = Face::XLO_YLO_ZLO;
174  else if (i==lo.x && j==lo.y && k==hi.z) face = Face::XLO_YLO_ZHI;
175  else if (i==lo.x && j==hi.y && k==lo.z) face = Face::XLO_YHI_ZLO;
176  else if (i==lo.x && j==hi.y && k==hi.z) face = Face::XLO_YHI_ZHI;
177  else if (i==hi.x && j==lo.y && k==lo.z) face = Face::XHI_YLO_ZLO;
178  else if (i==hi.x && j==lo.y && k==hi.z) face = Face::XHI_YLO_ZHI;
179  else if (i==hi.x && j==hi.y && k==lo.z) face = Face::XHI_YHI_ZLO;
180  else if (i==hi.x && j==hi.y && k==hi.z) face = Face::XHI_YHI_ZHI;
181 
182  else if (j==lo.y && k==lo.z) face = Face::YLO_ZLO;
183  else if (j==lo.y && k==hi.z) face = Face::YLO_ZHI;
184  else if (j==hi.y && k==lo.z) face = Face::YHI_ZLO;
185  else if (j==hi.y && k==hi.z) face = Face::YHI_ZHI;
186  else if (k==lo.z && i==lo.x) face = Face::ZLO_XLO;
187  else if (k==lo.z && i==hi.x) face = Face::ZLO_XHI;
188  else if (k==hi.z && i==lo.x) face = Face::ZHI_XLO;
189  else if (k==hi.z && i==hi.x) face = Face::ZHI_XHI;
190  else if (i==lo.x && j==lo.y) face = Face::XLO_YLO;
191  else if (i==lo.x && j==hi.y) face = Face::XLO_YHI;
192  else if (i==hi.x && j==lo.y) face = Face::XHI_YLO;
193  else if (i==hi.x && j==hi.y) face = Face::XHI_YHI;
194 
195  else if (i==lo.x) face = Face::XLO;
196  else if (i==hi.x) face = Face::XHI;
197  else if (j==lo.y) face = Face::YLO;
198  else if (j==hi.y) face = Face::YHI;
199  else if (k==lo.z) face = Face::ZLO;
200  else if (k==hi.z) face = Face::ZHI;
201 
202  #endif
203 
204  amrex::IndexType type = a_rhs->ixType();
205  if (!(face == Face::INT))
206  {
207  if (a_homogeneous && m_bc_type[face][dir] == Type::Displacement)
208  rhs(i,j,k,dir) = 0.0;
209  else
210  {
211  Set::Vector x = Set::Position(i, j, k, a_geom, type);
212  #if AMREX_SPACEDIM==1
213  rhs(i,j,k,dir) = (m_bc_func[face][dir])(x(0),0.0,0.0,m_time);
214  #elif AMREX_SPACEDIM==2
215  rhs(i,j,k,dir) = (m_bc_func[face][dir])(x(0),x(1),0.0,m_time);
216  #elif AMREX_SPACEDIM==3
217  rhs(i,j,k,dir) = (m_bc_func[face][dir])(x(0),x(1),x(2),m_time);
218  #endif
219  }
220  }
221  }
222  });
223  }
224  }
225 
226 
227 
228  virtual
229  std::array<Type,AMREX_SPACEDIM> getType (
230  const int &i, const int &j, [[maybe_unused]] const int &k,
231  const amrex::Box &domain) override
232  {
233  (void)i, (void)j, (void)k;
234  amrex::IntVect m(AMREX_D_DECL(i,j,k));
235  const amrex::Dim3 lo = amrex::lbound(domain), hi = amrex::ubound(domain);
236 
237  // Corners
238  #if AMREX_SPACEDIM == 2
239  if (m[0] == lo.x && m[1] == lo.y) return m_bc_type[Face::XLO_YLO];
240  if (m[0] == lo.x && m[1] == hi.y) return m_bc_type[Face::XLO_YHI];
241  if (m[0] == hi.x && m[1] == lo.y) return m_bc_type[Face::XHI_YLO];
242  if (m[0] == hi.x && m[1] == hi.y) return m_bc_type[Face::XHI_YHI];
243  if (m[0] == lo.x) return m_bc_type[Face::XLO];
244  if (m[0] == hi.x) return m_bc_type[Face::XHI];
245  if (m[1] == lo.y) return m_bc_type[Face::YLO];
246  if (m[1] == hi.y) return m_bc_type[Face::YHI];
247 
248  #elif AMREX_SPACEDIM == 3
249 
250  if (m[0] == lo.x && m[1] == lo.y && m[2] == lo.z) return m_bc_type[Face::XLO_YLO_ZLO];
251  if (m[0] == lo.x && m[1] == lo.y && m[2] == hi.z) return m_bc_type[Face::XLO_YLO_ZHI];
252  if (m[0] == lo.x && m[1] == hi.y && m[2] == lo.z) return m_bc_type[Face::XLO_YHI_ZLO];
253  if (m[0] == lo.x && m[1] == hi.y && m[2] == hi.z) return m_bc_type[Face::XLO_YHI_ZHI];
254  if (m[0] == hi.x && m[1] == lo.y && m[2] == lo.z) return m_bc_type[Face::XHI_YLO_ZLO];
255  if (m[0] == hi.x && m[1] == lo.y && m[2] == hi.z) return m_bc_type[Face::XHI_YLO_ZHI];
256  if (m[0] == hi.x && m[1] == hi.y && m[2] == lo.z) return m_bc_type[Face::XHI_YHI_ZLO];
257  if (m[0] == hi.x && m[1] == hi.y && m[2] == hi.z) return m_bc_type[Face::XHI_YHI_ZHI];
258 
259  if (m[1] == lo.y && m[2] == lo.z) return m_bc_type[Face::YLO_ZLO];
260  if (m[1] == lo.y && m[2] == hi.z) return m_bc_type[Face::YLO_ZHI];
261  if (m[1] == hi.y && m[2] == lo.z) return m_bc_type[Face::YHI_ZLO];
262  if (m[1] == hi.y && m[2] == hi.z) return m_bc_type[Face::YHI_ZHI];
263  if (m[2] == lo.z && m[0] == lo.x) return m_bc_type[Face::ZLO_XLO];
264  if (m[2] == lo.z && m[0] == hi.x) return m_bc_type[Face::ZLO_XHI];
265  if (m[2] == hi.z && m[0] == lo.x) return m_bc_type[Face::ZHI_XLO];
266  if (m[2] == hi.z && m[0] == hi.x) return m_bc_type[Face::ZHI_XHI];
267  if (m[0] == lo.x && m[1] == lo.y) return m_bc_type[Face::XLO_YLO];
268  if (m[0] == lo.x && m[1] == hi.y) return m_bc_type[Face::XLO_YHI];
269  if (m[0] == hi.x && m[1] == lo.y) return m_bc_type[Face::XHI_YLO];
270  if (m[0] == hi.x && m[1] == hi.y) return m_bc_type[Face::XHI_YHI];
271 
272  if (m[0] == lo.x) return m_bc_type[Face::XLO];
273  if (m[0] == hi.x) return m_bc_type[Face::XHI];
274  if (m[1] == lo.y) return m_bc_type[Face::YLO];
275  if (m[1] == hi.y) return m_bc_type[Face::YHI];
276  if (m[2] == lo.z) return m_bc_type[Face::ZLO];
277  if (m[2] == hi.z) return m_bc_type[Face::ZHI];
278 
279  #endif
280 
282  }
283 
284 
285  AMREX_FORCE_INLINE
287  const Set::Matrix &gradu,
288  const Set::Matrix &sigma,
289  const int &i, const int &j, const int &k,
290  const amrex::Box &domain) override
291  {
292  (void)i; (void)j; (void)k; // Suppress "unused variable" warnings
293  //Set::Vector f;
294 
295  const amrex::Dim3 lo= amrex::lbound(domain), hi = amrex::ubound(domain);
296 
297  amrex::IntVect m(AMREX_D_DECL(i,j,k));
298 
299  // Corners
300  #if AMREX_SPACEDIM == 2
301 
302  if (m[0] == lo.x && m[1] == lo.y) return set(m_bc_type[Face::XLO_YLO], u, gradu, sigma, Set::Vector(-SQRT2INV, -SQRT2INV));
303  if (m[0] == lo.x && m[1] == hi.y) return set(m_bc_type[Face::XLO_YHI], u, gradu, sigma, Set::Vector(-SQRT2INV, +SQRT2INV));
304  if (m[0] == hi.x && m[1] == lo.y) return set(m_bc_type[Face::XHI_YLO], u, gradu, sigma, Set::Vector(+SQRT2INV, -SQRT2INV));
305  if (m[0] == hi.x && m[1] == hi.y) return set(m_bc_type[Face::XHI_YHI], u, gradu, sigma, Set::Vector(+SQRT2INV, +SQRT2INV));
306 
307  if (m[0] == lo.x) return set(m_bc_type[Face::XLO], u, gradu, sigma, Set::Vector(-1, 0));
308  if (m[0] == hi.x) return set(m_bc_type[Face::XHI], u, gradu, sigma, Set::Vector(+1, 0));
309  if (m[1] == lo.y) return set(m_bc_type[Face::YLO], u, gradu, sigma, Set::Vector( 0,-1));
310  if (m[1] == hi.y) return set(m_bc_type[Face::YHI], u, gradu, sigma, Set::Vector( 0,+1));
311 
312  #elif AMREX_SPACEDIM == 3
313 
314  if (m[0] == lo.x && m[1] == lo.y && m[2] == lo.z) return set(m_bc_type[Face::XLO_YLO_ZLO], u, gradu, sigma, Set::Vector(-SQRT3INV,-SQRT3INV,-SQRT3INV));
315  if (m[0] == lo.x && m[1] == lo.y && m[2] == hi.z) return set(m_bc_type[Face::XLO_YLO_ZHI], u, gradu, sigma, Set::Vector(-SQRT3INV,-SQRT3INV,+SQRT3INV));
316  if (m[0] == lo.x && m[1] == hi.y && m[2] == lo.z) return set(m_bc_type[Face::XLO_YHI_ZLO], u, gradu, sigma, Set::Vector(-SQRT3INV,+SQRT3INV,-SQRT3INV));
317  if (m[0] == lo.x && m[1] == hi.y && m[2] == hi.z) return set(m_bc_type[Face::XLO_YHI_ZHI], u, gradu, sigma, Set::Vector(-SQRT3INV,+SQRT3INV,+SQRT3INV));
318  if (m[0] == hi.x && m[1] == lo.y && m[2] == lo.z) return set(m_bc_type[Face::XHI_YLO_ZLO], u, gradu, sigma, Set::Vector(+SQRT3INV,-SQRT3INV,-SQRT3INV));
319  if (m[0] == hi.x && m[1] == lo.y && m[2] == hi.z) return set(m_bc_type[Face::XHI_YLO_ZHI], u, gradu, sigma, Set::Vector(+SQRT3INV,-SQRT3INV,+SQRT3INV));
320  if (m[0] == hi.x && m[1] == hi.y && m[2] == lo.z) return set(m_bc_type[Face::XHI_YHI_ZLO], u, gradu, sigma, Set::Vector(+SQRT3INV,+SQRT3INV,-SQRT3INV));
321  if (m[0] == hi.x && m[1] == hi.y && m[2] == hi.z) return set(m_bc_type[Face::XHI_YHI_ZHI], u, gradu, sigma, Set::Vector(+SQRT3INV,+SQRT3INV,+SQRT3INV));
322 
323  if (m[1] == lo.y && m[2] == lo.z) return set(m_bc_type[Face::YLO_ZLO], u, gradu, sigma, Set::Vector(0, -SQRT2INV,-SQRT2INV));
324  if (m[1] == lo.y && m[2] == hi.z) return set(m_bc_type[Face::YLO_ZHI], u, gradu, sigma, Set::Vector(0, -SQRT2INV,+SQRT2INV));
325  if (m[1] == hi.y && m[2] == lo.z) return set(m_bc_type[Face::YHI_ZLO], u, gradu, sigma, Set::Vector(0, +SQRT2INV,-SQRT2INV));
326  if (m[1] == hi.y && m[2] == hi.z) return set(m_bc_type[Face::YHI_ZHI], u, gradu, sigma, Set::Vector(0, +SQRT2INV,+SQRT2INV));
327  if (m[2] == lo.z && m[0] == lo.x) return set(m_bc_type[Face::ZLO_XLO], u, gradu, sigma, Set::Vector(-SQRT2INV, 0, -SQRT2INV));
328  if (m[2] == lo.z && m[0] == hi.x) return set(m_bc_type[Face::ZLO_XHI], u, gradu, sigma, Set::Vector(+SQRT2INV, 0, -SQRT2INV));
329  if (m[2] == hi.z && m[0] == lo.x) return set(m_bc_type[Face::ZHI_XLO], u, gradu, sigma, Set::Vector(-SQRT2INV, 0, +SQRT2INV));
330  if (m[2] == hi.z && m[0] == hi.x) return set(m_bc_type[Face::ZHI_XHI], u, gradu, sigma, Set::Vector(+SQRT2INV, 0, +SQRT2INV));
331  if (m[0] == lo.x && m[1] == lo.y) return set(m_bc_type[Face::XLO_YLO], u, gradu, sigma, Set::Vector(-SQRT2INV, -SQRT2INV,0 ));
332  if (m[0] == lo.x && m[1] == hi.y) return set(m_bc_type[Face::XLO_YHI], u, gradu, sigma, Set::Vector(-SQRT2INV, +SQRT2INV,0 ));
333  if (m[0] == hi.x && m[1] == lo.y) return set(m_bc_type[Face::XHI_YLO], u, gradu, sigma, Set::Vector(+SQRT2INV, -SQRT2INV,0 ));
334  if (m[0] == hi.x && m[1] == hi.y) return set(m_bc_type[Face::XHI_YHI], u, gradu, sigma, Set::Vector(+SQRT2INV, +SQRT2INV,0 ));
335 
336  if (m[0] == lo.x) return set(m_bc_type[Face::XLO], u, gradu, sigma, Set::Vector(-1, 0, 0));
337  if (m[0] == hi.x) return set(m_bc_type[Face::XHI], u, gradu, sigma, Set::Vector(+1, 0, 0));
338  if (m[1] == lo.y) return set(m_bc_type[Face::YLO], u, gradu, sigma, Set::Vector( 0,-1, 0));
339  if (m[1] == hi.y) return set(m_bc_type[Face::YHI], u, gradu, sigma, Set::Vector( 0,+1, 0));
340  if (m[2] == lo.z) return set(m_bc_type[Face::ZLO], u, gradu, sigma, Set::Vector( 0, 0,-1));
341  if (m[2] == hi.z) return set(m_bc_type[Face::ZHI], u, gradu, sigma, Set::Vector( 0, 0,+1));
342 
343  #endif
344 
345  Util::Abort(INFO,"Boundary condition error");
346  return Set::Vector::Zero();
347  }
348 
349  AMREX_FORCE_INLINE
350  Set::Vector set(std::array<Type,AMREX_SPACEDIM> &bc_type,
351  const Set::Vector &u, const Set::Matrix &gradu, const Set::Matrix &sigma, Set::Vector n) const
352  {
353  Set::Vector f = Set::Vector::Zero();
354  for (int i = 0; i < AMREX_SPACEDIM; i++)
355  {
356  if (bc_type[i] == Type::Displacement)
357  f(i) = u(i);
358  else if (bc_type[i] == Type::Traction)
359  f(i) = (sigma*n)(i);
360  else if (bc_type[i] == Type::Neumann)
361  f(i) = (gradu*n)(i);
362  else if (bc_type[i] == Periodic)
363  continue;
364  }
365  return f;
366  }
367 
368 protected:
369 
370  #if AMREX_SPACEDIM==2
371  static const int m_nfaces = 8;
372  #elif AMREX_SPACEDIM==3
373  static const int m_nfaces = 26;
374  #endif
375 
376  std::array<std::array<Type, AMREX_SPACEDIM>, m_nfaces> m_bc_type;
377  //std::array<std::array<FunctionParserAD,AMREX_SPACEDIM>, m_nfaces> m_bc_func;
378  std::array<std::array<amrex::Parser,AMREX_SPACEDIM>, m_nfaces> m_bc_func_parser;
379  std::array<std::array<amrex::ParserExecutor<4>,AMREX_SPACEDIM>, m_nfaces> m_bc_func;
380 
381 public:
382 
383  //
384  // This is basically the same as :ref:`BC::Operator::Elastic::Constant` except that
385  // you can use dynamically compiled expressions in space and time to define values.
386  //
387  // Usage is exactly the same except that the "val" inputs can depend on x, y, z, and t.
388  //
389  static void Parse(Expression & value, IO::ParmParse & pp)
390  {
391  std::map<std::string, Type> bcmap;
392  bcmap["displacement"] = Type::Displacement;
393  bcmap["disp"] = Type::Displacement;
394  bcmap["neumann"] = Type::Neumann;
395  bcmap["traction"] = Type::Traction;
396  bcmap["trac"] = Type::Traction;
397  bcmap["periodic"] = Type::Periodic;
398 
399  std::vector<std::string> str;
400 
401  // TYPES
402 
403  #if AMREX_SPACEDIM==3
404  if (pp.contains("type.xloylozlo")) { pp_queryarr("type.xloylozlo",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XLO_YLO_ZLO][i] = bcmap[str[i]]; } // 3D Corner
405  if (pp.contains("type.xloylozhi")) { pp_queryarr("type.xloylozhi",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XLO_YLO_ZHI][i] = bcmap[str[i]]; } // 3D Corner
406  if (pp.contains("type.xloyhizlo")) { pp_queryarr("type.xloyhizlo",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XLO_YHI_ZLO][i] = bcmap[str[i]]; } // 3D Corner
407  if (pp.contains("type.xloyhizhi")) { pp_queryarr("type.xloyhizhi",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XLO_YHI_ZHI][i] = bcmap[str[i]]; } // 3D Corner
408  if (pp.contains("type.xhiylozlo")) { pp_queryarr("type.xhiylozlo",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XHI_YLO_ZLO][i] = bcmap[str[i]]; } // 3D Corner
409  if (pp.contains("type.xhiylozhi")) { pp_queryarr("type.xhiylozhi",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XHI_YLO_ZHI][i] = bcmap[str[i]]; } // 3D Corner
410  if (pp.contains("type.xhiyhizlo")) { pp_queryarr("type.xhiyhizlo",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XHI_YHI_ZLO][i] = bcmap[str[i]]; } // 3D Corner
411  if (pp.contains("type.xhiyhizhi")) { pp_queryarr("type.xhiyhizhi",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XHI_YHI_ZHI][i] = bcmap[str[i]]; } // 3D Corner
412  #endif
413 
414  #if AMREX_SPACEDIM==3
415  if (pp.contains("type.ylozlo")) { pp_queryarr("type.ylozlo",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::YLO_ZLO][i] = bcmap[str[i]]; } // 3D Edge
416  if (pp.contains("type.ylozhi")) { pp_queryarr("type.ylozhi",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::YLO_ZHI][i] = bcmap[str[i]]; } // 3D Edge
417  if (pp.contains("type.yhizlo")) { pp_queryarr("type.yhizlo",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::YHI_ZLO][i] = bcmap[str[i]]; } // 3D Edge
418  if (pp.contains("type.yhizhi")) { pp_queryarr("type.yhizhi",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::YHI_ZHI][i] = bcmap[str[i]]; } // 3D Edge
419  if (pp.contains("type.zloxlo")) { pp_queryarr("type.zloxlo",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::ZLO_XLO][i] = bcmap[str[i]]; } // 3D Edge
420  if (pp.contains("type.zloxhi")) { pp_queryarr("type.zloxhi",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::ZLO_XHI][i] = bcmap[str[i]]; } // 3D Edge
421  if (pp.contains("type.zhixlo")) { pp_queryarr("type.zhixlo",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::ZHI_XLO][i] = bcmap[str[i]]; } // 3D Edge
422  if (pp.contains("type.zhixhi")) { pp_queryarr("type.zhixhi",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::ZHI_XHI][i] = bcmap[str[i]]; } // 3D Edge
423  #endif
424  if (pp.contains("type.xloylo")) { pp_queryarr("type.xloylo",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XLO_YLO][i] = bcmap[str[i]]; } // 3D Edge / 2D Corner
425  if (pp.contains("type.xloyhi")) { pp_queryarr("type.xloyhi",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XLO_YHI][i] = bcmap[str[i]]; } // 3D Edge / 2D Corner
426  if (pp.contains("type.xhiylo")) { pp_queryarr("type.xhiylo",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XHI_YLO][i] = bcmap[str[i]]; } // 3D Edge / 2D Corner
427  if (pp.contains("type.xhiyhi")) { pp_queryarr("type.xhiyhi",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XHI_YHI][i] = bcmap[str[i]]; } // 3D Edge / 2D Corner
428 
429  if (pp.contains("type.xlo")) { pp_queryarr("type.xlo",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XLO][i] = bcmap[str[i]]; } // 3D Face 2D Edge
430  if (pp.contains("type.xhi")) { pp_queryarr("type.xhi",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XHI][i] = bcmap[str[i]]; } // 3D Face 2D Edge
431  if (pp.contains("type.ylo")) { pp_queryarr("type.ylo",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::YLO][i] = bcmap[str[i]]; } // 3D Face 2D Edge
432  if (pp.contains("type.yhi")) { pp_queryarr("type.yhi",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::YHI][i] = bcmap[str[i]]; } // 3D Face 2D Edge
433  #if AMREX_SPACEDIM==3
434  if (pp.contains("type.zlo")) { pp_queryarr("type.zlo",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::ZLO][i] = bcmap[str[i]]; } // 3D Face
435  if (pp.contains("type.zhi")) { pp_queryarr("type.zhi",str); for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::ZHI][i] = bcmap[str[i]]; } // 3D Face
436  #endif
437 
438  // VALS
439 
440  for (int face = 0; face != Face::INT; face++)
441  {
442  std::vector<std::string> val(AMREX_SPACEDIM,"0.0");
443  val.resize(AMREX_SPACEDIM,"0.0");
444  std::string querystr = std::string("val.") + std::string(facestr[face]);
445  pp_queryarr(querystr.c_str(),val);
446  for (int i = 0 ; i < AMREX_SPACEDIM; i++)
447  {
448  value.m_bc_func_parser[face][i] = amrex::Parser(val[i].c_str());
449  value.m_bc_func_parser[face][i].registerVariables({"x","y","z","t"});
450  value.m_bc_func[face][i] = value.m_bc_func_parser[face][i].compile<4>();
451  }
452  }
453 
454 #if AMREX_SPACEDIM==2
455  // We may wish to use an input file that has 3D BC inputs
456  // This will prevent the parser from complaining that there are unused inputs.
457  std::vector<std::string> ignore_in_2d = {
458  "zlo","zhi",
459  "zhixlo","zloxlo","zhixhi","zloxhi","ylozlo","ylozhi","yhizlo","yhizhi",
460  "xloylozlo","xloylozhi","xloyhizlo","xloyhizhi","xhiylozlo","xhiylozhi","xhiyhizlo","xhiyhizhi"};
461  for (unsigned int n = 0; n < ignore_in_2d.size(); n++)
462  {
463  std::string querystr = std::string("val.") + ignore_in_2d[n];
464  pp.ignore(querystr);
465  querystr = std::string("type.") + ignore_in_2d[n];
466  pp.ignore(querystr);
467  }
468 #endif
469 
470 
471  }
472  Expression(IO::ParmParse &pp, std::string name): Expression()
473  {pp_queryclass(name,*this);}
474 };
475 }
476 }
477 }
478 #endif
Model::Solid::gradu
@ gradu
Definition: Solid.H:26
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
Linear.H
Operator
Documentation for operator namespace.
Definition: Diagonal.cpp:13
BC::Operator::Elastic::Expression
Definition: Expression.H:17
BC::Operator::Elastic::Expression::Init
virtual void Init(amrex::FabArray< amrex::BaseFab< Set::Vector >> *a_rhs, const amrex::Geometry &a_geom, bool a_homogeneous=false) const override
Definition: Expression.H:52
SQRT2INV
#define SQRT2INV
Definition: Elastic.H:78
BC::Operator::Elastic::Expression::m_bc_func_parser
std::array< std::array< amrex::Parser, AMREX_SPACEDIM >, m_nfaces > m_bc_func_parser
Definition: Expression.H:378
BC::Operator::Elastic::Elastic::m_time
Set::Scalar m_time
Definition: Elastic.H:93
BC::Operator::Elastic::Expression::m_bc_type
std::array< std::array< Type, AMREX_SPACEDIM >, m_nfaces > m_bc_type
Definition: Expression.H:376
Set::Vector
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition: Base.H:20
ParmParse.H
Set::None
@ None
Definition: Base.H:197
BC::Operator::Elastic::Elastic::AMREX_D_DECL
@ AMREX_D_DECL
Definition: Elastic.H:44
BC::Operator::Elastic::Elastic::Periodic
@ Periodic
Definition: Elastic.H:24
SQRT3INV
#define SQRT3INV
Definition: Elastic.H:77
pp_queryclass
#define pp_queryclass(...)
Definition: ParmParse.H:105
Set::Matrix
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, AMREX_SPACEDIM > Matrix
Definition: Base.H:23
Set::Face
@ Face
Definition: Set.H:32
IO::ParmParse::contains
bool contains(std::string name)
Definition: ParmParse.H:151
BC::Operator::Elastic::Expression::Expression
Expression(IO::ParmParse &pp, std::string name)
Definition: Expression.H:472
BC::Operator::Elastic::Elastic
Definition: Elastic.H:18
BC::Operator::Elastic::Expression::getType
virtual std::array< Type, AMREX_SPACEDIM > getType(const int &i, const int &j, [[maybe_unused]] const int &k, const amrex::Box &domain) override
Definition: Expression.H:229
Util::Abort
void Abort(const char *msg)
Definition: Util.cpp:165
BC::Operator::Elastic::Expression::~Expression
~Expression()
Definition: Expression.H:46
BC
Collection of boundary condition (BC) objects.
Definition: BC.cpp:4
Elastic.H
IO::ParmParse
Definition: ParmParse.H:110
BC::Operator::Elastic::Elastic::Init
virtual void Init(amrex::MultiFab *a_rhs, const amrex::Geometry &a_geom, bool a_homogeneous=false) const =0
BC::Operator::Elastic::Expression::Init
virtual void Init(amrex::MultiFab *a_rhs, const amrex::Geometry &a_geom, bool a_homogeneous=false) const override
Definition: Expression.H:140
IO::ParmParse::ignore
void ignore(std::string name)
Definition: ParmParse.H:132
BC::Operator::Elastic::Expression::Parse
static void Parse(Expression &value, IO::ParmParse &pp)
Definition: Expression.H:389
BC::Operator::Elastic::Expression::Expression
Expression()
Definition: Expression.H:37
BC::Operator::Elastic::Expression::operator()
AMREX_FORCE_INLINE Set::Vector operator()(const Set::Vector &u, const Set::Matrix &gradu, const Set::Matrix &sigma, const int &i, const int &j, const int &k, const amrex::Box &domain) override
Definition: Expression.H:286
INFO
#define INFO
Definition: Util.H:20
BC::Operator::Elastic::Expression::m_bc_func
std::array< std::array< amrex::ParserExecutor< 4 >, AMREX_SPACEDIM >, m_nfaces > m_bc_func
Definition: Expression.H:379
BC::Operator::Elastic::Expression::set
AMREX_FORCE_INLINE Set::Vector set(std::array< Type, AMREX_SPACEDIM > &bc_type, const Set::Vector &u, const Set::Matrix &gradu, const Set::Matrix &sigma, Set::Vector n) const
Definition: Expression.H:350
pp_queryarr
#define pp_queryarr(...)
Definition: ParmParse.H:103