Alamo
Constant.H
Go to the documentation of this file.
1 //
2 // This BC defines boundary conditions that are constant with respect to space.
3 // However they may change in time using the :ref:`Numeric::Interpolator::Linear` method.
4 //
5 // In 2D, BCs are prescribed for each edge (4) and corner (4) on the boundary, for a total of 8 possible regions.
6 // In 3D, BCs are prescribed for each face (6), each edge (12), and each corner (8), for a total of 26 possible regions.
7 // Care must be taken to define BCs for edges/corners consistently with the faces/edges, or you will get poor
8 // convergence and inaccurate behavior.
9 // (See :ref:`BC::Operator::Elastic::TensionTest` for a reduced BC with streamlined load options.)
10 //
11 // To define a boundary condition, you must define both a "type" and a "value" in each direction.
12 // The type can be "displacement", "disp", "neumann", "traction", "trac".
13 // The value is the corresponding value.
14 // All BCs are, by default, dirichlet (displacement) with a value of zero.
15 //
16 
17 #ifndef BC_OPERATOR_ELASTIC_CONSTANT_H
18 #define BC_OPERATOR_ELASTIC_CONSTANT_H
19 
20 // #include "Operator/Elastic.H"
21 #include "IO/ParmParse.H"
24 
25 namespace BC
26 {
27 namespace Operator
28 {
29 namespace Elastic
30 {
31 class Constant : public Elastic
32 {
33 public:
34  static constexpr const char* name = "constant";
35 
37  {
38  // By default, all boundary conditions are displacement
39  for (int face = 0; face < m_nfaces; face++)
40  for (int direction = 0; direction < AMREX_SPACEDIM; direction++)
41  {
42  m_bc_type[face][direction] = Type::Displacement;
43  m_bc_val [face][direction] = 0.0;
44  }
45  };
46 
47  Constant(IO::ParmParse &pp, std::string name): Constant()
48  { pp_queryclass(name,*this); }
49 
50  ~Constant() {};
51 
52 
53  void
54  Set(const Face face,
55  const Direction direction,
56  const Type type,
57  const Set::Scalar value)
58  {
59  m_bc_type[face][direction] = type;
60  m_bc_val [face][direction] = value;
61  }
62 
63  void
64  Set(const Face /*face*/,
65  const Direction /*direction*/,
66  const Type /*type*/,
67  const Set::Scalar /*value*/,
68  amrex::Vector<amrex::MultiFab *> &/*a_rhs*/,
69  const amrex::Vector<amrex::Geometry> &/*a_geom*/)
70  {
71  Util::Abort(INFO,"This has been depricated");
72  }
73 
74  void
75  Set(const Face face,
76  const Direction direction,
77  const Type type,
78  const Set::Scalar value,
79  amrex::Vector<amrex::MultiFab> &a_rhs,
80  const amrex::Vector<amrex::Geometry> &a_geom)
81  {
82  amrex::Vector<amrex::MultiFab *> pa_rhs = amrex::GetVecOfPtrs(a_rhs);
83  Set(face,direction,type,value,pa_rhs,a_geom);
84  }
85 
86  void
87  Set(const Face face,
88  const Direction direction,
89  const Type type,
90  const Set::Scalar value,
91  amrex::Vector<std::unique_ptr<amrex::MultiFab> > &a_rhs,
92  const amrex::Vector<amrex::Geometry> &a_geom)
93  {
94  amrex::Vector<amrex::MultiFab *> pa_rhs = amrex::GetVecOfPtrs(a_rhs);
95  Set(face,direction,type,value,pa_rhs,a_geom);
96  }
97 
98  //using Elastic::Init;
99 
100  virtual void
101  Init(amrex::FabArray<amrex::BaseFab<Set::Vector>> *a_rhs,
102  const amrex::Geometry &a_geom,
103  bool a_homogeneous = false) const
104  {
105  amrex::Box domain(a_geom.Domain());
106  domain.convert(amrex::IntVect::TheNodeVector());
107  const amrex::Dim3 lo= amrex::lbound(domain), hi = amrex::ubound(domain);
108  for (amrex::MFIter mfi(*a_rhs, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
109  {
110  amrex::Box bx = mfi.tilebox();
111  bx.grow(2);
112  bx = bx & domain;
113  amrex::Array4<Set::Vector> const& rhs = a_rhs->array(mfi);
114  amrex::ParallelFor (bx,[=] AMREX_GPU_DEVICE(int i, int j, int k) {
115 
116  for (int dir = 0; dir < AMREX_SPACEDIM; dir++)
117  {
118  Face face = Face::INT;
119 
120  #if AMREX_SPACEDIM == 2
121 
122  if (i==lo.x && j==lo.y) face = Face::XLO_YLO;
123  else if (i==lo.x && j==hi.y) face = Face::XLO_YHI;
124  else if (i==hi.x && j==lo.y) face = Face::XHI_YLO;
125  else if (i==hi.x && j==hi.y) face = Face::XHI_YHI;
126 
127  else if (i==lo.x) face = Face::XLO;
128  else if (i==hi.x) face = Face::XHI;
129  else if (j==lo.y) face = Face::YLO;
130  else if (j==hi.y) face = Face::YHI;
131 
132  #elif AMREX_SPACEDIM == 3
133 
134  if (i==lo.x && j==lo.y && k==lo.z) face = Face::XLO_YLO_ZLO;
135  else if (i==lo.x && j==lo.y && k==hi.z) face = Face::XLO_YLO_ZHI;
136  else if (i==lo.x && j==hi.y && k==lo.z) face = Face::XLO_YHI_ZLO;
137  else if (i==lo.x && j==hi.y && k==hi.z) face = Face::XLO_YHI_ZHI;
138  else if (i==hi.x && j==lo.y && k==lo.z) face = Face::XHI_YLO_ZLO;
139  else if (i==hi.x && j==lo.y && k==hi.z) face = Face::XHI_YLO_ZHI;
140  else if (i==hi.x && j==hi.y && k==lo.z) face = Face::XHI_YHI_ZLO;
141  else if (i==hi.x && j==hi.y && k==hi.z) face = Face::XHI_YHI_ZHI;
142 
143  else if (j==lo.y && k==lo.z) face = Face::YLO_ZLO;
144  else if (j==lo.y && k==hi.z) face = Face::YLO_ZHI;
145  else if (j==hi.y && k==lo.z) face = Face::YHI_ZLO;
146  else if (j==hi.y && k==hi.z) face = Face::YHI_ZHI;
147  else if (k==lo.z && i==lo.x) face = Face::ZLO_XLO;
148  else if (k==lo.z && i==hi.x) face = Face::ZLO_XHI;
149  else if (k==hi.z && i==lo.x) face = Face::ZHI_XLO;
150  else if (k==hi.z && i==hi.x) face = Face::ZHI_XHI;
151  else if (i==lo.x && j==lo.y) face = Face::XLO_YLO;
152  else if (i==lo.x && j==hi.y) face = Face::XLO_YHI;
153  else if (i==hi.x && j==lo.y) face = Face::XHI_YLO;
154  else if (i==hi.x && j==hi.y) face = Face::XHI_YHI;
155 
156  else if (i==lo.x) face = Face::XLO;
157  else if (i==hi.x) face = Face::XHI;
158  else if (j==lo.y) face = Face::YLO;
159  else if (j==hi.y) face = Face::YHI;
160  else if (k==lo.z) face = Face::ZLO;
161  else if (k==hi.z) face = Face::ZHI;
162 
163  #endif
164 
165  if (!(face == Face::INT))
166  {
167  if (a_homogeneous && m_bc_type[face][dir] == Type::Displacement)
168  rhs(i,j,k)(dir) = 0.0;
169  else
170  rhs(i,j,k)(dir) = m_bc_val[face][dir](m_time);
171  }
172  }
173  });
174  }
175  }
176 
177  using Elastic::Init;
178  virtual void
179  Init(amrex::MultiFab * a_rhs,
180  const amrex::Geometry &a_geom,
181  bool a_homogeneous = false) const
182  {
183  amrex::Box domain(a_geom.Domain());
184  domain.convert(amrex::IntVect::TheNodeVector());
185  const amrex::Dim3 lo= amrex::lbound(domain), hi = amrex::ubound(domain);
186  for (amrex::MFIter mfi(*a_rhs, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
187  {
188  amrex::Box bx = mfi.tilebox();
189  bx.grow(2);
190  bx = bx & domain;
191  amrex::Array4<amrex::Real> const& rhs = a_rhs->array(mfi);
192  amrex::ParallelFor (bx,[=] AMREX_GPU_DEVICE(int i, int j, int k) {
193 
194  for (int dir = 0; dir < AMREX_SPACEDIM; dir++)
195  {
196  Face face = Face::INT;
197 
198  #if AMREX_SPACEDIM == 2
199 
200  if (i==lo.x && j==lo.y) face = Face::XLO_YLO;
201  else if (i==lo.x && j==hi.y) face = Face::XLO_YHI;
202  else if (i==hi.x && j==lo.y) face = Face::XHI_YLO;
203  else if (i==hi.x && j==hi.y) face = Face::XHI_YHI;
204 
205  else if (i==lo.x) face = Face::XLO;
206  else if (i==hi.x) face = Face::XHI;
207  else if (j==lo.y) face = Face::YLO;
208  else if (j==hi.y) face = Face::YHI;
209 
210  #elif AMREX_SPACEDIM == 3
211 
212  if (i==lo.x && j==lo.y && k==lo.z) face = Face::XLO_YLO_ZLO;
213  else if (i==lo.x && j==lo.y && k==hi.z) face = Face::XLO_YLO_ZHI;
214  else if (i==lo.x && j==hi.y && k==lo.z) face = Face::XLO_YHI_ZLO;
215  else if (i==lo.x && j==hi.y && k==hi.z) face = Face::XLO_YHI_ZHI;
216  else if (i==hi.x && j==lo.y && k==lo.z) face = Face::XHI_YLO_ZLO;
217  else if (i==hi.x && j==lo.y && k==hi.z) face = Face::XHI_YLO_ZHI;
218  else if (i==hi.x && j==hi.y && k==lo.z) face = Face::XHI_YHI_ZLO;
219  else if (i==hi.x && j==hi.y && k==hi.z) face = Face::XHI_YHI_ZHI;
220 
221  else if (j==lo.y && k==lo.z) face = Face::YLO_ZLO;
222  else if (j==lo.y && k==hi.z) face = Face::YLO_ZHI;
223  else if (j==hi.y && k==lo.z) face = Face::YHI_ZLO;
224  else if (j==hi.y && k==hi.z) face = Face::YHI_ZHI;
225  else if (k==lo.z && i==lo.x) face = Face::ZLO_XLO;
226  else if (k==lo.z && i==hi.x) face = Face::ZLO_XHI;
227  else if (k==hi.z && i==lo.x) face = Face::ZHI_XLO;
228  else if (k==hi.z && i==hi.x) face = Face::ZHI_XHI;
229  else if (i==lo.x && j==lo.y) face = Face::XLO_YLO;
230  else if (i==lo.x && j==hi.y) face = Face::XLO_YHI;
231  else if (i==hi.x && j==lo.y) face = Face::XHI_YLO;
232  else if (i==hi.x && j==hi.y) face = Face::XHI_YHI;
233 
234  else if (i==lo.x) face = Face::XLO;
235  else if (i==hi.x) face = Face::XHI;
236  else if (j==lo.y) face = Face::YLO;
237  else if (j==hi.y) face = Face::YHI;
238  else if (k==lo.z) face = Face::ZLO;
239  else if (k==hi.z) face = Face::ZHI;
240 
241  #endif
242 
243  if (!(face == Face::INT))
244  {
245  if (a_homogeneous && m_bc_type[face][dir] == Type::Displacement)
246  rhs(i,j,k,dir) = 0.0;
247  else
248  rhs(i,j,k,dir) = m_bc_val[face][dir](m_time);
249  }
250  }
251  });
252  }
253  }
254 
255  virtual
256  std::array<Type,AMREX_SPACEDIM> getType (
257  const int &i, const int &j, [[maybe_unused]] const int &k,
258  const amrex::Box &domain)
259  {
260  amrex::IntVect m(AMREX_D_DECL(i,j,k));
261  const amrex::Dim3 lo = amrex::lbound(domain), hi = amrex::ubound(domain);
262 
263  // Corners
264  #if AMREX_SPACEDIM == 2
265  if (m[0] == lo.x && m[1] == lo.y) return m_bc_type[Face::XLO_YLO];
266  if (m[0] == lo.x && m[1] == hi.y) return m_bc_type[Face::XLO_YHI];
267  if (m[0] == hi.x && m[1] == lo.y) return m_bc_type[Face::XHI_YLO];
268  if (m[0] == hi.x && m[1] == hi.y) return m_bc_type[Face::XHI_YHI];
269  if (m[0] == lo.x) return m_bc_type[Face::XLO];
270  if (m[0] == hi.x) return m_bc_type[Face::XHI];
271  if (m[1] == lo.y) return m_bc_type[Face::YLO];
272  if (m[1] == hi.y) return m_bc_type[Face::YHI];
273 
274  #elif AMREX_SPACEDIM == 3
275 
276  if (m[0] == lo.x && m[1] == lo.y && m[2] == lo.z) return m_bc_type[Face::XLO_YLO_ZLO];
277  if (m[0] == lo.x && m[1] == lo.y && m[2] == hi.z) return m_bc_type[Face::XLO_YLO_ZHI];
278  if (m[0] == lo.x && m[1] == hi.y && m[2] == lo.z) return m_bc_type[Face::XLO_YHI_ZLO];
279  if (m[0] == lo.x && m[1] == hi.y && m[2] == hi.z) return m_bc_type[Face::XLO_YHI_ZHI];
280  if (m[0] == hi.x && m[1] == lo.y && m[2] == lo.z) return m_bc_type[Face::XHI_YLO_ZLO];
281  if (m[0] == hi.x && m[1] == lo.y && m[2] == hi.z) return m_bc_type[Face::XHI_YLO_ZHI];
282  if (m[0] == hi.x && m[1] == hi.y && m[2] == lo.z) return m_bc_type[Face::XHI_YHI_ZLO];
283  if (m[0] == hi.x && m[1] == hi.y && m[2] == hi.z) return m_bc_type[Face::XHI_YHI_ZHI];
284 
285  if (m[1] == lo.y && m[2] == lo.z) return m_bc_type[Face::YLO_ZLO];
286  if (m[1] == lo.y && m[2] == hi.z) return m_bc_type[Face::YLO_ZHI];
287  if (m[1] == hi.y && m[2] == lo.z) return m_bc_type[Face::YHI_ZLO];
288  if (m[1] == hi.y && m[2] == hi.z) return m_bc_type[Face::YHI_ZHI];
289  if (m[2] == lo.z && m[0] == lo.x) return m_bc_type[Face::ZLO_XLO];
290  if (m[2] == lo.z && m[0] == hi.x) return m_bc_type[Face::ZLO_XHI];
291  if (m[2] == hi.z && m[0] == lo.x) return m_bc_type[Face::ZHI_XLO];
292  if (m[2] == hi.z && m[0] == hi.x) return m_bc_type[Face::ZHI_XHI];
293  if (m[0] == lo.x && m[1] == lo.y) return m_bc_type[Face::XLO_YLO];
294  if (m[0] == lo.x && m[1] == hi.y) return m_bc_type[Face::XLO_YHI];
295  if (m[0] == hi.x && m[1] == lo.y) return m_bc_type[Face::XHI_YLO];
296  if (m[0] == hi.x && m[1] == hi.y) return m_bc_type[Face::XHI_YHI];
297 
298  if (m[0] == lo.x) return m_bc_type[Face::XLO];
299  if (m[0] == hi.x) return m_bc_type[Face::XHI];
300  if (m[1] == lo.y) return m_bc_type[Face::YLO];
301  if (m[1] == hi.y) return m_bc_type[Face::YHI];
302  if (m[2] == lo.z) return m_bc_type[Face::ZLO];
303  if (m[2] == hi.z) return m_bc_type[Face::ZHI];
304 
305  #endif
306 
308  }
309 
310 
311  AMREX_FORCE_INLINE
313  const Set::Matrix &gradu,
314  const Set::Matrix &sigma,
315  const int &i, const int &j, const int &k,
316  const amrex::Box &domain)
317  {
318  (void)i; (void)j; (void)k; // Suppress "unused variable" warnings
319  //Set::Vector f;
320 
321  const amrex::Dim3 lo= amrex::lbound(domain), hi = amrex::ubound(domain);
322 
323  amrex::IntVect m(AMREX_D_DECL(i,j,k));
324 
325  // Corners
326  #if AMREX_SPACEDIM == 2
327 
328  if (m[0] == lo.x && m[1] == lo.y) return set(m_bc_type[Face::XLO_YLO], u, gradu, sigma, Set::Vector(-SQRT2INV, -SQRT2INV));
329  if (m[0] == lo.x && m[1] == hi.y) return set(m_bc_type[Face::XLO_YHI], u, gradu, sigma, Set::Vector(-SQRT2INV, +SQRT2INV));
330  if (m[0] == hi.x && m[1] == lo.y) return set(m_bc_type[Face::XHI_YLO], u, gradu, sigma, Set::Vector(+SQRT2INV, -SQRT2INV));
331  if (m[0] == hi.x && m[1] == hi.y) return set(m_bc_type[Face::XHI_YHI], u, gradu, sigma, Set::Vector(+SQRT2INV, +SQRT2INV));
332 
333  if (m[0] == lo.x) return set(m_bc_type[Face::XLO], u, gradu, sigma, Set::Vector(-1, 0));
334  if (m[0] == hi.x) return set(m_bc_type[Face::XHI], u, gradu, sigma, Set::Vector(+1, 0));
335  if (m[1] == lo.y) return set(m_bc_type[Face::YLO], u, gradu, sigma, Set::Vector( 0,-1));
336  if (m[1] == hi.y) return set(m_bc_type[Face::YHI], u, gradu, sigma, Set::Vector( 0,+1));
337 
338  #elif AMREX_SPACEDIM == 3
339 
340  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));
341  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));
342  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));
343  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));
344  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));
345  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));
346  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));
347  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));
348 
349  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));
350  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));
351  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));
352  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));
353  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));
354  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));
355  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));
356  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));
357  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 ));
358  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 ));
359  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 ));
360  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 ));
361 
362  if (m[0] == lo.x) return set(m_bc_type[Face::XLO], u, gradu, sigma, Set::Vector(-1, 0, 0));
363  if (m[0] == hi.x) return set(m_bc_type[Face::XHI], u, gradu, sigma, Set::Vector(+1, 0, 0));
364  if (m[1] == lo.y) return set(m_bc_type[Face::YLO], u, gradu, sigma, Set::Vector( 0,-1, 0));
365  if (m[1] == hi.y) return set(m_bc_type[Face::YHI], u, gradu, sigma, Set::Vector( 0,+1, 0));
366  if (m[2] == lo.z) return set(m_bc_type[Face::ZLO], u, gradu, sigma, Set::Vector( 0, 0,-1));
367  if (m[2] == hi.z) return set(m_bc_type[Face::ZHI], u, gradu, sigma, Set::Vector( 0, 0,+1));
368 
369  #endif
370 
371  Util::Abort(INFO,"Boundary condition error");
372  return Set::Vector::Zero();
373  }
374 
375  AMREX_FORCE_INLINE
376  Set::Vector set(std::array<Type,AMREX_SPACEDIM> &bc_type,
377  const Set::Vector &u, const Set::Matrix &gradu, const Set::Matrix &sigma, Set::Vector n) const
378  {
379  Set::Vector f = Set::Vector::Zero();
380  for (int i = 0; i < AMREX_SPACEDIM; i++)
381  {
382  if (bc_type[i] == Type::Displacement)
383  f(i) = u(i);
384  else if (bc_type[i] == Type::Traction)
385  f(i) = (sigma*n)(i);
386  else if (bc_type[i] == Type::Neumann)
387  f(i) = (gradu*n)(i);
388  else if (bc_type[i] == Periodic)
389  continue;
390  }
391  return f;
392  }
393 
394 protected:
395 
396  #if AMREX_SPACEDIM==2
397  static const int m_nfaces = 8;
398  #elif AMREX_SPACEDIM==3
399  static const int m_nfaces = 26;
400  #endif
401 
402  std::array<std::array<Type, AMREX_SPACEDIM>, m_nfaces> m_bc_type;
403  std::array<std::array<Numeric::Interpolator::Linear<Set::Scalar>,AMREX_SPACEDIM>, m_nfaces> m_bc_val;
404 
405 public:
406  static void Parse(Constant & value, IO::ParmParse & pp)
407  {
408  std::map<std::string, Type> bcmap;
409  bcmap["displacement"] = Type::Displacement;
410  bcmap["disp"] = Type::Displacement;
411  bcmap["neumann"] = Type::Neumann;
412  bcmap["traction"] = Type::Traction;
413  bcmap["trac"] = Type::Traction;
414  bcmap["periodic"] = Type::Periodic;
415 
416  std::vector<std::string> str;
417 
418  // TYPES
419 
420  #if AMREX_SPACEDIM==3
421  if (pp.contains("type.xloylozlo")) {
422  pp_queryarr("type.xloylozlo",str); // 3D Corner
423  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XLO_YLO_ZLO][i] = bcmap[str[i]]; }
424  if (pp.contains("type.xloylozhi")) {
425  pp_queryarr("type.xloylozhi",str); // 3D Corner
426  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XLO_YLO_ZHI][i] = bcmap[str[i]]; }
427  if (pp.contains("type.xloyhizlo")) {
428  pp_queryarr("type.xloyhizlo",str); // 3D Corner
429  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XLO_YHI_ZLO][i] = bcmap[str[i]]; }
430  if (pp.contains("type.xloyhizhi")) {
431  pp_queryarr("type.xloyhizhi",str); // 3D Corner
432  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XLO_YHI_ZHI][i] = bcmap[str[i]]; }
433  if (pp.contains("type.xhiylozlo")) {
434  pp_queryarr("type.xhiylozlo",str); // 3D Corner
435  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XHI_YLO_ZLO][i] = bcmap[str[i]]; }
436  if (pp.contains("type.xhiylozhi")) {
437  pp_queryarr("type.xhiylozhi",str); // 3D Corner
438  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XHI_YLO_ZHI][i] = bcmap[str[i]]; }
439  if (pp.contains("type.xhiyhizlo")) {
440  pp_queryarr("type.xhiyhizlo",str); // 3D Corner
441  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XHI_YHI_ZLO][i] = bcmap[str[i]]; }
442  if (pp.contains("type.xhiyhizhi")) {
443  pp_queryarr("type.xhiyhizhi",str); // 3D Corner
444  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XHI_YHI_ZHI][i] = bcmap[str[i]]; }
445  #else // ignore unused inputs
446  pp.ignore("type.xloylozlo"); pp.ignore("type.xloylozhi"); pp.ignore("type.xloyhizlo"); pp.ignore("type.xloyhizhi");
447  pp.ignore("type.xhiylozlo"); pp.ignore("type.xhiylozhi"); pp.ignore("type.xhiyhizlo"); pp.ignore("type.xhiyhizhi");
448  #endif
449 
450  #if AMREX_SPACEDIM==3
451  if (pp.contains("type.ylozlo")) {
452  pp_queryarr("type.ylozlo",str); // 3D Edge
453  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::YLO_ZLO][i] = bcmap[str[i]]; }
454  if (pp.contains("type.ylozhi")) {
455  pp_queryarr("type.ylozhi",str); // 3D Edge
456  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::YLO_ZHI][i] = bcmap[str[i]]; }
457  if (pp.contains("type.yhizlo")) {
458  pp_queryarr("type.yhizlo",str); // 3D Edge
459  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::YHI_ZLO][i] = bcmap[str[i]]; }
460  if (pp.contains("type.yhizhi")) {
461  pp_queryarr("type.yhizhi",str); // 3D Edge
462  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::YHI_ZHI][i] = bcmap[str[i]]; }
463  if (pp.contains("type.zloxlo")) {
464  pp_queryarr("type.zloxlo",str); // 3D Edge
465  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::ZLO_XLO][i] = bcmap[str[i]]; }
466  if (pp.contains("type.zloxhi")) {
467  pp_queryarr("type.zloxhi",str); // 3D Edge
468  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::ZLO_XHI][i] = bcmap[str[i]]; }
469  if (pp.contains("type.zhixlo")) {
470  pp_queryarr("type.zhixlo",str); // 3D Edge
471  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::ZHI_XLO][i] = bcmap[str[i]]; }
472  if (pp.contains("type.zhixhi")) {
473  pp_queryarr("type.zhixhi",str); // 3D Edge
474  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::ZHI_XHI][i] = bcmap[str[i]]; }
475  #else
476  pp.ignore("type.ylozlo"); pp.ignore("type.ylozhi"); pp.ignore("type.yhizlo"); pp.ignore("type.yhizhi");
477  pp.ignore("type.zloxlo"); pp.ignore("type.zloxhi"); pp.ignore("type.zhixlo"); pp.ignore("type.zhixhi");
478  #endif
479 
480  if (pp.contains("type.xloylo")) {
481  pp_queryarr("type.xloylo",str); // 3D Edge / 2D Corner
482  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XLO_YLO][i] = bcmap[str[i]]; }
483  if (pp.contains("type.xloyhi")) {
484  pp_queryarr("type.xloyhi",str); // 3D Edge / 2D Corner
485  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XLO_YHI][i] = bcmap[str[i]]; }
486  if (pp.contains("type.xhiylo")) {
487  pp_queryarr("type.xhiylo",str); // 3D Edge / 2D Corner
488  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XHI_YLO][i] = bcmap[str[i]]; }
489  if (pp.contains("type.xhiyhi")) {
490  pp_queryarr("type.xhiyhi",str); // 3D Edge / 2D Corner
491  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XHI_YHI][i] = bcmap[str[i]]; }
492 
493  if (pp.contains("type.xlo")) {
494  pp_queryarr("type.xlo",str); // 3D Face / 2D Edge
495  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XLO][i] = bcmap[str[i]]; }
496  if (pp.contains("type.xhi")) {
497  pp_queryarr("type.xhi",str); // 3D Face / 2D Edge
498  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::XHI][i] = bcmap[str[i]]; }
499  if (pp.contains("type.ylo")) {
500  pp_queryarr("type.ylo",str); // 3D Face / 2D Edge
501  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::YLO][i] = bcmap[str[i]]; }
502  if (pp.contains("type.yhi")) {
503  pp_queryarr("type.yhi",str); // 3D Face / 2D Edge
504  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::YHI][i] = bcmap[str[i]]; }
505  #if AMREX_SPACEDIM==3
506  if (pp.contains("type.zlo")) {
507  pp_queryarr("type.zlo",str); // 3D Face
508  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::ZLO][i] = bcmap[str[i]]; }
509  if (pp.contains("type.zhi")) {
510  pp_queryarr("type.zhi",str); // 3D Face
511  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_type[Face::ZHI][i] = bcmap[str[i]]; }
512  #else
513  pp.ignore("type.zlo"); pp.ignore("type.zhi");
514  #endif
515 
516  // VALS
517  //std::vector<Set::Scalar> val;
518  std::vector<std::string> val;
519 
520  #if AMREX_SPACEDIM==3
521  if (pp.contains("val.xloylozlo")) {
522  pp_queryarr("val.xloylozlo",val); // 3D Corner
523  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::XLO_YLO_ZLO][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
524  if (pp.contains("val.xloylozhi")) {
525  pp_queryarr("val.xloylozhi",val); // 3D Corner
526  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::XLO_YLO_ZHI][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
527  if (pp.contains("val.xloyhizlo")) {
528  pp_queryarr("val.xloyhizlo",val); // 3D Corner
529  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::XLO_YHI_ZLO][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
530  if (pp.contains("val.xloyhizhi")) {
531  pp_queryarr("val.xloyhizhi",val); // 3D Corner
532  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::XLO_YHI_ZHI][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
533  if (pp.contains("val.xhiylozlo")) {
534  pp_queryarr("val.xhiylozlo",val); // 3D Corner
535  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::XHI_YLO_ZLO][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
536  if (pp.contains("val.xhiylozhi")) {
537  pp_queryarr("val.xhiylozhi",val); // 3D Corner
538  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::XHI_YLO_ZHI][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
539  if (pp.contains("val.xhiyhizlo")) {
540  pp_queryarr("val.xhiyhizlo",val); // 3D Corner
541  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::XHI_YHI_ZLO][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
542  if (pp.contains("val.xhiyhizhi")) {
543  pp_queryarr("val.xhiyhizhi",val); // 3D Corner
544  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::XHI_YHI_ZHI][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
545  #else
546  pp.ignore("val.xloylozlo"); pp.ignore("val.xloylozhi"); pp.ignore("val.xloyhizlo"); pp.ignore("val.xloyhizhi");
547  pp.ignore("val.xhiylozlo"); pp.ignore("val.xhiylozhi"); pp.ignore("val.xhiyhizlo"); pp.ignore("val.xhiyhizhi");
548  #endif
549 
550  #if AMREX_SPACEDIM==3
551  if (pp.contains("val.ylozlo")) {
552  pp_queryarr("val.ylozlo",val); // 3D Edge
553  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::YLO_ZLO][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
554  if (pp.contains("val.ylozhi")) {
555  pp_queryarr("val.ylozhi",val); // 3D Edge
556  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::YLO_ZHI][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
557  if (pp.contains("val.yhizlo")) {
558  pp_queryarr("val.yhizlo",val); // 3D Edge
559  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::YHI_ZLO][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
560  if (pp.contains("val.yhizhi")) {
561  pp_queryarr("val.yhizhi",val); // 3D Edge
562  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::YHI_ZHI][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
563  if (pp.contains("val.zloxlo")) {
564  pp_queryarr("val.zloxlo",val); // 3D Edge
565  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::ZLO_XLO][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
566  if (pp.contains("val.zloxhi")) {
567  pp_queryarr("val.zloxhi",val); // 3D Edge
568  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::ZLO_XHI][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
569  if (pp.contains("val.zhixlo")) {
570  pp_queryarr("val.zhixlo",val); // 3D Edge
571  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::ZHI_XLO][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
572  if (pp.contains("val.zhixhi")) {
573  pp_queryarr("val.zhixhi",val); // 3D Edge
574  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::ZHI_XHI][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
575  #else
576  pp.ignore("val.ylozlo"); pp.ignore("val.ylozhi"); pp.ignore("val.yhizlo"); pp.ignore("val.yhizhi");
577  pp.ignore("val.zloxlo"); pp.ignore("val.zloxhi"); pp.ignore("val.zhixlo"); pp.ignore("val.zhixhi");
578  #endif
579  if (pp.contains("val.xloylo")) {
580  pp_queryarr("val.xloylo",val); // 3D Edge / 2D Corner
581  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::XLO_YLO][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
582  if (pp.contains("val.xloyhi")) {
583  pp_queryarr("val.xloyhi",val); // 3D Edge / 2D Corner
584  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::XLO_YHI][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
585  if (pp.contains("val.xhiylo")) {
586  pp_queryarr("val.xhiylo",val); // 3D Edge / 2D Corner
587  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::XHI_YLO][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
588  if (pp.contains("val.xhiyhi")) {
589  pp_queryarr("val.xhiyhi",val); // 3D Edge / 2D Corner
590  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::XHI_YHI][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
591 
592  if (pp.contains("val.xlo")) {
593  pp_queryarr("val.xlo",val); // 3D Face / 2D Edge
594  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::XLO][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
595  if (pp.contains("val.xhi")) {
596  pp_queryarr("val.xhi",val); // 3D Face / 2D Edge
597  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::XHI][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
598  if (pp.contains("val.ylo")) {
599  pp_queryarr("val.ylo",val); // 3D Face / 2D Edge
600  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::YLO][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
601  if (pp.contains("val.yhi")) {
602  pp_queryarr("val.yhi",val); // 3D Face / 2D Edge
603  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::YHI][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
604  #if AMREX_SPACEDIM==3
605  if (pp.contains("val.zlo")) {
606  pp_queryarr("val.zlo",val); // 3D Face
607  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::ZLO][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
608  if (pp.contains("val.zhi")) {
609  pp_queryarr("val.zhi",val); // 3D Face
610  for (int i = 0; i < AMREX_SPACEDIM; i++) value.m_bc_val[Face::ZHI][i] = Numeric::Interpolator::Linear<Set::Scalar>(val[i]); }
611  #else
612  pp.ignore("val.zlo"); pp.ignore("val.zhi");
613  #endif
614 
615  }
616 };
617 }
618 }
619 }
620 #endif
Model::Solid::gradu
@ gradu
Definition: Solid.H:26
BC::Operator::Elastic::Elastic::Direction
Direction
Definition: Elastic.H:44
BC::Operator::Elastic::Constant::name
static constexpr const char * name
Definition: Constant.H:34
BC::Operator::Elastic::Constant::m_bc_type
std::array< std::array< Type, AMREX_SPACEDIM >, m_nfaces > m_bc_type
Definition: Constant.H:402
BC::Operator::Elastic::Constant::Set
void Set(const Face, const Direction, const Type, const Set::Scalar, amrex::Vector< amrex::MultiFab * > &, const amrex::Vector< amrex::Geometry > &)
Definition: Constant.H:64
BC::Operator::Elastic::Constant::getType
virtual std::array< Type, AMREX_SPACEDIM > getType(const int &i, const int &j, [[maybe_unused]] const int &k, const amrex::Box &domain)
Definition: Constant.H:256
BC::Operator::Elastic::Constant::Init
virtual void Init(amrex::MultiFab *a_rhs, const amrex::Geometry &a_geom, bool a_homogeneous=false) const
Definition: Constant.H:179
Linear.H
Operator
Documentation for operator namespace.
Definition: Diagonal.cpp:13
BC::Operator::Elastic::Constant::Constant
Constant(IO::ParmParse &pp, std::string name)
Definition: Constant.H:47
SQRT2INV
#define SQRT2INV
Definition: Elastic.H:78
BC::Operator::Elastic::Elastic::m_time
Set::Scalar m_time
Definition: Elastic.H:93
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
BC::Operator::Elastic::Constant::Parse
static void Parse(Constant &value, IO::ParmParse &pp)
Definition: Constant.H:406
Set::Scalar
amrex::Real Scalar
Definition: Base.H:19
SQRT3INV
#define SQRT3INV
Definition: Elastic.H:77
pp_queryclass
#define pp_queryclass(...)
Definition: ParmParse.H:106
BC::Operator::Elastic::Elastic::Type
Type
Definition: Elastic.H:24
BC::Operator::Elastic::Constant::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)
Definition: Constant.H:312
Set::Matrix
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, AMREX_SPACEDIM > Matrix
Definition: Base.H:23
BC::Operator::Elastic::Constant::Set
void Set(const Face face, const Direction direction, const Type type, const Set::Scalar value, amrex::Vector< std::unique_ptr< amrex::MultiFab > > &a_rhs, const amrex::Vector< amrex::Geometry > &a_geom)
Definition: Constant.H:87
BC::Operator::Elastic::Constant::Init
virtual void Init(amrex::FabArray< amrex::BaseFab< Set::Vector >> *a_rhs, const amrex::Geometry &a_geom, bool a_homogeneous=false) const
Definition: Constant.H:101
Set::Face
@ Face
Definition: Set.H:32
BC::Operator::Elastic::Constant
Definition: Constant.H:31
BC::Operator::Elastic::Constant::Constant
Constant()
Definition: Constant.H:36
BC::Operator::Elastic::Constant::Set
void Set(const Face face, const Direction direction, const Type type, const Set::Scalar value, amrex::Vector< amrex::MultiFab > &a_rhs, const amrex::Vector< amrex::Geometry > &a_geom)
Definition: Constant.H:75
IO::ParmParse::contains
bool contains(std::string name)
Definition: ParmParse.H:153
BC::Operator::Elastic::Elastic
Definition: Elastic.H:18
Numeric::Interpolator::Linear< Set::Scalar >
BC::Operator::Elastic::Constant::m_bc_val
std::array< std::array< Numeric::Interpolator::Linear< Set::Scalar >, AMREX_SPACEDIM >, m_nfaces > m_bc_val
Definition: Constant.H:403
Util::Abort
void Abort(const char *msg)
Definition: Util.cpp:166
BC::Operator::Elastic::Constant::~Constant
~Constant()
Definition: Constant.H:50
BC::Operator::Elastic::Constant::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: Constant.H:376
BC
Collection of boundary condition (BC) objects.
Definition: BC.cpp:4
Elastic.H
IO::ParmParse
Definition: ParmParse.H:112
BC::Operator::Elastic::Elastic::Init
virtual void Init(amrex::MultiFab *a_rhs, const amrex::Geometry &a_geom, bool a_homogeneous=false) const =0
IO::ParmParse::ignore
void ignore(std::string name)
Definition: ParmParse.H:134
BC::Operator::Elastic::Constant::Set
void Set(const Face face, const Direction direction, const Type type, const Set::Scalar value)
Definition: Constant.H:54
INFO
#define INFO
Definition: Util.H:20
pp_queryarr
#define pp_queryarr(...)
Definition: ParmParse.H:103