ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_ThinBodyWallDist.cpp File Reference
#include <AMReX.H>
#include <AMReX_MultiFab.H>
Include dependency graph for ERF_ThinBodyWallDist.cpp:

Functions

void thinbody_wall_dist (std::unique_ptr< MultiFab > &wdist, Vector< IntVect > &xfaces, Vector< IntVect > &yfaces, Vector< IntVect > &zfaces, const Geometry &geomdata, std::unique_ptr< MultiFab > &z_phys_cc)
 

Function Documentation

◆ thinbody_wall_dist()

void thinbody_wall_dist ( std::unique_ptr< MultiFab > &  wdist,
Vector< IntVect > &  xfaces,
Vector< IntVect > &  yfaces,
Vector< IntVect > &  zfaces,
const Geometry &  geomdata,
std::unique_ptr< MultiFab > &  z_phys_cc 
)
18 {
19  BL_PROFILE("thinbody_wall_dist()");
20 
21  const Real* prob_lo = geomdata.ProbLo();
22  const Real* dx = geomdata.CellSize();
23 
24  const bool use_terrain = (z_phys_cc != nullptr);
25  if (use_terrain) {
26  Error("Thinbody wall dist calc not implemented for terrain yet");
27  }
28 
29  Gpu::DeviceVector<IntVect> xfaces_d(xfaces.size());
30  Gpu::DeviceVector<IntVect> yfaces_d(yfaces.size());
31  Gpu::DeviceVector<IntVect> zfaces_d(zfaces.size());
32  Gpu::copyAsync(Gpu::hostToDevice, xfaces.begin(), xfaces.end(), xfaces_d.begin());
33  Gpu::copyAsync(Gpu::hostToDevice, yfaces.begin(), yfaces.end(), yfaces_d.begin());
34  Gpu::copyAsync(Gpu::hostToDevice, zfaces.begin(), zfaces.end(), zfaces_d.begin());
35  auto const* xfaces_d_ptr = xfaces_d.data();
36  auto const* yfaces_d_ptr = yfaces_d.data();
37  auto const* zfaces_d_ptr = zfaces_d.data();
38 
39  for (MFIter mfi(*wdist); mfi.isValid(); ++mfi) {
40  const Box& bx = mfi.validbox();
41 
42  //const auto& z_cc = (use_terrain) ? z_phys_cc->const_array(mfi) : Array4<const Real>{};
43  auto wd_arr = wdist->array(mfi);
44 
45  if (!use_terrain) {
46  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
47  Real xr = prob_lo[0] + (i + 0.5) * dx[0];
48  Real yr = prob_lo[1] + (j + 0.5) * dx[1];
49  Real zr = prob_lo[2] + (k + 0.5) * dx[2];
50 
51  for (std::size_t iface=0; iface < xfaces_d_ptr->size(); ++iface) {
52  int ii = xfaces_d_ptr[iface][0];
53  int jj = xfaces_d_ptr[iface][1];
54  int kk = xfaces_d_ptr[iface][2];
55  Real xfc = prob_lo[0] + ii * dx[0];
56  Real yfc = prob_lo[1] + (jj+0.5) * dx[1];
57  Real zfc = prob_lo[2] + (kk+0.5) * dx[2];
58  Real y0 = prob_lo[1] + jj * dx[1];
59  Real y1 = prob_lo[1] + (jj+1)* dx[1];
60  Real z0 = prob_lo[2] + kk * dx[2];
61  Real z1 = prob_lo[2] + (kk+1)* dx[2];
62  Real wd2 = wd_arr(i, j, k) * wd_arr(i, j, k);
63  wd2 = min(wd2, (xfc-xr)*(xfc-xr) + (yfc-yr)*(yfc-yr) + (zfc-zr)*(zfc-zr));
64  wd2 = min(wd2, (xfc-xr)*(xfc-xr) + ( y0-yr)*( y0-yr) + ( z0-zr)*( z0-zr));
65  wd2 = min(wd2, (xfc-xr)*(xfc-xr) + ( y0-yr)*( y0-yr) + ( z1-zr)*( z1-zr));
66  wd2 = min(wd2, (xfc-xr)*(xfc-xr) + ( y1-yr)*( y1-yr) + ( z0-zr)*( z0-zr));
67  wd2 = min(wd2, (xfc-xr)*(xfc-xr) + ( y1-yr)*( y1-yr) + ( z1-zr)*( z1-zr));
68  wd_arr(i, j, k) = std::sqrt(wd2);
69  }
70 
71  for (std::size_t iface=0; iface < yfaces_d_ptr->size(); ++iface) {
72  int ii = yfaces_d_ptr[iface][0];
73  int jj = yfaces_d_ptr[iface][1];
74  int kk = yfaces_d_ptr[iface][2];
75  Real xfc = prob_lo[0] + (ii+0.5) * dx[0];
76  Real yfc = prob_lo[1] + jj * dx[1];
77  Real zfc = prob_lo[2] + (kk+0.5) * dx[2];
78  Real x0 = prob_lo[0] + ii * dx[0];
79  Real x1 = prob_lo[0] + (ii+1)* dx[0];
80  Real z0 = prob_lo[2] + kk * dx[2];
81  Real z1 = prob_lo[2] + (kk+1)* dx[2];
82  Real wd2 = wd_arr(i, j, k) * wd_arr(i, j, k);
83  wd2 = min(wd2, (xfc-xr)*(xfc-xr) + (yfc-yr)*(yfc-yr) + (zfc-zr)*(zfc-zr));
84  wd2 = min(wd2, ( x0-xr)*( x0-xr) + (yfc-yr)*(yfc-yr) + ( z0-zr)*( z0-zr));
85  wd2 = min(wd2, ( x0-xr)*( x0-xr) + (yfc-yr)*(yfc-yr) + ( z1-zr)*( z1-zr));
86  wd2 = min(wd2, ( x1-xr)*( x1-xr) + (yfc-yr)*(yfc-yr) + ( z0-zr)*( z0-zr));
87  wd2 = min(wd2, ( x1-xr)*( x1-xr) + (yfc-yr)*(yfc-yr) + ( z1-zr)*( z1-zr));
88  wd_arr(i, j, k) = std::sqrt(wd2);
89  }
90 
91  for (std::size_t iface=0; iface < zfaces_d_ptr->size(); ++iface) {
92  int ii = zfaces_d_ptr[iface][0];
93  int jj = zfaces_d_ptr[iface][1];
94  int kk = zfaces_d_ptr[iface][2];
95  Real xfc = prob_lo[0] + (ii+0.5) * dx[0];
96  Real yfc = prob_lo[1] + (jj+0.5) * dx[1];
97  Real zfc = prob_lo[2] + kk * dx[2];
98  Real x0 = prob_lo[0] + ii * dx[0];
99  Real x1 = prob_lo[0] + (ii+1)* dx[0];
100  Real y0 = prob_lo[1] + jj * dx[1];
101  Real y1 = prob_lo[1] + (jj+1)* dx[1];
102  Real wd2 = wd_arr(i, j, k) * wd_arr(i, j, k);
103  wd2 = min(wd2, (xfc-xr)*(xfc-xr) + (yfc-yr)*(yfc-yr) + (zfc-zr)*(zfc-zr));
104  wd2 = min(wd2, ( x0-xr)*( x0-xr) + ( y0-yr)*( y0-yr) + (zfc-zr)*(zfc-zr));
105  wd2 = min(wd2, ( x0-xr)*( x0-xr) + ( y0-yr)*( y0-yr) + (zfc-zr)*(zfc-zr));
106  wd2 = min(wd2, ( x1-xr)*( x1-xr) + ( y1-yr)*( y1-yr) + (zfc-zr)*(zfc-zr));
107  wd2 = min(wd2, ( x1-xr)*( x1-xr) + ( y1-yr)*( y1-yr) + (zfc-zr)*(zfc-zr));
108  wd_arr(i, j, k) = std::sqrt(wd2);
109  }
110  });
111  }
112  }
113 }

Referenced by ERF::MakeNewLevelFromScratch().

Here is the caller graph for this function: