REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_BoundaryConditions_zvel.cpp
Go to the documentation of this file.
1#include "AMReX_PhysBCFunct.H"
4
5using namespace amrex;
6
7/**
8 * @param[inout] dest_arr data on which to apply BCs
9 * @param[in ] bx box to update on
10 * @param[in ] domain domain box
11 * @param[in ] dxInv pm or pn
12 * @param[in ] mskr land-sea mask on rho-points
13 * @param[in ] calc_arr data to use in the RHS of calculations
14 * @param[in ] time current time
15 * @param[in ] bccomp index into both domain_bcs_type_bcr and bc_extdir_vals for icomp=0
16 */
17void REMORAPhysBCFunct::impose_zvel_bcs (const Array4<Real>& dest_arr, const Box& bx, const Box& domain,
18 const GpuArray<Real,AMREX_SPACEDIM> /*dxInv*/,const Array4<const Real>& /*mskr*/,
19 Real /*time*/, int bccomp)
20{
21 const auto& dom_lo = amrex::lbound(domain);
22 const auto& dom_hi = amrex::ubound(domain);
23
24 // Based on BCRec for the domain, we need to make BCRec for this Box
25 // bccomp is used as starting index for m_domain_bcs_type
26 // 0 is used as starting index for bcrs
27 int ncomp = 1;
28 Vector<BCRec> bcrs(ncomp);
29 amrex::setBC(bx, domain, bccomp, 0, ncomp, m_domain_bcs_type, bcrs);
30
31 // xlo: ori = 0
32 // ylo: ori = 1
33 // zlo: ori = 2
34 // xhi: ori = 3
35 // yhi: ori = 4
36 // zhi: ori = 5
37
38 //! Based on BCRec for the domain, we need to make BCRec for this Box
39 // bccomp is used as starting index for m_domain_bcs_type
40 // 0 is used as starting index for bcrs
41 amrex::setBC(bx, domain, bccomp, 0, ncomp, m_domain_bcs_type, bcrs);
42
43 amrex::Gpu::DeviceVector<BCRec> bcrs_d(ncomp);
44#ifdef AMREX_USE_GPU
45 Gpu::htod_memcpy_async(bcrs_d.data(), bcrs.data(), sizeof(BCRec)*ncomp);
46#else
47 std::memcpy(bcrs_d.data(), bcrs.data(), sizeof(BCRec)*ncomp);
48#endif
49 const amrex::BCRec* bc_ptr = bcrs_d.data();
50
51 GpuArray<GpuArray<Real, AMREX_SPACEDIM*2>,AMREX_SPACEDIM+NCONS+8> l_bc_extdir_vals_d;
52
53 for (int i = 0; i < ncomp; i++)
54 for (int ori = 0; ori < 2*AMREX_SPACEDIM; ori++)
55 l_bc_extdir_vals_d[i][ori] = m_bc_extdir_vals[bccomp+i][ori];
56
57 GeometryData const& geomdata = m_geom.data();
58 bool is_periodic_in_x = geomdata.isPeriodic(0);
59 bool is_periodic_in_y = geomdata.isPeriodic(1);
60
61 // First do all ext_dir bcs
62 if (!is_periodic_in_x or bccomp==BCVars::foextrap_bc)
63 {
64 Box bx_xlo(bx); bx_xlo.setBig (0,dom_lo.x-1);
65 Box bx_xhi(bx); bx_xhi.setSmall(0,dom_hi.x+1);
66 ParallelFor(
67 bx_xlo, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
68 int iflip = dom_lo.x - 1 - i;
69 if (bc_ptr[n].lo(0) == REMORABCType::ext_dir) {
70 dest_arr(i,j,k) = l_bc_extdir_vals_d[n][0];
71 } else if (bc_ptr[n].lo(0) == REMORABCType::foextrap || bc_ptr[n].lo(0) == REMORABCType::clamped) {
72 dest_arr(i,j,k) = dest_arr(dom_lo.x,j,k);
73 } else if (bc_ptr[n].lo(0) == REMORABCType::reflect_even) {
74 dest_arr(i,j,k) = dest_arr(iflip,j,k);
75 } else if (bc_ptr[n].lo(0) == REMORABCType::reflect_odd) {
76 dest_arr(i,j,k) = -dest_arr(iflip,j,k);
77 }
78 },
79 bx_xhi, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
80 int iflip = 2*dom_hi.x + 1 - i;
81 if (bc_ptr[n].hi(0) == REMORABCType::ext_dir) {
82 dest_arr(i,j,k) = l_bc_extdir_vals_d[n][3];
83 } else if (bc_ptr[n].hi(0) == REMORABCType::foextrap || bc_ptr[n].hi(0) == REMORABCType::clamped) {
84 dest_arr(i,j,k) = dest_arr(dom_hi.x,j,k);
85 } else if (bc_ptr[n].hi(0) == REMORABCType::reflect_even) {
86 dest_arr(i,j,k) = dest_arr(iflip,j,k);
87 } else if (bc_ptr[n].hi(0) == REMORABCType::reflect_odd) {
88 dest_arr(i,j,k) = -dest_arr(iflip,j,k);
89 }
90 }
91 );
92 }
93
94 // First do all ext_dir bcs
95 if (!is_periodic_in_y or bccomp==BCVars::foextrap_bc)
96 {
97 Box bx_ylo(bx); bx_ylo.setBig (1,dom_lo.y-1);
98 Box bx_yhi(bx); bx_yhi.setSmall(1,dom_hi.y+1);
99 ParallelFor(bx_ylo, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
100 int jflip = dom_lo.y - 1 - j;
101 if (bc_ptr[n].lo(1) == REMORABCType::ext_dir) {
102 dest_arr(i,j,k) = l_bc_extdir_vals_d[n][1];
103 } else if (bc_ptr[n].lo(1) == REMORABCType::foextrap || bc_ptr[n].lo(1) == REMORABCType::clamped) {
104 dest_arr(i,j,k) = dest_arr(i,dom_lo.y,k);
105 } else if (bc_ptr[n].lo(1) == REMORABCType::reflect_even) {
106 dest_arr(i,j,k) = dest_arr(i,jflip,k);
107 } else if (bc_ptr[n].lo(1) == REMORABCType::reflect_odd) {
108 dest_arr(i,j,k) = -dest_arr(i,jflip,k);
109 }
110 },
111 bx_yhi, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) {
112 int jflip = 2*dom_hi.y + 1 - j;
113 if (bc_ptr[n].hi(1) == REMORABCType::ext_dir) {
114 dest_arr(i,j,k) = l_bc_extdir_vals_d[n][4];
115 } else if (bc_ptr[n].hi(1) == REMORABCType::foextrap || bc_ptr[n].hi(1) == REMORABCType::clamped) {
116 dest_arr(i,j,k) = dest_arr(i,dom_hi.y,k);
117 } else if (bc_ptr[n].hi(1) == REMORABCType::reflect_even) {
118 dest_arr(i,j,k) = dest_arr(i,jflip,k);
119 } else if (bc_ptr[n].hi(1) == REMORABCType::reflect_odd) {
120 dest_arr(i,j,k) = -dest_arr(i,jflip,k);
121 }
122 });
123 }
124
125 // Never need to apply a boundary condition to the top
126
127 Gpu::streamSynchronize();
128}
#define NCONS
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NCONS+8 > m_bc_extdir_vals
amrex::Vector< amrex::BCRec > m_domain_bcs_type
amrex::Geometry m_geom
void impose_zvel_bcs(const amrex::Array4< amrex::Real > &dest_arr, const amrex::Box &bx, const amrex::Box &domain, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > dxInv, const amrex::Array4< const amrex::Real > &mskr, amrex::Real time, int bccomp)
apply z-velocity type boundary conditions