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