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