Block-Structured AMR Software Framework
AMReX_EB2_3D_C.H
Go to the documentation of this file.
1 #ifndef AMREX_EB2_3D_C_H_
2 #define AMREX_EB2_3D_C_H_
3 #include <AMReX_Config.H>
4 
5 namespace amrex::EB2 {
6 
8 void
9 amrex_eb2_build_types (Box const& tbx, Box const& bxg2,
10  Array4<Real const> const& s,
11  Array4<EBCellFlag> const& cell,
12  Array4<Type_t> const& fx,
13  Array4<Type_t> const& fy,
14  Array4<Type_t> const& fz,
15  Array4<Type_t> const& ex,
16  Array4<Type_t> const& ey,
17  Array4<Type_t> const& ez)
18 {
19  auto lo = amrex::max_lbound(tbx, bxg2);
20  auto hi = amrex::min_ubound(tbx, bxg2);
21  amrex::Loop(lo, hi,
22  [=] (int i, int j, int k) noexcept
23  {
24  if ( s(i,j ,k ) < 0.0_rt && s(i+1,j ,k ) < 0.0_rt
25  && s(i,j+1,k ) < 0.0_rt && s(i+1,j+1,k ) < 0.0_rt
26  && s(i,j ,k+1) < 0.0_rt && s(i+1,j ,k+1) < 0.0_rt
27  && s(i,j+1,k+1) < 0.0_rt && s(i+1,j+1,k+1) < 0.0_rt)
28  {
29  cell(i,j,k).setRegular();
30  }
31  else if (s(i,j ,k ) >= 0.0_rt && s(i+1,j ,k ) >= 0.0_rt
32  && s(i,j+1,k ) >= 0.0_rt && s(i+1,j+1,k ) >= 0.0_rt
33  && s(i,j ,k+1) >= 0.0_rt && s(i+1,j ,k+1) >= 0.0_rt
34  && s(i,j+1,k+1) >= 0.0_rt && s(i+1,j+1,k+1) >= 0.0_rt)
35  {
36  cell(i,j,k).setCovered();
37  }
38  else
39  {
40  cell(i,j,k).setSingleValued();
41  }
42  });
43 
44  // x-face
45  Box b = amrex::surroundingNodes(bxg2,0);
46  lo = amrex::max_lbound(tbx, b);
47  hi = amrex::min_ubound(tbx, b);
48  amrex::Loop(lo, hi,
49  [=] (int i, int j, int k) noexcept
50  {
51  if ( s(i,j,k ) < 0.0_rt && s(i,j+1,k ) < 0.0_rt
52  && s(i,j,k+1) < 0.0_rt && s(i,j+1,k+1) < 0.0_rt )
53  {
54  fx(i,j,k) = Type::regular;
55  }
56  else if (s(i,j,k ) >= 0.0_rt && s(i,j+1,k ) >= 0.0_rt
57  && s(i,j,k+1) >= 0.0_rt && s(i,j+1,k+1) >= 0.0_rt )
58  {
59  fx(i,j,k) = Type::covered;
60  }
61  else
62  {
63  fx(i,j,k) = Type::irregular;
64  }
65  });
66 
67  // y-face
68  b = amrex::surroundingNodes(bxg2,1);
69  lo = amrex::max_lbound(tbx, b);
70  hi = amrex::min_ubound(tbx, b);
71  amrex::Loop(lo, hi,
72  [=] (int i, int j, int k) noexcept
73  {
74  if ( s(i,j,k ) < 0.0_rt && s(i+1,j,k ) < 0.0_rt
75  && s(i,j,k+1) < 0.0_rt && s(i+1,j,k+1) < 0.0_rt )
76  {
77  fy(i,j,k) = Type::regular;
78  }
79  else if (s(i,j,k ) >= 0.0_rt && s(i+1,j,k ) >= 0.0_rt
80  && s(i,j,k+1) >= 0.0_rt && s(i+1,j,k+1) >= 0.0_rt )
81  {
82  fy(i,j,k) = Type::covered;
83  }
84  else
85  {
86  fy(i,j,k) = Type::irregular;
87  }
88  });
89 
90  // z-face
91  b = amrex::surroundingNodes(bxg2,2);
92  lo = amrex::max_lbound(tbx, b);
93  hi = amrex::min_ubound(tbx, b);
94  amrex::Loop(lo, hi,
95  [=] (int i, int j, int k) noexcept
96  {
97  if ( s(i,j ,k) < 0.0_rt && s(i+1,j ,k) < 0.0_rt
98  && s(i,j+1,k) < 0.0_rt && s(i+1,j+1,k) < 0.0_rt)
99  {
100  fz(i,j,k) = Type::regular;
101  }
102  else if (s(i,j ,k) >= 0.0_rt && s(i+1,j ,k) >= 0.0_rt
103  && s(i,j+1,k) >= 0.0_rt && s(i+1,j+1,k) >= 0.0_rt)
104  {
105  fz(i,j,k) = Type::covered;
106  }
107  else
108  {
109  fz(i,j,k) = Type::irregular;
110  }
111  });
112 
113  // x-edge
114  b = amrex::convert(bxg2,IntVect(0,1,1));
115  lo = amrex::max_lbound(tbx, b);
116  hi = amrex::min_ubound(tbx, b);
117  amrex::Loop(lo, hi,
118  [=] (int i, int j, int k) noexcept
119  {
120  if (s(i,j,k) < 0.0_rt && s(i+1,j,k) < 0.0_rt) {
121  ex(i,j,k) = Type::regular;
122  } else if (s(i,j,k) >= 0.0_rt && s(i+1,j,k) >= 0.0_rt) {
123  ex(i,j,k) = Type::covered;
124  } else {
125  ex(i,j,k) = Type::irregular;
126  }
127  });
128 
129  // y-edge
130  b = amrex::convert(bxg2,IntVect(1,0,1));
131  lo = amrex::max_lbound(tbx, b);
132  hi = amrex::min_ubound(tbx, b);
133  amrex::Loop(lo, hi,
134  [=] (int i, int j, int k) noexcept
135  {
136  if (s(i,j,k) < 0.0_rt && s(i,j+1,k) < 0.0_rt) {
137  ey(i,j,k) = Type::regular;
138  } else if (s(i,j,k) >= 0.0_rt && s(i,j+1,k) >= 0.0_rt) {
139  ey(i,j,k) = Type::covered;
140  } else {
141  ey(i,j,k) = Type::irregular;
142  }
143  });
144 
145  // z-edge
146  b = amrex::convert(bxg2,IntVect(1,1,0));
147  lo = amrex::max_lbound(tbx, b);
148  hi = amrex::min_ubound(tbx, b);
149  amrex::Loop(lo, hi,
150  [=] (int i, int j, int k) noexcept
151  {
152  if (s(i,j,k) < 0.0_rt && s(i,j,k+1) < 0.0_rt) {
153  ez(i,j,k) = Type::regular;
154  } else if (s(i,j,k) >= 0.0_rt && s(i,j,k+1) >= 0.0_rt) {
155  ez(i,j,k) = Type::covered;
156  } else {
157  ez(i,j,k) = Type::irregular;
158  }
159  });
160 }
161 
162 namespace detail {
164  int num_cuts (Real a, Real b) noexcept {
165  return (a >= 0.0_rt && b < 0.0_rt) || (b >= 0.0_rt && a < 0.0_rt);
166  }
167 }
168 
170 int check_mvmc (int i, int j, int k, Array4<Real const> const& fine)
171 {
172  using detail::num_cuts;
173 
174  int ierr = 0;
175 
176  i *= 2;
177  j *= 2;
178  k *= 2;
179 
180  // x-edges
181  int nx00 = num_cuts(fine(i,j,k),fine(i+1,j,k)) + num_cuts(fine(i+1,j,k),fine(i+2,j,k));
182  int nx10 = num_cuts(fine(i,j+2,k),fine(i+1,j+2,k)) + num_cuts(fine(i+1,j+2,k),fine(i+2,j+2,k));
183  int nx01 = num_cuts(fine(i,j,k+2),fine(i+1,j,k+2)) + num_cuts(fine(i+1,j,k+2),fine(i+2,j,k+2));
184  int nx11 = num_cuts(fine(i,j+2,k+2),fine(i+1,j+2,k+2)) + num_cuts(fine(i+1,j+2,k+2),fine(i+2,j+2,k+2));
185 
186  // y-edges
187  int ny00 = num_cuts(fine(i,j,k),fine(i,j+1,k)) + num_cuts(fine(i,j+1,k),fine(i,j+2,k));
188  int ny10 = num_cuts(fine(i+2,j,k),fine(i+2,j+1,k)) + num_cuts(fine(i+2,j+1,k),fine(i+2,j+2,k));
189  int ny01 = num_cuts(fine(i,j,k+2),fine(i,j+1,k+2)) + num_cuts(fine(i,j+1,k+2),fine(i,j+2,k+2));
190  int ny11 = num_cuts(fine(i+2,j,k+2),fine(i+2,j+1,k+2)) + num_cuts(fine(i+2,j+1,k+2),fine(i+2,j+2,k+2));
191 
192  // z-edges
193  int nz00 = num_cuts(fine(i,j,k),fine(i,j,k+1)) + num_cuts(fine(i,j,k+1),fine(i,j,k+2));
194  int nz10 = num_cuts(fine(i+2,j,k),fine(i+2,j,k+1)) + num_cuts(fine(i+2,j,k+1),fine(i+2,j,k+2));
195  int nz01 = num_cuts(fine(i,j+2,k),fine(i,j+2,k+1)) + num_cuts(fine(i,j+2,k+1),fine(i,j+2,k+2));
196  int nz11 = num_cuts(fine(i+2,j+2,k),fine(i+2,j+2,k+1)) + num_cuts(fine(i+2,j+2,k+1),fine(i+2,j+2,k+2));
197 
198  // x-faces
199  int nxm = -1;
200  int n = ny00 + ny01 + nz00 + nz01;
201  if (n == 0) {
202  nxm = 0;
203  } else if (n == 2) {
204  nxm = 1;
205  } else {
206  ierr = 1;
207  }
208 
209  int nxp = -1;
210  n = ny10 + ny11 + nz10 + nz11;
211  if (n == 0) {
212  nxp = 0;
213  } else if (n == 2) {
214  nxp = 1;
215  } else {
216  ierr = 1;
217  }
218 
219  // y-faces
220  int nym = -1;
221  n = nx00 + nx01 + nz00 + nz10;
222  if (n == 0) {
223  nym = 0;
224  } else if (n == 2) {
225  nym = 1;
226  } else {
227  ierr = 1;
228  }
229 
230  int nyp = -1;
231  n = nx10 + nx11 + nz01 + nz11;
232  if (n == 0) {
233  nyp = 0;
234  } else if (n == 2) {
235  nyp = 1;
236  } else {
237  ierr = 1;
238  }
239 
240  // z-faces
241  int nzm = -1;
242  n = nx00 + nx10 + ny00 + ny10;
243  if (n == 0) {
244  nzm = 0;
245  } else if (n == 2) {
246  nzm = 1;
247  } else {
248  ierr = 1;
249  }
250 
251  int nzp = -1;
252  n = nx01 + nx11 + ny01 + ny11;
253  if (n == 0) {
254  nzp = 0;
255  } else if (n == 2) {
256  nzp = 1;
257  } else {
258  ierr = 1;
259  }
260 
261  if (nxm == 1 && nym == 1 && nzm == 1 && nxp == 1 && nyp == 1 && nzp == 1) {
262  n = (fine(i ,j ,k ) < 0.0_rt) + (fine(i+2,j ,k ) < 0.0_rt) +
263  (fine(i ,j+2,k ) < 0.0_rt) + (fine(i+2,j+2,k ) < 0.0_rt) +
264  (fine(i ,j ,k+2) < 0.0_rt) + (fine(i+2,j ,k+2) < 0.0_rt) +
265  (fine(i ,j+2,k+2) < 0.0_rt) + (fine(i+2,j+2,k+2) < 0.0_rt);
266  if (n == 2 || n == 6) {
267  ierr = 1;
268  } else if (n != 4) {
269  ierr = 1;
270  amrex::Abort("amrex::check_mvmc: how did this happen? nopen != 4");
271  }
272  }
273 
274  return ierr;
275 }
276 
277 namespace detail {
279 Real coarsen_edge_cent (Real f1, Real f2)
280 {
281  if (f1 == Real(1.0) && f2 == Real(1.0)) {
282  return Real(1.0);
283  } else if (f1 == Real(-1.0) && f2 == Real(-1.0)) {
284  return Real(-1.0);
285  } else {
286  if (f1 == Real(-1.0)) {
287  f1 = Real(0.0);
288  } else if (f1 == Real(1.0)) {
289  f1 = Real(-0.25);
290  } else {
291  f1 = Real(0.5)*f1 - Real(0.25);
292  }
293  if (f2 == Real(-1.0)) {
294  f2 = Real(0.0);
295  } else if (f2 == Real(1.0)) {
296  f2 = Real(0.25);
297  } else {
298  f2 = Real(0.5)*f2 + Real(0.25);
299  }
300  Real r = (f2*f2-f1*f1)/(f2-f1+Real(1.e-30));
301  return amrex::min(Real(0.5),amrex::max(Real(-0.5),r));
302  }
303 }
304 }
305 
307 int coarsen_from_fine (int i, int j, int k, Box const& bx, int ngrow,
308  Array4<Real> const& cvol, Array4<Real> const& ccent,
309  Array4<Real> const& cba, Array4<Real> const& cbc,
310  Array4<Real> const& cbn, Array4<Real> const& capx,
311  Array4<Real> const& capy, Array4<Real> const& capz,
312  Array4<Real> const& cfcx, Array4<Real> const& cfcy,
313  Array4<Real> const& cfcz, Array4<Real> const& cecx,
314  Array4<Real> const& cecy, Array4<Real> const& cecz,
315  Array4<EBCellFlag> const& cflag,
316  Array4<Real const> const& fvol, Array4<Real const> const& fcent,
317  Array4<Real const> const& fba, Array4<Real const> const& fbc,
318  Array4<Real const> const& fbn, Array4<Real const> const& fapx,
319  Array4<Real const> const& fapy, Array4<Real const> const& fapz,
320  Array4<Real const> const& ffcx, Array4<Real const> const& ffcy,
321  Array4<Real const> const& ffcz, Array4<Real const> const& fecx,
322  Array4<Real const> const& fecy, Array4<Real const> const& fecz,
323  Array4<EBCellFlag const> const& fflag)
324 {
326 
327  const Box& gbx = amrex::grow(bx,ngrow);
328  const Box& xbx = amrex::surroundingNodes(bx,0); // face boxes
329  const Box& ybx = amrex::surroundingNodes(bx,1);
330  const Box& zbx = amrex::surroundingNodes(bx,2);
331  const Box& xgbx = amrex::surroundingNodes(gbx,0);
332  const Box& ygbx = amrex::surroundingNodes(gbx,1);
333  const Box& zgbx = amrex::surroundingNodes(gbx,2);
334  const Box& exbx = amrex::convert(bx,IntVect(0,1,1)); // edge boxes
335  const Box& eybx = amrex::convert(bx,IntVect(1,0,1));
336  const Box& ezbx = amrex::convert(bx,IntVect(1,1,0));
337  const Box& exgbx = amrex::convert(gbx,IntVect(0,1,1));
338  const Box& eygbx = amrex::convert(gbx,IntVect(1,0,1));
339  const Box& ezgbx = amrex::convert(gbx,IntVect(1,1,0));
340 
341  int ierr = 0;
342  int ii = i*2;
343  int jj = j*2;
344  int kk = k*2;
345 
346  if (bx.contains(i,j,k))
347  {
348  if (fflag(ii,jj ,kk ).isRegular() && fflag(ii+1,jj ,kk ).isRegular() &&
349  fflag(ii,jj+1,kk ).isRegular() && fflag(ii+1,jj+1,kk ).isRegular() &&
350  fflag(ii,jj ,kk+1).isRegular() && fflag(ii+1,jj ,kk+1).isRegular() &&
351  fflag(ii,jj+1,kk+1).isRegular() && fflag(ii+1,jj+1,kk+1).isRegular())
352  {
353  cflag(i,j,k).setRegular();
354  cvol(i,j,k) = 1.0_rt;
355  ccent(i,j,k,0) = 0.0_rt;
356  ccent(i,j,k,1) = 0.0_rt;
357  ccent(i,j,k,2) = 0.0_rt;
358  cba(i,j,k) = 0.0_rt;
359  cbc(i,j,k,0) = -1.0_rt;
360  cbc(i,j,k,1) = -1.0_rt;
361  cbc(i,j,k,2) = -1.0_rt;
362  cbn(i,j,k,0) = 0.0_rt;
363  cbn(i,j,k,1) = 0.0_rt;
364  cbn(i,j,k,2) = 0.0_rt;
365  }
366  else if (fflag(ii,jj ,kk ).isCovered() && fflag(ii+1,jj ,kk ).isCovered() &&
367  fflag(ii,jj+1,kk ).isCovered() && fflag(ii+1,jj+1,kk ).isCovered() &&
368  fflag(ii,jj ,kk+1).isCovered() && fflag(ii+1,jj ,kk+1).isCovered() &&
369  fflag(ii,jj+1,kk+1).isCovered() && fflag(ii+1,jj+1,kk+1).isCovered())
370  {
371  cflag(i,j,k).setCovered();
372  cvol(i,j,k) = 0.0_rt;
373  ccent(i,j,k,0) = 0.0_rt;
374  ccent(i,j,k,1) = 0.0_rt;
375  ccent(i,j,k,2) = 0.0_rt;
376  cba(i,j,k) = 0.0_rt;
377  cbc(i,j,k,0) = -1.0_rt;
378  cbc(i,j,k,1) = -1.0_rt;
379  cbc(i,j,k,2) = -1.0_rt;
380  cbn(i,j,k,0) = 0.0_rt;
381  cbn(i,j,k,1) = 0.0_rt;
382  cbn(i,j,k,2) = 0.0_rt;
383  }
384  else
385  {
386  cflag(i,j,k).setSingleValued();
387 
388  cvol(i,j,k) = 0.125_rt*(fvol(ii ,jj ,kk ) + fvol(ii+1,jj ,kk ) +
389  fvol(ii ,jj+1,kk ) + fvol(ii+1,jj+1,kk ) +
390  fvol(ii ,jj ,kk+1) + fvol(ii+1,jj ,kk+1) +
391  fvol(ii ,jj+1,kk+1) + fvol(ii+1,jj+1,kk+1));
392  Real cvolinv = 1.0_rt / cvol(i,j,k);
393 
394  ccent(i,j,k,0) = 0.125_rt * cvolinv *
395  (fvol(ii ,jj ,kk )*(0.5_rt*fcent(ii ,jj ,kk ,0)-0.25_rt) +
396  fvol(ii+1,jj ,kk )*(0.5_rt*fcent(ii+1,jj ,kk ,0)+0.25_rt) +
397  fvol(ii ,jj+1,kk )*(0.5_rt*fcent(ii ,jj+1,kk ,0)-0.25_rt) +
398  fvol(ii+1,jj+1,kk )*(0.5_rt*fcent(ii+1,jj+1,kk ,0)+0.25_rt) +
399  fvol(ii ,jj ,kk+1)*(0.5_rt*fcent(ii ,jj ,kk+1,0)-0.25_rt) +
400  fvol(ii+1,jj ,kk+1)*(0.5_rt*fcent(ii+1,jj ,kk+1,0)+0.25_rt) +
401  fvol(ii ,jj+1,kk+1)*(0.5_rt*fcent(ii ,jj+1,kk+1,0)-0.25_rt) +
402  fvol(ii+1,jj+1,kk+1)*(0.5_rt*fcent(ii+1,jj+1,kk+1,0)+0.25_rt));
403  ccent(i,j,k,1) = 0.125_rt * cvolinv *
404  (fvol(ii ,jj ,kk )*(0.5_rt*fcent(ii ,jj ,kk ,1)-0.25_rt) +
405  fvol(ii+1,jj ,kk )*(0.5_rt*fcent(ii+1,jj ,kk ,1)-0.25_rt) +
406  fvol(ii ,jj+1,kk )*(0.5_rt*fcent(ii ,jj+1,kk ,1)+0.25_rt) +
407  fvol(ii+1,jj+1,kk )*(0.5_rt*fcent(ii+1,jj+1,kk ,1)+0.25_rt) +
408  fvol(ii ,jj ,kk+1)*(0.5_rt*fcent(ii ,jj ,kk+1,1)-0.25_rt) +
409  fvol(ii+1,jj ,kk+1)*(0.5_rt*fcent(ii+1,jj ,kk+1,1)-0.25_rt) +
410  fvol(ii ,jj+1,kk+1)*(0.5_rt*fcent(ii ,jj+1,kk+1,1)+0.25_rt) +
411  fvol(ii+1,jj+1,kk+1)*(0.5_rt*fcent(ii+1,jj+1,kk+1,1)+0.25_rt));
412  ccent(i,j,k,2) = 0.125_rt * cvolinv *
413  (fvol(ii ,jj ,kk )*(0.5_rt*fcent(ii ,jj ,kk ,2)-0.25_rt) +
414  fvol(ii+1,jj ,kk )*(0.5_rt*fcent(ii+1,jj ,kk ,2)-0.25_rt) +
415  fvol(ii ,jj+1,kk )*(0.5_rt*fcent(ii ,jj+1,kk ,2)-0.25_rt) +
416  fvol(ii+1,jj+1,kk )*(0.5_rt*fcent(ii+1,jj+1,kk ,2)-0.25_rt) +
417  fvol(ii ,jj ,kk+1)*(0.5_rt*fcent(ii ,jj ,kk+1,2)+0.25_rt) +
418  fvol(ii+1,jj ,kk+1)*(0.5_rt*fcent(ii+1,jj ,kk+1,2)+0.25_rt) +
419  fvol(ii ,jj+1,kk+1)*(0.5_rt*fcent(ii ,jj+1,kk+1,2)+0.25_rt) +
420  fvol(ii+1,jj+1,kk+1)*(0.5_rt*fcent(ii+1,jj+1,kk+1,2)+0.25_rt));
421 
422  cba(i,j,k) = 0.25_rt*(fba(ii ,jj ,kk ) + fba(ii+1,jj ,kk ) +
423  fba(ii ,jj+1,kk ) + fba(ii+1,jj+1,kk ) +
424  fba(ii ,jj ,kk+1) + fba(ii+1,jj ,kk+1) +
425  fba(ii ,jj+1,kk+1) + fba(ii+1,jj+1,kk+1));
426  Real cbainv = 1.0_rt / cba(i,j,k);
427 
428  cbc(i,j,k,0) = 0.25_rt * cbainv *
429  ( fba(ii ,jj ,kk )*(0.5_rt*fbc(ii ,jj ,kk ,0)-0.25_rt)
430  + fba(ii+1,jj ,kk )*(0.5_rt*fbc(ii+1,jj ,kk ,0)+0.25_rt)
431  + fba(ii ,jj+1,kk )*(0.5_rt*fbc(ii ,jj+1,kk ,0)-0.25_rt)
432  + fba(ii+1,jj+1,kk )*(0.5_rt*fbc(ii+1,jj+1,kk ,0)+0.25_rt)
433  + fba(ii ,jj ,kk+1)*(0.5_rt*fbc(ii ,jj ,kk+1,0)-0.25_rt)
434  + fba(ii+1,jj ,kk+1)*(0.5_rt*fbc(ii+1,jj ,kk+1,0)+0.25_rt)
435  + fba(ii ,jj+1,kk+1)*(0.5_rt*fbc(ii ,jj+1,kk+1,0)-0.25_rt)
436  + fba(ii+1,jj+1,kk+1)*(0.5_rt*fbc(ii+1,jj+1,kk+1,0)+0.25_rt) );
437  cbc(i,j,k,1) = 0.25_rt * cbainv *
438  ( fba(ii ,jj ,kk )*(0.5_rt*fbc(ii ,jj ,kk ,1)-0.25_rt)
439  + fba(ii+1,jj ,kk )*(0.5_rt*fbc(ii+1,jj ,kk ,1)-0.25_rt)
440  + fba(ii ,jj+1,kk )*(0.5_rt*fbc(ii ,jj+1,kk ,1)+0.25_rt)
441  + fba(ii+1,jj+1,kk )*(0.5_rt*fbc(ii+1,jj+1,kk ,1)+0.25_rt)
442  + fba(ii ,jj ,kk+1)*(0.5_rt*fbc(ii ,jj ,kk+1,1)-0.25_rt)
443  + fba(ii+1,jj ,kk+1)*(0.5_rt*fbc(ii+1,jj ,kk+1,1)-0.25_rt)
444  + fba(ii ,jj+1,kk+1)*(0.5_rt*fbc(ii ,jj+1,kk+1,1)+0.25_rt)
445  + fba(ii+1,jj+1,kk+1)*(0.5_rt*fbc(ii+1,jj+1,kk+1,1)+0.25_rt) );
446  cbc(i,j,k,2) = 0.25_rt * cbainv *
447  ( fba(ii ,jj ,kk )*(0.5_rt*fbc(ii ,jj ,kk ,2)-0.25_rt)
448  + fba(ii+1,jj ,kk )*(0.5_rt*fbc(ii+1,jj ,kk ,2)-0.25_rt)
449  + fba(ii ,jj+1,kk )*(0.5_rt*fbc(ii ,jj+1,kk ,2)-0.25_rt)
450  + fba(ii+1,jj+1,kk )*(0.5_rt*fbc(ii+1,jj+1,kk ,2)-0.25_rt)
451  + fba(ii ,jj ,kk+1)*(0.5_rt*fbc(ii ,jj ,kk+1,2)+0.25_rt)
452  + fba(ii+1,jj ,kk+1)*(0.5_rt*fbc(ii+1,jj ,kk+1,2)+0.25_rt)
453  + fba(ii ,jj+1,kk+1)*(0.5_rt*fbc(ii ,jj+1,kk+1,2)+0.25_rt)
454  + fba(ii+1,jj+1,kk+1)*(0.5_rt*fbc(ii+1,jj+1,kk+1,2)+0.25_rt) );
455 
456  Real nx = fbn(ii ,jj ,kk ,0)*fba(ii ,jj ,kk )
457  + fbn(ii+1,jj ,kk ,0)*fba(ii+1,jj ,kk )
458  + fbn(ii ,jj+1,kk ,0)*fba(ii ,jj+1,kk )
459  + fbn(ii+1,jj+1,kk ,0)*fba(ii+1,jj+1,kk )
460  + fbn(ii ,jj ,kk+1,0)*fba(ii ,jj ,kk+1)
461  + fbn(ii+1,jj ,kk+1,0)*fba(ii+1,jj ,kk+1)
462  + fbn(ii ,jj+1,kk+1,0)*fba(ii ,jj+1,kk+1)
463  + fbn(ii+1,jj+1,kk+1,0)*fba(ii+1,jj+1,kk+1);
464  Real ny = fbn(ii ,jj ,kk ,1)*fba(ii ,jj ,kk )
465  + fbn(ii+1,jj ,kk ,1)*fba(ii+1,jj ,kk )
466  + fbn(ii ,jj+1,kk ,1)*fba(ii ,jj+1,kk )
467  + fbn(ii+1,jj+1,kk ,1)*fba(ii+1,jj+1,kk )
468  + fbn(ii ,jj ,kk+1,1)*fba(ii ,jj ,kk+1)
469  + fbn(ii+1,jj ,kk+1,1)*fba(ii+1,jj ,kk+1)
470  + fbn(ii ,jj+1,kk+1,1)*fba(ii ,jj+1,kk+1)
471  + fbn(ii+1,jj+1,kk+1,1)*fba(ii+1,jj+1,kk+1);
472  Real nz = fbn(ii ,jj ,kk ,2)*fba(ii ,jj ,kk )
473  + fbn(ii+1,jj ,kk ,2)*fba(ii+1,jj ,kk )
474  + fbn(ii ,jj+1,kk ,2)*fba(ii ,jj+1,kk )
475  + fbn(ii+1,jj+1,kk ,2)*fba(ii+1,jj+1,kk )
476  + fbn(ii ,jj ,kk+1,2)*fba(ii ,jj ,kk+1)
477  + fbn(ii+1,jj ,kk+1,2)*fba(ii+1,jj ,kk+1)
478  + fbn(ii ,jj+1,kk+1,2)*fba(ii ,jj+1,kk+1)
479  + fbn(ii+1,jj+1,kk+1,2)*fba(ii+1,jj+1,kk+1);
480  Real nfac = 1.0_rt / std::sqrt(nx*nx+ny*ny+nz*nz+1.e-30_rt);
481  cbn(i,j,k,0) = nx*nfac;
482  cbn(i,j,k,1) = ny*nfac;
483  cbn(i,j,k,2) = nz*nfac;
484  ierr = (nx == 0.0_rt && ny == 0.0_rt && nz == 0.0_rt)
485  // we must check the enclosing surface to make sure the coarse cell does not
486  // fully contains the geometry object.
487  || ( fapx(ii,jj ,kk )==1.0_rt && fapx(ii+2,jj ,kk )==1.0_rt
488  && fapx(ii,jj+1,kk )==1.0_rt && fapx(ii+2,jj+1,kk )==1.0_rt
489  && fapx(ii,jj ,kk+1)==1.0_rt && fapx(ii+2,jj ,kk+1)==1.0_rt
490  && fapx(ii,jj+1,kk+1)==1.0_rt && fapx(ii+2,jj+1,kk+1)==1.0_rt
491  && fapy(ii,jj ,kk )==1.0_rt && fapy(ii+1,jj ,kk )==1.0_rt
492  && fapy(ii,jj+2,kk )==1.0_rt && fapy(ii+1,jj+2,kk )==1.0_rt
493  && fapy(ii,jj ,kk+1)==1.0_rt && fapy(ii+1,jj ,kk+1)==1.0_rt
494  && fapy(ii,jj+2,kk+1)==1.0_rt && fapy(ii+1,jj+2,kk+1)==1.0_rt
495  && fapz(ii,jj ,kk )==1.0_rt && fapz(ii+1,jj ,kk )==1.0_rt
496  && fapz(ii,jj+1,kk )==1.0_rt && fapz(ii+1,jj+1,kk )==1.0_rt
497  && fapz(ii,jj ,kk+2)==1.0_rt && fapz(ii+1,jj ,kk+2)==1.0_rt
498  && fapz(ii,jj+1,kk+2)==1.0_rt && fapz(ii+1,jj+1,kk+2)==1.0_rt);
499 
500  }
501  }
502  else if (gbx.contains(i,j,k))
503  {
504  cvol(i,j,k) = 1.0_rt;
505  ccent(i,j,k,0) = 0.0_rt;
506  ccent(i,j,k,1) = 0.0_rt;
507  ccent(i,j,k,2) = 0.0_rt;
508  cba(i,j,k) = 0.0_rt;
509  cbc(i,j,k,0) = -1.0_rt;
510  cbc(i,j,k,1) = -1.0_rt;
511  cbc(i,j,k,2) = -1.0_rt;
512  cbn(i,j,k,0) = 0.0_rt;
513  cbn(i,j,k,1) = 0.0_rt;
514  cbn(i,j,k,2) = 0.0_rt;
515  }
516 
517  if (xbx.contains(i,j,k))
518  {
519  capx(i,j,k) = 0.25_rt*(fapx(ii,jj ,kk ) + fapx(ii,jj+1,kk ) +
520  fapx(ii,jj ,kk+1) + fapx(ii,jj+1,kk+1));
521  if (capx(i,j,k) != 0.0_rt) {
522  Real apinv = 1.0_rt / capx(i,j,k);
523  cfcx(i,j,k,0) = 0.25_rt * apinv *
524  (fapx(ii,jj ,kk )*(0.5_rt*ffcx(ii,jj ,kk ,0)-0.25_rt) +
525  fapx(ii,jj+1,kk )*(0.5_rt*ffcx(ii,jj+1,kk ,0)+0.25_rt) +
526  fapx(ii,jj ,kk+1)*(0.5_rt*ffcx(ii,jj ,kk+1,0)-0.25_rt) +
527  fapx(ii,jj+1,kk+1)*(0.5_rt*ffcx(ii,jj+1,kk+1,0)+0.25_rt) );
528  cfcx(i,j,k,1) = 0.25_rt * apinv *
529  (fapx(ii,jj ,kk )*(0.5_rt*ffcx(ii,jj ,kk ,1)-0.25_rt) +
530  fapx(ii,jj+1,kk )*(0.5_rt*ffcx(ii,jj+1,kk ,1)-0.25_rt) +
531  fapx(ii,jj ,kk+1)*(0.5_rt*ffcx(ii,jj ,kk+1,1)+0.25_rt) +
532  fapx(ii,jj+1,kk+1)*(0.5_rt*ffcx(ii,jj+1,kk+1,1)+0.25_rt) );
533  } else {
534  cfcx(i,j,k,0) = 0.0_rt;
535  cfcx(i,j,k,1) = 0.0_rt;
536  }
537  }
538  else if (xgbx.contains(i,j,k))
539  {
540  capx(i,j,k) = 1.0_rt;
541  cfcx(i,j,k,0) = 0.0_rt;
542  cfcx(i,j,k,1) = 0.0_rt;
543  }
544 
545  if (ybx.contains(i,j,k))
546  {
547  capy(i,j,k) = 0.25_rt*(fapy(ii ,jj,kk ) + fapy(ii+1,jj,kk ) +
548  fapy(ii ,jj,kk+1) + fapy(ii+1,jj,kk+1));
549  if (capy(i,j,k) != 0.0_rt) {
550  Real apinv = 1.0_rt / capy(i,j,k);
551  cfcy(i,j,k,0) = 0.25_rt * apinv *
552  (fapy(ii ,jj,kk )*(0.5_rt*ffcy(ii ,jj,kk ,0)-0.25_rt) +
553  fapy(ii+1,jj,kk )*(0.5_rt*ffcy(ii+1,jj,kk ,0)+0.25_rt) +
554  fapy(ii ,jj,kk+1)*(0.5_rt*ffcy(ii ,jj,kk+1,0)-0.25_rt) +
555  fapy(ii+1,jj,kk+1)*(0.5_rt*ffcy(ii+1,jj,kk+1,0)+0.25_rt) );
556  cfcy(i,j,k,1) = 0.25_rt * apinv *
557  (fapy(ii ,jj,kk )*(0.5_rt*ffcy(ii ,jj,kk ,1)-0.25_rt) +
558  fapy(ii+1,jj,kk )*(0.5_rt*ffcy(ii+1,jj,kk ,1)-0.25_rt) +
559  fapy(ii ,jj,kk+1)*(0.5_rt*ffcy(ii ,jj,kk+1,1)+0.25_rt) +
560  fapy(ii+1,jj,kk+1)*(0.5_rt*ffcy(ii+1,jj,kk+1,1)+0.25_rt) );
561  } else {
562  cfcy(i,j,k,0) = 0.0_rt;
563  cfcy(i,j,k,1) = 0.0_rt;
564  }
565  }
566  else if (ygbx.contains(i,j,k))
567  {
568  capy(i,j,k) = 1.0_rt;
569  cfcy(i,j,k,0) = 0.0_rt;
570  cfcy(i,j,k,1) = 0.0_rt;
571  }
572 
573  if (zbx.contains(i,j,k))
574  {
575  capz(i,j,k) = 0.25_rt * (fapz(ii ,jj ,kk) + fapz(ii+1,jj ,kk) +
576  fapz(ii ,jj+1,kk) + fapz(ii+1,jj+1,kk));
577  if (capz(i,j,k) != 0.0_rt) {
578  Real apinv = 1.0_rt / capz(i,j,k);
579  cfcz(i,j,k,0) = 0.25_rt * apinv *
580  (fapz(ii ,jj ,kk)*(0.5_rt*ffcz(ii ,jj ,kk,0)-0.25_rt) +
581  fapz(ii+1,jj ,kk)*(0.5_rt*ffcz(ii+1,jj ,kk,0)+0.25_rt) +
582  fapz(ii ,jj+1,kk)*(0.5_rt*ffcz(ii ,jj+1,kk,0)-0.25_rt) +
583  fapz(ii+1,jj+1,kk)*(0.5_rt*ffcz(ii+1,jj+1,kk,0)+0.25_rt) );
584  cfcz(i,j,k,1) = 0.25_rt * apinv *
585  (fapz(ii ,jj ,kk)*(0.5_rt*ffcz(ii ,jj ,kk,1)-0.25_rt) +
586  fapz(ii+1,jj ,kk)*(0.5_rt*ffcz(ii+1,jj ,kk,1)-0.25_rt) +
587  fapz(ii ,jj+1,kk)*(0.5_rt*ffcz(ii ,jj+1,kk,1)+0.25_rt) +
588  fapz(ii+1,jj+1,kk)*(0.5_rt*ffcz(ii+1,jj+1,kk,1)+0.25_rt) );
589  } else {
590  cfcz(i,j,k,0) = 0.0_rt;
591  cfcz(i,j,k,1) = 0.0_rt;
592  }
593  }
594  else if (zgbx.contains(i,j,k))
595  {
596  capz(i,j,k) = 1.0_rt;
597  cfcz(i,j,k,0) = 0.0_rt;
598  cfcz(i,j,k,1) = 0.0_rt;
599  }
600 
601  if (exbx.contains(i,j,k))
602  {
603  cecx(i,j,k) = coarsen_edge_cent(fecx(ii,jj,kk), fecx(ii+1,jj,kk));
604  }
605  else if (exgbx.contains(i,j,k))
606  {
607  cecx(i,j,k) = 1.0_rt;
608  }
609 
610  if (eybx.contains(i,j,k))
611  {
612  cecy(i,j,k) = coarsen_edge_cent(fecy(ii,jj,kk), fecy(ii,jj+1,kk));
613  }
614  else if (eygbx.contains(i,j,k))
615  {
616  cecy(i,j,k) = 1.0_rt;
617  }
618 
619  if (ezbx.contains(i,j,k))
620  {
621  cecz(i,j,k) = coarsen_edge_cent(fecz(ii,jj,kk), fecz(ii,jj,kk+1));
622  }
623  else if (ezgbx.contains(i,j,k))
624  {
625  cecz(i,j,k) = 1.0_rt;
626  }
627 
628  return ierr;
629 }
630 
632 void build_cellflag_from_ap (int i, int j, int k, Array4<EBCellFlag> const& cflag,
633  Array4<Real const> const& apx, Array4<Real const> const& apy,
634  Array4<Real const> const& apz)
635 {
636  auto flg = cflag(i,j,k);
637  flg.setDisconnected();
638 
639  if (!flg.isCovered())
640  {
641  flg.setConnected(0,0,0);
642 
643  if (apx(i ,j,k) != 0.0_rt) { flg.setConnected(-1, 0, 0); }
644  if (apx(i+1,j,k) != 0.0_rt) { flg.setConnected( 1, 0, 0); }
645  if (apy(i,j ,k) != 0.0_rt) { flg.setConnected( 0, -1, 0); }
646  if (apy(i,j+1,k) != 0.0_rt) { flg.setConnected( 0, 1, 0); }
647  if (apz(i,j,k ) != 0.0_rt) { flg.setConnected( 0, 0, -1); }
648  if (apz(i,j,k+1) != 0.0_rt) { flg.setConnected( 0, 0, 1); }
649 
650  if ( (apx(i,j,k) != 0.0_rt && apy(i-1,j,k) != 0.0_rt) ||
651  (apy(i,j,k) != 0.0_rt && apx(i,j-1,k) != 0.0_rt) )
652  {
653  flg.setConnected(-1, -1, 0);
654  if (apz(i-1,j-1,k ) != 0.0_rt) { flg.setConnected(-1,-1,-1); }
655  if (apz(i-1,j-1,k+1) != 0.0_rt) { flg.setConnected(-1,-1, 1); }
656  }
657 
658  if ( (apx(i+1,j,k) != 0.0_rt && apy(i+1,j ,k) != 0.0_rt) ||
659  (apy(i ,j,k) != 0.0_rt && apx(i+1,j-1,k) != 0.0_rt) )
660  {
661  flg.setConnected(1, -1, 0);
662  if (apz(i+1,j-1,k ) != 0.0_rt) { flg.setConnected(1,-1,-1); }
663  if (apz(i+1,j-1,k+1) != 0.0_rt) { flg.setConnected(1,-1, 1); }
664  }
665 
666  if ( (apx(i,j ,k) != 0.0_rt && apy(i-1,j+1,k) != 0.0_rt) ||
667  (apy(i,j+1,k) != 0.0_rt && apx(i ,j+1,k) != 0.0_rt) )
668  {
669  flg.setConnected(-1, 1, 0);
670  if (apz(i-1,j+1,k ) != 0.0_rt) { flg.setConnected(-1, 1,-1); }
671  if (apz(i-1,j+1,k+1) != 0.0_rt) { flg.setConnected(-1, 1, 1); }
672  }
673 
674  if ( (apx(i+1,j ,k) != 0.0_rt && apy(i+1,j+1,k) != 0.0_rt) ||
675  (apy(i ,j+1,k) != 0.0_rt && apx(i+1,j+1,k) != 0.0_rt) )
676  {
677  flg.setConnected(1, 1, 0);
678  if (apz(i+1,j+1,k ) != 0.0_rt) { flg.setConnected(1, 1,-1); }
679  if (apz(i+1,j+1,k+1) != 0.0_rt) { flg.setConnected(1, 1, 1); }
680  }
681 
682  if ( (apx(i,j,k) != 0.0_rt && apz(i-1,j,k ) != 0.0_rt) ||
683  (apz(i,j,k) != 0.0_rt && apx(i ,j,k-1) != 0.0_rt) )
684  {
685  flg.setConnected(-1, 0, -1);
686  if (apy(i-1,j ,k-1) != 0.0_rt) { flg.setConnected(-1,-1,-1); }
687  if (apy(i-1,j+1,k-1) != 0.0_rt) { flg.setConnected(-1, 1,-1); }
688  }
689 
690  if ( (apx(i+1,j,k) != 0.0_rt && apz(i+1,j,k ) != 0.0_rt) ||
691  (apz(i ,j,k) != 0.0_rt && apx(i+1,j,k-1) != 0.0_rt) )
692  {
693  flg.setConnected(1, 0, -1);
694  if (apy(i+1,j ,k-1) != 0.0_rt) { flg.setConnected(1,-1,-1); }
695  if (apy(i+1,j+1,k-1) != 0.0_rt) { flg.setConnected(1, 1,-1); }
696  }
697 
698  if ( (apx(i,j,k ) != 0.0_rt && apz(i-1,j,k+1) != 0.0_rt) ||
699  (apz(i,j,k+1) != 0.0_rt && apx(i ,j,k+1) != 0.0_rt) )
700  {
701  flg.setConnected(-1, 0, 1);
702  if (apy(i-1,j ,k+1) != 0.0_rt) { flg.setConnected(-1,-1, 1); }
703  if (apy(i-1,j+1,k+1) != 0.0_rt) { flg.setConnected(-1, 1, 1); }
704  }
705 
706  if ( (apx(i+1,j,k ) != 0.0_rt && apz(i+1,j,k+1) != 0.0_rt) ||
707  (apz(i ,j,k+1) != 0.0_rt && apx(i+1,j,k+1) != 0.0_rt) )
708  {
709  flg.setConnected(1, 0, 1);
710  if (apy(i+1,j ,k+1) != 0.0_rt) { flg.setConnected(1,-1, 1); }
711  if (apy(i+1,j+1,k+1) != 0.0_rt) { flg.setConnected(1, 1, 1); }
712  }
713 
714  if ( (apy(i,j,k) != 0.0_rt && apz(i,j-1,k ) != 0.0_rt) ||
715  (apz(i,j,k) != 0.0_rt && apy(i,j ,k-1) != 0.0_rt) )
716  {
717  flg.setConnected(0, -1, -1);
718  if (apx(i ,j-1,k-1) != 0.0_rt) { flg.setConnected(-1,-1,-1); }
719  if (apx(i+1,j-1,k-1) != 0.0_rt) { flg.setConnected( 1,-1,-1); }
720  }
721 
722  if ( (apy(i,j+1,k) != 0.0_rt && apz(i,j+1,k ) != 0.0_rt) ||
723  (apz(i,j ,k) != 0.0_rt && apy(i,j+1,k-1) != 0.0_rt) )
724  {
725  flg.setConnected(0, 1, -1);
726  if (apx(i ,j+1,k-1) != 0.0_rt) { flg.setConnected(-1, 1,-1); }
727  if (apx(i+1,j+1,k-1) != 0.0_rt) { flg.setConnected( 1, 1,-1); }
728  }
729 
730  if ( (apy(i,j,k ) != 0.0_rt && apz(i,j-1,k+1) != 0.0_rt) ||
731  (apz(i,j,k+1) != 0.0_rt && apy(i,j ,k+1) != 0.0_rt) )
732  {
733  flg.setConnected(0, -1, 1);
734  if (apx(i ,j-1,k+1) != 0.0_rt) { flg.setConnected(-1,-1, 1); }
735  if (apx(i+1,j-1,k+1) != 0.0_rt) { flg.setConnected( 1,-1, 1); }
736  }
737 
738  if ( (apy(i,j+1,k ) != 0.0_rt && apz(i,j+1,k+1) != 0.0_rt) ||
739  (apz(i,j ,k+1) != 0.0_rt && apy(i,j+1,k+1) != 0.0_rt) )
740  {
741  flg.setConnected(0, 1, 1);
742  if (apx(i ,j+1,k+1) != 0.0_rt) { flg.setConnected(-1, 1, 1); }
743  if (apx(i+1,j+1,k+1) != 0.0_rt) { flg.setConnected( 1, 1, 1); }
744  }
745  }
746 
747  cflag(i,j,k) = flg;
748 }
749 
750 }
751 
752 #endif
#define AMREX_FORCE_INLINE
Definition: AMReX_Extension.H:119
#define AMREX_GPU_HOST_DEVICE
Definition: AMReX_GpuQualifiers.H:20
Array4< Real > fine
Definition: AMReX_InterpFaceRegister.cpp:90
AMREX_GPU_HOST_DEVICE bool contains(const IntVectND< dim > &p) const noexcept
Returns true if argument is contained within BoxND.
Definition: AMReX_Box.H:204
static constexpr Type_t regular
Definition: AMReX_EB2_Graph.H:38
static constexpr Type_t covered
Definition: AMReX_EB2_Graph.H:39
static constexpr Type_t irregular
Definition: AMReX_EB2_Graph.H:40
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int num_cuts(Real a, Real b) noexcept
Definition: AMReX_EB2_2D_C.H:71
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Real coarsen_edge_cent(Real f1, Real f2)
Definition: AMReX_EB2_3D_C.H:279
Definition: AMReX_FabArrayBase.H:32
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int coarsen_from_fine(int i, int j, Box const &bx, int ngrow, Array4< Real > const &cvol, Array4< Real > const &ccent, Array4< Real > const &cba, Array4< Real > const &cbc, Array4< Real > const &cbn, Array4< Real > const &capx, Array4< Real > const &capy, Array4< Real > const &cfcx, Array4< Real > const &cfcy, Array4< Real > const &cecx, Array4< Real > const &cecy, Array4< EBCellFlag > const &cflag, Array4< Real const > const &fvol, Array4< Real const > const &fcent, Array4< Real const > const &fba, Array4< Real const > const &fbc, Array4< Real const > const &fbn, Array4< Real const > const &fapx, Array4< Real const > const &fapy, Array4< Real const > const &ffcx, Array4< Real const > const &ffcy, Array4< Real const > const &fecx, Array4< Real const > const &fecy, Array4< EBCellFlag const > const &fflag)
Definition: AMReX_EB2_2D_C.H:96
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void build_cellflag_from_ap(int i, int j, Array4< EBCellFlag > const &cflag, Array4< Real const > const &apx, Array4< Real const > const &apy)
Definition: AMReX_EB2_2D_C.H:268
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void amrex_eb2_build_types(Box const &tbx, Box const &bxg2, Array4< Real const > const &s, Array4< EBCellFlag > const &cell, Array4< Type_t > const &fx, Array4< Type_t > const &fy)
Definition: AMReX_EB2_2D_C.H:9
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int check_mvmc(int i, int j, int, Array4< Real const > const &fine)
Definition: AMReX_EB2_2D_C.H:77
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Dim3 min_ubound(BoxND< dim > const &b1, BoxND< dim > const &b2) noexcept
Definition: AMReX_Box.H:1928
AMREX_GPU_HOST_DEVICE constexpr AMREX_FORCE_INLINE const T & max(const T &a, const T &b) noexcept
Definition: AMReX_Algorithm.H:35
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > grow(const BoxND< dim > &b, int i) noexcept
Grow BoxND in all directions by given amount.
Definition: AMReX_Box.H:1211
AMREX_GPU_HOST_DEVICE constexpr AMREX_FORCE_INLINE const T & min(const T &a, const T &b) noexcept
Definition: AMReX_Algorithm.H:21
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Returns a BoxND with different type.
Definition: AMReX_Box.H:1435
IntVectND< AMREX_SPACEDIM > IntVect
Definition: AMReX_BaseFwd.H:30
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > surroundingNodes(const BoxND< dim > &b, int dir) noexcept
Returns a BoxND with NODE based coordinates in direction dir that encloses BoxND b....
Definition: AMReX_Box.H:1399
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition: AMReX.cpp:225
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Dim3 max_lbound(BoxND< dim > const &b1, BoxND< dim > const &b2) noexcept
Definition: AMReX_Box.H:1909
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE GpuComplex< T > sqrt(const GpuComplex< T > &a_z) noexcept
Return the square root of a complex number.
Definition: AMReX_GpuComplex.H:373
AMREX_GPU_HOST_DEVICE AMREX_ATTRIBUTE_FLATTEN_FOR void Loop(Dim3 lo, Dim3 hi, F const &f) noexcept
Definition: AMReX_Loop.H:125
Definition: AMReX_FabArrayCommI.H:841
Definition: AMReX_Array4.H:61